/[MITgcm]/MITgcm/pkg/ctrl/ctrl_getobcss.F
ViewVC logotype

Contents of /MITgcm/pkg/ctrl/ctrl_getobcss.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.10 - (show annotations) (download)
Wed Mar 9 16:17:14 2011 UTC (13 years, 3 months ago) by mlosch
Branch: MAIN
Changes since 1.9: +2 -2 lines
fix a bug in the declaration of tmpfldxz (spotted by Holly)

1 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_getobcss.F,v 1.9 2011/03/07 09:24:10 mlosch Exp $
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 _RL tmpfldxz (1-olx:snx+olx,nr,nsx,nsy)
74
75 logical doglobalread
76 logical ladinit
77
78 character*(80) fnameobcss
79
80 cgg( Variables for splitting barotropic/baroclinic vels.
81 _RL vbaro
82 _RL vtop
83 cgg)
84
85 c == external functions ==
86
87 integer ilnblnk
88 external ilnblnk
89
90
91 c == end of interface ==
92
93 jtlo = mybylo(mythid)
94 jthi = mybyhi(mythid)
95 itlo = mybxlo(mythid)
96 ithi = mybxhi(mythid)
97 jmin = 1-oly
98 jmax = sny+oly
99 imin = 1-olx
100 imax = snx+olx
101 jp1 = 1
102
103 cgg( Initialize variables for balancing volume flux.
104 vbaro = 0.d0
105 vtop = 0.d0
106 cgg)
107
108 c-- Now, read the control vector.
109 doglobalread = .false.
110 ladinit = .false.
111
112 if (optimcycle .ge. 0) then
113 ilobcss=ilnblnk( xx_obcss_file )
114 write(fnameobcss(1:80),'(2a,i10.10)')
115 & xx_obcss_file(1:ilobcss), '.', optimcycle
116 endif
117
118 c-- Get the counters, flags, and the interpolation factor.
119 call ctrl_get_gen_rec(
120 I xx_obcssstartdate, xx_obcssperiod,
121 O obcssfac, obcssfirst, obcsschanged,
122 O obcsscount0,obcsscount1,
123 I mytime, myiter, mythid )
124
125 do iobcs = 1,nobcs
126 if ( obcssfirst ) then
127 call active_read_xz( fnameobcss, tmpfldxz,
128 & (obcsscount0-1)*nobcs+iobcs,
129 & doglobalread, ladinit, optimcycle,
130 & mythid, xx_obcss_dummy )
131
132 #ifdef ALLOW_CTRL_OBCS_BALANCE
133
134 if ( optimcycle .gt. 0) then
135 if (iobcs .eq. 3) then
136 cgg Special attention is needed for the normal velocity.
137 cgg For the north, this is the v velocity, iobcs = 4.
138 cgg This is done on a columnwise basis here.
139 do bj = jtlo,jthi
140 do bi = itlo, ithi
141 do i = imin,imax
142 j = OB_Js(I,bi,bj)
143
144 cgg The barotropic velocity is stored in the level 1.
145 vbaro = tmpfldxz(i,1,bi,bj)
146 tmpfldxz(i,1,bi,bj) = 0.d0
147 vtop = 0.d0
148
149 do k = 1,Nr
150 cgg If cells are not full, this should be modified with hFac.
151 cgg
152 cgg The xx field (tmpfldxz) does not contain the velocity at the
153 cgg surface level. This velocity is not independent; it must
154 cgg exactly balance the volume flux, since we are dealing with
155 cgg the baroclinic velocity structure..
156 vtop = tmpfldxz(i,k,bi,bj)*
157 & maskS(i,j+jp1,k,bi,bj) * delR(k) + vtop
158 cgg Add the barotropic velocity component.
159 if (maskS(i,j+jp1,k,bi,bj) .ne. 0.) then
160 tmpfldxz(i,k,bi,bj) = tmpfldxz(i,k,bi,bj)+ vbaro
161 endif
162 enddo
163 cgg Compute the baroclinic velocity at level 1. Should balance flux.
164 tmpfldxz(i,1,bi,bj) = tmpfldxz(i,1,bi,bj)
165 & - vtop / delR(1)
166 enddo
167 enddo
168 enddo
169 endif
170
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 i = imin,imax
178 j = OB_Js(I,bi,bj)
179
180 cgg The barotropic velocity is stored in the level 1.
181 vbaro = tmpfldxz(i,1,bi,bj)
182 tmpfldxz(i,1,bi,bj) = 0.d0
183 vtop = 0.d0
184
185 do k = 1,Nr
186 cgg If cells are not full, this should be modified with hFac.
187 cgg
188 cgg The xx field (tmpfldxz) does not contain the velocity at the
189 cgg surface level. This velocity is not independent; it must
190 cgg exactly balance the volume flux, since we are dealing with
191 cgg the baroclinic velocity structure..
192 vtop = tmpfldxz(i,k,bi,bj)*
193 & maskW(i,j,k,bi,bj) * delR(k) + vtop
194 cgg Add the barotropic velocity component.
195 if (maskW(i,j,k,bi,bj) .ne. 0.) then
196 tmpfldxz(i,k,bi,bj) = tmpfldxz(i,k,bi,bj)+ vbaro
197 endif
198 enddo
199 cgg Compute the baroclinic velocity at level 1. Should balance flux.
200 tmpfldxz(i,1,bi,bj) = tmpfldxz(i,1,bi,bj)
201 & - vtop / delR(1)
202 enddo
203 enddo
204 enddo
205 endif
206 endif
207
208 #endif /* ALLOW_CTRL_OBCS_BALANCE */
209
210 do bj = jtlo,jthi
211 do bi = itlo,ithi
212 do k = 1,nr
213 do i = imin,imax
214 xx_obcss1(i,k,bi,bj,iobcs) = tmpfldxz (i,k,bi,bj)
215 cgg & * maskxz (i,k,bi,bj)
216 enddo
217 enddo
218 enddo
219 enddo
220 endif
221
222 if ( (obcssfirst) .or. (obcsschanged)) then
223
224 do bj = jtlo,jthi
225 do bi = itlo,ithi
226 do k = 1,nr
227 do i = imin,imax
228 xx_obcss0(i,k,bi,bj,iobcs) = xx_obcss1(i,k,bi,bj,iobcs)
229 tmpfldxz (i,k,bi,bj) = 0. _d 0
230 enddo
231 enddo
232 enddo
233 enddo
234
235 call active_read_xz( fnameobcss, tmpfldxz,
236 & (obcsscount1-1)*nobcs+iobcs,
237 & doglobalread, ladinit, optimcycle,
238 & mythid, xx_obcss_dummy )
239
240 #ifdef ALLOW_CTRL_OBCS_BALANCE
241
242 if ( optimcycle .gt. 0) then
243 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 i = imin,imax
250 j = OB_Js(I,bi,bj)
251
252 cgg The barotropic velocity is stored in the level 1.
253 vbaro = tmpfldxz(i,1,bi,bj)
254 tmpfldxz(i,1,bi,bj) = 0.d0
255 vtop = 0.d0
256
257 do k = 1,Nr
258 cgg If cells are not full, this should be modified with hFac.
259 cgg
260 cgg The xx field (tmpfldxz) does not contain the velocity at the
261 cgg surface level. This velocity is not independent; it must
262 cgg exactly balance the volume flux, since we are dealing with
263 cgg the baroclinic velocity structure..
264 vtop = tmpfldxz(i,k,bi,bj)*
265 & maskS(i,j+jp1,k,bi,bj) * delR(k) + vtop
266 cgg Add the barotropic velocity component.
267 if (maskS(i,j+jp1,k,bi,bj) .ne. 0.) then
268 tmpfldxz(i,k,bi,bj) = tmpfldxz(i,k,bi,bj)+ vbaro
269 endif
270 enddo
271 cgg Compute the baroclinic velocity at level 1. Should balance flux.
272 tmpfldxz(i,1,bi,bj) = tmpfldxz(i,1,bi,bj)
273 & - vtop / delR(1)
274 enddo
275 enddo
276 enddo
277 endif
278
279 if (iobcs .eq. 4) then
280 cgg Special attention is needed for the normal velocity.
281 cgg For the north, this is the v velocity, iobcs = 4.
282 cgg This is done on a columnwise basis here.
283 do bj = jtlo,jthi
284 do bi = itlo, ithi
285 do i = imin,imax
286 j = OB_Js(I,bi,bj)
287
288 cgg The barotropic velocity is stored in the level 1.
289 vbaro = tmpfldxz(i,1,bi,bj)
290 tmpfldxz(i,1,bi,bj) = 0.d0
291 vtop = 0.d0
292
293 do k = 1,Nr
294 cgg If cells are not full, this should be modified with hFac.
295 cgg
296 cgg The xx field (tmpfldxz) does not contain the velocity at the
297 cgg surface level. This velocity is not independent; it must
298 cgg exactly balance the volume flux, since we are dealing with
299 cgg the baroclinic velocity structure..
300 vtop = tmpfldxz(i,k,bi,bj)*
301 & maskW(i,j,k,bi,bj) * delR(k) + vtop
302 cgg Add the barotropic velocity component.
303 if (maskW(i,j,k,bi,bj) .ne. 0.) then
304 tmpfldxz(i,k,bi,bj) = tmpfldxz(i,k,bi,bj)+ vbaro
305 endif
306 enddo
307 cgg Compute the baroclinic velocity at level 1. Should balance flux.
308 tmpfldxz(i,1,bi,bj) = tmpfldxz(i,1,bi,bj)
309 & - vtop / delR(1)
310 enddo
311 enddo
312 enddo
313 endif
314 endif
315
316 #endif /* ALLOW_CTRL_OBCS_BALANCE */
317
318 do bj = jtlo,jthi
319 do bi = itlo,ithi
320 do k = 1,nr
321 do i = imin,imax
322 xx_obcss1 (i,k,bi,bj,iobcs) = tmpfldxz (i,k,bi,bj)
323 cgg & * maskxz (i,k,bi,bj)
324 enddo
325 enddo
326 enddo
327 enddo
328 endif
329
330 c-- Add control to model variable.
331 do bj = jtlo,jthi
332 do bi = itlo,ithi
333 c-- Calculate mask for tracer cells (0 => land, 1 => water).
334 do k = 1,nr
335 do i = 1,snx
336 j = OB_Js(I,bi,bj)
337 if (iobcs .EQ. 1) then
338 OBSt(i,k,bi,bj) = OBSt (i,k,bi,bj)
339 & + obcssfac *xx_obcss0(i,k,bi,bj,iobcs)
340 & + (1. _d 0 - obcssfac)*xx_obcss1(i,k,bi,bj,iobcs)
341 OBSt(i,k,bi,bj) = OBSt(i,k,bi,bj)
342 & *maskS(i,j+jp1,k,bi,bj)
343 else if (iobcs .EQ. 2) then
344 OBSs(i,k,bi,bj) = OBSs (i,k,bi,bj)
345 & + obcssfac *xx_obcss0(i,k,bi,bj,iobcs)
346 & + (1. _d 0 - obcssfac)*xx_obcss1(i,k,bi,bj,iobcs)
347 OBSs(i,k,bi,bj) = OBSs(i,k,bi,bj)
348 & *maskS(i,j+jp1,k,bi,bj)
349 else if (iobcs .EQ. 4) then
350 OBSu(i,k,bi,bj) = OBSu (i,k,bi,bj)
351 & + obcssfac *xx_obcss0(i,k,bi,bj,iobcs)
352 & + (1. _d 0 - obcssfac)*xx_obcss1(i,k,bi,bj,iobcs)
353 OBSu(i,k,bi,bj) = OBSu(i,k,bi,bj)
354 & *maskW(i,j,k,bi,bj)
355 else if (iobcs .EQ. 3) then
356 OBSv(i,k,bi,bj) = OBSv (i,k,bi,bj)
357 & + obcssfac *xx_obcss0(i,k,bi,bj,iobcs)
358 & + (1. _d 0 - obcssfac)*xx_obcss1(i,k,bi,bj,iobcs)
359 OBSv(i,k,bi,bj) = OBSv(i,k,bi,bj)
360 & *maskS(i,j+jp1,k,bi,bj)
361 endif
362 enddo
363 enddo
364 enddo
365 enddo
366
367 C-- End over iobcs loop
368 enddo
369
370 #else /* ALLOW_OBCSS_CONTROL undefined */
371
372 c == routine arguments ==
373
374 _RL mytime
375 integer myiter
376 integer mythid
377
378 c-- CPP flag ALLOW_OBCSS_CONTROL undefined.
379
380 #endif /* ALLOW_OBCSS_CONTROL */
381
382 end
383

  ViewVC Help
Powered by ViewVC 1.1.22