1 |
C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_getobcsn.F,v 1.16 2012/09/18 20:21:23 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "CTRL_OPTIONS.h" |
5 |
#ifdef ALLOW_OBCS |
6 |
# include "OBCS_OPTIONS.h" |
7 |
#endif |
8 |
|
9 |
subroutine ctrl_getobcsn( |
10 |
I mytime, |
11 |
I myiter, |
12 |
I mythid |
13 |
& ) |
14 |
|
15 |
c ================================================================== |
16 |
c SUBROUTINE ctrl_getobcsn |
17 |
c ================================================================== |
18 |
c |
19 |
c o Get northern obc of the control vector and add it |
20 |
c to dyn. fields |
21 |
c |
22 |
c started: heimbach@mit.edu, 29-Aug-2001 |
23 |
c |
24 |
c modified: gebbie@mit.edu, 18-Mar-2003 |
25 |
c ================================================================== |
26 |
c SUBROUTINE ctrl_getobcsn |
27 |
c ================================================================== |
28 |
|
29 |
implicit none |
30 |
|
31 |
c == global variables == |
32 |
#ifdef ALLOW_OBCSN_CONTROL |
33 |
#include "EEPARAMS.h" |
34 |
#include "SIZE.h" |
35 |
#include "PARAMS.h" |
36 |
#include "GRID.h" |
37 |
c#include "OBCS_PARAMS.h" |
38 |
#include "OBCS_GRID.h" |
39 |
#include "OBCS_FIELDS.h" |
40 |
#include "CTRL_SIZE.h" |
41 |
#include "ctrl.h" |
42 |
#include "ctrl_dummy.h" |
43 |
#include "CTRL_OBCS.h" |
44 |
#include "optim.h" |
45 |
#endif /* ALLOW_OBCSN_CONTROL */ |
46 |
|
47 |
c == routine arguments == |
48 |
_RL mytime |
49 |
integer myiter |
50 |
integer mythid |
51 |
|
52 |
#ifdef ALLOW_OBCSN_CONTROL |
53 |
c == local variables == |
54 |
|
55 |
integer bi,bj |
56 |
integer i,j,k |
57 |
integer itlo,ithi |
58 |
integer jtlo,jthi |
59 |
integer jmin,jmax |
60 |
integer imin,imax |
61 |
integer ilobcsn |
62 |
integer iobcs |
63 |
integer jp1 |
64 |
|
65 |
_RL dummy |
66 |
_RL obcsnfac |
67 |
logical obcsnfirst |
68 |
logical obcsnchanged |
69 |
integer obcsncount0 |
70 |
integer obcsncount1 |
71 |
|
72 |
cgg _RL maskxz (1-olx:snx+olx,nr,nsx,nsy) |
73 |
_RL tmpfldxz (1-olx:snx+olx,nr,nsx,nsy) |
74 |
|
75 |
logical doglobalread |
76 |
logical ladinit |
77 |
|
78 |
character*(80) fnameobcsn |
79 |
|
80 |
#ifdef ALLOW_OBCS_CONTROL_MODES |
81 |
integer nk,nz |
82 |
_RL tmpz (nr,nsx,nsy) |
83 |
_RL stmp |
84 |
#endif |
85 |
|
86 |
c == external functions == |
87 |
|
88 |
integer ilnblnk |
89 |
external ilnblnk |
90 |
|
91 |
|
92 |
c == end of interface == |
93 |
|
94 |
jtlo = mybylo(mythid) |
95 |
jthi = mybyhi(mythid) |
96 |
itlo = mybxlo(mythid) |
97 |
ithi = mybxhi(mythid) |
98 |
cgg jmin = 1-oly |
99 |
cgg jmax = sny+oly |
100 |
cgg imin = 1-olx |
101 |
cgg imax = snx+olx |
102 |
|
103 |
jmin = 1 |
104 |
jmax = sny |
105 |
imin = 1 |
106 |
imax = snx |
107 |
jp1 = 0 |
108 |
|
109 |
c-- Now, read the control vector. |
110 |
doglobalread = .false. |
111 |
ladinit = .false. |
112 |
|
113 |
if (optimcycle .ge. 0) then |
114 |
ilobcsn=ilnblnk( xx_obcsn_file ) |
115 |
write(fnameobcsn(1:80),'(2a,i10.10)') |
116 |
& xx_obcsn_file(1:ilobcsn), '.', optimcycle |
117 |
endif |
118 |
|
119 |
c-- Get the counters, flags, and the interpolation factor. |
120 |
call ctrl_get_gen_rec( |
121 |
I xx_obcsnstartdate, xx_obcsnperiod, |
122 |
O obcsnfac, obcsnfirst, obcsnchanged, |
123 |
O obcsncount0,obcsncount1, |
124 |
I mytime, myiter, mythid ) |
125 |
|
126 |
do iobcs = 1,nobcs |
127 |
if ( obcsnfirst ) then |
128 |
call active_read_xz( fnameobcsn, tmpfldxz, |
129 |
& (obcsncount0-1)*nobcs+iobcs, |
130 |
& doglobalread, ladinit, optimcycle, |
131 |
& mythid, xx_obcsn_dummy ) |
132 |
|
133 |
do bj = jtlo,jthi |
134 |
do bi = itlo,ithi |
135 |
#ifdef ALLOW_OBCS_CONTROL_MODES |
136 |
if (iobcs .gt. 2) then |
137 |
do i = imin,imax |
138 |
j = OB_Jn(i,bi,bj) |
139 |
IF ( j.EQ.OB_indexNone ) j = 1 |
140 |
cih Determine number of open vertical layers. |
141 |
nz = 0 |
142 |
do k = 1,Nr |
143 |
if (iobcs .eq. 3) then |
144 |
nz = nz + maskS(i,j+jp1,k,bi,bj) |
145 |
else |
146 |
nz = nz + maskW(i,j,k,bi,bj) |
147 |
endif |
148 |
end do |
149 |
cih Compute absolute velocities from the barotropic-baroclinic modes. |
150 |
do k = 1,Nr |
151 |
if (k.le.nz) then |
152 |
stmp = 0. |
153 |
do nk = 1,nz |
154 |
stmp = stmp + |
155 |
& modesv(k,nk,nz)*tmpfldxz(i,nk,bi,bj) |
156 |
end do |
157 |
tmpz(k,bi,bj) = stmp |
158 |
else |
159 |
tmpz(k,bi,bj) = 0. |
160 |
end if |
161 |
end do |
162 |
do k = 1,Nr |
163 |
if (iobcs .eq. 3) then |
164 |
tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj) |
165 |
& *recip_hFacS(i,j+jp1,k,bi,bj) |
166 |
else |
167 |
tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj) |
168 |
& *recip_hFacW(i,j,k,bi,bj) |
169 |
endif |
170 |
end do |
171 |
enddo |
172 |
endif |
173 |
#endif |
174 |
do k = 1,nr |
175 |
do i = imin,imax |
176 |
xx_obcsn1(i,k,bi,bj,iobcs) = tmpfldxz (i,k,bi,bj) |
177 |
cgg & * maskxz (i,k,bi,bj) |
178 |
enddo |
179 |
enddo |
180 |
enddo |
181 |
enddo |
182 |
endif |
183 |
|
184 |
if ( (obcsnfirst) .or. (obcsnchanged)) then |
185 |
|
186 |
do bj = jtlo,jthi |
187 |
do bi = itlo,ithi |
188 |
do k = 1,nr |
189 |
do i = imin,imax |
190 |
xx_obcsn0(i,k,bi,bj,iobcs) = xx_obcsn1(i,k,bi,bj,iobcs) |
191 |
tmpfldxz (i,k,bi,bj) = 0. _d 0 |
192 |
enddo |
193 |
enddo |
194 |
enddo |
195 |
enddo |
196 |
|
197 |
call active_read_xz( fnameobcsn, tmpfldxz, |
198 |
& (obcsncount1-1)*nobcs+iobcs, |
199 |
& doglobalread, ladinit, optimcycle, |
200 |
& mythid, xx_obcsn_dummy ) |
201 |
|
202 |
do bj = jtlo,jthi |
203 |
do bi = itlo,ithi |
204 |
#ifdef ALLOW_OBCS_CONTROL_MODES |
205 |
if (iobcs .gt. 2) then |
206 |
do i = imin,imax |
207 |
j = OB_Jn(i,bi,bj) |
208 |
IF ( j.EQ.OB_indexNone ) j = 1 |
209 |
cih Determine number of open vertical layers. |
210 |
nz = 0 |
211 |
do k = 1,Nr |
212 |
if (iobcs .eq. 3) then |
213 |
nz = nz + maskS(i,j+jp1,k,bi,bj) |
214 |
else |
215 |
nz = nz + maskW(i,j,k,bi,bj) |
216 |
endif |
217 |
end do |
218 |
cih Compute absolute velocities from the barotropic-baroclinic modes. |
219 |
do k = 1,Nr |
220 |
if (k.le.nz) then |
221 |
stmp = 0. |
222 |
do nk = 1,nz |
223 |
stmp = stmp + |
224 |
& modesv(k,nk,nz)*tmpfldxz(i,nk,bi,bj) |
225 |
end do |
226 |
tmpz(k,bi,bj) = stmp |
227 |
else |
228 |
tmpz(k,bi,bj) = 0. |
229 |
end if |
230 |
end do |
231 |
do k = 1,Nr |
232 |
if (iobcs .eq. 3) then |
233 |
tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj) |
234 |
& *recip_hFacS(i,j+jp1,k,bi,bj) |
235 |
else |
236 |
tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj) |
237 |
& *recip_hFacW(i,j,k,bi,bj) |
238 |
endif |
239 |
end do |
240 |
enddo |
241 |
endif |
242 |
#endif |
243 |
do k = 1,nr |
244 |
do i = imin,imax |
245 |
xx_obcsn1 (i,k,bi,bj,iobcs) = tmpfldxz (i,k,bi,bj) |
246 |
cgg & * maskxz (i,k,bi,bj) |
247 |
enddo |
248 |
enddo |
249 |
enddo |
250 |
enddo |
251 |
|
252 |
endif |
253 |
|
254 |
c-- Add control to model variable. |
255 |
do bj = jtlo,jthi |
256 |
do bi = itlo,ithi |
257 |
c-- Calculate mask for tracer cells (0 => land, 1 => water). |
258 |
do k = 1,nr |
259 |
do i = 1,snx |
260 |
j = OB_Jn(I,bi,bj) |
261 |
IF ( j.EQ.OB_indexNone ) j = 1 |
262 |
if (iobcs .EQ. 1) then |
263 |
OBNt(i,k,bi,bj) = OBNt (i,k,bi,bj) |
264 |
& + obcsnfac *xx_obcsn0(i,k,bi,bj,iobcs) |
265 |
& + (1. _d 0 - obcsnfac)*xx_obcsn1(i,k,bi,bj,iobcs) |
266 |
OBNt(i,k,bi,bj) = OBNt(i,k,bi,bj) |
267 |
& *maskS(i,j+jp1,k,bi,bj) |
268 |
else if (iobcs .EQ. 2) then |
269 |
OBNs(i,k,bi,bj) = OBNs (i,k,bi,bj) |
270 |
& + obcsnfac *xx_obcsn0(i,k,bi,bj,iobcs) |
271 |
& + (1. _d 0 - obcsnfac)*xx_obcsn1(i,k,bi,bj,iobcs) |
272 |
OBNs(i,k,bi,bj) = OBNs(i,k,bi,bj) |
273 |
& *maskS(i,j+jp1,k,bi,bj) |
274 |
else if (iobcs .EQ. 4) then |
275 |
OBNu(i,k,bi,bj) = OBNu (i,k,bi,bj) |
276 |
& + obcsnfac *xx_obcsn0(i,k,bi,bj,iobcs) |
277 |
& + (1. _d 0 - obcsnfac)*xx_obcsn1(i,k,bi,bj,iobcs) |
278 |
OBNu(i,k,bi,bj) = OBNu(i,k,bi,bj) |
279 |
& *maskW(i,j,k,bi,bj) |
280 |
else if (iobcs .EQ. 3) then |
281 |
OBNv(i,k,bi,bj) = OBNv (i,k,bi,bj) |
282 |
& + obcsnfac *xx_obcsn0(i,k,bi,bj,iobcs) |
283 |
& + (1. _d 0 - obcsnfac)*xx_obcsn1(i,k,bi,bj,iobcs) |
284 |
OBNv(i,k,bi,bj) = OBNv(i,k,bi,bj) |
285 |
& *maskS(i,j+jp1,k,bi,bj) |
286 |
endif |
287 |
enddo |
288 |
enddo |
289 |
enddo |
290 |
enddo |
291 |
|
292 |
C-- End over iobcs loop |
293 |
enddo |
294 |
|
295 |
#endif /* ALLOW_OBCSN_CONTROL */ |
296 |
|
297 |
return |
298 |
end |