/[MITgcm]/MITgcm_contrib/nesting_sannino/code_nest_merged/obcs_calc.F
ViewVC logotype

Contents of /MITgcm_contrib/nesting_sannino/code_nest_merged/obcs_calc.F

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


Revision 1.5 - (show annotations) (download)
Mon Dec 6 14:34:33 2010 UTC (15 years, 1 month ago) by jmc
Branch: MAIN
CVS Tags: HEAD
Changes since 1.4: +1 -1 lines
FILE REMOVED
nesting specific modification have been merged into main code:
 remove the last 2 remaining S/R.

1 C $Header: /u/gcmpack/MITgcm_contrib/nesting_sannino/code_nest_merged/obcs_calc.F,v 1.4 2010/11/28 01:56:55 jmc Exp $
2 C $Name: $
3
4 #include "OBCS_OPTIONS.h"
5
6 SUBROUTINE OBCS_CALC( futureTime, futureIter,
7 & uVel, vVel, wVel, theta, salt,
8 & myThid )
9 C *==========================================================*
10 C | SUBROUTINE OBCS_CALC |
11 C | o Calculate future boundary data at open boundaries |
12 C | at time = futureTime |
13 C *==========================================================*
14 C | |
15 C *==========================================================*
16 IMPLICIT NONE
17
18 C === Global variables ===
19 #include "SIZE.h"
20 #include "EEPARAMS.h"
21 #include "PARAMS.h"
22 #include "GRID.h"
23 #include "OBCS.h"
24 #ifdef ALLOW_PTRACERS
25 # include "PTRACERS_SIZE.h"
26 # include "PTRACERS_PARAMS.h"
27 # include "PTRACERS_FIELDS.h"
28 # include "OBCS_PTRACERS.h"
29 #endif /* ALLOW_PTRACERS */
30 #ifdef ALLOW_NEST_CHILD
31 # include "NEST_CHILD.h"
32 #endif /* ALLOW_NEST_CHILD */
33
34 C == Routine arguments ==
35 INTEGER futureIter
36 _RL futureTime
37 _RL uVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
38 _RL vVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
39 _RL wVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
40 _RL theta(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
41 _RL salt (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
42 INTEGER myThid
43
44 #ifdef ALLOW_OBCS
45
46 C == Local variables ==
47 C bi, bj :: tile indices
48 C I,J,K :: loop indices
49 C I_obc, J_obc :: local index of open boundary
50 C msgBuf :: Informational/error message buffer
51 INTEGER bi, bj
52 INTEGER I, J, K, I_obc, J_obc
53 #ifdef ALLOW_OBCS_BALANCE
54 _RL Tr_T, Ar_T, Ar
55 #endif /* ALLOW_OBCS_BALANCE */
56 #ifdef ALLOW_PTRACERS
57 CHARACTER*(MAX_LEN_MBUF) msgBuf
58 INTEGER iTracer
59 #endif /* ALLOW_PTRACERS */
60
61 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
62
63 #ifdef ALLOW_DEBUG
64 IF (debugMode) CALL DEBUG_ENTER('OBCS_CALC',myThid)
65 #endif
66
67 DO bj=myByLo(myThid),myByHi(myThid)
68 DO bi=myBxLo(myThid),myBxHi(myThid)
69
70 #ifdef ALLOW_NEST_CHILD
71 IF ( useNEST_CHILD ) THEN
72 IF ( PASSI.LT.2 ) THEN
73 CALL NEST_CHILD_RECV ( myThid )
74 ENDIF
75 ENDIF
76 #endif /* ALLOW_NEST_CHILD */
77
78 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
79
80 #ifdef ALLOW_OBCS_EAST
81 C Eastern OB
82 #ifdef ALLOW_DEBUG
83 IF (debugMode) CALL DEBUG_MSG('OBCS_CALC: East',myThid)
84 #endif
85 IF (useOrlanskiEast) THEN
86 #ifdef ALLOW_ORLANSKI
87 CALL ORLANSKI_EAST(
88 & bi, bj, futureTime,
89 & uVel, vVel, wVel, theta, salt,
90 & myThid )
91 #endif
92 #ifdef ALLOW_NEST_CHILD
93 ELSEIF ( useNEST_CHILD ) THEN
94 DO k=1,Nr
95 DO j=1-Oly,sNy+Oly
96 I_obc=OB_Ie(j,bi,bj)
97 IF (I_obc.ne.0) THEN
98 OBEu(j,k,bi,bj)= U_F1(j,k,2)
99 OBEv(j,k,bi,bj)= V_F1(j,k,2)
100 OBEt(j,k,bi,bj)= T_F1(j,k,2)
101 OBEs(j,k,bi,bj)= S_F1(j,k,2)
102 #ifdef NONLIN_FRSURF
103 OBEeta(j,bi,bj)= ETA_F1(j,1,2)
104 #endif
105 ENDIF
106 ENDDO
107 ENDDO
108 #endif /* ALLOW_NEST_CHILD */
109 ELSE
110 DO K=1,Nr
111 DO J=1-Oly,sNy+Oly
112 I_obc=OB_Ie(J,bi,bj)
113 IF (I_obc.ne.0) THEN
114 OBEu(J,K,bi,bj)=0.
115 OBEv(J,K,bi,bj)=0.
116 OBEt(J,K,bi,bj)=tRef(K)
117 OBEs(J,K,bi,bj)=sRef(K)
118 #ifdef ALLOW_NONHYDROSTATIC
119 OBEw(J,K,bi,bj)=0.
120 #endif
121 #ifdef NONLIN_FRSURF
122 OBEeta(J,bi,bj)=0.
123 #endif
124 ENDIF
125 ENDDO
126 ENDDO
127 ENDIF
128 #endif /* ALLOW_OBCS_EAST */
129
130 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
131
132 #ifdef ALLOW_OBCS_WEST
133 C Western OB
134 #ifdef ALLOW_DEBUG
135 IF (debugMode) CALL DEBUG_MSG('OBCS_CALC: West',myThid)
136 #endif
137 IF (useOrlanskiWest) THEN
138 #ifdef ALLOW_ORLANSKI
139 CALL ORLANSKI_WEST(
140 & bi, bj, futureTime,
141 & uVel, vVel, wVel, theta, salt,
142 & myThid )
143 #endif
144 #ifdef ALLOW_NEST_CHILD
145 ELSEIF ( useNEST_CHILD ) THEN
146 DO k=1,Nr
147 DO j=1-Oly,sNy+Oly
148 I_obc=OB_Iw(j,bi,bj)
149 IF (I_obc.ne.0) THEN
150 OBWu(j,k,bi,bj)= U_F1(j,k,1)
151 OBWv(j,k,bi,bj)= V_F1(j,k,1)
152 OBWt(j,k,bi,bj)= T_F1(j,k,1)
153 OBWs(j,k,bi,bj)= S_F1(j,k,1)
154 #ifdef NONLIN_FRSURF
155 OBWeta(j,bi,bj)= ETA_F1(j,1,1)
156 #endif
157 ENDIF
158 ENDDO
159 ENDDO
160 #endif /* ALLOW_NEST_CHILD */
161 ELSE
162 DO K=1,Nr
163 DO J=1-Oly,sNy+Oly
164 I_obc=OB_Iw(J,bi,bj)
165 IF (I_obc.ne.0) THEN
166 OBWu(J,K,bi,bj)=0.
167 OBWv(J,K,bi,bj)=0.
168 OBWt(J,K,bi,bj)=tRef(K)
169 OBWs(J,K,bi,bj)=sRef(K)
170 #ifdef ALLOW_NONHYDROSTATIC
171 OBWw(J,K,bi,bj)=0.
172 #endif
173 #ifdef NONLIN_FRSURF
174 OBWeta(J,bi,bj)=0.
175 #endif
176 ENDIF
177 ENDDO
178 ENDDO
179 ENDIF
180 #endif /* ALLOW_OBCS_WEST */
181
182 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
183
184 #ifdef ALLOW_OBCS_NORTH
185 C Northern OB
186 #ifdef ALLOW_DEBUG
187 IF (debugMode) CALL DEBUG_MSG('OBCS_CALC: North',myThid)
188 #endif
189 IF (useOrlanskiNorth) THEN
190 #ifdef ALLOW_ORLANSKI
191 CALL ORLANSKI_NORTH(
192 & bi, bj, futureTime,
193 & uVel, vVel, wVel, theta, salt,
194 & myThid )
195 #endif
196 ELSE
197 DO K=1,Nr
198 DO I=1-Olx,sNx+Olx
199 J_obc=OB_Jn(I,bi,bj)
200 IF (J_obc.ne.0) THEN
201 OBNv(I,K,bi,bj)=0.
202 OBNu(I,K,bi,bj)=0.
203 OBNt(I,K,bi,bj)=tRef(K)
204 OBNs(I,K,bi,bj)=sRef(K)
205 #ifdef ALLOW_NONHYDROSTATIC
206 OBNw(I,K,bi,bj)=0.
207 #endif
208 #ifdef NONLIN_FRSURF
209 OBNeta(I,bi,bj)=0.
210 #endif
211 ENDIF
212 ENDDO
213 ENDDO
214 ENDIF
215 #endif /* ALLOW_OBCS_NORTH */
216
217 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
218
219 #ifdef ALLOW_OBCS_SOUTH
220 C Southern OB
221 #ifdef ALLOW_DEBUG
222 IF (debugMode) CALL DEBUG_MSG('OBCS_CALC: South',myThid)
223 #endif
224 IF (useOrlanskiSouth) THEN
225 #ifdef ALLOW_ORLANSKI
226 CALL ORLANSKI_SOUTH(
227 & bi, bj, futureTime,
228 & uVel, vVel, wVel, theta, salt,
229 & myThid )
230 #endif
231 ELSE
232 DO K=1,Nr
233 DO I=1-Olx,sNx+Olx
234 J_obc=OB_Js(I,bi,bj)
235 IF (J_obc.ne.0) THEN
236 OBSu(I,K,bi,bj)=0.
237 OBSv(I,K,bi,bj)=0.
238 OBSt(I,K,bi,bj)=tRef(K)
239 OBSs(I,K,bi,bj)=sRef(K)
240 #ifdef ALLOW_NONHYDROSTATIC
241 OBSw(I,K,bi,bj)=0.
242 #endif
243 #ifdef NONLIN_FRSURF
244 OBSeta(I,bi,bj)=0.
245 #endif
246 ENDIF
247 ENDDO
248 ENDDO
249 ENDIF
250 #endif /* ALLOW_OBCS_SOUTH */
251
252 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
253
254 #ifdef ALLOW_PTRACERS
255 IF ( usePTRACERS ) THEN
256 C
257 C Calculate some default open boundary conditions for passive tracers:
258 C The default is a homogeneous v.Neumann conditions, that is, the
259 C tracer gradient across the open boundary is nearly zero;
260 C only nearly, because the boundary conditions are applied throughout
261 C the time step during which the interior field does change; therefore
262 C we have to use the values from the previous time step here. If you
263 C really want exact v.Neumann conditions, you have to modify
264 C obcs_apply_ptracer directly.
265 C
266 # ifdef ALLOW_OBCS_EAST
267 C Eastern OB
268 # ifdef ALLOW_DEBUG
269 IF (debugMode)
270 & CALL DEBUG_MSG('OBCS_CALC: East, pTracers',myThid)
271 # endif
272 IF (useOrlanskiEast) THEN
273 WRITE(msgBuf,'(A)')
274 & 'OBCS_CALC: ERROR: useOrlanskiEast Rad OBC with'
275 CALL PRINT_ERROR( msgBuf, myThid )
276 WRITE(msgBuf,'(A)')
277 & 'OBCS_CALC: ERROR: pTracers not yet implemented'
278 CALL PRINT_ERROR( msgBuf, myThid )
279 STOP 'ABNORMAL END: S/R OBCS_CALC'
280 ELSE
281 DO iTracer=1,PTRACERS_numInUse
282 DO K=1,Nr
283 DO J=1-Oly,sNy+Oly
284 I_obc=OB_Ie(J,bi,bj)
285 IF (I_obc.ne.0) THEN
286 OBEptr(J,K,bi,bj,iTracer) =
287 & pTracer(I_obc-1,J,K,bi,bj,iTracer)
288 & *_maskW(I_obc,J,K,bi,bj)
289 ENDIF
290 ENDDO
291 ENDDO
292 ENDDO
293 ENDIF
294 # endif /* ALLOW_OBCS_EAST */
295
296 C ------------------------------------------------------------------------------
297
298 # ifdef ALLOW_OBCS_WEST
299 C Western OB
300 # ifdef ALLOW_DEBUG
301 IF (debugMode)
302 & CALL DEBUG_MSG('OBCS_CALC: West, pTracers',myThid)
303 # endif
304 IF (useOrlanskiWest) THEN
305 WRITE(msgBuf,'(A)')
306 & 'OBCS_CALC: ERROR: useOrlanskiWest Rad OBC with'
307 CALL PRINT_ERROR( msgBuf, myThid )
308 WRITE(msgBuf,'(A)')
309 & 'OBCS_CALC: ERROR: pTracers not yet implemented'
310 CALL PRINT_ERROR( msgBuf, myThid )
311 STOP 'ABNORMAL END: S/R OBCS_CALC'
312 ELSE
313 DO iTracer=1,PTRACERS_numInUse
314 DO K=1,Nr
315 DO J=1-Oly,sNy+Oly
316 I_obc=OB_Iw(J,bi,bj)
317 IF (I_obc.ne.0) THEN
318 OBWptr(J,K,bi,bj,iTracer) =
319 & pTracer(I_obc+1,J,K,bi,bj,iTracer)
320 & *_maskW(I_obc+1,J,K,bi,bj)
321 ENDIF
322 ENDDO
323 ENDDO
324 ENDDO
325 ENDIF
326 # endif /* ALLOW_OBCS_WEST */
327
328 C ------------------------------------------------------------------------------
329
330 # ifdef ALLOW_OBCS_NORTH
331 C Northern OB
332 # ifdef ALLOW_DEBUG
333 IF (debugMode)
334 & CALL DEBUG_MSG('OBCS_CALC: North, pTracers',myThid)
335 # endif
336 IF (useOrlanskiNorth) THEN
337 WRITE(msgBuf,'(A)')
338 & 'OBCS_CALC: ERROR: useOrlanskiNorth Rad OBC with'
339 CALL PRINT_ERROR( msgBuf, myThid )
340 WRITE(msgBuf,'(A)')
341 & 'OBCS_CALC: ERROR: pTracers not yet implemented'
342 CALL PRINT_ERROR( msgBuf, myThid )
343 STOP 'ABNORMAL END: S/R OBCS_CALC'
344 ELSE
345 DO iTracer=1,PTRACERS_numInUse
346 DO K=1,Nr
347 DO I=1-Olx,sNx+Olx
348 J_obc=OB_Jn(I,bi,bj)
349 IF (J_obc.ne.0) THEN
350 OBNptr(I,K,bi,bj,iTracer) =
351 & pTracer(I,J_obc-1,K,bi,bj,iTracer)
352 & *_maskS(I,J_obc,K,bi,bj)
353 ENDIF
354 ENDDO
355 ENDDO
356 ENDDO
357 ENDIF
358 # endif /* ALLOW_OBCS_NORTH */
359
360 C ------------------------------------------------------------------------------
361
362 # ifdef ALLOW_OBCS_SOUTH
363 C Southern OB
364 # ifdef ALLOW_DEBUG
365 IF (debugMode)
366 & CALL DEBUG_MSG('OBCS_CALC: South, pTracers',myThid)
367 #endif
368 IF (useOrlanskiSouth) THEN
369 WRITE(msgBuf,'(A)')
370 & 'OBCS_CALC: ERROR: useOrlanskiSouth Rad OBC with'
371 CALL PRINT_ERROR( msgBuf, myThid )
372 WRITE(msgBuf,'(A)')
373 & 'OBCS_CALC: ERROR: pTracers not yet implemented'
374 CALL PRINT_ERROR( msgBuf, myThid )
375 STOP 'ABNORMAL END: S/R OBCS_CALC'
376 ELSE
377 DO iTracer=1,PTRACERS_numInUse
378 DO K=1,Nr
379 DO I=1-Olx,sNx+Olx
380 J_obc=OB_Js(I,bi,bj)
381 IF (J_obc.ne.0) THEN
382 OBSptr(I,K,bi,bj,iTracer) =
383 & pTracer(I,J_obc+1,K,bi,bj,iTracer)
384 & *_maskS(I,J_obc+1,K,bi,bj)
385 ENDIF
386 ENDDO
387 ENDDO
388 ENDDO
389 ENDIF
390 # endif /* ALLOW_OBCS_SOUTH */
391 C end if (usePTracers)
392 ENDIF
393 #endif /* ALLOW_PTRACERS */
394
395 C-- end bi,bj loops.
396 ENDDO
397 ENDDO
398
399 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
400
401 #ifdef ALLOW_OBCS_PRESCRIBE
402 IF (useOBCSprescribe) THEN
403 C-- Calculate future values on open boundaries
404 #ifdef ALLOW_DEBUG
405 IF (debugMode) CALL DEBUG_CALL('OBCS_PRESCRIBE_READ',myThid)
406 #endif
407 CALL OBCS_PRESCRIBE_READ( futureTime, futureIter, myThid )
408 ENDIF
409 #endif /* ALLOW_OBCS_PRESCRIBE */
410
411 C ------------------------------------------------------------------------------
412
413 #ifndef ALLOW_AUTODIFF_TAMC
414 C I do not think that we want this to be differentiated for now
415 #ifdef ALLOW_OBCS_STEVENS
416 C The Stevens (1990) boundary conditions come after reading data from
417 C files, because they use the data to compute a mix of simplified
418 C Orlanski and prescribed boundary conditions
419 IF (useStevensNorth.OR.useStevensSouth.OR.
420 & useStevensEast.OR.useStevensWest) THEN
421 #ifdef ALLOW_DEBUG
422 IF (debugMode) CALL DEBUG_CALL('OBCS_CALC_STEVENS',myThid)
423 #endif
424 CALL OBCS_CALC_STEVENS( futureTime, futureIter, myThid )
425 ENDIF
426 #endif /* ALLOW_OBCS_STEVENS */
427 #endif /* ndef ALLOW_AUTODIFF_TAMC */
428
429 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
430
431 #ifdef ALLOW_OBCS_BALANCE
432 IF ( useOBCSbalance) THEN
433 #ifdef ALLOW_DEBUG
434 IF (debugMode) CALL DEBUG_MSG('useOBCSbalance=.TRUE.',myThid)
435 #endif
436
437 #ifdef ALLOW_OBCS_EAST
438 Tr_T = 0. _d 0
439 Ar_T = 0. _d 0
440 DO bj=myByLo(myThid),myByHi(myThid)
441 DO bi=myBxLo(myThid),myBxHi(myThid)
442 DO K=1,Nr
443 DO J=1,sNy
444 I_obc=OB_Ie(J,bi,bj)
445 IF (I_obc.ne.0) THEN
446 Ar = drF(k)*hFacC(I_obc,j,k,bi,bj)*dyG(I_obc,j,bi,bj)
447 Ar_T = Ar_T + Ar
448 Tr_T = Tr_T + Ar * OBEu(J,K,bi,bj)
449 ENDIF
450 ENDDO
451 ENDDO
452 ENDDO
453 ENDDO
454 _GLOBAL_SUM_RL( Ar_T , myThid )
455 IF ( Ar_T .GT. 0. _d 0 ) THEN
456 _GLOBAL_SUM_RL( Tr_T , myThid )
457 Tr_T = (0. _d 0 - Tr_T)/Ar_T
458 DO bj=myByLo(myThid),myByHi(myThid)
459 DO bi=myBxLo(myThid),myBxHi(myThid)
460 DO K=1,Nr
461 DO J=1-Oly,sNy+Oly
462 I_obc=OB_Ie(J,bi,bj)
463 IF (I_obc.ne.0) THEN
464 OBEu(J,K,bi,bj) = OBEu(J,K,bi,bj) + Tr_T
465 c OBEv(J,K,bi,bj) = 0.
466 ENDIF
467 ENDDO
468 ENDDO
469 ENDDO
470 ENDDO
471 ENDIF
472 #endif
473
474 #ifdef ALLOW_OBCS_WEST
475 Tr_T = 0. _d 0
476 Ar_T = 0. _d 0
477 DO bj=myByLo(myThid),myByHi(myThid)
478 DO bi=myBxLo(myThid),myBxHi(myThid)
479 DO K=1,Nr
480 DO J=1,sNy
481 I_obc=OB_Iw(J,bi,bj)
482 IF (I_obc.ne.0) THEN
483 Ar = drF(k)*hFacC(I_obc,j,k,bi,bj)*dyG(I_obc,j,bi,bj)
484 Ar_T = Ar_T + Ar
485 Tr_T = Tr_T + Ar * OBWu(J,K,bi,bj)
486 ENDIF
487 ENDDO
488 ENDDO
489 ENDDO
490 ENDDO
491 _GLOBAL_SUM_RL( Ar_T , myThid )
492 IF ( Ar_T .GT. 0. _d 0 ) THEN
493 _GLOBAL_SUM_RL( Tr_T , myThid )
494 Tr_T = (0. _d 0 - Tr_T)/Ar_T
495 DO bj=myByLo(myThid),myByHi(myThid)
496 DO bi=myBxLo(myThid),myBxHi(myThid)
497 DO K=1,Nr
498 DO J=1-Oly,sNy+Oly
499 I_obc=OB_Iw(J,bi,bj)
500 IF (I_obc.ne.0) THEN
501 OBWu(J,K,bi,bj) = OBWu(J,K,bi,bj) + Tr_T
502 c OBWv(J,K,bi,bj) = 0.
503 ENDIF
504 ENDDO
505 ENDDO
506 ENDDO
507 ENDDO
508 ENDIF
509 #endif
510
511 #ifdef ALLOW_OBCS_NORTH
512 Tr_T = 0. _d 0
513 Ar_T = 0. _d 0
514 DO bj=myByLo(myThid),myByHi(myThid)
515 DO bi=myBxLo(myThid),myBxHi(myThid)
516 DO K=1,Nr
517 DO I=1,sNx
518 J_obc=OB_Jn(I,bi,bj)
519 IF (J_obc.ne.0) THEN
520 Ar = drF(k)*hFacC(i,J_obc,k,bi,bj)*dxG(i,J_obc,bi,bj)
521 Ar_T = Ar_T + Ar
522 Tr_T = Tr_T + Ar * OBNv(I,K,bi,bj)
523 ENDIF
524 ENDDO
525 ENDDO
526 ENDDO
527 ENDDO
528 _GLOBAL_SUM_RL( Ar_T , myThid )
529 IF ( Ar_T .GT. 0. _d 0 ) THEN
530 _GLOBAL_SUM_RL( Tr_T , myThid )
531 Tr_T = (0. _d 0 - Tr_T)/Ar_T
532 DO bj=myByLo(myThid),myByHi(myThid)
533 DO bi=myBxLo(myThid),myBxHi(myThid)
534 DO K=1,Nr
535 DO I=1-Olx,sNx+Olx
536 J_obc=OB_Jn(I,bi,bj)
537 IF (J_obc.ne.0) THEN
538 c OBNu(I,K,bi,bj) = 0.
539 OBNv(I,K,bi,bj) = OBNv(I,K,bi,bj) + Tr_T
540 ENDIF
541 ENDDO
542 ENDDO
543 ENDDO
544 ENDDO
545 ENDIF
546 #endif
547
548 #ifdef ALLOW_OBCS_SOUTH
549 Tr_T = 0. _d 0
550 Ar_T = 0. _d 0
551 DO bj=myByLo(myThid),myByHi(myThid)
552 DO bi=myBxLo(myThid),myBxHi(myThid)
553 DO K=1,Nr
554 DO I=1,sNx
555 J_obc=OB_Js(I,bi,bj)
556 IF (J_obc.ne.0) THEN
557 Ar = drF(k)*hFacC(i,J_obc,k,bi,bj)*dxG(i,J_obc,bi,bj)
558 Ar_T = Ar_T + Ar
559 Tr_T = Tr_T + Ar * OBSv(I,K,bi,bj)
560 ENDIF
561 ENDDO
562 ENDDO
563 ENDDO
564 ENDDO
565 _GLOBAL_SUM_RL( Ar_T , myThid )
566 IF ( Ar_T .GT. 0. _d 0 ) THEN
567 _GLOBAL_SUM_RL( Tr_T , myThid )
568 Tr_T = (0. _d 0 - Tr_T)/Ar_T
569 DO bj=myByLo(myThid),myByHi(myThid)
570 DO bi=myBxLo(myThid),myBxHi(myThid)
571 DO K=1,Nr
572 DO I=1-Olx,sNx+Olx
573 J_obc=OB_Js(I,bi,bj)
574 IF (J_obc.ne.0) THEN
575 c OBSu(I,K,bi,bj) = 0.
576 OBSv(I,K,bi,bj) = OBSv(I,K,bi,bj) + Tr_T
577 ENDIF
578 ENDDO
579 ENDDO
580 ENDDO
581 ENDDO
582 ENDIF
583 #endif
584
585 ENDIF
586 #endif /* ALLOW_OBCS_BALANCE */
587
588 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
589
590 #ifdef ALLOW_DEBUG
591 IF (debugMode) CALL DEBUG_LEAVE('OBCS_CALC',myThid)
592 #endif
593 #endif /* ALLOW_OBCS */
594
595 RETURN
596 END

  ViewVC Help
Powered by ViewVC 1.1.22