1 |
mmazloff |
1.11 |
C $Header: /u/gcmpack/MITgcm_contrib/SOSE/BoxAdj/code_ad/ctrl_getobcse.F,v 1.6 2011/03/18 18:43:05 mmazloff Exp $ |
2 |
jmc |
1.7 |
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_getobcse( |
11 |
|
|
I mytime, |
12 |
|
|
I myiter, |
13 |
|
|
I mythid |
14 |
|
|
& ) |
15 |
|
|
|
16 |
|
|
c ================================================================== |
17 |
|
|
c SUBROUTINE ctrl_getobcse |
18 |
|
|
c ================================================================== |
19 |
|
|
c |
20 |
|
|
c o Get eastern 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 ================================================================== |
26 |
|
|
c SUBROUTINE ctrl_getobcse |
27 |
|
|
c ================================================================== |
28 |
|
|
|
29 |
|
|
implicit none |
30 |
|
|
|
31 |
|
|
#ifdef ALLOW_OBCSE_CONTROL |
32 |
|
|
|
33 |
|
|
c == global variables == |
34 |
|
|
|
35 |
|
|
#include "EEPARAMS.h" |
36 |
|
|
#include "SIZE.h" |
37 |
|
|
#include "PARAMS.h" |
38 |
|
|
#include "GRID.h" |
39 |
|
|
#include "OBCS.h" |
40 |
|
|
|
41 |
|
|
#include "ctrl.h" |
42 |
|
|
#include "ctrl_dummy.h" |
43 |
|
|
#include "optim.h" |
44 |
|
|
|
45 |
|
|
c == routine arguments == |
46 |
|
|
|
47 |
|
|
_RL mytime |
48 |
|
|
integer myiter |
49 |
|
|
integer mythid |
50 |
|
|
|
51 |
|
|
c == local variables == |
52 |
|
|
|
53 |
|
|
integer bi,bj |
54 |
|
|
integer i,j,k |
55 |
|
|
integer itlo,ithi |
56 |
|
|
integer jtlo,jthi |
57 |
|
|
integer jmin,jmax |
58 |
|
|
integer imin,imax |
59 |
|
|
integer ilobcse |
60 |
|
|
integer iobcs |
61 |
|
|
|
62 |
|
|
_RL dummy |
63 |
|
|
_RL obcsefac |
64 |
|
|
logical obcsefirst |
65 |
|
|
logical obcsechanged |
66 |
|
|
integer obcsecount0 |
67 |
|
|
integer obcsecount1 |
68 |
|
|
integer ip1 |
69 |
|
|
|
70 |
|
|
cgg _RL maskyz (1-oly:sny+oly,nr,nsx,nsy) |
71 |
mlosch |
1.9 |
_RL tmpfldyz (1-oly:sny+oly,nr,nsx,nsy) |
72 |
heimbach |
1.2 |
|
73 |
|
|
logical doglobalread |
74 |
|
|
logical ladinit |
75 |
|
|
|
76 |
|
|
character*(80) fnameobcse |
77 |
|
|
|
78 |
mmazloff |
1.11 |
#ifdef ALLOW_OBCS_CONTROL_MODES |
79 |
|
|
integer nk,nz |
80 |
|
|
_RL tmpz (nr,nsx,nsy) |
81 |
|
|
_RL stmp |
82 |
|
|
#endif |
83 |
heimbach |
1.2 |
|
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 = 0 |
101 |
|
|
|
102 |
|
|
c-- Now, read the control vector. |
103 |
|
|
doglobalread = .false. |
104 |
|
|
ladinit = .false. |
105 |
|
|
|
106 |
|
|
if (optimcycle .ge. 0) then |
107 |
mlosch |
1.10 |
ilobcse=ilnblnk( xx_obcse_file ) |
108 |
|
|
write(fnameobcse(1:80),'(2a,i10.10)') |
109 |
|
|
& xx_obcse_file(1:ilobcse), '.', optimcycle |
110 |
heimbach |
1.2 |
endif |
111 |
|
|
|
112 |
|
|
c-- Get the counters, flags, and the interpolation factor. |
113 |
|
|
call ctrl_get_gen_rec( |
114 |
|
|
I xx_obcsestartdate, xx_obcseperiod, |
115 |
|
|
O obcsefac, obcsefirst, obcsechanged, |
116 |
|
|
O obcsecount0,obcsecount1, |
117 |
|
|
I mytime, myiter, mythid ) |
118 |
|
|
|
119 |
|
|
do iobcs = 1,nobcs |
120 |
|
|
|
121 |
mlosch |
1.10 |
if ( obcsefirst ) then |
122 |
|
|
call active_read_yz( fnameobcse, tmpfldyz, |
123 |
|
|
& (obcsecount0-1)*nobcs+iobcs, |
124 |
|
|
& doglobalread, ladinit, optimcycle, |
125 |
|
|
& mythid, xx_obcse_dummy ) |
126 |
heimbach |
1.2 |
|
127 |
mlosch |
1.10 |
do bj = jtlo,jthi |
128 |
|
|
do bi = itlo,ithi |
129 |
mmazloff |
1.11 |
#ifdef ALLOW_OBCS_CONTROL_MODES |
130 |
|
|
if (iobcs .gt. 2) then |
131 |
|
|
do j = jmin,jmax |
132 |
|
|
i = OB_Ie(j,bi,bj) |
133 |
|
|
cih Determine number of open vertical layers. |
134 |
|
|
nz = 0 |
135 |
|
|
do k = 1,Nr |
136 |
|
|
if (iobcs .eq. 3) then |
137 |
|
|
nz = nz + maskW(i+ip1,j,k,bi,bj) |
138 |
|
|
else |
139 |
|
|
nz = nz + maskS(i,j,k,bi,bj) |
140 |
|
|
endif |
141 |
|
|
end do |
142 |
|
|
cih Compute absolute velocities from the barotropic-baroclinic modes. |
143 |
|
|
do k = 1,Nr |
144 |
|
|
if (k.le.nz) then |
145 |
|
|
stmp = 0. |
146 |
|
|
do nk = 1,nz |
147 |
|
|
stmp = stmp + |
148 |
|
|
& modesv(k,nk,nz)*tmpfldyz(j,nk,bi,bj) |
149 |
|
|
end do |
150 |
|
|
tmpz(k,bi,bj) = stmp |
151 |
|
|
else |
152 |
|
|
tmpz(k,bi,bj) = 0. |
153 |
|
|
end if |
154 |
|
|
enddo |
155 |
|
|
do k = 1,Nr |
156 |
|
|
if (iobcs .eq. 3) then |
157 |
|
|
tmpfldyz(j,k,bi,bj) = tmpz(k,bi,bj) |
158 |
|
|
& *recip_hFacW(i+ip1,j,k,bi,bj) |
159 |
|
|
else |
160 |
|
|
tmpfldyz(j,k,bi,bj) = tmpz(k,bi,bj) |
161 |
|
|
& *recip_hFacS(i,j,k,bi,bj) |
162 |
|
|
endif |
163 |
|
|
end do |
164 |
|
|
enddo |
165 |
|
|
endif |
166 |
|
|
#endif |
167 |
mlosch |
1.10 |
do k = 1,nr |
168 |
|
|
do j = jmin,jmax |
169 |
|
|
xx_obcse1(j,k,bi,bj,iobcs) = tmpfldyz (j,k,bi,bj) |
170 |
|
|
cgg & * maskyz (j,k,bi,bj) |
171 |
mlosch |
1.8 |
enddo |
172 |
heimbach |
1.2 |
enddo |
173 |
mlosch |
1.10 |
enddo |
174 |
|
|
enddo |
175 |
|
|
endif |
176 |
|
|
|
177 |
|
|
if ( (obcsefirst) .or. (obcsechanged)) then |
178 |
|
|
|
179 |
|
|
do bj = jtlo,jthi |
180 |
|
|
do bi = itlo,ithi |
181 |
|
|
do j = jmin,jmax |
182 |
|
|
do k = 1,nr |
183 |
|
|
xx_obcse0(j,k,bi,bj,iobcs) = xx_obcse1(j,k,bi,bj,iobcs) |
184 |
|
|
tmpfldyz (j,k,bi,bj) = 0. _d 0 |
185 |
|
|
enddo |
186 |
heimbach |
1.2 |
enddo |
187 |
mlosch |
1.10 |
enddo |
188 |
|
|
enddo |
189 |
|
|
|
190 |
|
|
call active_read_yz( fnameobcse, tmpfldyz, |
191 |
|
|
& (obcsecount1-1)*nobcs+iobcs, |
192 |
|
|
& doglobalread, ladinit, optimcycle, |
193 |
|
|
& mythid, xx_obcse_dummy ) |
194 |
heimbach |
1.2 |
|
195 |
|
|
do bj = jtlo,jthi |
196 |
mlosch |
1.9 |
do bi = itlo,ithi |
197 |
mmazloff |
1.11 |
#ifdef ALLOW_OBCS_CONTROL_MODES |
198 |
|
|
if (iobcs .gt. 2) then |
199 |
|
|
do j = jmin,jmax |
200 |
|
|
i = OB_Ie(j,bi,bj) |
201 |
|
|
cih Determine number of open vertical layers. |
202 |
|
|
nz = 0 |
203 |
|
|
do k = 1,Nr |
204 |
|
|
if (iobcs .eq. 3) then |
205 |
|
|
nz = nz + maskW(i+ip1,j,k,bi,bj) |
206 |
|
|
else |
207 |
|
|
nz = nz + maskS(i,j,k,bi,bj) |
208 |
|
|
endif |
209 |
|
|
end do |
210 |
|
|
cih Compute absolute velocities from the barotropic-baroclinic modes. |
211 |
|
|
do k = 1,Nr |
212 |
|
|
if (k.le.nz) then |
213 |
|
|
stmp = 0. |
214 |
|
|
do nk = 1,nz |
215 |
|
|
stmp = stmp + |
216 |
|
|
& modesv(k,nk,nz)*tmpfldyz(j,nk,bi,bj) |
217 |
|
|
end do |
218 |
|
|
tmpz(k,bi,bj) = stmp |
219 |
|
|
else |
220 |
|
|
tmpz(k,bi,bj) = 0. |
221 |
|
|
endif |
222 |
|
|
enddo |
223 |
|
|
do k = 1,Nr |
224 |
|
|
if (iobcs .eq. 3) then |
225 |
|
|
tmpfldyz(j,k,bi,bj) = tmpz(k,bi,bj) |
226 |
|
|
& *recip_hFacW(i+ip1,j,k,bi,bj) |
227 |
|
|
else |
228 |
|
|
tmpfldyz(j,k,bi,bj) = tmpz(k,bi,bj) |
229 |
|
|
& *recip_hFacS(i,j,k,bi,bj) |
230 |
|
|
endif |
231 |
|
|
end do |
232 |
|
|
enddo |
233 |
|
|
endif |
234 |
|
|
#endif |
235 |
mlosch |
1.9 |
do k = 1,nr |
236 |
mlosch |
1.10 |
do j = jmin,jmax |
237 |
|
|
xx_obcse1 (j,k,bi,bj,iobcs) = tmpfldyz (j,k,bi,bj) |
238 |
|
|
cgg & * maskyz (j,k,bi,bj) |
239 |
heimbach |
1.2 |
enddo |
240 |
mlosch |
1.9 |
enddo |
241 |
|
|
enddo |
242 |
heimbach |
1.2 |
enddo |
243 |
mlosch |
1.10 |
endif |
244 |
|
|
|
245 |
|
|
c-- Add control to model variable. |
246 |
|
|
do bj = jtlo,jthi |
247 |
|
|
do bi = itlo,ithi |
248 |
|
|
c-- Calculate mask for tracer cells (0 => land, 1 => water). |
249 |
|
|
do k = 1,nr |
250 |
|
|
do j = 1,sny |
251 |
|
|
i = OB_Ie(j,bi,bj) |
252 |
|
|
if (iobcs .EQ. 1) then |
253 |
|
|
OBEt(j,k,bi,bj) = OBEt(j,k,bi,bj) |
254 |
|
|
& + obcsefac *xx_obcse0(j,k,bi,bj,iobcs) |
255 |
|
|
& + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs) |
256 |
|
|
OBEt(j,k,bi,bj) = OBEt(j,k,bi,bj) |
257 |
|
|
& *maskW(i+ip1,j,k,bi,bj) |
258 |
|
|
else if (iobcs .EQ. 2) then |
259 |
|
|
OBEs(j,k,bi,bj) = OBEs(j,k,bi,bj) |
260 |
|
|
& + obcsefac *xx_obcse0(j,k,bi,bj,iobcs) |
261 |
|
|
& + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs) |
262 |
|
|
OBEs(j,k,bi,bj) = OBEs(j,k,bi,bj) |
263 |
|
|
& *maskW(i+ip1,j,k,bi,bj) |
264 |
|
|
else if (iobcs .EQ. 3) then |
265 |
|
|
OBEu(j,k,bi,bj) = OBEu(j,k,bi,bj) |
266 |
|
|
& + obcsefac *xx_obcse0(j,k,bi,bj,iobcs) |
267 |
|
|
& + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs) |
268 |
|
|
OBEu(j,k,bi,bj) = OBEu(j,k,bi,bj) |
269 |
|
|
& *maskW(i+ip1,j,k,bi,bj) |
270 |
|
|
else if (iobcs .EQ. 4) then |
271 |
|
|
OBEv(j,k,bi,bj) = OBEv(j,k,bi,bj) |
272 |
|
|
& + obcsefac *xx_obcse0(j,k,bi,bj,iobcs) |
273 |
|
|
& + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs) |
274 |
|
|
OBEv(j,k,bi,bj) = OBEv(j,k,bi,bj) |
275 |
|
|
& *maskS(i,j,k,bi,bj) |
276 |
|
|
endif |
277 |
|
|
enddo |
278 |
|
|
enddo |
279 |
|
|
enddo |
280 |
|
|
enddo |
281 |
mlosch |
1.9 |
|
282 |
heimbach |
1.2 |
C-- End over iobcs loop |
283 |
|
|
enddo |
284 |
|
|
|
285 |
|
|
#else /* ALLOW_OBCSE_CONTROL undefined */ |
286 |
|
|
|
287 |
|
|
c == routine arguments == |
288 |
|
|
|
289 |
|
|
_RL mytime |
290 |
|
|
integer myiter |
291 |
|
|
integer mythid |
292 |
|
|
|
293 |
|
|
c-- CPP flag ALLOW_OBCSE_CONTROL undefined. |
294 |
|
|
|
295 |
|
|
#endif /* ALLOW_OBCSE_CONTROL */ |
296 |
|
|
|
297 |
|
|
end |
298 |
|
|
|