1 |
mlosch |
1.10 |
C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_getobcsw.F,v 1.9 2011/01/19 08:42:06 mlosch Exp $ |
2 |
jmc |
1.8 |
C $Name: $ |
3 |
heimbach |
1.2 |
|
4 |
|
|
#include "CTRL_CPPOPTIONS.h" |
5 |
|
|
#ifdef ALLOW_OBCS |
6 |
|
|
# include "OBCS_OPTIONS.h" |
7 |
|
|
#endif |
8 |
|
|
|
9 |
|
|
|
10 |
|
|
subroutine ctrl_getobcsw( |
11 |
|
|
I mytime, |
12 |
|
|
I myiter, |
13 |
|
|
I mythid |
14 |
|
|
& ) |
15 |
|
|
|
16 |
|
|
c ================================================================== |
17 |
|
|
c SUBROUTINE ctrl_getobcsw |
18 |
|
|
c ================================================================== |
19 |
|
|
c |
20 |
|
|
c o Get western obc of the control vector and add it |
21 |
|
|
c to dyn. fields |
22 |
|
|
c |
23 |
|
|
c started: heimbach@mit.edu, 29-Aug-2001 |
24 |
|
|
c |
25 |
|
|
c modified: gebbie@mit.edu, 18-Mar-2003 |
26 |
|
|
c ================================================================== |
27 |
|
|
c SUBROUTINE ctrl_getobcsw |
28 |
|
|
c ================================================================== |
29 |
|
|
|
30 |
|
|
implicit none |
31 |
|
|
|
32 |
|
|
#ifdef ALLOW_OBCSW_CONTROL |
33 |
|
|
|
34 |
|
|
c == global variables == |
35 |
|
|
|
36 |
|
|
#include "EEPARAMS.h" |
37 |
|
|
#include "SIZE.h" |
38 |
|
|
#include "PARAMS.h" |
39 |
|
|
#include "GRID.h" |
40 |
|
|
#include "OBCS.h" |
41 |
|
|
|
42 |
|
|
#include "ctrl.h" |
43 |
|
|
#include "ctrl_dummy.h" |
44 |
|
|
#include "optim.h" |
45 |
|
|
|
46 |
|
|
c == routine arguments == |
47 |
|
|
|
48 |
|
|
_RL mytime |
49 |
|
|
integer myiter |
50 |
|
|
integer mythid |
51 |
|
|
|
52 |
|
|
c == local variables == |
53 |
|
|
|
54 |
|
|
integer bi,bj |
55 |
|
|
integer i,j,k |
56 |
|
|
integer itlo,ithi |
57 |
|
|
integer jtlo,jthi |
58 |
|
|
integer jmin,jmax |
59 |
|
|
integer imin,imax |
60 |
|
|
integer ilobcsw |
61 |
|
|
integer iobcs |
62 |
|
|
integer ip1 |
63 |
|
|
|
64 |
|
|
_RL dummy |
65 |
|
|
_RL obcswfac |
66 |
|
|
logical obcswfirst |
67 |
|
|
logical obcswchanged |
68 |
|
|
integer obcswcount0 |
69 |
|
|
integer obcswcount1 |
70 |
|
|
|
71 |
|
|
cgg _RL maskyz (1-oly:sny+oly,nr,nsx,nsy) |
72 |
mlosch |
1.10 |
_RL tmpfldyz (1-oly:sny+oly,nr,nsx,nsy) |
73 |
heimbach |
1.2 |
|
74 |
|
|
logical doglobalread |
75 |
|
|
logical ladinit |
76 |
|
|
|
77 |
|
|
character*(80) fnameobcsw |
78 |
|
|
|
79 |
|
|
cgg( Variables for splitting barotropic/baroclinic vels. |
80 |
|
|
_RL ubaro |
81 |
|
|
_RL utop |
82 |
|
|
cgg) |
83 |
|
|
|
84 |
|
|
c == external functions == |
85 |
|
|
|
86 |
|
|
integer ilnblnk |
87 |
|
|
external ilnblnk |
88 |
|
|
|
89 |
|
|
|
90 |
|
|
c == end of interface == |
91 |
|
|
|
92 |
|
|
jtlo = mybylo(mythid) |
93 |
|
|
jthi = mybyhi(mythid) |
94 |
|
|
itlo = mybxlo(mythid) |
95 |
|
|
ithi = mybxhi(mythid) |
96 |
|
|
jmin = 1-oly |
97 |
|
|
jmax = sny+oly |
98 |
|
|
imin = 1-olx |
99 |
|
|
imax = snx+olx |
100 |
|
|
ip1 = 1 |
101 |
|
|
|
102 |
|
|
cgg( Initialize variables for balancing volume flux. |
103 |
|
|
ubaro = 0.d0 |
104 |
|
|
utop = 0.d0 |
105 |
|
|
cgg) |
106 |
|
|
|
107 |
|
|
c-- Now, read the control vector. |
108 |
|
|
doglobalread = .false. |
109 |
|
|
ladinit = .false. |
110 |
|
|
|
111 |
|
|
if (optimcycle .ge. 0) then |
112 |
|
|
ilobcsw=ilnblnk( xx_obcsw_file ) |
113 |
jmc |
1.8 |
write(fnameobcsw(1:80),'(2a,i10.10)') |
114 |
heimbach |
1.2 |
& xx_obcsw_file(1:ilobcsw), '.', optimcycle |
115 |
|
|
endif |
116 |
|
|
|
117 |
|
|
c-- Get the counters, flags, and the interpolation factor. |
118 |
|
|
call ctrl_get_gen_rec( |
119 |
|
|
I xx_obcswstartdate, xx_obcswperiod, |
120 |
|
|
O obcswfac, obcswfirst, obcswchanged, |
121 |
|
|
O obcswcount0,obcswcount1, |
122 |
|
|
I mytime, myiter, mythid ) |
123 |
|
|
|
124 |
mlosch |
1.10 |
CML print *,'ml-getobcs: ',myIter,obcswfirst,obcswchanged, |
125 |
|
|
CML & obcswcount0,obcswcount1,obcswfac |
126 |
heimbach |
1.2 |
do iobcs = 1,nobcs |
127 |
|
|
if ( obcswfirst ) then |
128 |
jmc |
1.8 |
call active_read_yz( fnameobcsw, tmpfldyz, |
129 |
heimbach |
1.2 |
& (obcswcount0-1)*nobcs+iobcs, |
130 |
|
|
& doglobalread, ladinit, optimcycle, |
131 |
|
|
& mythid, xx_obcsw_dummy ) |
132 |
|
|
|
133 |
mlosch |
1.6 |
#ifdef ALLOW_CTRL_OBCS_BALANCE |
134 |
|
|
|
135 |
jmc |
1.8 |
if ( optimcycle .gt. 0) then |
136 |
heimbach |
1.2 |
if (iobcs .eq. 3) then |
137 |
|
|
cgg Special attention is needed for the normal velocity. |
138 |
|
|
cgg For the north, this is the v velocity, iobcs = 4. |
139 |
|
|
cgg This is done on a columnwise basis here. |
140 |
|
|
do bj = jtlo,jthi |
141 |
|
|
do bi = itlo, ithi |
142 |
|
|
do j = jmin,jmax |
143 |
|
|
i = OB_Iw(J,bi,bj) |
144 |
|
|
|
145 |
|
|
cgg The barotropic velocity is stored in the level 1. |
146 |
|
|
ubaro = tmpfldyz(j,1,bi,bj) |
147 |
|
|
tmpfldyz(j,1,bi,bj) = 0.d0 |
148 |
|
|
utop = 0.d0 |
149 |
|
|
|
150 |
|
|
do k = 1,Nr |
151 |
|
|
cgg If cells are not full, this should be modified with hFac. |
152 |
jmc |
1.8 |
cgg |
153 |
mlosch |
1.10 |
cgg The xx field (tmpfldyz) does not contain the velocity at the |
154 |
jmc |
1.8 |
cgg surface level. This velocity is not independent; it must |
155 |
heimbach |
1.2 |
cgg exactly balance the volume flux, since we are dealing with |
156 |
jmc |
1.8 |
cgg the baroclinic velocity structure.. |
157 |
heimbach |
1.2 |
utop = tmpfldyz(j,k,bi,bj)* |
158 |
heimbach |
1.5 |
& maskW(i+ip1,j,k,bi,bj) * delR(k) + utop |
159 |
jmc |
1.8 |
cgg Add the barotropic velocity component. |
160 |
heimbach |
1.2 |
if (maskW(i+ip1,j,k,bi,bj) .ne. 0.) then |
161 |
|
|
tmpfldyz(j,k,bi,bj) = tmpfldyz(j,k,bi,bj)+ ubaro |
162 |
|
|
endif |
163 |
|
|
enddo |
164 |
|
|
cgg Compute the baroclinic velocity at level 1. Should balance flux. |
165 |
jmc |
1.8 |
tmpfldyz(j,1,bi,bj) = tmpfldyz(j,1,bi,bj) |
166 |
heimbach |
1.5 |
& - utop / delR(1) |
167 |
heimbach |
1.2 |
enddo |
168 |
|
|
enddo |
169 |
|
|
enddo |
170 |
|
|
endif |
171 |
|
|
if (iobcs .eq. 4) then |
172 |
|
|
cgg Special attention is needed for the normal velocity. |
173 |
|
|
cgg For the north, this is the v velocity, iobcs = 4. |
174 |
|
|
cgg This is done on a columnwise basis here. |
175 |
|
|
do bj = jtlo,jthi |
176 |
|
|
do bi = itlo, ithi |
177 |
|
|
do j = jmin,jmax |
178 |
|
|
i = OB_Iw(J,bi,bj) |
179 |
|
|
|
180 |
|
|
cgg The barotropic velocity is stored in the level 1. |
181 |
|
|
ubaro = tmpfldyz(j,1,bi,bj) |
182 |
|
|
tmpfldyz(j,1,bi,bj) = 0.d0 |
183 |
|
|
utop = 0.d0 |
184 |
|
|
|
185 |
|
|
do k = 1,Nr |
186 |
|
|
cgg If cells are not full, this should be modified with hFac. |
187 |
jmc |
1.8 |
cgg |
188 |
mlosch |
1.10 |
cgg The xx field (tmpfldyz) does not contain the velocity at the |
189 |
jmc |
1.8 |
cgg surface level. This velocity is not independent; it must |
190 |
heimbach |
1.2 |
cgg exactly balance the volume flux, since we are dealing with |
191 |
jmc |
1.8 |
cgg the baroclinic velocity structure.. |
192 |
heimbach |
1.2 |
utop = tmpfldyz(j,k,bi,bj)* |
193 |
heimbach |
1.5 |
& maskS(i,j,k,bi,bj) * delR(k) + utop |
194 |
jmc |
1.8 |
cgg Add the barotropic velocity component. |
195 |
heimbach |
1.2 |
if (maskS(i,j,k,bi,bj) .ne. 0.) then |
196 |
|
|
tmpfldyz(j,k,bi,bj) = tmpfldyz(j,k,bi,bj)+ ubaro |
197 |
|
|
endif |
198 |
|
|
enddo |
199 |
|
|
cgg Compute the baroclinic velocity at level 1. Should balance flux. |
200 |
jmc |
1.8 |
tmpfldyz(j,1,bi,bj) = tmpfldyz(j,1,bi,bj) |
201 |
heimbach |
1.5 |
& - utop / delR(1) |
202 |
heimbach |
1.2 |
enddo |
203 |
|
|
enddo |
204 |
|
|
enddo |
205 |
|
|
endif |
206 |
|
|
endif |
207 |
|
|
|
208 |
mlosch |
1.6 |
#endif /* ALLOW_CTRL_OBCS_BALANCE */ |
209 |
|
|
|
210 |
heimbach |
1.2 |
do bj = jtlo,jthi |
211 |
|
|
do bi = itlo,ithi |
212 |
|
|
do k = 1,nr |
213 |
|
|
do j = jmin,jmax |
214 |
|
|
xx_obcsw1(j,k,bi,bj,iobcs) = tmpfldyz (j,k,bi,bj) |
215 |
|
|
cgg & * maskyz (j,k,bi,bj) |
216 |
|
|
enddo |
217 |
|
|
enddo |
218 |
|
|
enddo |
219 |
|
|
enddo |
220 |
|
|
endif |
221 |
|
|
|
222 |
|
|
if ( (obcswfirst) .or. (obcswchanged)) then |
223 |
jmc |
1.8 |
|
224 |
heimbach |
1.2 |
do bj = jtlo,jthi |
225 |
mlosch |
1.9 |
do bi = itlo,ithi |
226 |
|
|
do k = 1,nr |
227 |
|
|
do j = jmin,jmax |
228 |
|
|
xx_obcsw0(j,k,bi,bj,iobcs) = xx_obcsw1(j,k,bi,bj,iobcs) |
229 |
|
|
tmpfldyz (j,k,bi,bj) = 0. _d 0 |
230 |
|
|
enddo |
231 |
heimbach |
1.2 |
enddo |
232 |
mlosch |
1.9 |
enddo |
233 |
heimbach |
1.2 |
enddo |
234 |
|
|
|
235 |
jmc |
1.8 |
call active_read_yz( fnameobcsw, tmpfldyz, |
236 |
heimbach |
1.2 |
& (obcswcount1-1)*nobcs+iobcs, |
237 |
|
|
& doglobalread, ladinit, optimcycle, |
238 |
|
|
& mythid, xx_obcsw_dummy ) |
239 |
|
|
|
240 |
mlosch |
1.6 |
#ifdef ALLOW_CTRL_OBCS_BALANCE |
241 |
|
|
|
242 |
jmc |
1.8 |
if ( optimcycle .gt. 0) then |
243 |
heimbach |
1.2 |
if (iobcs .eq. 3) then |
244 |
|
|
cgg Special attention is needed for the normal velocity. |
245 |
|
|
cgg For the north, this is the v velocity, iobcs = 4. |
246 |
|
|
cgg This is done on a columnwise basis here. |
247 |
|
|
do bj = jtlo,jthi |
248 |
|
|
do bi = itlo, ithi |
249 |
|
|
do j = jmin,jmax |
250 |
|
|
i = OB_Iw(J,bi,bj) |
251 |
|
|
|
252 |
|
|
cgg The barotropic velocity is stored in the level 1. |
253 |
|
|
ubaro = tmpfldyz(j,1,bi,bj) |
254 |
|
|
tmpfldyz(j,1,bi,bj) = 0.d0 |
255 |
|
|
utop = 0.d0 |
256 |
|
|
|
257 |
|
|
do k = 1,Nr |
258 |
|
|
cgg If cells are not full, this should be modified with hFac. |
259 |
jmc |
1.8 |
cgg |
260 |
mlosch |
1.10 |
cgg The xx field (tmpfldyz) does not contain the velocity at the |
261 |
jmc |
1.8 |
cgg surface level. This velocity is not independent; it must |
262 |
heimbach |
1.2 |
cgg exactly balance the volume flux, since we are dealing with |
263 |
jmc |
1.8 |
cgg the baroclinic velocity structure.. |
264 |
heimbach |
1.2 |
utop = tmpfldyz(j,k,bi,bj)* |
265 |
heimbach |
1.5 |
& maskW(i+ip1,j,k,bi,bj) * delR(k) + utop |
266 |
jmc |
1.8 |
cgg Add the barotropic velocity component. |
267 |
heimbach |
1.2 |
if (maskW(i+ip1,j,k,bi,bj) .ne. 0.) then |
268 |
|
|
tmpfldyz(j,k,bi,bj) = tmpfldyz(j,k,bi,bj)+ ubaro |
269 |
|
|
endif |
270 |
|
|
enddo |
271 |
|
|
cgg Compute the baroclinic velocity at level 1. Should balance flux. |
272 |
jmc |
1.8 |
tmpfldyz(j,1,bi,bj) = tmpfldyz(j,1,bi,bj) |
273 |
heimbach |
1.5 |
& - utop / delR(1) |
274 |
heimbach |
1.2 |
enddo |
275 |
|
|
enddo |
276 |
|
|
enddo |
277 |
|
|
endif |
278 |
|
|
if (iobcs .eq. 4) then |
279 |
|
|
cgg Special attention is needed for the normal velocity. |
280 |
|
|
cgg For the north, this is the v velocity, iobcs = 4. |
281 |
|
|
cgg This is done on a columnwise basis here. |
282 |
|
|
do bj = jtlo,jthi |
283 |
|
|
do bi = itlo, ithi |
284 |
|
|
do j = jmin,jmax |
285 |
|
|
i = OB_Iw(J,bi,bj) |
286 |
|
|
|
287 |
|
|
cgg The barotropic velocity is stored in the level 1. |
288 |
|
|
ubaro = tmpfldyz(j,1,bi,bj) |
289 |
|
|
tmpfldyz(j,1,bi,bj) = 0.d0 |
290 |
|
|
utop = 0.d0 |
291 |
|
|
|
292 |
|
|
do k = 1,Nr |
293 |
|
|
cgg If cells are not full, this should be modified with hFac. |
294 |
jmc |
1.8 |
cgg |
295 |
mlosch |
1.10 |
cgg The xx field (tmpfldyz) does not contain the velocity at the |
296 |
jmc |
1.8 |
cgg surface level. This velocity is not independent; it must |
297 |
heimbach |
1.2 |
cgg exactly balance the volume flux, since we are dealing with |
298 |
jmc |
1.8 |
cgg the baroclinic velocity structure.. |
299 |
heimbach |
1.2 |
utop = tmpfldyz(j,k,bi,bj)* |
300 |
heimbach |
1.5 |
& maskS(i,j,k,bi,bj) * delR(k) + utop |
301 |
jmc |
1.8 |
cgg Add the barotropic velocity component. |
302 |
heimbach |
1.2 |
if (maskS(i,j,k,bi,bj) .ne. 0.) then |
303 |
|
|
tmpfldyz(j,k,bi,bj) = tmpfldyz(j,k,bi,bj)+ ubaro |
304 |
|
|
endif |
305 |
|
|
enddo |
306 |
|
|
cgg Compute the baroclinic velocity at level 1. Should balance flux. |
307 |
jmc |
1.8 |
tmpfldyz(j,1,bi,bj) = tmpfldyz(j,1,bi,bj) |
308 |
heimbach |
1.5 |
& - utop / delR(1) |
309 |
heimbach |
1.2 |
enddo |
310 |
|
|
enddo |
311 |
|
|
enddo |
312 |
|
|
endif |
313 |
|
|
endif |
314 |
|
|
|
315 |
mlosch |
1.6 |
#endif /* ALLOW_CTRL_OBCS_BALANCE */ |
316 |
|
|
|
317 |
heimbach |
1.2 |
do bj = jtlo,jthi |
318 |
|
|
do bi = itlo,ithi |
319 |
|
|
do k = 1,nr |
320 |
|
|
do j = jmin,jmax |
321 |
|
|
xx_obcsw1 (j,k,bi,bj,iobcs) = tmpfldyz (j,k,bi,bj) |
322 |
|
|
cgg & * maskyz (j,k,bi,bj) |
323 |
|
|
enddo |
324 |
|
|
enddo |
325 |
|
|
enddo |
326 |
|
|
enddo |
327 |
|
|
endif |
328 |
|
|
|
329 |
mlosch |
1.10 |
c-- Add control to model variable. |
330 |
heimbach |
1.2 |
do bj = jtlo, jthi |
331 |
mlosch |
1.10 |
do bi = itlo, ithi |
332 |
|
|
c-- Calculate mask for tracer cells (0 => land, 1 => water). |
333 |
|
|
do k = 1,nr |
334 |
|
|
do j = 1,sny |
335 |
|
|
i = OB_Iw(j,bi,bj) |
336 |
|
|
if (iobcs .EQ. 1) then |
337 |
|
|
OBWt(j,k,bi,bj) = OBWt (j,k,bi,bj) |
338 |
|
|
& + obcswfac *xx_obcsw0(j,k,bi,bj,iobcs) |
339 |
|
|
& + (1. _d 0 - obcswfac)*xx_obcsw1(j,k,bi,bj,iobcs) |
340 |
|
|
OBWt(j,k,bi,bj) = OBWt(j,k,bi,bj) |
341 |
|
|
& *maskW(i+ip1,j,k,bi,bj) |
342 |
|
|
else if (iobcs .EQ. 2) then |
343 |
|
|
OBWs(j,k,bi,bj) = OBWs (j,k,bi,bj) |
344 |
|
|
& + obcswfac *xx_obcsw0(j,k,bi,bj,iobcs) |
345 |
|
|
& + (1. _d 0 - obcswfac)*xx_obcsw1(j,k,bi,bj,iobcs) |
346 |
|
|
OBWs(j,k,bi,bj) = OBWs(j,k,bi,bj) |
347 |
|
|
& *maskW(i+ip1,j,k,bi,bj) |
348 |
|
|
else if (iobcs .EQ. 3) then |
349 |
|
|
OBWu(j,k,bi,bj) = OBWu (j,k,bi,bj) |
350 |
|
|
& + obcswfac *xx_obcsw0(j,k,bi,bj,iobcs) |
351 |
|
|
& + (1. _d 0 - obcswfac)*xx_obcsw1(j,k,bi,bj,iobcs) |
352 |
|
|
OBWu(j,k,bi,bj) = OBWu(j,k,bi,bj) |
353 |
|
|
& *maskW(i+ip1,j,k,bi,bj) |
354 |
|
|
else if (iobcs .EQ. 4) then |
355 |
|
|
OBWv(j,k,bi,bj) = OBWv (j,k,bi,bj) |
356 |
|
|
& + obcswfac *xx_obcsw0(j,k,bi,bj,iobcs) |
357 |
|
|
& + (1. _d 0 - obcswfac)*xx_obcsw1(j,k,bi,bj,iobcs) |
358 |
|
|
OBWv(j,k,bi,bj) = OBWv(j,k,bi,bj) |
359 |
|
|
& *maskS(i,j,k,bi,bj) |
360 |
|
|
endif |
361 |
heimbach |
1.2 |
enddo |
362 |
mlosch |
1.10 |
enddo |
363 |
|
|
enddo |
364 |
heimbach |
1.2 |
enddo |
365 |
|
|
|
366 |
|
|
C-- End over iobcs loop |
367 |
|
|
enddo |
368 |
|
|
|
369 |
|
|
#else /* ALLOW_OBCSW_CONTROL undefined */ |
370 |
|
|
|
371 |
|
|
c == routine arguments == |
372 |
|
|
|
373 |
|
|
_RL mytime |
374 |
|
|
integer myiter |
375 |
|
|
integer mythid |
376 |
|
|
|
377 |
|
|
c-- CPP flag ALLOW_OBCSW_CONTROL undefined. |
378 |
|
|
|
379 |
|
|
#endif /* ALLOW_OBCSW_CONTROL */ |
380 |
|
|
|
381 |
|
|
end |
382 |
|
|
|