/[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.63 - (hide annotations) (download)
Fri Dec 27 15:52:17 2013 UTC (10 years, 5 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64s, checkpoint64u, checkpoint64t
Changes since 1.62: +6 -16 lines
add fZon & fMer as output argument of S/R GAD_CALC_RHS

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

  ViewVC Help
Powered by ViewVC 1.1.22