/[MITgcm]/MITgcm_contrib/sannino/GRID_Refinemet/code/obcs_calc.F
ViewVC logotype

Annotation of /MITgcm_contrib/sannino/GRID_Refinemet/code/obcs_calc.F

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


Revision 1.1 - (hide annotations) (download)
Thu Jul 20 21:08:14 2006 UTC (19 years, 1 month ago) by sannino
Branch: MAIN
CVS Tags: HEAD
o Adding OASIS package
o Adding grid refinement package

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

  ViewVC Help
Powered by ViewVC 1.1.22