/[MITgcm]/MITgcm/model/src/calc_gt.F
ViewVC logotype

Diff of /MITgcm/model/src/calc_gt.F

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

revision 1.13 by adcroft, Mon Jun 22 15:26:25 1998 UTC revision 1.27 by heimbach, Mon Nov 13 16:32:57 2000 UTC
# Line 1  Line 1 
1  C $Header$  C $Header$
2    
3  #include "CPP_EEOPTIONS.h"  #include "CPP_OPTIONS.h"
4    
5  CStartOfInterFace  CStartOfInterFace
6        SUBROUTINE CALC_GT(        SUBROUTINE CALC_GT(
7       I           bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,       I           bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
8       I           xA,yA,uTrans,vTrans,wTrans,maskup,maskC,       I           xA,yA,uTrans,vTrans,rTrans,maskup,maskC,
9       I           K13,K23,KappaZT,KapGM,       I           KappaRT,
10       U           af,df,fZon,fMer,fVerT,       U           af,df,fZon,fMer,fVerT,
11       I           myThid )       I           myCurrentTime, myThid )
12  C     /==========================================================\  C     /==========================================================\
13  C     | SUBROUTINE CALC_GT                                       |  C     | SUBROUTINE CALC_GT                                       |
14  C     | o Calculate the temperature tendency terms.              |  C     | o Calculate the temperature tendency terms.              |
# Line 43  C     == GLobal variables == Line 43  C     == GLobal variables ==
43  #include "PARAMS.h"  #include "PARAMS.h"
44  #include "GRID.h"  #include "GRID.h"
45  #include "FFIELDS.h"  #include "FFIELDS.h"
46    c #include "GM_ARRAYS.h"
47    
48    
49  C     == Routine arguments ==  C     == Routine arguments ==
50  C     fZon    - Work array for flux of temperature in the east-west  C     fZon    - Work array for flux of temperature in the east-west
# Line 57  C     xA      - Tracer cell face area no Line 59  C     xA      - Tracer cell face area no
59  C     yA      - Tracer cell face area normal to X  C     yA      - Tracer cell face area normal to X
60  C     uTrans  - Zonal volume transport through cell face  C     uTrans  - Zonal volume transport through cell face
61  C     vTrans  - Meridional volume transport through cell face  C     vTrans  - Meridional volume transport through cell face
62  C     wTrans  - Vertical volume transport through cell face  C     rTrans  - Vertical volume transport through cell face
63  C     af      - Advective flux component work array  C     af      - Advective flux component work array
64  C     df      - Diffusive flux component work array  C     df      - Diffusive flux component work array
65  C     bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation  C     bi, bj, iMin, iMax, jMin, jMax - Range of points for which calculation
# Line 70  C     myThid - Instance number for this Line 72  C     myThid - Instance number for this
72        _RS yA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS yA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
73        _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
74        _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
75        _RL wTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
76        _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
77        _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS maskC (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
78        _RL K13   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)        _RL KappaRT(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
       _RL K23   (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)  
       _RL KappaZT(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nz)  
       _RL KapGM (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
79        _RL af    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL af    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
80        _RL df    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL df    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
81        INTEGER k,kUp,kDown,kM1        INTEGER k,kUp,kDown,kM1
82        INTEGER bi,bj,iMin,iMax,jMin,jMax        INTEGER bi,bj,iMin,iMax,jMin,jMax
83        INTEGER myThid        INTEGER myThid
84          _RL     myCurrentTime
85  CEndOfInterface  CEndOfInterface
86    
87  C     == Local variables ==  C     == Local variables ==
# Line 91  C     I, J, K - Loop counters Line 91  C     I, J, K - Loop counters
91        _RL afFacT, dfFacT        _RL afFacT, dfFacT
92        _RL dTdx(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL dTdx(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
93        _RL dTdy(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL dTdy(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94          _RL df4   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95    
96    #ifdef ALLOW_AUTODIFF_TAMC
97    C--   only the kUp part of fverT is set in this subroutine
98    C--   the kDown is still required
99    
100          fVerT(1,1,kDown) = fVerT(1,1,kDown)
101          DO j=1-OLy,sNy+OLy
102           DO i=1-OLx,sNx+OLx
103            fZon(i,j)      = 0.0
104            fMer(i,j)      = 0.0
105            fVerT(i,j,kUp) = 0.0
106           ENDDO
107          ENDDO
108    #endif
109    
110        afFacT = 1. _d 0        afFacT = 1. _d 0
111        dfFacT = 1. _d 0        dfFacT = 1. _d 0
# Line 98  C     I, J, K - Loop counters Line 113  C     I, J, K - Loop counters
113    
114  C---  Calculate advective and diffusive fluxes between cells.  C---  Calculate advective and diffusive fluxes between cells.
115    
116  C--   Zonal flux (fZon is at west face of "theta" cell)  #ifdef INCLUDE_T_DIFFUSION_CODE
117  C     Advective component of zonal flux  C     o Zonal tracer gradient
118        DO j=jMin,jMax        DO j=1-Oly,sNy+Oly
119         DO i=iMin,iMax         DO i=1-Olx+1,sNx+Olx
120          af(i,j) =          dTdx(i,j) = _recip_dxC(i,j,bi,bj)*
121       &   uTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i-1,j,k,bi,bj))*0.5 _d 0       &  (theta(i,j,k,bi,bj)-theta(i-1,j,k,bi,bj))
122         ENDDO         ENDDO
123        ENDDO        ENDDO
124  C     Zonal tracer gradient  C     o Meridional tracer gradient
125        DO j=jMin,jMax        DO j=1-Oly+1,sNy+Oly
126         DO i=iMin,iMax         DO i=1-Olx,sNx+Olx
127          dTdx(i,j) = _rdxC(i,j,bi,bj)*          dTdy(i,j) = _recip_dyC(i,j,bi,bj)*
128       &  (theta(i,j,k,bi,bj)-theta(i-1,j,k,bi,bj))       &  (theta(i,j,k,bi,bj)-theta(i,j-1,k,bi,bj))
129         ENDDO         ENDDO
130        ENDDO        ENDDO
131  C     Diffusive component of zonal flux  
132    C--   del^2 of T, needed for bi-harmonic (del^4) term
133          IF (diffK4T .NE. 0.) THEN
134           DO j=1-Oly+1,sNy+Oly-1
135            DO i=1-Olx+1,sNx+Olx-1
136             df4(i,j)= _recip_hFacC(i,j,k,bi,bj)
137         &             *recip_drF(k)/_rA(i,j,bi,bj)
138         &            *(
139         &             +( xA(i+1,j)*dTdx(i+1,j)-xA(i,j)*dTdx(i,j) )
140         &             +( yA(i,j+1)*dTdy(i,j+1)-yA(i,j)*dTdy(i,j) )
141         &             )
142            ENDDO
143           ENDDO
144          ENDIF
145    #endif
146    
147    C--   Zonal flux (fZon is at west face of "theta" cell)
148    #ifdef INCLUDE_T_ADVECTION_CODE
149    C     o Advective component of zonal flux
150        DO j=jMin,jMax        DO j=jMin,jMax
151         DO i=iMin,iMax         DO i=iMin,iMax
152          df(i,j) = -(diffKhT+0.5*(KapGM(i,j)+KapGM(i-1,j)))*          af(i,j) =
153       &            xA(i,j)*dTdx(i,j)       &   uTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i-1,j,k,bi,bj))*0.5 _d 0
154         ENDDO         ENDDO
155        ENDDO        ENDDO
156  C     Net zonal flux  #endif /* INCLUDE_T_ADVECTION_CODE */
157    #ifdef INCLUDE_T_DIFFUSION_CODE
158    C     o Diffusive component of zonal flux
159          DO j=jMin,jMax
160           DO i=iMin,iMax
161            df(i,j) = -diffKhT*xA(i,j)*dTdx(i,j)
162           ENDDO
163          ENDDO
164    #ifdef ALLOW_GMREDI
165          IF (useGMRedi) CALL GMREDI_XTRANSPORT(
166         I     iMin,iMax,jMin,jMax,bi,bj,K,
167         I     xA,theta,
168         U     df,
169         I     myThid)
170    #endif
171    C     o Add the bi-harmonic contribution
172          IF (diffK4T .NE. 0.) THEN
173           DO j=jMin,jMax
174            DO i=iMin,iMax
175             df(i,j) = df(i,j) + xA(i,j)*
176         &    diffK4T*(df4(i,j)-df4(i-1,j))*_recip_dxC(i,j,bi,bj)
177            ENDDO
178           ENDDO
179          ENDIF
180    #endif /* INCLUDE_T_DIFFUSION_CODE */
181    C     o Net zonal flux
182        DO j=jMin,jMax        DO j=jMin,jMax
183         DO i=iMin,iMax         DO i=iMin,iMax
184          fZon(i,j) = afFacT*af(i,j) + dfFacT*df(i,j)          fZon(i,j) = 0.
185         & _ADT( + afFacT*af(i,j) )
186         & _LPT( + dfFacT*df(i,j) )
187         ENDDO         ENDDO
188        ENDDO        ENDDO
189    
190  C--   Meridional flux (fMer is at south face of "theta" cell)  C--   Meridional flux (fMer is at south face of "theta" cell)
191  C     Advective component of meridional flux  #ifdef INCLUDE_T_ADVECTION_CODE
192    C     o Advective component of meridional flux
193        DO j=jMin,jMax        DO j=jMin,jMax
194         DO i=iMin,iMax         DO i=iMin,iMax
 C       Advective component of meridional flux  
195          af(i,j) =          af(i,j) =
196       &   vTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i,j-1,k,bi,bj))*0.5 _d 0       &   vTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i,j-1,k,bi,bj))*0.5 _d 0
197         ENDDO         ENDDO
198        ENDDO        ENDDO
199  C     Zonal tracer gradient  #endif /* INCLUDE_T_ADVECTION_CODE */
200        DO j=jMin,jMax  #ifdef INCLUDE_T_DIFFUSION_CODE
201         DO i=iMin,iMax  C     o Diffusive component of meridional flux
202          dTdy(i,j) = _rdyC(i,j,bi,bj)*        DO j=jMin,jMax
203       &  (theta(i,j,k,bi,bj)-theta(i,j-1,k,bi,bj))         DO i=iMin,iMax
204         ENDDO          df(i,j) = -diffKhT*yA(i,j)*dTdy(i,j)
205        ENDDO         ENDDO
206  C     Diffusive component of meridional flux        ENDDO
207        DO j=jMin,jMax  #ifdef ALLOW_GMREDI
208         DO i=iMin,iMax        IF (useGMRedi) CALL GMREDI_YTRANSPORT(
209          df(i,j) = -(diffKhT+0.5*(KapGM(i,j)+KapGM(i,j-1)))*       I     iMin,iMax,jMin,jMax,bi,bj,K,
210       &            yA(i,j)*dTdy(i,j)       I     yA,theta,
211         U     df,
212         I     myThid)
213    #endif
214    C     o Add the bi-harmonic contribution
215          IF (diffK4T .NE. 0.) THEN
216           DO j=jMin,jMax
217            DO i=iMin,iMax
218             df(i,j) = df(i,j) + yA(i,j)*
219         &    diffK4T*(df4(i,j)-df4(i,j-1))*_recip_dyC(i,j,bi,bj)
220            ENDDO
221         ENDDO         ENDDO
222        ENDDO        ENDIF
223  C     Net meridional flux  #endif /* INCLUDE_T_DIFFUSION_CODE */
224    C     o Net meridional flux
225        DO j=jMin,jMax        DO j=jMin,jMax
226         DO i=iMin,iMax         DO i=iMin,iMax
227          fMer(i,j) = afFacT*af(i,j) + dfFacT*df(i,j)          fMer(i,j) = 0.
228         & _ADT( + afFacT*af(i,j) )
229         & _LPT( + dfFacT*df(i,j) )
230         ENDDO         ENDDO
231        ENDDO        ENDDO
232    
233  C--   Interpolate terms for Redi/GM scheme  #ifdef INCLUDE_T_DIFFUSION_CODE
234    C--   Terms that diffusion tensor projects onto z
235        DO j=jMin,jMax        DO j=jMin,jMax
236         DO i=iMin,iMax         DO i=iMin,iMax
237          dTdx(i,j) = 0.5*(          dTdx(i,j) = 0.5*(
238       &   +0.5*(_maskW(i+1,j,k,bi,bj)*_rdxC(i+1,j,bi,bj)*       &   +0.5*(_maskW(i+1,j,k,bi,bj)
239         &         *_recip_dxC(i+1,j,bi,bj)*
240       &           (theta(i+1,j,k,bi,bj)-theta(i,j,k,bi,bj))       &           (theta(i+1,j,k,bi,bj)-theta(i,j,k,bi,bj))
241       &        +_maskW(i,j,k,bi,bj)*_rdxC(i,j,bi,bj)*       &        +_maskW(i,j,k,bi,bj)
242         &         *_recip_dxC(i,j,bi,bj)*
243       &           (theta(i,j,k,bi,bj)-theta(i-1,j,k,bi,bj)))       &           (theta(i,j,k,bi,bj)-theta(i-1,j,k,bi,bj)))
244       &   +0.5*(_maskW(i+1,j,km1,bi,bj)*_rdxC(i+1,j,bi,bj)*       &   +0.5*(_maskW(i+1,j,km1,bi,bj)
245         &         *_recip_dxC(i+1,j,bi,bj)*
246       &           (theta(i+1,j,km1,bi,bj)-theta(i,j,km1,bi,bj))       &           (theta(i+1,j,km1,bi,bj)-theta(i,j,km1,bi,bj))
247       &        +_maskW(i,j,km1,bi,bj)*_rdxC(i,j,bi,bj)*       &        +_maskW(i,j,km1,bi,bj)
248         &         *_recip_dxC(i,j,bi,bj)*
249       &           (theta(i,j,km1,bi,bj)-theta(i-1,j,km1,bi,bj)))       &           (theta(i,j,km1,bi,bj)-theta(i-1,j,km1,bi,bj)))
250       &       )       &       )
251         ENDDO         ENDDO
# Line 175  C--   Interpolate terms for Redi/GM sche Line 253  C--   Interpolate terms for Redi/GM sche
253        DO j=jMin,jMax        DO j=jMin,jMax
254         DO i=iMin,iMax         DO i=iMin,iMax
255          dTdy(i,j) = 0.5*(          dTdy(i,j) = 0.5*(
256       &   +0.5*(_maskS(i,j,k,bi,bj)*_rdyC(i,j,bi,bj)*       &   +0.5*(_maskS(i,j,k,bi,bj)
257         &         *_recip_dyC(i,j,bi,bj)*
258       &           (theta(i,j,k,bi,bj)-theta(i,j-1,k,bi,bj))       &           (theta(i,j,k,bi,bj)-theta(i,j-1,k,bi,bj))
259       &        +_maskS(i,j+1,k,bi,bj)*_rdyC(i,j+1,bi,bj)*       &        +_maskS(i,j+1,k,bi,bj)
260         &         *_recip_dyC(i,j+1,bi,bj)*
261       &           (theta(i,j+1,k,bi,bj)-theta(i,j,k,bi,bj)))       &           (theta(i,j+1,k,bi,bj)-theta(i,j,k,bi,bj)))
262       &   +0.5*(_maskS(i,j,km1,bi,bj)*_rdyC(i,j,bi,bj)*       &   +0.5*(_maskS(i,j,km1,bi,bj)
263         &         *_recip_dyC(i,j,bi,bj)*
264       &           (theta(i,j,km1,bi,bj)-theta(i,j-1,km1,bi,bj))       &           (theta(i,j,km1,bi,bj)-theta(i,j-1,km1,bi,bj))
265       &        +_maskS(i,j+1,km1,bi,bj)*_rdyC(i,j+1,bi,bj)*       &        +_maskS(i,j+1,km1,bi,bj)
266         &         *_recip_dyC(i,j+1,bi,bj)*
267       &           (theta(i,j+1,km1,bi,bj)-theta(i,j,km1,bi,bj)))       &           (theta(i,j+1,km1,bi,bj)-theta(i,j,km1,bi,bj)))
268       &       )       &       )
269         ENDDO         ENDDO
270        ENDDO        ENDDO
271    #endif /* INCLUDE_T_DIFFUSION_CODE */
272    
273  C--   Vertical flux (fVerT) above  C--   Vertical flux ( fVerT(,,kUp) is at upper face of "theta" cell )
274  C     Advective component of vertical flux  #ifdef INCLUDE_T_ADVECTION_CODE
275    C     o Advective component of vertical flux
276  C     Note: For K=1 then KM1=1 this gives a barZ(T) = T  C     Note: For K=1 then KM1=1 this gives a barZ(T) = T
277  C     (this plays the role of the free-surface correction)  C     (this plays the role of the free-surface correction)
278        DO j=jMin,jMax        DO j=jMin,jMax
279         DO i=iMin,iMax         DO i=iMin,iMax
280          af(i,j) =          af(i,j) =
281       &   wTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i,j,kM1,bi,bj))*0.5 _d 0       &   rTrans(i,j)*(theta(i,j,k,bi,bj)+theta(i,j,kM1,bi,bj))*0.5 _d 0
282         ENDDO         ENDDO
283        ENDDO        ENDDO
284  C     Diffusive component of vertical flux  #endif /* INCLUDE_T_ADVECTION_CODE */
285  C     Note: For K=1 then KM1=1 this gives a dT/dz = 0 upper  #ifdef INCLUDE_T_DIFFUSION_CODE
286    C     o Diffusive component of vertical flux
287    C     Note: For K=1 then KM1=1 and this gives a dT/dr = 0 upper
288  C           boundary condition.  C           boundary condition.
289        DO j=jMin,jMax        IF (implicitDiffusion) THEN
290         DO i=iMin,iMax         DO j=jMin,jMax
291          df(i,j) = _zA(i,j,bi,bj)*(          DO i=iMin,iMax
292       &   -KapGM(i,j)*K13(i,j,k)*dTdx(i,j)           df(i,j) = 0.
293       &   -KapGM(i,j)*K23(i,j,k)*dTdy(i,j)          ENDDO
      &   )  
