1 |
C $Header: $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "CTRL_CPPOPTIONS.h" |
5 |
#ifdef ALLOW_OBCS |
6 |
# include "OBCS_OPTIONS.h" |
7 |
#endif |
8 |
|
9 |
|
10 |
subroutine ctrl_getobcss( |
11 |
I mytime, |
12 |
I myiter, |
13 |
I mythid |
14 |
& ) |
15 |
|
16 |
c ================================================================== |
17 |
c SUBROUTINE ctrl_getobcss |
18 |
c ================================================================== |
19 |
c |
20 |
c o Get southern 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 new flags: gebbie@mit.edu, 25 Jan 2003. |
26 |
c |
27 |
c ================================================================== |
28 |
c SUBROUTINE ctrl_getobcss |
29 |
c ================================================================== |
30 |
|
31 |
implicit none |
32 |
|
33 |
#ifdef ALLOW_OBCSS_CONTROL |
34 |
|
35 |
c == global variables == |
36 |
|
37 |
#include "EEPARAMS.h" |
38 |
#include "SIZE.h" |
39 |
#include "PARAMS.h" |
40 |
#include "GRID.h" |
41 |
#include "OBCS.h" |
42 |
|
43 |
#include "ctrl.h" |
44 |
#include "ctrl_dummy.h" |
45 |
#include "optim.h" |
46 |
|
47 |
c == routine arguments == |
48 |
|
49 |
_RL mytime |
50 |
integer myiter |
51 |
integer mythid |
52 |
|
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 ilobcss |
62 |
integer iobcs |
63 |
|
64 |
_RL dummy |
65 |
_RL obcssfac |
66 |
logical obcssfirst |
67 |
logical obcsschanged |
68 |
integer obcsscount0 |
69 |
integer obcsscount1 |
70 |
integer jp1 |
71 |
|
72 |
cgg _RL maskxz (1-olx:snx+olx,nr,nsx,nsy) |
73 |
|
74 |
logical doglobalread |
75 |
logical ladinit |
76 |
|
77 |
character*(80) fnameobcss |
78 |
|
79 |
cgg( Variables for splitting barotropic/baroclinic vels. |
80 |
_RL vbaro |
81 |
_RL vtop |
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 |
jp1 = 1 |
101 |
|
102 |
cgg( Initialize variables for balancing volume flux. |
103 |
vbaro = 0.d0 |
104 |
vtop = 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 |
ilobcss=ilnblnk( xx_obcss_file ) |
113 |
write(fnameobcss(1:80),'(2a,i10.10)') |
114 |
& xx_obcss_file(1:ilobcss), '.', optimcycle |
115 |
endif |
116 |
|
117 |
c-- Get the counters, flags, and the interpolation factor. |
118 |
call ctrl_get_gen_rec( |
119 |
I xx_obcssstartdate, xx_obcssperiod, |
120 |
O obcssfac, obcssfirst, obcsschanged, |
121 |
O obcsscount0,obcsscount1, |
122 |
I mytime, myiter, mythid ) |
123 |
|
124 |
do iobcs = 1,nobcs |
125 |
if ( obcssfirst ) then |
126 |
call active_read_xz( fnameobcss, tmpfldxz, |
127 |
& (obcsscount0-1)*nobcs+iobcs, |
128 |
& doglobalread, ladinit, optimcycle, |
129 |
& mythid, xx_obcss_dummy ) |
130 |
|
131 |
#ifdef ALLOW_CTRL_OBCS_BALANCE |
132 |
|
133 |
if ( optimcycle .gt. 0) then |
134 |
if (iobcs .eq. 3) then |
135 |
cgg Special attention is needed for the normal velocity. |
136 |
cgg For the north, this is the v velocity, iobcs = 4. |
137 |
cgg This is done on a columnwise basis here. |
138 |
do bj = jtlo,jthi |
139 |
do bi = itlo, ithi |
140 |
do i = imin,imax |
141 |
j = OB_Js(I,bi,bj) |
142 |
|
143 |
cgg The barotropic velocity is stored in the level 1. |
144 |
vbaro = tmpfldxz(i,1,bi,bj) |
145 |
tmpfldxz(i,1,bi,bj) = 0.d0 |
146 |
vtop = 0.d0 |
147 |
|
148 |
do k = 1,Nr |
149 |
cgg If cells are not full, this should be modified with hFac. |
150 |
cgg |
151 |
cgg The xx field (tmpfldxz) does not contain the velocity at the |
152 |
cgg surface level. This velocity is not independent; it must |
153 |
cgg exactly balance the volume flux, since we are dealing with |
154 |
cgg the baroclinic velocity structure.. |
155 |
vtop = tmpfldxz(i,k,bi,bj)* |
156 |
& maskS(i,j+jp1,k,bi,bj) * delR(k) + vtop |
157 |
cgg Add the barotropic velocity component. |
158 |
if (maskS(i,j+jp1,k,bi,bj) .ne. 0.) then |
159 |
tmpfldxz(i,k,bi,bj) = tmpfldxz(i,k,bi,bj)+ vbaro |
160 |
endif |
161 |
enddo |
162 |
cgg Compute the baroclinic velocity at level 1. Should balance flux. |
163 |
tmpfldxz(i,1,bi,bj) = tmpfldxz(i,1,bi,bj) |
164 |
& - vtop / delR(1) |
165 |
enddo |
166 |
enddo |
167 |
enddo |
168 |
endif |
169 |
|
170 |
if (iobcs .eq. 4) then |
171 |
cgg Special attention is needed for the normal velocity. |
172 |
cgg For the north, this is the v velocity, iobcs = 4. |
173 |
cgg This is done on a columnwise basis here. |
174 |
do bj = jtlo,jthi |
175 |
do bi = itlo, ithi |
176 |
do i = imin,imax |
177 |
j = OB_Js(I,bi,bj) |
178 |
|
179 |
cgg The barotropic velocity is stored in the level 1. |
180 |
vbaro = tmpfldxz(i,1,bi,bj) |
181 |
tmpfldxz(i,1,bi,bj) = 0.d0 |
182 |
vtop = 0.d0 |
183 |
|
184 |
do k = 1,Nr |
185 |
cgg If cells are not full, this should be modified with hFac. |
186 |
cgg |
187 |
cgg The xx field (tmpfldxz) does not contain the velocity at the |
188 |
cgg surface level. This velocity is not independent; it must |
189 |
cgg exactly balance the volume flux, since we are dealing with |
190 |
cgg the baroclinic velocity structure.. |
191 |
vtop = tmpfldxz(i,k,bi,bj)* |
192 |
& maskW(i,j,k,bi,bj) * delR(k) + vtop |
193 |
cgg Add the barotropic velocity component. |
194 |
if (maskW(i,j,k,bi,bj) .ne. 0.) then |
195 |
tmpfldxz(i,k,bi,bj) = tmpfldxz(i,k,bi,bj)+ vbaro |
196 |
endif |
197 |
enddo |
198 |
cgg Compute the baroclinic velocity at level 1. Should balance flux. |
199 |
tmpfldxz(i,1,bi,bj) = tmpfldxz(i,1,bi,bj) |
200 |
& - vtop / delR(1) |
201 |
enddo |
202 |
enddo |
203 |
enddo |
204 |
endif |
205 |
endif |
206 |
|
207 |
#endif /* ALLOW_CTRL_OBCS_BALANCE */ |
208 |
|
209 |
do bj = jtlo,jthi |
210 |
do bi = itlo,ithi |
211 |
do k = 1,nr |
212 |
do i = imin,imax |
213 |
xx_obcss1(i,k,bi,bj,iobcs) = tmpfldxz (i,k,bi,bj) |
214 |
cgg & * maskxz (i,k,bi,bj) |
215 |
enddo |
216 |
enddo |
217 |
enddo |
218 |
enddo |
219 |
endif |
220 |
|
221 |
if ( (obcssfirst) .or. (obcsschanged)) then |
222 |
|
223 |
do bj = jtlo,jthi |
224 |
do bi = itlo,ithi |
225 |
do k = 1,nr |
226 |
do i = imin,imax |
227 |
tmpfldxz(i,k,bi,bj) = xx_obcss1(i,k,bi,bj,iobcs) |
228 |
enddo |
229 |
enddo |
230 |
enddo |
231 |
enddo |
232 |
|
233 |
call exf_swapffields_xz( tmpfldxz2, tmpfldxz, mythid) |
234 |
|
235 |
do bj = jtlo,jthi |
236 |
do bi = itlo,ithi |
237 |
do k = 1,nr |
238 |
do i = imin,imax |
239 |
xx_obcss0(i,k,bi,bj,iobcs) = tmpfldxz2(i,k,bi,bj) |
240 |
enddo |
241 |
enddo |
242 |
enddo |
243 |
enddo |
244 |
|
245 |
call active_read_xz( fnameobcss, tmpfldxz, |
246 |
& (obcsscount1-1)*nobcs+iobcs, |
247 |
& doglobalread, ladinit, optimcycle, |
248 |
& mythid, xx_obcss_dummy ) |
249 |
|
250 |
#ifdef ALLOW_CTRL_OBCS_BALANCE |
251 |
|
252 |
if ( optimcycle .gt. 0) then |
253 |
if (iobcs .eq. 3) then |
254 |
cgg Special attention is needed for the normal velocity. |
255 |
cgg For the north, this is the v velocity, iobcs = 4. |
256 |
cgg This is done on a columnwise basis here. |
257 |
do bj = jtlo,jthi |
258 |
do bi = itlo, ithi |
259 |
do i = imin,imax |
260 |
j = OB_Js(I,bi,bj) |
261 |
|
262 |
cgg The barotropic velocity is stored in the level 1. |
263 |
vbaro = tmpfldxz(i,1,bi,bj) |
264 |
tmpfldxz(i,1,bi,bj) = 0.d0 |
265 |
vtop = 0.d0 |
266 |
|
267 |
do k = 1,Nr |
268 |
cgg If cells are not full, this should be modified with hFac. |
269 |
cgg |
270 |
cgg The xx field (tmpfldxz) does not contain the velocity at the |
271 |
cgg surface level. This velocity is not independent; it must |
272 |
cgg exactly balance the volume flux, since we are dealing with |
273 |
cgg the baroclinic velocity structure.. |
274 |
vtop = tmpfldxz(i,k,bi,bj)* |
275 |
& maskS(i,j+jp1,k,bi,bj) * delR(k) + vtop |
276 |
cgg Add the barotropic velocity component. |
277 |
if (maskS(i,j+jp1,k,bi,bj) .ne. 0.) then |
278 |
tmpfldxz(i,k,bi,bj) = tmpfldxz(i,k,bi,bj)+ vbaro |
279 |
endif |
280 |
enddo |
281 |
cgg Compute the baroclinic velocity at level 1. Should balance flux. |
282 |
tmpfldxz(i,1,bi,bj) = tmpfldxz(i,1,bi,bj) |
283 |
& - vtop / delR(1) |
284 |
enddo |
285 |
enddo |
286 |
enddo |
287 |
endif |
288 |
|
289 |
if (iobcs .eq. 4) then |
290 |
cgg Special attention is needed for the normal velocity. |
291 |
cgg For the north, this is the v velocity, iobcs = 4. |
292 |
cgg This is done on a columnwise basis here. |
293 |
do bj = jtlo,jthi |
294 |
do bi = itlo, ithi |
295 |
do i = imin,imax |
296 |
j = OB_Js(I,bi,bj) |
297 |
|
298 |
cgg The barotropic velocity is stored in the level 1. |
299 |
vbaro = tmpfldxz(i,1,bi,bj) |
300 |
tmpfldxz(i,1,bi,bj) = 0.d0 |
301 |
vtop = 0.d0 |
302 |
|
303 |
do k = 1,Nr |
304 |
cgg If cells are not full, this should be modified with hFac. |
305 |
cgg |
306 |
cgg The xx field (tmpfldxz) does not contain the velocity at the |
307 |
cgg surface level. This velocity is not independent; it must |
308 |
cgg exactly balance the volume flux, since we are dealing with |
309 |
cgg the baroclinic velocity structure.. |
310 |
vtop = tmpfldxz(i,k,bi,bj)* |
311 |
& maskW(i,j,k,bi,bj) * delR(k) + vtop |
312 |
cgg Add the barotropic velocity component. |
313 |
if (maskW(i,j,k,bi,bj) .ne. 0.) then |
314 |
tmpfldxz(i,k,bi,bj) = tmpfldxz(i,k,bi,bj)+ vbaro |
315 |
endif |
316 |
enddo |
317 |
cgg Compute the baroclinic velocity at level 1. Should balance flux. |
318 |
tmpfldxz(i,1,bi,bj) = tmpfldxz(i,1,bi,bj) |
319 |
& - vtop / delR(1) |
320 |
enddo |
321 |
enddo |
322 |
enddo |
323 |
endif |
324 |
endif |
325 |
|
326 |
#endif /* ALLOW_CTRL_OBCS_BALANCE */ |
327 |
|
328 |
do bj = jtlo,jthi |
329 |
do bi = itlo,ithi |
330 |
do k = 1,nr |
331 |
do i = imin,imax |
332 |
xx_obcss1 (i,k,bi,bj,iobcs) = tmpfldxz (i,k,bi,bj) |
333 |
cgg & * maskxz (i,k,bi,bj) |
334 |
enddo |
335 |
enddo |
336 |
enddo |
337 |
enddo |
338 |
endif |
339 |
|
340 |
c-- Add control to model variable. |
341 |
do bj = jtlo,jthi |
342 |
do bi = itlo,ithi |
343 |
c-- Calculate mask for tracer cells (0 => land, 1 => water). |
344 |
do k = 1,nr |
345 |
do i = 1,snx |
346 |
j = OB_Js(I,bi,bj) |
347 |
if (iobcs .EQ. 1) then |
348 |
OBSt(i,k,bi,bj) = OBSt (i,k,bi,bj) |
349 |
& + obcssfac *xx_obcss0(i,k,bi,bj,iobcs) |
350 |
& + (1. _d 0 - obcssfac)*xx_obcss1(i,k,bi,bj,iobcs) |
351 |
OBSt(i,k,bi,bj) = OBSt(i,k,bi,bj) |
352 |
& *maskS(i,j+jp1,k,bi,bj) |
353 |
else if (iobcs .EQ. 2) then |
354 |
OBSs(i,k,bi,bj) = OBSs (i,k,bi,bj) |
355 |
& + obcssfac *xx_obcss0(i,k,bi,bj,iobcs) |
356 |
& + (1. _d 0 - obcssfac)*xx_obcss1(i,k,bi,bj,iobcs) |
357 |
OBSs(i,k,bi,bj) = OBSs(i,k,bi,bj) |
358 |
& *maskS(i,j+jp1,k,bi,bj) |
359 |
else if (iobcs .EQ. 4) then |
360 |
OBSu(i,k,bi,bj) = OBSu (i,k,bi,bj) |
361 |
& + obcssfac *xx_obcss0(i,k,bi,bj,iobcs) |
362 |
& + (1. _d 0 - obcssfac)*xx_obcss1(i,k,bi,bj,iobcs) |
363 |
OBSu(i,k,bi,bj) = OBSu(i,k,bi,bj) |
364 |
& *maskW(i,j,k,bi,bj) |
365 |
else if (iobcs .EQ. 3) then |
366 |
OBSv(i,k,bi,bj) = OBSv (i,k,bi,bj) |
367 |
& + obcssfac *xx_obcss0(i,k,bi,bj,iobcs) |
368 |
& + (1. _d 0 - obcssfac)*xx_obcss1(i,k,bi,bj,iobcs) |
369 |
OBSv(i,k,bi,bj) = OBSv(i,k,bi,bj) |
370 |
& *maskS(i,j+jp1,k,bi,bj) |
371 |
endif |
372 |
enddo |
373 |
enddo |
374 |
enddo |
375 |
enddo |
376 |
|
377 |
C-- End over iobcs loop |
378 |
enddo |
379 |
|
380 |
#else /* ALLOW_OBCSS_CONTROL undefined */ |
381 |
|
382 |
c == routine arguments == |
383 |
|
384 |
_RL mytime |
385 |
integer myiter |
386 |
integer mythid |
387 |
|
388 |
c-- CPP flag ALLOW_OBCSS_CONTROL undefined. |
389 |
|
390 |
#endif /* ALLOW_OBCSS_CONTROL */ |
391 |
|
392 |
end |
393 |
|