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

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

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


Revision 1.64 - (hide annotations) (download)
Fri Apr 4 20:29:08 2014 UTC (10 years, 2 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64w, checkpoint64v
Changes since 1.63: +5 -10 lines
- Replace ALLOW_AUTODIFF_TAMC by ALLOW_AUTODIFF (except for tape/storage
  which are specific to TAF/TAMC).

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

  ViewVC Help
Powered by ViewVC 1.1.22