294         ENDDO         ENDDO
295        ENDDO        ELSE
       IF (.NOT.implicitDiffusion) THEN  
296         DO j=jMin,jMax         DO j=jMin,jMax
297          DO i=iMin,iMax          DO i=iMin,iMax
298           df(i,j) = df(i,j) + _zA(i,j,bi,bj)*(           df(i,j) = - _rA(i,j,bi,bj)*(
299       &    -KappaZT(i,j,k)*rdzC(k)       &    KappaRT(i,j,k)*recip_drC(k)
300       &    *(theta(i,j,kM1,bi,bj)-theta(i,j,k,bi,bj))       &    *(theta(i,j,kM1,bi,bj)-theta(i,j,k,bi,bj))*rkFac
301       &    )       &    )
302          ENDDO          ENDDO
303         ENDDO         ENDDO
304        ENDIF        ENDIF
305  C     Net vertical flux  #endif /* INCLUDE_T_DIFFUSION_CODE */
306    
307    #ifdef ALLOW_GMREDI
308          IF (useGMRedi) CALL GMREDI_RTRANSPORT(
309         I     iMin,iMax,jMin,jMax,bi,bj,K,
310         I     maskUp,theta,
311         U     df,
312         I     myThid)
313    #endif
314    
315    #ifdef ALLOW_KPP
316    C--   Add non local KPP transport term (ghat) to diffusive T flux.
317          IF (useKPP) CALL KPP_TRANSPORT_T(
318         I     iMin,iMax,jMin,jMax,bi,bj,k,km1,
319         I     maskC,KappaRT,
320         U     df )
321    #endif
322    
323    C     o Net vertical flux
324        DO j=jMin,jMax        DO j=jMin,jMax
325         DO i=iMin,iMax         DO i=iMin,iMax
326          fVerT(i,j,kUp) = ( afFacT*af(i,j)+  dfFacT*df(i,j) )*maskUp(i,j)          fVerT(i,j,kUp) = 0.
327         & _ADT( +afFacT*af(i,j)*maskUp(i,j) )
328         & _LPT( +dfFacT*df(i,j)*maskUp(i,j) )
329         ENDDO         ENDDO
330        ENDDO        ENDDO
331    #ifdef INCLUDE_T_ADVECTION_CODE
332        IF ( TOP_LAYER ) THEN        IF ( TOP_LAYER ) THEN
333         DO j=jMin,jMax         DO j=jMin,jMax
334          DO i=iMin,iMax          DO i=iMin,iMax
# Line 231  C     Net vertical flux Line 336  C     Net vertical flux
336          ENDDO          ENDDO
337         ENDDO         ENDDO
338        ENDIF        ENDIF
339    #endif /* INCLUDE_T_ADVECTION_CODE */
340    
341  C--   Tendency is minus divergence of the fluxes.  C--   Tendency is minus divergence of the fluxes.
342  C     Note. Tendency terms will only be correct for range  C     Note. Tendency terms will only be correct for range
# Line 240  C           they are not algorithmically Line 346  C           they are not algorithmically
346  C           are not used.  C           are not used.
347        DO j=jMin,jMax        DO j=jMin,jMax
348         DO i=iMin,iMax         DO i=iMin,iMax
349  C    &   -_rhFacC(i,j,k,bi,bj)*rdzF(k)*_rdxF(i,j,bi,bj)*_rdyF(i,j,bi,bj)  #define _recip_VolT1(i,j,k,bi,bj) _recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
350  C    &   -_rhFacC(i,j,k,bi,bj)*rdzF(k)/_zA(i,j,bi,bj)  #define _recip_VolT2(i,j,k,bi,bj) /_rA(i,j,bi,bj)
 C #define _rVolT(i,j,k,bi,bj) _rhFacC(i,j,k,bi,bj)*rdzF(k)*_rdxF(i,j,bi,bj)*_rdyF(i,j,bi,bj)  
 #define _rVolT(i,j,k,bi,bj) _rhFacC(i,j,k,bi,bj)*rdzF(k)/_zA(i,j,bi,bj)  
