/[MITgcm]/MITgcm/pkg/generic_advdiff/gad_calc_rhs.F
ViewVC logotype

Diff of /MITgcm/pkg/generic_advdiff/gad_calc_rhs.F

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

revision 1.17 by jmc, Sun Jan 26 21:08:36 2003 UTC revision 1.69 by jmc, Thu Aug 14 16:46:36 2014 UTC
# Line 7  CBOP Line 7  CBOP
7  C !ROUTINE: GAD_CALC_RHS  C !ROUTINE: GAD_CALC_RHS
8    
9  C !INTERFACE: ==========================================================  C !INTERFACE: ==========================================================
10        SUBROUTINE GAD_CALC_RHS(        SUBROUTINE GAD_CALC_RHS(
11       I           bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,       I           bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
12       I           xA,yA,uTrans,vTrans,rTrans,maskUp,       I           xA, yA, maskUp, uFld, vFld, wFld,
13       I           diffKh, diffK4, KappaRT, Tracer,       I           uTrans, vTrans, rTrans, rTransKp1,
14       I           tracerIdentity, advectionScheme, calcAdvection,       I           diffKh, diffK4, KappaR, diffKr4, TracerN, TracAB,
15         I           deltaTLev, trIdentity,
16         I           advectionScheme, vertAdvecScheme,
17         I           calcAdvection, implicitAdvection, applyAB_onTracer,
18         I           trUseDiffKr4, trUseGMRedi, trUseKPP,
19         O           fZon, fMer,
20       U           fVerT, gTracer,       U           fVerT, gTracer,
21       I           myThid )       I           myTime, myIter, myThid )
22    
23  C !DESCRIPTION:  C !DESCRIPTION:
24  C Calculates the tendancy of a tracer due to advection and diffusion.  C Calculates the tendency of a tracer due to advection and diffusion.
25  C It calculates the fluxes in each direction indepentently and then  C It calculates the fluxes in each direction indepentently and then
26  C sets the tendancy to the divergence of these fluxes. The advective  C sets the tendency to the divergence of these fluxes. The advective
27  C fluxes are only calculated here when using the linear advection schemes  C fluxes are only calculated here when using the linear advection schemes
28  C otherwise only the diffusive and parameterized fluxes are calculated.  C otherwise only the diffusive and parameterized fluxes are calculated.
29  C  C
# Line 27  C \begin{equation*} Line 32  C \begin{equation*}
32  C {\bf F} = {\bf F}_{adv} + {\bf F}_{diff} +{\bf F}_{GM} + {\bf F}_{KPP}  C {\bf F} = {\bf F}_{adv} + {\bf F}_{diff} +{\bf F}_{GM} + {\bf F}_{KPP}
33  C \end{equation*}  C \end{equation*}
34  C  C
35  C The tendancy is the divergence of the fluxes:  C The tendency is the divergence of the fluxes:
36  C \begin{equation*}  C \begin{equation*}
37  C G_\theta = G_\theta + \nabla \cdot {\bf F}  C G_\theta = G_\theta + \nabla \cdot {\bf F}
38  C \end{equation*}  C \end{equation*}
39  C  C
40  C The tendancy is assumed to contain data on entry.  C The tendency is assumed to contain data on entry.
41    
42  C !USES: ===============================================================  C !USES: ===============================================================
43        IMPLICIT NONE        IMPLICIT NONE
# Line 40  C !USES: =============================== Line 45  C !USES: ===============================
45  #include "EEPARAMS.h"  #include "EEPARAMS.h"
46  #include "PARAMS.h"  #include "PARAMS.h"
47  #include "GRID.h"  #include "GRID.h"
 #include "DYNVARS.h"  
48  #include "SURFACE.h"  #include "SURFACE.h"
49  #include "GAD.h"  #include "GAD.h"
50    #ifdef ALLOW_AUTODIFF
51  #ifdef ALLOW_AUTODIFF_TAMC  # include "AUTODIFF_PARAMS.h"
52  #include "tamc.h"  #endif /* ALLOW_AUTODIFF */
 #include "tamc_keys.h"  
 #endif /* ALLOW_AUTODIFF_TAMC */  
