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

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

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


Revision 1.5 - (show annotations) (download)
Fri Dec 3 00:48:57 2004 UTC (19 years, 5 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint57m_post, checkpoint57g_pre, checkpoint57s_post, checkpoint58b_post, checkpoint57b_post, checkpoint57g_post, checkpoint57y_post, checkpoint57h_post, checkpoint57r_post, checkpoint57d_post, checkpoint57i_post, checkpoint59, checkpoint58, checkpoint57, checkpoint58f_post, checkpoint57n_post, checkpoint58d_post, checkpoint58a_post, checkpoint57z_post, checkpoint58y_post, checkpoint58t_post, checkpoint58m_post, checkpoint57l_post, checkpoint57t_post, checkpoint57v_post, checkpoint57f_post, checkpoint57a_post, checkpoint57h_pre, checkpoint58w_post, checkpoint57y_pre, checkpoint58o_post, checkpoint57c_post, checkpoint58p_post, checkpoint58q_post, checkpoint58e_post, checkpoint57c_pre, checkpoint58r_post, checkpoint58n_post, checkpoint57e_post, checkpoint59a, checkpoint57p_post, checkpint57u_post, checkpoint57q_post, eckpoint57e_pre, checkpoint58k_post, checkpoint58v_post, checkpoint58l_post, checkpoint57h_done, checkpoint57j_post, checkpoint57f_pre, checkpoint58g_post, checkpoint58x_post, checkpoint58h_post, checkpoint56c_post, checkpoint58j_post, checkpoint57a_pre, checkpoint57o_post, checkpoint57k_post, checkpoint57w_post, checkpoint58i_post, checkpoint57x_post, checkpoint58c_post, checkpoint58u_post, checkpoint58s_post
Changes since 1.4: +8 -0 lines
o OBCS as control variables
  - update ad_diff.list
  - remove balance of obcs controls from default
  - fix index bug nobcs in ctrl_init
  - fix dummy fields filen in ctrl_pack
  - add dummy weights for obcs

1
2 #include "CTRL_CPPOPTIONS.h"
3 #ifdef ALLOW_OBCS
4 # include "OBCS_OPTIONS.h"
5 #endif
6
7
8 subroutine ctrl_getobcsn(
9 I mytime,
10 I myiter,
11 I mythid
12 & )
13
14 c ==================================================================
15 c SUBROUTINE ctrl_getobcsn
16 c ==================================================================
17 c
18 c o Get northern obc of the control vector and add it
19 c to dyn. fields
20 c
21 c started: heimbach@mit.edu, 29-Aug-2001
22 c
23 c modified: gebbie@mit.edu, 18-Mar-2003
24 c ==================================================================
25 c SUBROUTINE ctrl_getobcsn
26 c ==================================================================
27
28 implicit none
29
30 #ifdef ALLOW_OBCSN_CONTROL
31
32 c == global variables ==
33
34 #include "EEPARAMS.h"
35 #include "SIZE.h"
36 #include "PARAMS.h"
37 #include "GRID.h"
38 #include "OBCS.h"
39
40 #include "ctrl.h"
41 #include "ctrl_dummy.h"
42 #include "optim.h"
43
44 c == routine arguments ==
45
46 _RL mytime
47 integer myiter
48 integer mythid
49
50 c == local variables ==
51
52 integer bi,bj
53 integer i,j,k
54 integer itlo,ithi
55 integer jtlo,jthi
56 integer jmin,jmax
57 integer imin,imax
58 integer ilobcsn
59 integer iobcs
60 integer jp1
61
62 _RL dummy
63 _RL obcsnfac
64 logical obcsnfirst
65 logical obcsnchanged
66 integer obcsncount0
67 integer obcsncount1
68
69 cgg _RL maskxz (1-olx:snx+olx,nr,nsx,nsy)
70
71 logical doglobalread
72 logical ladinit
73
74 character*(80) fnameobcsn
75
76 cgg( Variables for splitting barotropic/baroclinic vels.
77 _RL vbaro
78 _RL vtop
79
80 c == external functions ==
81
82 integer ilnblnk
83 external ilnblnk
84
85
86 c == end of interface ==
87
88 jtlo = mybylo(mythid)
89 jthi = mybyhi(mythid)
90 itlo = mybxlo(mythid)
91 ithi = mybxhi(mythid)
92 cgg jmin = 1-oly
93 cgg jmax = sny+oly
94 cgg imin = 1-olx
95 cgg imax = snx+olx
96
97 jmin = 1
98 jmax = sny
99 imin = 1
100 imax = snx
101 jp1 = 0
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 ilobcsn=ilnblnk( xx_obcsn_file )
114 write(fnameobcsn(1:80),'(2a,i10.10)')
115 & xx_obcsn_file(1:ilobcsn), '.', optimcycle
116 endif
117
118 c-- Get the counters, flags, and the interpolation factor.
119 call ctrl_get_gen_rec(
120 I xx_obcsnstartdate, xx_obcsnperiod,
121 O obcsnfac, obcsnfirst, obcsnchanged,
122 O obcsncount0,obcsncount1,
123 I mytime, myiter, mythid )
124
125 do iobcs = 1,nobcs
126 if ( obcsnfirst ) then
127 call active_read_xz_loc( fnameobcsn, tmpfldxz,
128 & (obcsncount0-1)*nobcs+iobcs,
129 & doglobalread, ladinit, optimcycle,
130 & mythid, xx_obcsn_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
143 cgg The barotropic velocity is stored in the level 1.
144 vbaro = tmpfldxz(i,1,bi,bj)
145 cgg Except for the special point which balances barotropic vol.flux.
146 cgg Special column in the NW corner.
147 j = OB_Jn(I,bi,bj)
148 if (ob_iw(j,bi,bj).eq.(i-1).and.
149 & ob_iw(j,bi,bj).ne. 0) then
150 print*,'Apply shiftvel1 @ i,j'
151 print*,shiftvel(1),i,j
152 vbaro = shiftvel(1)
153 endif
154 tmpfldxz(i,1,bi,bj) = 0.d0
155 vtop = 0.d0
156
157 do k = 1,Nr
158 cgg If cells are not full, this should be modified with hFac.
159 cgg
160 cgg The xx field (tmpfldxz) does not contain the velocity at the
161 cgg surface level. This velocity is not independent; it must
162 cgg exactly balance the volume flux, since we are dealing with
163 cgg the baroclinic velocity structure..
164 vtop = tmpfldxz(i,k,bi,bj)*
165 & maskS(i,j+jp1,k,bi,bj) * delR(k) + vtop
166 cgg Add the barotropic velocity component.
167 if (maskS(i,j+jp1,k,bi,bj) .ne. 0.) then
168 tmpfldxz(i,k,bi,bj) = tmpfldxz(i,k,bi,bj)+ vbaro
169 endif
170 enddo
171 cgg Compute the baroclinic velocity at level 1. Should balance flux.
172 tmpfldxz(i,1,bi,bj) = tmpfldxz(i,1,bi,bj)
173 & - vtop / delR(1)
174 enddo
175 enddo
176 enddo
177 endif
178
179 if (iobcs .eq. 4) then
180 cgg Special attention is needed for the normal velocity.
181 cgg For the north, this is the v velocity, iobcs = 4.
182 cgg This is done on a columnwise basis here.
183 do bj = jtlo,jthi
184 do bi = itlo, ithi
185 do i = imin,imax
186
187 cgg The barotropic velocity is stored in the level 1.
188 vbaro = tmpfldxz(i,1,bi,bj)
189 cgg Except for the special point which balances barotropic vol.flux.
190 cgg Special column in the NW corner.
191 j = OB_Jn(I,bi,bj)
192 tmpfldxz(i,1,bi,bj) = 0.d0
193 vtop = 0.d0
194
195 do k = 1,Nr
196 cgg If cells are not full, this should be modified with hFac.
197 cgg
198 cgg The xx field (tmpfldxz) does not contain the velocity at the
199 cgg surface level. This velocity is not independent; it must
200 cgg exactly balance the volume flux, since we are dealing with
201 cgg the baroclinic velocity structure..
202 vtop = tmpfldxz(i,k,bi,bj)*
203 & maskW(i,j,k,bi,bj) * delR(k) + vtop
204 cgg Add the barotropic velocity component.
205 if (maskW(i,j,k,bi,bj) .ne. 0.) then
206 tmpfldxz(i,k,bi,bj) = tmpfldxz(i,k,bi,bj)+ vbaro
207 endif
208 enddo
209 cgg Compute the baroclinic velocity at level 1. Should balance flux.
210 tmpfldxz(i,1,bi,bj) = tmpfldxz(i,1,bi,bj)
211 & - vtop / delR(1)
212 enddo
213 enddo
214 enddo
215 endif
216
217 endif
218
219 #endif /* ALLOW_CTRL_OBCS_BALANCE */
220
221 do bj = jtlo,jthi
222 do bi = itlo,ithi
223 do k = 1,nr
224 do i = imin,imax
225 xx_obcsn1(i,k,bi,bj,iobcs) = tmpfldxz (i,k,bi,bj)
226 cgg & * maskxz (i,k,bi,bj)
227 enddo
228 enddo
229 enddo
230 enddo
231 endif
232
233 if ( (obcsnfirst) .or. (obcsnchanged)) then
234
235 do bj = jtlo,jthi
236 do bi = itlo,ithi
237 do k = 1,nr
238 do i = imin,imax
239 tmpfldxz(i,k,bi,bj) = xx_obcsn1(i,k,bi,bj,iobcs)
240 enddo
241 enddo
242 enddo
243 enddo
244
245 call exf_swapffields_xz( tmpfldxz2, tmpfldxz, mythid)
246
247 do bj = jtlo,jthi
248 do bi = itlo,ithi
249 do k = 1,nr
250 do i = imin,imax
251 xx_obcsn0(i,k,bi,bj,iobcs) = tmpfldxz2(i,k,bi,bj)
252 enddo
253 enddo
254 enddo
255 enddo
256
257 call active_read_xz_loc( fnameobcsn, tmpfldxz,
258 & (obcsncount1-1)*nobcs+iobcs,
259 & doglobalread, ladinit, optimcycle,
260 & mythid, xx_obcsn_dummy )
261
262 #ifdef ALLOW_CTRL_OBCS_BALANCE
263
264 if ( optimcycle .gt. 0) then
265 if (iobcs .eq. 3) then
266 cgg Special attention is needed for the normal velocity.
267 cgg For the north, this is the v velocity, iobcs = 3.
268 cgg This is done on a columnwise basis here.
269 do bj = jtlo,jthi
270 do bi = itlo, ithi
271 do i = imin,imax
272
273 cgg The barotropic velocity is stored in the level 1.
274 vbaro = tmpfldxz(i,1,bi,bj)
275 cgg Except for the special point which balances barotropic vol.flux.
276 cgg Special column in the NW corner.
277 j = OB_Jn(I,bi,bj)
278 if (ob_iw(j,bi,bj).eq.(i-1).and.
279 & ob_iw(j,bi,bj).ne. 0) then
280 print*,'correct vbaro',vbaro
281 print*,'Apply shiftvel2 @ i,j'
282 print*,shiftvel(2),i,j
283 vbaro = shiftvel(2)
284 endif
285 tmpfldxz(i,1,bi,bj) = 0.d0
286 vtop = 0.d0
287
288 do k = 1,Nr
289 cgg If cells are not full, this should be modified with hFac.
290 cgg
291 cgg The xx field (tmpfldxz) does not contain the velocity at the
292 cgg surface level. This velocity is not independent; it must
293 cgg exactly balance the volume flux, since we are dealing with
294 cgg the baroclinic velocity structure..
295 vtop = tmpfldxz(i,k,bi,bj)*
296 & maskS(i,j+jp1,k,bi,bj) * delR(k) + vtop
297 cgg Add the barotropic velocity component.
298 if (maskS(i,j+jp1,k,bi,bj) .ne. 0.) then
299 tmpfldxz(i,k,bi,bj) = tmpfldxz(i,k,bi,bj)+ vbaro
300 endif
301 enddo
302 cgg Compute the baroclinic velocity at level 1. Should balance flux.
303 tmpfldxz(i,1,bi,bj) = tmpfldxz(i,1,bi,bj)
304 & - vtop / delR(1)
305 enddo
306 enddo
307 enddo
308 endif
309 if (iobcs .eq. 4) then
310 cgg Special attention is needed for the normal velocity.
311 cgg For the north, this is the v velocity, iobcs = 3.
312 cgg This is done on a columnwise basis here.
313 do bj = jtlo,jthi
314 do bi = itlo, ithi
315 do i = imin,imax
316
317 cgg The barotropic velocity is stored in the level 1.
318 vbaro = tmpfldxz(i,1,bi,bj)
319 cgg Except for the special point which balances barotropic vol.flux.
320 cgg Special column in the NW corner.
321 j = OB_Jn(I,bi,bj)
322 tmpfldxz(i,1,bi,bj) = 0.d0
323 vtop = 0.d0
324
325 do k = 1,Nr
326 cgg If cells are not full, this should be modified with hFac.
327 cgg
328 cgg The xx field (tmpfldxz) does not contain the velocity at the
329 cgg surface level. This velocity is not independent; it must
330 cgg exactly balance the volume flux, since we are dealing with
331 cgg the baroclinic velocity structure..
332 vtop = tmpfldxz(i,k,bi,bj)*
333 & maskW(i,j,k,bi,bj) * delR(k) + vtop
334 cgg Add the barotropic velocity component.
335 if (maskW(i,j,k,bi,bj) .ne. 0.) then
336 tmpfldxz(i,k,bi,bj) = tmpfldxz(i,k,bi,bj)+ vbaro
337 endif
338 enddo
339 cgg Compute the baroclinic velocity at level 1. Should balance flux.
340 tmpfldxz(i,1,bi,bj) = tmpfldxz(i,1,bi,bj)
341 & - vtop / delR(1)
342 enddo
343 enddo
344 enddo
345 endif
346 endif
347
348 #endif /* ALLOW_CTRL_OBCS_BALANCE */
349
350 do bj = jtlo,jthi
351 do bi = itlo,ithi
352 do k = 1,nr
353 do i = imin,imax
354 xx_obcsn1 (i,k,bi,bj,iobcs) = tmpfldxz (i,k,bi,bj)
355 cgg & * maskxz (i,k,bi,bj)
356 enddo
357 enddo
358 enddo
359 enddo
360
361 endif
362
363 c-- Add control to model variable.
364 do bj = jtlo,jthi
365 do bi = itlo,ithi
366 c-- Calculate mask for tracer cells (0 => land, 1 => water).
367 do k = 1,nr
368 do i = 1,snx
369 j = OB_Jn(I,bi,bj)
370 if (iobcs .EQ. 1) then
371 OBNt(i,k,bi,bj) = OBNt (i,k,bi,bj)
372 & + obcsnfac *xx_obcsn0(i,k,bi,bj,iobcs)
373 & + (1. _d 0 - obcsnfac)*xx_obcsn1(i,k,bi,bj,iobcs)
374 OBNt(i,k,bi,bj) = OBNt(i,k,bi,bj)
375 & *maskS(i,j+jp1,k,bi,bj)
376 cgg & *maskxz(i,k,bi,bj)
377 else if (iobcs .EQ. 2) then
378 OBNs(i,k,bi,bj) = OBNs (i,k,bi,bj)
379 & + obcsnfac *xx_obcsn0(i,k,bi,bj,iobcs)
380 & + (1. _d 0 - obcsnfac)*xx_obcsn1(i,k,bi,bj,iobcs)
381 OBNs(i,k,bi,bj) = OBNs(i,k,bi,bj)
382 & *maskS(i,j+jp1,k,bi,bj)
383 cgg & *maskxz(i,k,bi,bj)
384 else if (iobcs .EQ. 4) then
385 OBNu(i,k,bi,bj) = OBNu (i,k,bi,bj)
386 & + obcsnfac *xx_obcsn0(i,k,bi,bj,iobcs)
387 & + (1. _d 0 - obcsnfac)*xx_obcsn1(i,k,bi,bj,iobcs)
388 OBNu(i,k,bi,bj) = OBNu(i,k,bi,bj)
389 & *maskW(i,j,k,bi,bj)
390 cgg & *maskxz(i,k,bi,bj)
391 else if (iobcs .EQ. 3) then
392 OBNv(i,k,bi,bj) = OBNv (i,k,bi,bj)
393 & + obcsnfac *xx_obcsn0(i,k,bi,bj,iobcs)
394 & + (1. _d 0 - obcsnfac)*xx_obcsn1(i,k,bi,bj,iobcs)
395 OBNv(i,k,bi,bj) = OBNv(i,k,bi,bj)
396 & *maskS(i,j+jp1,k,bi,bj)
397 cgg & *maskxz(i,k,bi,bj)
398 endif
399 enddo
400 enddo
401 enddo
402 enddo
403
404 C-- End over iobcs loop
405 enddo
406
407 #else /* ALLOW_OBCSN_CONTROL undefined */
408
409 c == routine arguments ==
410
411 _RL mytime
412 integer myiter
413 integer mythid
414
415 c-- CPP flag ALLOW_OBCSN_CONTROL undefined.
416
417 #endif /* ALLOW_OBCSN_CONTROL */
418
419 end
420
421
422
423
424

  ViewVC Help
Powered by ViewVC 1.1.22