351          gT(i,j,k,bi,bj)=          gT(i,j,k,bi,bj)=
352       &   -_rVolT(i,j,k,bi,bj)       &   -_recip_VolT1(i,j,k,bi,bj)
353         &    _recip_VolT2(i,j,k,bi,bj)
354       &   *(       &   *(
355       &    +( fZon(i+1,j)-fZon(i,j) )       &    +( fZon(i+1,j)-fZon(i,j) )
356       &    +( fMer(i,j+1)-fMer(i,j) )       &    +( fMer(i,j+1)-fMer(i,j) )
357       &    +( fVerT(i,j,kUp)-fVerT(i,j,kDown) )       &    +( fVerT(i,j,kUp)-fVerT(i,j,kDown) )*rkFac
358       &    )       &    )
359         ENDDO         ENDDO
360        ENDDO        ENDDO
361    
362    #ifdef INCLUDE_T_FORCING_CODE
363  C--   External thermal forcing term(s)  C--   External thermal forcing term(s)
364  C     o Surface relaxation term        CALL EXTERNAL_FORCING_T(
365        IF ( TOP_LAYER ) THEN       I     iMin,iMax,jMin,jMax,bi,bj,k,
366         DO j=jMin,jMax       I     maskC,
367          DO i=iMin,iMax       I     myCurrentTime,myThid)
368           gT(i,j,k,bi,bj)=gT(i,j,k,bi,bj)  #endif /*  INCLUDE_T_FORCING_CODE */
369       &  +maskC(i,j)*(  
370       &   -lambdaThetaClimRelax*(theta(i,j,k,bi,bj)-SST(i,j,bi,bj))  #ifdef INCLUDE_LAT_CIRC_FFT_FILTER_CODE
371       &   -Qnet(i,j,bi,bj) )  C--   Zonal FFT filter of tendency
372          ENDDO        CALL FILTER_LATCIRCS_FFT_APPLY(
373         ENDDO       U     gT,
374        ENDIF       I     1, sNy, k, k, bi, bj, 1, myThid)
375    #endif /* INCLUDE_LAT_CIRC_FFT_FILTER_CODE */
376    
377        RETURN        RETURN
378        END        END

Legend:
Removed from v.1.13  
changed lines
  Added in v.1.27

  ViewVC Help
Powered by ViewVC 1.1.22