53    
54  C !INPUT PARAMETERS: ===================================================  C !INPUT PARAMETERS: ===================================================
55  C  bi,bj                :: tile indices  C bi,bj            :: tile indices
56  C  iMin,iMax,jMin,jMax  :: loop range for called routines  C iMin,iMax        :: loop range for called routines
57  C  kup                  :: index into 2 1/2D array, toggles between 1 and 2  C jMin,jMax        :: loop range for called routines
58  C  kdown                :: index into 2 1/2D array, toggles between 2 and 1  C k                :: vertical index
59  C  kp1                  :: =k+1 for k<Nr, =Nr for k=Nr  C kM1              :: =k-1 for k>1, =1 for k=1
60  C  xA,yA                :: areas of X and Y face of tracer cells  C kUp              :: index into 2 1/2D array, toggles between 1|2
61  C  uTrans,vTrans,rTrans :: 2-D arrays of volume transports at U,V and W points  C kDown            :: index into 2 1/2D array, toggles between 2|1
62  C  maskUp               :: 2-D array for mask at W points  C xA,yA            :: areas of X and Y face of tracer cells
63  C  diffKh               :: horizontal diffusion coefficient  C maskUp           :: 2-D array for mask at W points
64  C  diffK4               :: bi-harmonic diffusion coefficient  C uFld,vFld,wFld   :: Local copy of velocity field (3 components)
65  C  KappaRT              :: 3-D array for vertical diffusion coefficient  C uTrans,vTrans    :: 2-D arrays of volume transports at U,V points
66  C  Tracer               :: tracer field  C rTrans           :: 2-D arrays of volume transports at W points
67  C  tracerIdentity       :: identifier for the tracer (required only for KPP)  C rTransKp1        :: 2-D array of volume trans at W pts, interf k+1
68  C  advectionScheme      :: advection scheme to use  C diffKh           :: horizontal diffusion coefficient
69  C  calcAdvection        :: =False if Advec terms computed with multiDim scheme  C diffK4           :: horizontal bi-harmonic diffusion coefficient
70  C  myThid               :: thread number  C KappaR           :: 2-D array for vertical diffusion coefficient, interf k
71    C diffKr4          :: 1-D array for vertical bi-harmonic diffusion coefficient
72    C TracerN          :: tracer field @ time-step n (Note: only used
73    C                     if applying AB on tracer field rather than on tendency gTr)
74    C TracAB           :: current tracer field (@ time-step n if applying AB on gTr
75    C                     or extrapolated fwd in time to n+1/2 if applying AB on Tr)
76    C trIdentity       :: tracer identifier (required for KPP,GM)
77    C advectionScheme  :: advection scheme to use (Horizontal plane)
78    C vertAdvecScheme  :: advection scheme to use (Vertical direction)
79    C calcAdvection    :: =False if Advec computed with multiDim scheme
80    C implicitAdvection:: =True if vertical Advec computed implicitly
81    C applyAB_onTracer :: apply Adams-Bashforth on Tracer (rather than on gTr)
82    C trUseDiffKr4     :: true if this tracer uses vertical bi-harmonic diffusion
83    C trUseGMRedi      :: true if this tracer uses GM-Redi
84    C trUseKPP         :: true if this tracer uses KPP
85    C myTime           :: current time
86    C myIter           :: iteration number
87    C myThid           :: thread number
88        INTEGER bi,bj,iMin,iMax,jMin,jMax        INTEGER bi,bj,iMin,iMax,jMin,jMax
89        INTEGER k,kUp,kDown,kM1        INTEGER k,kUp,kDown,kM1
90        _RS xA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS xA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
91        _RS yA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RS yA    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
92          _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
93          _RL uFld  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94          _RL vFld  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95          _RL wFld  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96        _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97        _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98        _RL rTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
99        _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL rTransKp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
100        _RL diffKh, diffK4        _RL diffKh, diffK4
101        _RL KappaRT(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL KappaR(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
102        _RL Tracer(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)        _RL diffKr4(Nr)
103        INTEGER tracerIdentity        _RL TracerN(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
104        INTEGER advectionScheme        _RL TracAB (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
105          _RL deltaTLev(Nr)
106          INTEGER trIdentity
107          INTEGER advectionScheme, vertAdvecScheme
108        LOGICAL calcAdvection        LOGICAL calcAdvection
109        INTEGER myThid        LOGICAL implicitAdvection, applyAB_onTracer
110          LOGICAL trUseDiffKr4, trUseGMRedi, trUseKPP
111          _RL     myTime
112          INTEGER myIter, myThid
113    
114  C !OUTPUT PARAMETERS: ==================================================  C !OUTPUT PARAMETERS: ==================================================
115  C  gTracer              :: tendancy array  C gTracer          :: tendency array
116  C  fVerT                :: 2 1/2D arrays for vertical advective flux  C fZon             :: zonal flux
117        _RL gTracer(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)  C fMer             :: meridional flux
118    C fVerT            :: 2 1/2D arrays for vertical advective flux
119          _RL gTracer(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
120          _RL fZon  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
121          _RL fMer  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
122        _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)        _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
123    
124    C !FUNCTIONS:       ====================================================
125    #ifdef ALLOW_DIAGNOSTICS
126          CHARACTER*4 GAD_DIAG_SUFX
127          EXTERNAL    GAD_DIAG_SUFX
128    #endif /* ALLOW_DIAGNOSTICS */
129    
130  C !LOCAL VARIABLES: ====================================================  C !LOCAL VARIABLES: ====================================================
131  C  i,j                  :: loop indices  C i,j              :: loop indices
132  C  df4                  :: used for storing del^2 T for bi-harmonic term  C df4              :: used for storing del^2 T for bi-harmonic term
133  C  fZon                 :: zonal flux  C af               :: advective flux
134  C  fmer                 :: meridional flux  C df               :: diffusive flux
135  C  af                   :: advective flux  C localT           :: local copy of tracer field
136  C  df                   :: diffusive flux  C locABT           :: local copy of (AB-extrapolated) tracer field
 C  localT               :: local copy of tracer field  
137        INTEGER i,j        INTEGER i,j
138          _RS maskLocW(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
139          _RS maskLocS(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
140        _RL df4   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL df4   (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
       _RL fZon  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
       _RL fMer  (1-OLx:sNx+OLx,1-OLy:sNy+OLy)  
141        _RL af    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL af    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
142        _RL df    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL df    (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
143        _RL localT(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL localT(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
144          _RL locABT(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
145          _RL advFac, rAdvFac
146    #ifdef GAD_SMOLARKIEWICZ_HACK
147          _RL outFlux, trac, fac, gTrFac
148    #endif
149    #ifdef ALLOW_DIAGNOSTICS
150          CHARACTER*8 diagName
151          CHARACTER*4 diagSufx
152    #endif
153  CEOP  CEOP
154    
155  #ifdef ALLOW_AUTODIFF_TAMC  #ifdef ALLOW_AUTODIFF
156  C--   only the kUp part of fverT is set in this subroutine  C--   only the kUp part of fverT is set in this subroutine
157  C--   the kDown is still required  C--   the kDown is still required
158        fVerT(1,1,kDown) = fVerT(1,1,kDown)        fVerT(1,1,kDown) = fVerT(1,1,kDown)
159  #endif  #endif
160    
161    #ifdef ALLOW_DIAGNOSTICS
162    C--   Set diagnostic suffix for the current tracer
163          IF ( useDiagnostics ) THEN
164            diagSufx = GAD_DIAG_SUFX( trIdentity, myThid )
165          ENDIF
166    #endif
167    
168          advFac  = 0. _d 0
169          IF (calcAdvection) advFac = 1. _d 0
170          rAdvFac = rkSign*advFac
171          IF (implicitAdvection) rAdvFac = rkSign
172    
173        DO j=1-OLy,sNy+OLy        DO j=1-OLy,sNy+OLy
174         DO i=1-OLx,sNx+OLx         DO i=1-OLx,sNx+OLx
175          fZon(i,j)      = 0. _d 0          fZon(i,j)      = 0. _d 0
# Line 118  C--   the kDown is still required Line 177  C--   the kDown is still required
177          fVerT(i,j,kUp) = 0. _d 0          fVerT(i,j,kUp) = 0. _d 0
178          df(i,j)        = 0. _d 0          df(i,j)        = 0. _d 0
179          df4(i,j)       = 0. _d 0          df4(i,j)       = 0. _d 0
         localT(i,j)    = 0. _d 0  
180         ENDDO         ENDDO
181        ENDDO        ENDDO
182    
183  C--   Make local copy of tracer array  C--   Make local copy of tracer array
184        DO j=1-OLy,sNy+OLy        IF ( applyAB_onTracer ) THEN
185         DO i=1-OLx,sNx+OLx          DO j=1-OLy,sNy+OLy
186          localT(i,j)=tracer(i,j,k,bi,bj)           DO i=1-OLx,sNx+OLx
187         ENDDO            localT(i,j)=TracerN(i,j,k,bi,bj)
188        ENDDO            locABT(i,j)= TracAB(i,j,k,bi,bj)
189             ENDDO
190  C--   Unless we have already calculated the advection terms we initialize          ENDDO
191  C     the tendency to zero.        ELSE
192        IF (calcAdvection) THEN          DO j=1-OLy,sNy+OLy
193         DO j=1-Oly,sNy+Oly           DO i=1-OLx,sNx+OLx
194          DO i=1-Olx,sNx+Olx            localT(i,j)=TracerN(i,j,k,bi,bj)
195           gTracer(i,j,k,bi,bj)=0. _d 0            locABT(i,j)=TracerN(i,j,k,bi,bj)
196             ENDDO
197          ENDDO          ENDDO
        ENDDO  
198        ENDIF        ENDIF
199    
200  C--   Pre-calculate del^2 T if bi-harmonic coefficient is non-zero  C--   Pre-calculate del^2 T if bi-harmonic coefficient is non-zero
# Line 147  C--   Pre-calculate del^2 T if bi-harmon Line 205  C--   Pre-calculate del^2 T if bi-harmon
205        ENDIF        ENDIF
206    
207  C--   Initialize net flux in X direction  C--   Initialize net flux in X direction
208        DO j=1-Oly,sNy+Oly        DO j=1-OLy,sNy+OLy
209         DO i=1-Olx,sNx+Olx         DO i=1-OLx,sNx+OLx
210          fZon(i,j) = 0. _d 0          fZon(i,j) = 0. _d 0
211         ENDDO         ENDDO
212        ENDDO        ENDDO
213    
214  C-    Advective flux in X  C-    Advective flux in X
215        IF (calcAdvection) THEN        IF (calcAdvection) THEN
216        IF (advectionScheme.EQ.ENUM_CENTERED_2ND) THEN          IF (advectionScheme.EQ.ENUM_CENTERED_2ND) THEN
217         CALL GAD_C2_ADV_X(bi,bj,k,uTrans,localT,af,myThid)             CALL GAD_C2_ADV_X( bi,bj,k, uTrans, locABT, af, myThid )
218        ELSEIF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN          ELSEIF ( advectionScheme.EQ.ENUM_UPWIND_1RST
219         CALL GAD_FLUXLIMIT_ADV_X(       &          .OR. advectionScheme.EQ.ENUM_DST2 ) THEN
220       &      bi,bj,k,deltaTtracer,uTrans,uVel,localT,af,myThid)             CALL GAD_DST2U1_ADV_X( bi,bj,k, advectionScheme, .TRUE.,
221        ELSEIF (advectionScheme.EQ.ENUM_UPWIND_3RD ) THEN       I              deltaTLev(k), uTrans, uFld, locABT,
222         CALL GAD_U3_ADV_X(bi,bj,k,uTrans,localT,af,myThid)       O              af, myThid )
223        ELSEIF (advectionScheme.EQ.ENUM_CENTERED_4TH) THEN          ELSE
224         CALL GAD_C4_ADV_X(bi,bj,k,uTrans,localT,af,myThid)           DO j=1-OLy,sNy+OLy
225        ELSEIF (advectionScheme.EQ.ENUM_DST3 ) THEN            DO i=1-OLx,sNx+OLx
226         CALL GAD_DST3_ADV_X(  #ifdef ALLOW_OBCS
227       &       bi,bj,k,deltaTtracer,uTrans,uVel,localT,af,myThid)             maskLocW(i,j) = _maskW(i,j,k,bi,bj)*maskInW(i,j,bi,bj)
228        ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN  #else /* ALLOW_OBCS */
229         CALL GAD_DST3FL_ADV_X(             maskLocW(i,j) = _maskW(i,j,k,bi,bj)
230       &       bi,bj,k,deltaTtracer,uTrans,uVel,localT,af,myThid)  #endif /* ALLOW_OBCS */
231        ELSE            ENDDO
232         STOP 'GAD_CALC_RHS: Bad advectionScheme (X)'           ENDDO
233        ENDIF           IF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN
234        DO j=1-Oly,sNy+Oly             CALL GAD_FLUXLIMIT_ADV_X( bi,bj,k, .TRUE., deltaTLev(k),
235         DO i=1-Olx,sNx+Olx       I              uTrans, uFld, maskLocW, locABT,
236          fZon(i,j) = fZon(i,j) + af(i,j)       O              af, myThid )
237         ENDDO           ELSEIF (advectionScheme.EQ.ENUM_UPWIND_3RD ) THEN
238        ENDDO             CALL GAD_U3_ADV_X( bi,bj,k, uTrans, maskLocW, locABT,
239         O              af, myThid )
240             ELSEIF (advectionScheme.EQ.ENUM_CENTERED_4TH) THEN
241               CALL GAD_C4_ADV_X( bi,bj,k, uTrans, maskLocW, locABT,
242         O              af, myThid )
243    #ifdef ALLOW_AUTODIFF
244             ELSEIF( advectionScheme.EQ.ENUM_DST3 .OR.
245         &          (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT .AND. inAdMode)
246         &         ) THEN
247    cph This block is to trick the adjoint:
248    cph If inAdExact=.FALSE., we want to use DST3
249    cph with limiters in forward, but without limiters in reverse.
250    #else /* ALLOW_AUTODIFF */
251             ELSEIF( advectionScheme.EQ.ENUM_DST3 ) THEN
252    #endif /* ALLOW_AUTODIFF */
253               CALL GAD_DST3_ADV_X( bi,bj,k, .TRUE., deltaTLev(k),
254         I              uTrans, uFld, maskLocW, locABT,
255         O              af, myThid )
256             ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
257               CALL GAD_DST3FL_ADV_X( bi,bj,k, .TRUE., deltaTLev(k),
258         I              uTrans, uFld, maskLocW, locABT,
259         O              af, myThid )
260    #ifndef ALLOW_AUTODIFF
261             ELSEIF (advectionScheme.EQ.ENUM_OS7MP ) THEN
262               CALL GAD_OS7MP_ADV_X( bi,bj,k, .TRUE., deltaTLev(k),
263         I              uTrans, uFld, maskLocW, locABT,
264         O              af, myThid )
265    #endif
266             ELSE
267              STOP 'GAD_CALC_RHS: Bad advectionScheme (X)'
268             ENDIF
269            ENDIF
270    #ifdef ALLOW_OBCS
271            IF ( useOBCS ) THEN
272    C-      replace advective flux with 1st order upwind scheme estimate
273              CALL OBCS_U1_ADV_TRACER( .TRUE., trIdentity, bi, bj, k,
274         I                             maskW(1-OLx,1-OLy,k,bi,bj),
275         I                             uTrans, locABT,
276         U                             af, myThid )
277            ENDIF
278    #endif /* ALLOW_OBCS */
279            DO j=1-OLy,sNy+OLy
280             DO i=1-OLx,sNx+OLx
281              fZon(i,j) = fZon(i,j) + af(i,j)
282             ENDDO
283            ENDDO
284    #ifdef ALLOW_DIAGNOSTICS
285            IF ( useDiagnostics ) THEN
286              diagName = 'ADVx'//diagSufx
287              CALL DIAGNOSTICS_FILL( af, diagName, k,1, 2,bi,bj, myThid )
288            ENDIF
289    #endif
290        ENDIF        ENDIF
291    
292  C-    Diffusive flux in X  C-    Diffusive flux in X
293        IF (diffKh.NE.0.) THEN        IF (diffKh.NE.0.) THEN
294         CALL GAD_DIFF_X(bi,bj,k,xA,diffKh,localT,df,myThid)         CALL GAD_DIFF_X(bi,bj,k,xA,diffKh,localT,df,myThid)
295        ELSE        ELSE
296         DO j=1-Oly,sNy+Oly         DO j=1-OLy,sNy+OLy
297          DO i=1-Olx,sNx+Olx          DO i=1-OLx,sNx+OLx
298           df(i,j) = 0. _d 0           df(i,j) = 0. _d 0
299          ENDDO          ENDDO
300         ENDDO         ENDDO
301        ENDIF        ENDIF
302    
303    C-    Add bi-harmonic diffusive flux in X
304          IF (diffK4 .NE. 0.) THEN
305           CALL GAD_BIHARM_X(bi,bj,k,xA,df4,diffK4,df,myThid)
306          ENDIF
307    
308  #ifdef ALLOW_GMREDI  #ifdef ALLOW_GMREDI
309  C-    GM/Redi flux in X  C-    GM/Redi flux in X
310        IF (useGMRedi) THEN        IF ( trUseGMRedi ) THEN
311  C *note* should update GMREDI_XTRANSPORT to use localT and set df  *aja*  C *note* should update GMREDI_XTRANSPORT to set df  *aja*
312          CALL GMREDI_XTRANSPORT(          CALL GMREDI_XTRANSPORT(
313       I     iMin,iMax,jMin,jMax,bi,bj,K,       I         iMin,iMax,jMin,jMax,bi,bj,k,
314       I     xA,Tracer,tracerIdentity,       I         xA,TracerN,trIdentity,
315       U     df,       U         df,
316       I     myThid)       I         myThid)
317        ENDIF        ENDIF
318  #endif  #endif
319        DO j=1-Oly,sNy+Oly  C     anelastic: advect.fluxes are scaled by rhoFac but hor.diff. flx are not
320         DO i=1-Olx,sNx+Olx        DO j=1-OLy,sNy+OLy
321          fZon(i,j) = fZon(i,j) + df(i,j)         DO i=1-OLx,sNx+OLx
322            fZon(i,j) = fZon(i,j) + df(i,j)*rhoFacC(k)
323         ENDDO         ENDDO
324        ENDDO        ENDDO
325    
326  C-    Bi-harmonic duffusive flux in X  #ifdef ALLOW_DIAGNOSTICS
327        IF (diffK4 .NE. 0.) THEN  C-    Diagnostics of Tracer flux in X dir (mainly Diffusive term),
328         CALL GAD_BIHARM_X(bi,bj,k,xA,df4,diffK4,df,myThid)  C       excluding advective terms:
329         DO j=1-Oly,sNy+Oly        IF ( useDiagnostics .AND.
330          DO i=1-Olx,sNx+Olx       &    (diffKh.NE.0. .OR. diffK4 .NE.0. .OR. trUseGMRedi) ) THEN
331           fZon(i,j) = fZon(i,j) + df(i,j)            diagName = 'DFxE'//diagSufx
332          ENDDO            CALL DIAGNOSTICS_FILL( df, diagName, k,1, 2,bi,bj, myThid )
333         ENDDO  #ifdef ALLOW_LAYERS
334              IF ( useLayers ) THEN
335               CALL LAYERS_FILL_DFX( df, trIdentity, k, 1, 2,bi,bj, myThid )
336              ENDIF
337    #endif /* ALLOW_LAYERS */
338        ENDIF        ENDIF
339    #endif
340    
341  C--   Initialize net flux in Y direction  C--   Initialize net flux in Y direction
342        DO j=1-Oly,sNy+Oly        DO j=1-OLy,sNy+OLy
343         DO i=1-Olx,sNx+Olx         DO i=1-OLx,sNx+OLx
344          fMer(i,j) = 0. _d 0          fMer(i,j) = 0. _d 0
345         ENDDO         ENDDO
346        ENDDO        ENDDO
347    
348  C-    Advective flux in Y  C-    Advective flux in Y
349        IF (calcAdvection) THEN        IF (calcAdvection) THEN
350        IF (advectionScheme.EQ.ENUM_CENTERED_2ND) THEN          IF (advectionScheme.EQ.ENUM_CENTERED_2ND) THEN
351         CALL GAD_C2_ADV_Y(bi,bj,k,vTrans,localT,af,myThid)             CALL GAD_C2_ADV_Y( bi,bj,k, vTrans, locABT, af, myThid )
352        ELSEIF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN          ELSEIF ( advectionScheme.EQ.ENUM_UPWIND_1RST
353         CALL GAD_FLUXLIMIT_ADV_Y(       &          .OR. advectionScheme.EQ.ENUM_DST2 ) THEN
354       &       bi,bj,k,deltaTtracer,vTrans,vVel,localT,af,myThid)             CALL GAD_DST2U1_ADV_Y( bi,bj,k, advectionScheme, .TRUE.,
355        ELSEIF (advectionScheme.EQ.ENUM_UPWIND_3RD ) THEN       I              deltaTLev(k), vTrans, vFld, locABT,
356         CALL GAD_U3_ADV_Y(bi,bj,k,vTrans,localT,af,myThid)       O              af, myThid )
357        ELSEIF (advectionScheme.EQ.ENUM_CENTERED_4TH) THEN          ELSE
358         CALL GAD_C4_ADV_Y(bi,bj,k,vTrans,localT,af,myThid)           DO j=1-OLy,sNy+OLy
359        ELSEIF (advectionScheme.EQ.ENUM_DST3 ) THEN            DO i=1-OLx,sNx+OLx
360         CALL GAD_DST3_ADV_Y(  #ifdef ALLOW_OBCS
361       &       bi,bj,k,deltaTtracer,vTrans,vVel,localT,af,myThid)             maskLocS(i,j) = _maskS(i,j,k,bi,bj)*maskInS(i,j,bi,bj)
362        ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN  #else /* ALLOW_OBCS */
363         CALL GAD_DST3FL_ADV_Y(             maskLocS(i,j) = _maskS(i,j,k,bi,bj)
364       &       bi,bj,k,deltaTtracer,vTrans,vVel,localT,af,myThid)  #endif /* ALLOW_OBCS */
365        ELSE            ENDDO
366         STOP 'GAD_CALC_RHS: Bad advectionScheme (Y)'           ENDDO
367        ENDIF           IF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN
368        DO j=1-Oly,sNy+Oly             CALL GAD_FLUXLIMIT_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k),
369         DO i=1-Olx,sNx+Olx       I              vTrans, vFld, maskLocS, locABT,
370          fMer(i,j) = fMer(i,j) + af(i,j)       O              af, myThid )
371         ENDDO           ELSEIF (advectionScheme.EQ.ENUM_UPWIND_3RD ) THEN
372        ENDDO             CALL GAD_U3_ADV_Y( bi,bj,k, vTrans, maskLocS, locABT,
373         O              af, myThid )
374             ELSEIF (advectionScheme.EQ.ENUM_CENTERED_4TH) THEN
375               CALL GAD_C4_ADV_Y( bi,bj,k, vTrans, maskLocS, locABT,
376         O              af, myThid )
377    #ifdef ALLOW_AUTODIFF
378             ELSEIF( advectionScheme.EQ.ENUM_DST3 .OR.
379         &          (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT .AND. inAdMode)
380         &         ) THEN
381    cph This block is to trick the adjoint:
382    cph If inAdExact=.FALSE., we want to use DST3
383    cph with limiters in forward, but without limiters in reverse.
384    #else /* ALLOW_AUTODIFF */
385             ELSEIF( advectionScheme.EQ.ENUM_DST3 ) THEN
386    #endif /* ALLOW_AUTODIFF */
387               CALL GAD_DST3_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k),
388         I              vTrans, vFld, maskLocS, locABT,
389         O              af, myThid )
390             ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
391               CALL GAD_DST3FL_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k),
392         I              vTrans, vFld, maskLocS, locABT,
393         O              af, myThid )
394    #ifndef ALLOW_AUTODIFF
395             ELSEIF (advectionScheme.EQ.ENUM_OS7MP ) THEN
396               CALL GAD_OS7MP_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k),
397         I              vTrans, vFld, maskLocS, locABT,
398         O              af, myThid )
399    #endif
400             ELSE
401               STOP 'GAD_CALC_RHS: Bad advectionScheme (Y)'
402             ENDIF
403            ENDIF
404    #ifdef ALLOW_OBCS
405            IF ( useOBCS ) THEN
406    C-      replace advective flux with 1st order upwind scheme estimate
407              CALL OBCS_U1_ADV_TRACER( .FALSE., trIdentity, bi, bj, k,
408         I                             maskS(1-OLx,1-OLy,k,bi,bj),
409         I                             vTrans, locABT,
410         U                             af, myThid )
411            ENDIF
412    #endif /* ALLOW_OBCS */
413            DO j=1-OLy,sNy+OLy
414             DO i=1-OLx,sNx+OLx
415              fMer(i,j) = fMer(i,j) + af(i,j)
416             ENDDO
417            ENDDO
418    #ifdef ALLOW_DIAGNOSTICS
419            IF ( useDiagnostics ) THEN
420              diagName = 'ADVy'//diagSufx
421              CALL DIAGNOSTICS_FILL( af, diagName, k,1, 2,bi,bj, myThid )
422            ENDIF
423    #endif
424        ENDIF        ENDIF
425    
426  C-    Diffusive flux in Y  C-    Diffusive flux in Y
427        IF (diffKh.NE.0.) THEN        IF (diffKh.NE.0.) THEN
428         CALL GAD_DIFF_Y(bi,bj,k,yA,diffKh,localT,df,myThid)         CALL GAD_DIFF_Y(bi,bj,k,yA,diffKh,localT,df,myThid)
429        ELSE        ELSE
430         DO j=1-Oly,sNy+Oly         DO j=1-OLy,sNy+OLy
431          DO i=1-Olx,sNx+Olx          DO i=1-OLx,sNx+OLx
432           df(i,j) = 0. _d 0           df(i,j) = 0. _d 0
433          ENDDO          ENDDO
434         ENDDO         ENDDO
435        ENDIF        ENDIF
436    
437    C-    Add bi-harmonic flux in Y
438          IF (diffK4 .NE. 0.) THEN
439           CALL GAD_BIHARM_Y(bi,bj,k,yA,df4,diffK4,df,myThid)
440          ENDIF
441    
442  #ifdef ALLOW_GMREDI  #ifdef ALLOW_GMREDI
443  C-    GM/Redi flux in Y  C-    GM/Redi flux in Y
444        IF (useGMRedi) THEN        IF ( trUseGMRedi ) THEN
445  C *note* should update GMREDI_YTRANSPORT to use localT and set df  *aja*  C *note* should update GMREDI_YTRANSPORT to set df  *aja*
446         CALL GMREDI_YTRANSPORT(          CALL GMREDI_YTRANSPORT(
447       I     iMin,iMax,jMin,jMax,bi,bj,K,       I         iMin,iMax,jMin,jMax,bi,bj,k,
448       I     yA,Tracer,tracerIdentity,       I         yA,TracerN,trIdentity,
449       U     df,       U         df,
450       I     myThid)       I         myThid)
451        ENDIF        ENDIF
452  #endif  #endif
453        DO j=1-Oly,sNy+Oly  C     anelastic: advect.fluxes are scaled by rhoFac but hor.diff. flx are not
454         DO i=1-Olx,sNx+Olx        DO j=1-OLy,sNy+OLy
455          fMer(i,j) = fMer(i,j) + df(i,j)         DO i=1-OLx,sNx+OLx
456            fMer(i,j) = fMer(i,j) + df(i,j)*rhoFacC(k)
457         ENDDO         ENDDO
458        ENDDO        ENDDO
459    
460  C-    Bi-harmonic flux in Y  #ifdef ALLOW_DIAGNOSTICS
461        IF (diffK4 .NE. 0.) THEN  C-    Diagnostics of Tracer flux in Y dir (mainly Diffusive terms),
462         CALL GAD_BIHARM_Y(bi,bj,k,yA,df4,diffK4,df,myThid)  C       excluding advective terms:
463         DO j=1-Oly,sNy+Oly        IF ( useDiagnostics .AND.
464          DO i=1-Olx,sNx+Olx       &    (diffKh.NE.0. .OR. diffK4 .NE.0. .OR. trUseGMRedi) ) THEN
465           fMer(i,j) = fMer(i,j) + df(i,j)            diagName = 'DFyE'//diagSufx
466          ENDDO            CALL DIAGNOSTICS_FILL( df, diagName, k,1, 2,bi,bj, myThid )
467         ENDDO  #ifdef ALLOW_LAYERS
468        ENDIF            IF ( useLayers ) THEN
469               CALL LAYERS_FILL_DFY( df, trIdentity,k, 1, 2,bi,bj, myThid )
470  #ifdef NONLIN_FRSURF            ENDIF
471  C--   Compute vertical flux fVerT(kDown) at interface k+1 (between k & k+1):  #endif /* ALLOW_LAYERS */
       IF ( calcAdvection .AND. K.EQ.Nr .AND.  
      &     useRealFreshWaterFlux .AND.  
      &     buoyancyRelation .EQ. 'OCEANICP' ) THEN    
        DO j=1-Oly,sNy+Oly  
         DO i=1-Olx,sNx+Olx  
          fVerT(i,j,kDown) = convertEmP2rUnit*PmEpR(i,j,bi,bj)  
      &     *rA(i,j,bi,bj)*maskC(i,j,k,bi,bj)*Tracer(i,j,k,bi,bj)  
         ENDDO  
        ENDDO  
472        ENDIF        ENDIF
473  #endif /* NONLIN_FRSURF */  #endif
474    
475  C--   Compute vertical flux fVerT(kUp) at interface k (between k-1 & k):  C--   Compute vertical flux fVerT(kUp) at interface k (between k-1 & k):
476  C-    Advective flux in R  C-    Advective flux in R
477        IF (calcAdvection) THEN  #ifdef ALLOW_AIM
478  C     Note: wVel needs to be masked  C- a hack to prevent Water-Vapor vert.transport into the stratospheric level Nr
479        IF (K.GE.2) THEN        IF (calcAdvection .AND. .NOT.implicitAdvection .AND. k.GE.2 .AND.
480  C-    Compute vertical advective flux in the interior:       &     (.NOT.useAIM .OR. trIdentity.NE.GAD_SALINITY .OR. k.LT.Nr)
481         IF (advectionScheme.EQ.ENUM_CENTERED_2ND) THEN       &   ) THEN
482          CALL GAD_C2_ADV_R(bi,bj,k,rTrans,tracer,af,myThid)  #else
483         ELSEIF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN        IF (calcAdvection .AND. .NOT.implicitAdvection .AND. k.GE.2) THEN
484          CALL GAD_FLUXLIMIT_ADV_R(  #endif
485       &       bi,bj,k,deltaTtracer,rTrans,wVel,tracer,af,myThid)         IF ( applyAB_onTracer ) THEN
486         ELSEIF (advectionScheme.EQ.ENUM_UPWIND_3RD ) THEN  C-    Compute vertical advective flux in the interior using TracAB:
487          CALL GAD_U3_ADV_R(bi,bj,k,rTrans,tracer,af,myThid)          IF (vertAdvecScheme.EQ.ENUM_CENTERED_2ND) THEN
488         ELSEIF (advectionScheme.EQ.ENUM_CENTERED_4TH) THEN             CALL GAD_C2_ADV_R( bi,bj,k, rTrans, TracAB, af, myThid )
489          CALL GAD_C4_ADV_R(bi,bj,k,rTrans,tracer,af,myThid)          ELSEIF ( vertAdvecScheme.EQ.ENUM_UPWIND_1RST
490         ELSEIF (advectionScheme.EQ.ENUM_DST3 ) THEN       &          .OR. vertAdvecScheme.EQ.ENUM_DST2 ) THEN
491          CALL GAD_DST3_ADV_R(             CALL GAD_DST2U1_ADV_R( bi,bj,k,vertAdvecScheme,deltaTLev(k),
492       &       bi,bj,k,deltaTtracer,rTrans,wVel,tracer,af,myThid)       I              rTrans, wFld, TracAB(1-OLx,1-OLy,1,bi,bj),
493         ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN       O              af, myThid )
494          CALL GAD_DST3FL_ADV_R(          ELSEIF (vertAdvecScheme.EQ.ENUM_FLUX_LIMIT) THEN
495       &       bi,bj,k,deltaTtracer,rTrans,wVel,tracer,af,myThid)             CALL GAD_FLUXLIMIT_ADV_R( bi,bj,k, deltaTLev(k),
496         I              rTrans, wFld, TracAB(1-OLx,1-OLy,1,bi,bj),
497         O              af, myThid )
498            ELSEIF (vertAdvecScheme.EQ.ENUM_UPWIND_3RD ) THEN
499               CALL GAD_U3_ADV_R( bi,bj,k, rTrans, TracAB, af, myThid )
500            ELSEIF (vertAdvecScheme.EQ.ENUM_CENTERED_4TH) THEN
501               CALL GAD_C4_ADV_R( bi,bj,k, rTrans, TracAB, af, myThid )
502    #ifdef ALLOW_AUTODIFF
503            ELSEIF( vertAdvecScheme.EQ.ENUM_DST3 .OR.
504         &         (vertAdvecScheme.EQ.ENUM_DST3_FLUX_LIMIT .AND. inAdMode)
505         &        ) THEN
506    cph This block is to trick the adjoint:
507    cph If inAdExact=.FALSE., we want to use DST3
508    cph with limiters in forward, but without limiters in reverse.
509    #else /* ALLOW_AUTODIFF */
510            ELSEIF( vertAdvecScheme.EQ.ENUM_DST3 ) THEN
511    #endif /* ALLOW_AUTODIFF */
512               CALL GAD_DST3_ADV_R( bi,bj,k, deltaTLev(k),
513         I              rTrans, wFld, TracAB(1-OLx,1-OLy,1,bi,bj),
514         O              af, myThid )
515            ELSEIF (vertAdvecScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
516               CALL GAD_DST3FL_ADV_R( bi,bj,k, deltaTLev(k),
517         I              rTrans, wFld, TracAB(1-OLx,1-OLy,1,bi,bj),
518         O              af, myThid )
519    #ifndef ALLOW_AUTODIFF
520            ELSEIF (vertAdvecScheme.EQ.ENUM_OS7MP ) THEN
521               CALL GAD_OS7MP_ADV_R( bi,bj,k, deltaTLev(k),
522         I              rTrans, wFld, TracAB(1-OLx,1-OLy,1,bi,bj),
523         O              af, myThid )
524    #endif
525            ELSE
526              STOP 'GAD_CALC_RHS: Bad vertAdvecScheme (R)'
527            ENDIF
528         ELSE         ELSE
529          STOP 'GAD_CALC_RHS: Bad advectionScheme (R)'  C-    Compute vertical advective flux in the interior using TracerN:
530            IF (vertAdvecScheme.EQ.ENUM_CENTERED_2ND) THEN
531               CALL GAD_C2_ADV_R( bi,bj,k, rTrans, TracerN, af, myThid )
532            ELSEIF ( vertAdvecScheme.EQ.ENUM_UPWIND_1RST
533         &          .OR. vertAdvecScheme.EQ.ENUM_DST2 ) THEN
534               CALL GAD_DST2U1_ADV_R( bi,bj,k,vertAdvecScheme,deltaTLev(k),
535         I              rTrans, wFld, TracerN(1-OLx,1-OLy,1,bi,bj),
536         O              af, myThid )
537            ELSEIF (vertAdvecScheme.EQ.ENUM_FLUX_LIMIT) THEN
538               CALL GAD_FLUXLIMIT_ADV_R( bi,bj,k, deltaTLev(k),
539         I              rTrans, wFld, TracerN(1-OLx,1-OLy,1,bi,bj),
540         O              af, myThid )
541            ELSEIF (vertAdvecScheme.EQ.ENUM_UPWIND_3RD ) THEN
542               CALL GAD_U3_ADV_R( bi,bj,k, rTrans, TracerN, af, myThid )
543            ELSEIF (vertAdvecScheme.EQ.ENUM_CENTERED_4TH) THEN
544               CALL GAD_C4_ADV_R( bi,bj,k, rTrans, TracerN, af, myThid )
545    #ifdef ALLOW_AUTODIFF
546            ELSEIF( vertAdvecScheme.EQ.ENUM_DST3 .OR.
547         &         (vertAdvecScheme.EQ.ENUM_DST3_FLUX_LIMIT .AND. inAdMode)
548         &        ) THEN
549    cph This block is to trick the adjoint:
550    cph If inAdExact=.FALSE., we want to use DST3
551    cph with limiters in forward, but without limiters in reverse.
552    #else /* ALLOW_AUTODIFF */
553            ELSEIF( vertAdvecScheme.EQ.ENUM_DST3 ) THEN
554    #endif /* ALLOW_AUTODIFF */
555               CALL GAD_DST3_ADV_R( bi,bj,k, deltaTLev(k),
556         I              rTrans, wFld, TracerN(1-OLx,1-OLy,1,bi,bj),
557         O              af, myThid )
558            ELSEIF (vertAdvecScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
559               CALL GAD_DST3FL_ADV_R( bi,bj,k, deltaTLev(k),
560         I              rTrans, wFld, TracerN(1-OLx,1-OLy,1,bi,bj),
561         O              af, myThid )
562    #ifndef ALLOW_AUTODIFF
563            ELSEIF (vertAdvecScheme.EQ.ENUM_OS7MP ) THEN
564               CALL GAD_OS7MP_ADV_R( bi,bj,k, deltaTLev(k),
565         I              rTrans, wFld, TracerN(1-OLx,1-OLy,1,bi,bj),
566         O              af, myThid )
567    #endif
568            ELSE
569              STOP 'GAD_CALC_RHS: Bad vertAdvecScheme (R)'
570            ENDIF
571         ENDIF         ENDIF
572  C-    Surface "correction" term at k>1 :  C-     add the advective flux to fVerT
573         DO j=1-Oly,sNy+Oly         DO j=1-OLy,sNy+OLy
574          DO i=1-Olx,sNx+Olx          DO i=1-OLx,sNx+OLx
575           af(i,j) = af(i,j)            fVerT(i,j,kUp) = fVerT(i,j,kUp) + af(i,j)*maskInC(i,j,bi,bj)
      &           + (maskC(i,j,k,bi,bj)-maskC(i,j,k-1,bi,bj))*  
      &             rTrans(i,j)*Tracer(i,j,k,bi,bj)  
576          ENDDO          ENDDO
        ENDDO  
       ELSE  
 C-    Surface "correction" term at k=1 :  
        DO j=1-Oly,sNy+Oly  
         DO i=1-Olx,sNx+Olx  
          af(i,j) = rTrans(i,j)*Tracer(i,j,k,bi,bj)  
         ENDDO  
        ENDDO  
       ENDIF  
 C-    add the advective flux to fVerT  
       DO j=1-Oly,sNy+Oly  
        DO i=1-Olx,sNx+Olx  
         fVerT(i,j,kUp) = fVerT(i,j,kUp) + af(i,j)  
577         ENDDO         ENDDO
578        ENDDO  #ifdef ALLOW_DIAGNOSTICS
579           IF ( useDiagnostics ) THEN
580              diagName = 'ADVr'//diagSufx
581              CALL DIAGNOSTICS_FILL( af, diagName, k,1, 2,bi,bj, myThid )
582    C- note: needs to explicitly increment the counter since DIAGNOSTICS_FILL
583    C        does it only if k=1 (never the case here)
584              IF ( k.EQ.2 ) CALL DIAGNOSTICS_COUNT(diagName,bi,bj,myThid)
585           ENDIF
586    #endif
587        ENDIF        ENDIF
588    
589  C-    Diffusive flux in R  C-    Diffusive flux in R
590  C     Note: For K=1 then KM1=1 and this gives a dT/dr = 0 upper  C     Note: For K=1 then KM1=1 and this gives a dT/dr = 0 upper
591  C           boundary condition.  C           boundary condition.
592        IF (implicitDiffusion) THEN        IF (implicitDiffusion) THEN
593         DO j=1-Oly,sNy+Oly         DO j=1-OLy,sNy+OLy
594          DO i=1-Olx,sNx+Olx          DO i=1-OLx,sNx+OLx
595           df(i,j) = 0. _d 0           df(i,j) = 0. _d 0
596          ENDDO          ENDDO
597         ENDDO         ENDDO
598        ELSE        ELSE
599         CALL GAD_DIFF_R(bi,bj,k,KappaRT,tracer,df,myThid)          CALL GAD_DIFF_R(bi,bj,k,KappaR,TracerN,df,myThid)
600          ENDIF
601    
602          IF ( trUseDiffKr4 ) THEN
603            CALL GAD_BIHARM_R( bi,bj,k, diffKr4, TracerN, df, myThid )
604        ENDIF        ENDIF
605    
606  #ifdef ALLOW_GMREDI  #ifdef ALLOW_GMREDI
607  C-    GM/Redi flux in R  C-    GM/Redi flux in R
608        IF (useGMRedi) THEN        IF ( trUseGMRedi ) THEN
609  C *note* should update GMREDI_RTRANSPORT to set df  *aja*  C *note* should update GMREDI_RTRANSPORT to set df  *aja*
610         CALL GMREDI_RTRANSPORT(          CALL GMREDI_RTRANSPORT(
611       I     iMin,iMax,jMin,jMax,bi,bj,K,       I         iMin,iMax,jMin,jMax,bi,bj,k,
612       I     Tracer,tracerIdentity,       I         TracerN,trIdentity,
613       U     df,       U         df,
614       I     myThid)       I         myThid)
615        ENDIF        ENDIF
616  #endif  #endif
617    
618        DO j=1-Oly,sNy+Oly        DO j=1-OLy,sNy+OLy
619         DO i=1-Olx,sNx+Olx         DO i=1-OLx,sNx+OLx
620          fVerT(i,j,kUp) = fVerT(i,j,kUp) + df(i,j)*maskUp(i,j)          fVerT(i,j,kUp) = fVerT(i,j,kUp) + df(i,j)*maskUp(i,j)
621         ENDDO         ENDDO
622        ENDDO        ENDDO
623    
624    #ifdef ALLOW_DIAGNOSTICS
625    C-    Diagnostics of Tracer flux in R dir (mainly Diffusive terms),
626    C       Explicit terms only & excluding advective terms:
627          IF ( useDiagnostics .AND.
628         &    (.NOT.implicitDiffusion .OR. trUseGMRedi) ) THEN
629              diagName = 'DFrE'//diagSufx
630              CALL DIAGNOSTICS_FILL( df, diagName, k,1, 2,bi,bj, myThid )
631    #ifdef ALLOW_LAYERS
632              IF ( useLayers ) THEN
633               CALL LAYERS_FILL_DFR( df, trIdentity, k, 1, 2,bi,bj, myThid )
634              ENDIF
635    #endif /* ALLOW_LAYERS */
636          ENDIF
637    #endif
638    
639  #ifdef ALLOW_KPP  #ifdef ALLOW_KPP
640  C-    Add non local KPP transport term (ghat) to diffusive T flux.  C-    Set non local KPP transport term (ghat):
641        IF (useKPP) THEN        IF ( trUseKPP .AND. k.GE.2 ) THEN
642         DO j=1-Oly,sNy+Oly         DO j=1-OLy,sNy+OLy
643          DO i=1-Olx,sNx+Olx          DO i=1-OLx,sNx+OLx
644           df(i,j) = 0. _d 0           df(i,j) = 0. _d 0
645          ENDDO          ENDDO
646         ENDDO         ENDDO
647         IF (tracerIdentity.EQ.GAD_TEMPERATURE) THEN         IF (trIdentity.EQ.GAD_TEMPERATURE) THEN
 C *note* should update KPP_TRANSPORT_T to set df  *aja*  
648          CALL KPP_TRANSPORT_T(          CALL KPP_TRANSPORT_T(
649       I     iMin,iMax,jMin,jMax,bi,bj,k,km1,       I           iMin,iMax,jMin,jMax,bi,bj,k,km1,
650       I     KappaRT,       O           df,
651       U     df )       I           myTime, myIter, myThid )
652         ELSEIF (tracerIdentity.EQ.GAD_SALINITY) THEN         ELSEIF (trIdentity.EQ.GAD_SALINITY) THEN
653          CALL KPP_TRANSPORT_S(          CALL KPP_TRANSPORT_S(
654       I     iMin,iMax,jMin,jMax,bi,bj,k,km1,       I           iMin,iMax,jMin,jMax,bi,bj,k,km1,
655       I     KappaRT,       O           df,
656       U     df )       I           myTime, myIter, myThid )
657    #ifdef ALLOW_PTRACERS
658           ELSEIF (trIdentity .GE. GAD_TR1) THEN
659            CALL KPP_TRANSPORT_PTR(
660         I           iMin,iMax,jMin,jMax,bi,bj,k,km1,
661         I           trIdentity-GAD_TR1+1,
662         O           df,
663         I           myTime, myIter, myThid )
664    #endif
665         ELSE         ELSE
666          STOP 'GAD_CALC_RHS: Ooops'          WRITE(errorMessageUnit,*)
667         &    'tracer identity =', trIdentity, ' is not valid => STOP'
668            STOP 'ABNORMAL END: S/R GAD_CALC_RHS: invalid tracer identity'
669           ENDIF
670           DO j=1-OLy,sNy+OLy
671            DO i=1-OLx,sNx+OLx
672             fVerT(i,j,kUp) = fVerT(i,j,kUp)
673         &                  + df(i,j)*maskUp(i,j)*rhoFacF(k)
674            ENDDO
675           ENDDO
676    #ifdef ALLOW_DIAGNOSTICS
677    C-    Diagnostics of Non-Local Tracer (vertical) flux
678           IF ( useDiagnostics ) THEN
679             diagName = 'KPPg'//diagSufx
680             CALL DIAGNOSTICS_FILL( df, diagName, k,1, 2,bi,bj, myThid )
681    C- note: needs to explicitly increment the counter since DIAGNOSTICS_FILL
682    C        does it only if k=1 (never the case here)
683             IF ( k.EQ.2 ) CALL DIAGNOSTICS_COUNT(diagName,bi,bj,myThid)
684    #ifdef ALLOW_LAYERS
685             IF ( useLayers ) THEN
686              CALL LAYERS_FILL_DFR( df, trIdentity, k, 1, 2,bi,bj, myThid )
687             ENDIF
688    #endif /* ALLOW_LAYERS */
689         ENDIF         ENDIF
690         DO j=1-Oly,sNy+Oly  #endif
691          DO i=1-Olx,sNx+Olx        ENDIF
692           fVerT(i,j,kUp) = fVerT(i,j,kUp) + df(i,j)*maskUp(i,j)  #endif /* ALLOW_KPP */
693    
694    #ifdef GAD_SMOLARKIEWICZ_HACK
695    coj   Hack to make redi (and everything else in this s/r) positive
696    coj   (see Smolarkiewicz MWR 1989 and Bott MWR 1989).
697    coj   Only works if 'down' is k+1 and k loop in thermodynamics is k=Nr,1,-1
698    coj
699    coj   Apply to all tracers except temperature
700          IF ( trIdentity.NE.GAD_TEMPERATURE .AND.
701         &     trIdentity.NE.GAD_SALINITY ) THEN
702           DO j=1-OLy,sNy+OLy-1
703            DO i=1-OLx,sNx+OLx-1
704    coj   Add outgoing fluxes
705             outFlux=deltaTLev(k)*
706         &    _recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
707         &   *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)*recip_rhoFacC(k)
708         &    *( MAX(0. _d 0,fZon(i+1,j)) + MAX(0. _d 0,-fZon(i,j))
709         &      +MAX(0. _d 0,fMer(i,j+1)) + MAX(0. _d 0,-fMer(i,j))
710         &      +MAX(0. _d 0,fVerT(i,j,kDown)*rkSign)
711         &      +MAX(0. _d 0,-fVerT(i,j,kUp)*rkSign)
712         &     )
713             trac = localT(i,j)
714    coj   If they would reduce tracer by a fraction of more than
715    coj   SmolarkiewiczMaxFrac, scale them down
716             IF (outFlux.GT.0. _d 0 .AND.
717         &       outFlux.GT.SmolarkiewiczMaxFrac*trac) THEN
718    coj   If tracer is already negative, scale flux to zero
719               fac = MAX(0. _d 0,SmolarkiewiczMaxFrac*trac/outFlux)
720    
721               IF (fZon(i+1,j).GT.0. _d 0) fZon(i+1,j)=fac*fZon(i+1,j)
722               IF (-fZon(i,j) .GT.0. _d 0) fZon(i,j)  =fac*fZon(i,j)
723               IF (fMer(i,j+1).GT.0. _d 0) fMer(i,j+1)=fac*fMer(i,j+1)
724               IF (-fMer(i,j) .GT.0. _d 0) fMer(i,j)  =fac*fMer(i,j)
725               IF (-fVerT(i,j,kUp)*rkSign .GT.0. _d 0)
726         &       fVerT(i,j,kUp)=fac*fVerT(i,j,kUp)
727    
728               IF (k.LT.Nr .AND. fVerT(i,j,kDown)*rkSign.GT.0. _d 0) THEN
729    coj   Down flux is special: it has already been applied in lower layer,
730    coj   so we have to readjust this.
731    coj   Note: for k+1, gTracer is now the updated tracer, not the tendency!
732    coj   thus it has an extra factor deltaTLev(k+1)
733                 gTrFac=deltaTLev(k+1)
734    coj   Other factors that have been applied to gTracer since the last call:
735    #ifdef NONLIN_FRSURF
736                 IF (nonlinFreeSurf.GT.0) THEN
737                  IF (select_rStar.GT.0) THEN
738    #ifndef DISABLE_RSTAR_CODE
739                    gTrFac = gTrFac/rStarExpC(i,j,bi,bj)
740    #endif /* DISABLE_RSTAR_CODE */
741                  ENDIF
742                 ENDIF
743    #endif /* NONLIN_FRSURF */
744    coj   Now: undo down flux, ...
745                 gTracer(i,j,k+1) = gTracer(i,j,k+1)
746         &        +gTrFac
747         &         *_recip_hFacC(i,j,k+1,bi,bj)*recip_drF(k+1)
748         &         *recip_rA(i,j,bi,bj)*recip_deepFac2C(k+1)
749         &         *recip_rhoFacC(k+1)
750         &         *( -fVerT(i,j,kDown)*rkSign )
751    coj   ... scale ...
752                 fVerT(i,j,kDown)=fac*fVerT(i,j,kDown)
753    coj   ... and reapply
754                 gTracer(i,j,k+1) = gTracer(i,j,k+1)
755         &        +gTrFac
756         &         *_recip_hFacC(i,j,k+1,bi,bj)*recip_drF(k+1)
757         &         *recip_rA(i,j,bi,bj)*recip_deepFac2C(k+1)
758         &         *recip_rhoFacC(k+1)
759         &         *( fVerT(i,j,kDown)*rkSign )
760               ENDIF
761    
762             ENDIF
763          ENDDO          ENDDO
764         ENDDO         ENDDO
765        ENDIF        ENDIF
766  #endif  #endif
767    
768  C--   Divergence of fluxes  C--   Divergence of fluxes
769        DO j=1-Oly,sNy+Oly-1  C     Anelastic: scale vertical fluxes by rhoFac and leave Horizontal fluxes unchanged
770         DO i=1-Olx,sNx+Olx-1  C     for Stevens OBC: keep only vertical diffusive contribution on boundaries
771          gTracer(i,j,k,bi,bj)=gTracer(i,j,k,bi,bj)        DO j=1-OLy,sNy+OLy-1
772           DO i=1-OLx,sNx+OLx-1
773            gTracer(i,j,k) = gTracer(i,j,k)
774       &   -_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)       &   -_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
775       &    *recip_rA(i,j,bi,bj)       &   *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)*recip_rhoFacC(k)
776       &    *(       &   *( (fZon(i+1,j)-fZon(i,j))*maskInC(i,j,bi,bj)
777       &    +( fZon(i+1,j)-fZon(i,j) )       &     +(fMer(i,j+1)-fMer(i,j))*maskInC(i,j,bi,bj)
778       &    +( fMer(i,j+1)-fMer(i,j) )       &     +(fVerT(i,j,kDown)-fVerT(i,j,kUp))*rkSign
779       &    +( fVerT(i,j,kUp)-fVerT(i,j,kDown) )*rkFac       &     -localT(i,j)*( (uTrans(i+1,j)-uTrans(i,j))*advFac
780         &                   +(vTrans(i,j+1)-vTrans(i,j))*advFac
781         &                   +(rTransKp1(i,j)-rTrans(i,j))*rAdvFac
782         &                  )*maskInC(i,j,bi,bj)
783       &    )       &    )
784         ENDDO         ENDDO
785        ENDDO        ENDDO
786    
787  #ifdef NONLIN_FRSURF  #ifdef ALLOW_DEBUG
788  C-- account for 3.D divergence of the flow in rStar coordinate:        IF ( debugLevel .GE. debLevC
789        IF (calcAdvection .AND. select_rStar.GT.0) THEN       &   .AND. trIdentity.EQ.GAD_TEMPERATURE
790         DO j=1-Oly,sNy+Oly-1       &   .AND. k.EQ.2 .AND. myIter.EQ.1+nIter0
791          DO i=1-Olx,sNx+Olx-1       &   .AND. nPx.EQ.1 .AND. nPy.EQ.1
792           gTracer(i,j,k,bi,bj) = gTracer(i,j,k,bi,bj)       &   .AND. useCubedSphereExchange ) THEN
793       &     - (rStarExpC(i,j,bi,bj) - 1. _d 0)/deltaTfreesurf          CALL DEBUG_CS_CORNER_UV( ' fZon,fMer from GAD_CALC_RHS',
794       &       *tracer(i,j,k,bi,bj)*maskC(i,j,k,bi,bj)       &             fZon,fMer, k, standardMessageUnit,bi,bj,myThid )
         ENDDO  
        ENDDO  
       ENDIF  
       IF (calcAdvection .AND. select_rStar.LT.0) THEN  
        DO j=1-Oly,sNy+Oly-1  
         DO i=1-Olx,sNx+Olx-1  
          gTracer(i,j,k,bi,bj) = gTracer(i,j,k,bi,bj)  
      &     - rStarDhCDt(i,j,bi,bj)  
      &       *tracer(i,j,k,bi,bj)*maskC(i,j,k,bi,bj)  
         ENDDO  
        ENDDO  
795        ENDIF        ENDIF
796  #endif /* NONLIN_FRSURF */  #endif /* ALLOW_DEBUG */
         
797    
798        RETURN        RETURN
799        END        END

Legend:
Removed from v.1.17  
changed lines
  Added in v.1.69

  ViewVC Help
Powered by ViewVC 1.1.22