/[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.8 - (show annotations) (download)
Wed Jan 19 08:42:06 2011 UTC (13 years, 3 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint62t, checkpoint62s, checkpoint62r
Changes since 1.7: +8 -19 lines
simplify code so that it does not need to use S/R exf_swapfields

1 C $Header: /u/gcmpack/MITgcm/pkg/ctrl/ctrl_getobcss.F,v 1.7 2007/10/09 00:00:00 jmc 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
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 xx_obcss0(i,k,bi,bj,iobcs) = xx_obcss1(i,k,bi,bj,iobcs)
228 tmpfldxz (i,k,bi,bj) = 0. _d 0
229 enddo
230 enddo
231 enddo
232 enddo
233
234 call active_read_xz( fnameobcss, tmpfldxz,
235 & (obcsscount1-1)*nobcs+iobcs,
236 & doglobalread, ladinit, optimcycle,
237 & mythid, xx_obcss_dummy )
238
239 #ifdef ALLOW_CTRL_OBCS_BALANCE
240
241 if ( optimcycle .gt. 0) then
242 if (iobcs .eq. 3) then
243 cgg Special attention is needed for the normal velocity.
244 cgg For the north, this is the v velocity, iobcs = 4.
245 cgg This is done on a columnwise basis here.
246 do bj = jtlo,jthi
247 do bi = itlo, ithi
248 do i = imin,imax
249 j = OB_Js(I,bi,bj)
250
251 cgg The barotropic velocity is stored in the level 1.
252 vbaro = tmpfldxz(i,1,bi,bj)
253 tmpfldxz(i,1,bi,bj) = 0.d0
254 vtop = 0.d0
255
256 do k = 1,Nr
257 cgg If cells are not full, this should be modified with hFac.
258 cgg
259 cgg The xx field (tmpfldxz) does not contain the velocity at the
260 cgg surface level. This velocity is not independent; it must
261 cgg exactly balance the volume flux, since we are dealing with
262 cgg the baroclinic velocity structure..
263 vtop = tmpfldxz(i,k,bi,bj)*
264 & maskS(i,j+jp1,k,bi,bj) * delR(k) + vtop
265 cgg Add the barotropic velocity component.
266 if (maskS(i,j+jp1,k,bi,bj) .ne. 0.) then
267 tmpfldxz(i,k,bi,bj) = tmpfldxz(i,k,bi,bj)+ vbaro
268 endif
269 enddo
270 cgg Compute the baroclinic velocity at level 1. Should balance flux.
271 tmpfldxz(i,1,bi,bj) = tmpfldxz(i,1,bi,bj)
272 & - vtop / delR(1)
273 enddo
274 enddo
275 enddo
276 endif
277
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 i = imin,imax
285 j = OB_Js(I,bi,bj)
286
287 cgg The barotropic velocity is stored in the level 1.
288 vbaro = tmpfldxz(i,1,bi,bj)
289 tmpfldxz(i,1,bi,bj) = 0.d0
290 vtop = 0.d0
291
292 do k = 1,Nr
293 cgg If cells are not full, this should be modified with hFac.
294 cgg
295 cgg The xx field (tmpfldxz) does not contain the velocity at the
296 cgg surface level. This velocity is not independent; it must
297 cgg exactly balance the volume flux, since we are dealing with
298 cgg the baroclinic velocity structure..
299 vtop = tmpfldxz(i,k,bi,bj)*
300 & maskW(i,j,k,bi,bj) * delR(k) + vtop
301 cgg Add the barotropic velocity component.
302 if (maskW(i,j,k,bi,bj) .ne. 0.) then
303 tmpfldxz(i,k,bi,bj) = tmpfldxz(i,k,bi,bj)+ vbaro
304 endif
305 enddo
306 cgg Compute the baroclinic velocity at level 1. Should balance flux.
307 tmpfldxz(i,1,bi,bj) = tmpfldxz(i,1,bi,bj)
308 & - vtop / delR(1)
309 enddo
310 enddo
311 enddo
312 endif
313 endif
314
315 #endif /* ALLOW_CTRL_OBCS_BALANCE */
316
317 do bj = jtlo,jthi
318 do bi = itlo,ithi
319 do k = 1,nr
320 do i = imin,imax
321 xx_obcss1 (i,k,bi,bj,iobcs) = tmpfldxz (i,k,bi,bj)
322 cgg & * maskxz (i,k,bi,bj)
323 enddo
324 enddo
325 enddo
326 enddo
327 endif
328
329 c-- Add control to model variable.
330 do bj = jtlo,jthi
331 do bi = itlo,ithi
332 c-- Calculate mask for tracer cells (0 => land, 1 => water).
333 do k = 1,nr
334 do i = 1,snx
335 j = OB_Js(I,bi,bj)
336 if (iobcs .EQ. 1) then
337 OBSt(i,k,bi,bj) = OBSt (i,k,bi,bj)
338 & + obcssfac *xx_obcss0(i,k,bi,bj,iobcs)
339 & + (1. _d 0 - obcssfac)*xx_obcss1(i,k,bi,bj,iobcs)
340 OBSt(i,k,bi,bj) = OBSt(i,k,bi,bj)
341 & *maskS(i,j+jp1,k,bi,bj)
342 else if (iobcs .EQ. 2) then
343 OBSs(i,k,bi,bj) = OBSs (i,k,bi,bj)
344 & + obcssfac *xx_obcss0(i,k,bi,bj,iobcs)
345 & + (1. _d 0 - obcssfac)*xx_obcss1(i,k,bi,bj,iobcs)
346 OBSs(i,k,bi,bj) = OBSs(i,k,bi,bj)
347 & *maskS(i,j+jp1,k,bi,bj)
348 else if (iobcs .EQ. 4) then
349 OBSu(i,k,bi,bj) = OBSu (i,k,bi,bj)
350 & + obcssfac *xx_obcss0(i,k,bi,bj,iobcs)
351 & + (1. _d 0 - obcssfac)*xx_obcss1(i,k,bi,bj,iobcs)
352 OBSu(i,k,bi,bj) = OBSu(i,k,bi,bj)
353 & *maskW(i,j,k,bi,bj)
354 else if (iobcs .EQ. 3) then
355 OBSv(i,k,bi,bj) = OBSv (i,k,bi,bj)
356 & + obcssfac *xx_obcss0(i,k,bi,bj,iobcs)
357 & + (1. _d 0 - obcssfac)*xx_obcss1(i,k,bi,bj,iobcs)
358 OBSv(i,k,bi,bj) = OBSv(i,k,bi,bj)
359 & *maskS(i,j+jp1,k,bi,bj)
360 endif
361 enddo
362 enddo
363 enddo
364 enddo
365
366 C-- End over iobcs loop
367 enddo
368
369 #else /* ALLOW_OBCSS_CONTROL undefined */
370
371 c == routine arguments ==
372
373 _RL mytime
374 integer myiter
375 integer mythid
376
377 c-- CPP flag ALLOW_OBCSS_CONTROL undefined.
378
379 #endif /* ALLOW_OBCSS_CONTROL */
380
381 end
382

  ViewVC Help
Powered by ViewVC 1.1.22