1 |
jmc |
1.69 |
C $Header: /u/gcmpack/MITgcm/pkg/generic_advdiff/gad_calc_rhs.F,v 1.68 2014/07/22 12:05:23 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 |
jmc |
1.69 |
_RL gTracer(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
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 |
jmc |
1.68 |
localT(i,j)=TracerN(i,j,k,bi,bj) |
195 |
|
|
locABT(i,j)=TracerN(i,j,k,bi,bj) |
196 |
jmc |
1.40 |
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.68 |
CALL GMREDI_XTRANSPORT( |
313 |
jmc |
1.40 |
I iMin,iMax,jMin,jMax,bi,bj,k, |
314 |
jmc |
1.60 |
I xA,TracerN,trIdentity, |
315 |
jmc |
1.40 |
U df, |
316 |
|
|
I myThid) |
317 |
adcroft |
1.1 |
ENDIF |
318 |
|
|
#endif |
319 |
jmc |
1.43 |
C anelastic: advect.fluxes are scaled by rhoFac but hor.diff. flx are not |
320 |
jmc |
1.60 |
DO j=1-OLy,sNy+OLy |
321 |
|
|
DO i=1-OLx,sNx+OLx |
322 |
jmc |
1.43 |
fZon(i,j) = fZon(i,j) + df(i,j)*rhoFacC(k) |
323 |
adcroft |
1.1 |
ENDDO |
324 |
|
|
ENDDO |
325 |
|
|
|
326 |
jmc |
1.32 |
#ifdef ALLOW_DIAGNOSTICS |
327 |
|
|
C- Diagnostics of Tracer flux in X dir (mainly Diffusive term), |
328 |
|
|
C excluding advective terms: |
329 |
|
|
IF ( useDiagnostics .AND. |
330 |
jmc |
1.48 |
& (diffKh.NE.0. .OR. diffK4 .NE.0. .OR. trUseGMRedi) ) THEN |
331 |
jmc |
1.42 |
diagName = 'DFxE'//diagSufx |
332 |
jmc |
1.60 |
CALL DIAGNOSTICS_FILL( df, diagName, k,1, 2,bi,bj, myThid ) |
333 |
rpa |
1.66 |
#ifdef ALLOW_LAYERS |
334 |
rpa |
1.67 |
IF ( useLayers ) THEN |
335 |
|
|
CALL LAYERS_FILL_DFX( df, trIdentity, k, 1, 2,bi,bj, myThid ) |
336 |
|
|
ENDIF |
337 |
rpa |
1.66 |
#endif /* ALLOW_LAYERS */ |
338 |
adcroft |
1.1 |
ENDIF |
339 |
jmc |
1.32 |
#endif |
340 |
adcroft |
1.1 |
|
341 |
|
|
C-- Initialize net flux in Y direction |
342 |
jmc |
1.60 |
DO j=1-OLy,sNy+OLy |
343 |
|
|
DO i=1-OLx,sNx+OLx |
344 |
heimbach |
1.12 |
fMer(i,j) = 0. _d 0 |
345 |
adcroft |
1.1 |
ENDDO |
346 |
|
|
ENDDO |
347 |
|
|
|
348 |
|
|
C- Advective flux in Y |
349 |
jmc |
1.14 |
IF (calcAdvection) THEN |
350 |
jmc |
1.32 |
IF (advectionScheme.EQ.ENUM_CENTERED_2ND) THEN |
351 |
jmc |
1.56 |
CALL GAD_C2_ADV_Y( bi,bj,k, vTrans, locABT, af, myThid ) |
352 |
jmc |
1.41 |
ELSEIF ( advectionScheme.EQ.ENUM_UPWIND_1RST |
353 |
jmc |
1.37 |
& .OR. advectionScheme.EQ.ENUM_DST2 ) THEN |
354 |
jmc |
1.56 |
CALL GAD_DST2U1_ADV_Y( bi,bj,k, advectionScheme, .TRUE., |
355 |
|
|
I deltaTLev(k), vTrans, vFld, locABT, |
356 |
|
|
O af, myThid ) |
357 |
|
|
ELSE |
358 |
|
|
DO j=1-OLy,sNy+OLy |
359 |
|
|
DO i=1-OLx,sNx+OLx |
360 |
|
|
#ifdef ALLOW_OBCS |
361 |
|
|
maskLocS(i,j) = _maskS(i,j,k,bi,bj)*maskInS(i,j,bi,bj) |
362 |
|
|
#else /* ALLOW_OBCS */ |
363 |
|
|
maskLocS(i,j) = _maskS(i,j,k,bi,bj) |
364 |
|
|
#endif /* ALLOW_OBCS */ |
365 |
|
|
ENDDO |
366 |
|
|
ENDDO |
367 |
|
|
IF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN |
368 |
|
|
CALL GAD_FLUXLIMIT_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k), |
369 |
|
|
I vTrans, vFld, maskLocS, locABT, |
370 |
|
|
O af, myThid ) |
371 |
|
|
ELSEIF (advectionScheme.EQ.ENUM_UPWIND_3RD ) THEN |
372 |
|
|
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 |
jmc |
1.62 |
#ifdef ALLOW_AUTODIFF |
378 |
|
|
ELSEIF( advectionScheme.EQ.ENUM_DST3 .OR. |
379 |
|
|
& (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT .AND. inAdMode) |
380 |
|
|
& ) THEN |
381 |
heimbach |
1.35 |
cph This block is to trick the adjoint: |
382 |
jmc |
1.62 |
cph If inAdExact=.FALSE., we want to use DST3 |
383 |
heimbach |
1.35 |
cph with limiters in forward, but without limiters in reverse. |
384 |
jmc |
1.62 |
#else /* ALLOW_AUTODIFF */ |
385 |
|
|
ELSEIF( advectionScheme.EQ.ENUM_DST3 ) THEN |
386 |
|
|
#endif /* ALLOW_AUTODIFF */ |
387 |
jmc |
1.56 |
CALL GAD_DST3_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k), |
388 |
|
|
I vTrans, vFld, maskLocS, locABT, |
389 |
|
|
O af, myThid ) |
390 |
jmc |
1.62 |
ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN |
391 |
jmc |
1.56 |
CALL GAD_DST3FL_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k), |
392 |
|
|
I vTrans, vFld, maskLocS, locABT, |
393 |
|
|
O af, myThid ) |
394 |
jmc |
1.64 |
#ifndef ALLOW_AUTODIFF |
395 |
jmc |
1.56 |
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 |
heimbach |
1.35 |
ELSE |
401 |
jmc |
1.56 |
STOP 'GAD_CALC_RHS: Bad advectionScheme (Y)' |
402 |
heimbach |
1.35 |
ENDIF |
403 |
jmc |
1.32 |
ENDIF |
404 |
jmc |
1.60 |
#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 |
jmc |
1.32 |
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 |
jmc |
1.60 |
CALL DIAGNOSTICS_FILL( af, diagName, k,1, 2,bi,bj, myThid ) |
422 |
jmc |
1.32 |
ENDIF |
423 |
|
|
#endif |
424 |
adcroft |
1.8 |
ENDIF |
425 |
adcroft |
1.1 |
|
426 |
|
|
C- Diffusive flux in Y |
427 |
|
|
IF (diffKh.NE.0.) THEN |
428 |
|
|
CALL GAD_DIFF_Y(bi,bj,k,yA,diffKh,localT,df,myThid) |
429 |
|
|
ELSE |
430 |
jmc |
1.60 |
DO j=1-OLy,sNy+OLy |
431 |
|
|
DO i=1-OLx,sNx+OLx |
432 |
heimbach |
1.12 |
df(i,j) = 0. _d 0 |
433 |
adcroft |
1.1 |
ENDDO |
434 |
|
|
ENDDO |
435 |
|
|
ENDIF |
436 |
|
|
|
437 |
jmc |
1.32 |
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 |
adcroft |
1.1 |
#ifdef ALLOW_GMREDI |
443 |
|
|
C- GM/Redi flux in Y |
444 |
jmc |
1.48 |
IF ( trUseGMRedi ) THEN |
445 |
jmc |
1.38 |
C *note* should update GMREDI_YTRANSPORT to set df *aja* |
446 |
jmc |
1.68 |
CALL GMREDI_YTRANSPORT( |
447 |
jmc |
1.40 |
I iMin,iMax,jMin,jMax,bi,bj,k, |
448 |
jmc |
1.60 |
I yA,TracerN,trIdentity, |
449 |
jmc |
1.40 |
U df, |
450 |
|
|
I myThid) |
451 |
adcroft |
1.1 |
ENDIF |
452 |
|
|
#endif |
453 |
jmc |
1.43 |
C anelastic: advect.fluxes are scaled by rhoFac but hor.diff. flx are not |
454 |
jmc |
1.60 |
DO j=1-OLy,sNy+OLy |
455 |
|
|
DO i=1-OLx,sNx+OLx |
456 |
jmc |
1.43 |
fMer(i,j) = fMer(i,j) + df(i,j)*rhoFacC(k) |
457 |
adcroft |
1.1 |
ENDDO |
458 |
|
|
ENDDO |
459 |
|
|
|
460 |
jmc |
1.32 |
#ifdef ALLOW_DIAGNOSTICS |
461 |
|
|
C- Diagnostics of Tracer flux in Y dir (mainly Diffusive terms), |
462 |
|
|
C excluding advective terms: |
463 |
|
|
IF ( useDiagnostics .AND. |
464 |
jmc |
1.48 |
& (diffKh.NE.0. .OR. diffK4 .NE.0. .OR. trUseGMRedi) ) THEN |
465 |
jmc |
1.42 |
diagName = 'DFyE'//diagSufx |
466 |
jmc |
1.60 |
CALL DIAGNOSTICS_FILL( df, diagName, k,1, 2,bi,bj, myThid ) |
467 |
rpa |
1.66 |
#ifdef ALLOW_LAYERS |
468 |
rpa |
1.67 |
IF ( useLayers ) THEN |
469 |
|
|
CALL LAYERS_FILL_DFY( df, trIdentity,k, 1, 2,bi,bj, myThid ) |
470 |
|
|
ENDIF |
471 |
rpa |
1.66 |
#endif /* ALLOW_LAYERS */ |
472 |
adcroft |
1.1 |
ENDIF |
473 |
jmc |
1.32 |
#endif |
474 |
adcroft |
1.1 |
|
475 |
jmc |
1.16 |
C-- Compute vertical flux fVerT(kUp) at interface k (between k-1 & k): |
476 |
adcroft |
1.1 |
C- Advective flux in R |
477 |
jmc |
1.25 |
#ifdef ALLOW_AIM |
478 |
|
|
C- a hack to prevent Water-Vapor vert.transport into the stratospheric level Nr |
479 |
jmc |
1.40 |
IF (calcAdvection .AND. .NOT.implicitAdvection .AND. k.GE.2 .AND. |
480 |
jmc |
1.60 |
& (.NOT.useAIM .OR. trIdentity.NE.GAD_SALINITY .OR. k.LT.Nr) |
481 |
jmc |
1.25 |
& ) THEN |
482 |
|
|
#else |
483 |
jmc |
1.40 |
IF (calcAdvection .AND. .NOT.implicitAdvection .AND. k.GE.2) THEN |
484 |
jmc |
1.25 |
#endif |
485 |
jmc |
1.68 |
IF ( applyAB_onTracer ) THEN |
486 |
|
|
C- Compute vertical advective flux in the interior using TracAB: |
487 |
jmc |
1.32 |
IF (vertAdvecScheme.EQ.ENUM_CENTERED_2ND) THEN |
488 |
jmc |
1.60 |
CALL GAD_C2_ADV_R( bi,bj,k, rTrans, TracAB, af, myThid ) |
489 |
jmc |
1.41 |
ELSEIF ( vertAdvecScheme.EQ.ENUM_UPWIND_1RST |
490 |
jmc |
1.37 |
& .OR. vertAdvecScheme.EQ.ENUM_DST2 ) THEN |
491 |
jmc |
1.60 |
CALL GAD_DST2U1_ADV_R( bi,bj,k,vertAdvecScheme,deltaTLev(k), |
492 |
|
|
I rTrans, wFld, TracAB(1-OLx,1-OLy,1,bi,bj), |
493 |
|
|
O af, myThid ) |
494 |
jmc |
1.32 |
ELSEIF (vertAdvecScheme.EQ.ENUM_FLUX_LIMIT) THEN |
495 |
jmc |
1.60 |
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 |
jmc |
1.32 |
ELSEIF (vertAdvecScheme.EQ.ENUM_UPWIND_3RD ) THEN |
499 |
jmc |
1.60 |
CALL GAD_U3_ADV_R( bi,bj,k, rTrans, TracAB, af, myThid ) |
500 |
jmc |
1.32 |
ELSEIF (vertAdvecScheme.EQ.ENUM_CENTERED_4TH) THEN |
501 |
jmc |
1.60 |
CALL GAD_C4_ADV_R( bi,bj,k, rTrans, TracAB, af, myThid ) |
502 |
jmc |
1.62 |
#ifdef ALLOW_AUTODIFF |
503 |
|
|
ELSEIF( vertAdvecScheme.EQ.ENUM_DST3 .OR. |
504 |
|
|
& (vertAdvecScheme.EQ.ENUM_DST3_FLUX_LIMIT .AND. inAdMode) |
505 |
|
|
& ) THEN |
506 |
heimbach |
1.35 |
cph This block is to trick the adjoint: |
507 |
jmc |
1.62 |
cph If inAdExact=.FALSE., we want to use DST3 |
508 |
heimbach |
1.35 |
cph with limiters in forward, but without limiters in reverse. |
509 |
jmc |
1.62 |
#else /* ALLOW_AUTODIFF */ |
510 |
|
|
ELSEIF( vertAdvecScheme.EQ.ENUM_DST3 ) THEN |
511 |
|
|
#endif /* ALLOW_AUTODIFF */ |
512 |
jmc |
1.60 |
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 |
jmc |
1.62 |
ELSEIF (vertAdvecScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN |
516 |
jmc |
1.60 |
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 |
jmc |
1.64 |
#ifndef ALLOW_AUTODIFF |
520 |
adcroft |
1.44 |
ELSEIF (vertAdvecScheme.EQ.ENUM_OS7MP ) THEN |
521 |
jmc |
1.60 |
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 |
heimbach |
1.55 |
#endif |
525 |
jmc |
1.32 |
ELSE |
526 |
|
|
STOP 'GAD_CALC_RHS: Bad vertAdvecScheme (R)' |
527 |
|
|
ENDIF |
528 |
jmc |
1.68 |
ELSE |
529 |
|
|
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 |
572 |
jmc |
1.23 |
C- add the advective flux to fVerT |
573 |
jmc |
1.68 |
DO j=1-OLy,sNy+OLy |
574 |
|
|
DO i=1-OLx,sNx+OLx |
575 |
mlosch |
1.59 |
fVerT(i,j,kUp) = fVerT(i,j,kUp) + af(i,j)*maskInC(i,j,bi,bj) |
576 |
jmc |
1.2 |
ENDDO |
577 |
jmc |
1.68 |
ENDDO |
578 |
jmc |
1.32 |
#ifdef ALLOW_DIAGNOSTICS |
579 |
jmc |
1.68 |
IF ( useDiagnostics ) THEN |
580 |
jmc |
1.32 |
diagName = 'ADVr'//diagSufx |
581 |
jmc |
1.60 |
CALL DIAGNOSTICS_FILL( af, diagName, k,1, 2,bi,bj, myThid ) |
582 |
jmc |
1.34 |
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 |
jmc |
1.68 |
ENDIF |
586 |
jmc |
1.32 |
#endif |
587 |
adcroft |
1.8 |
ENDIF |
588 |
adcroft |
1.1 |
|
589 |
|
|
C- Diffusive flux in R |
590 |
|
|
C Note: For K=1 then KM1=1 and this gives a dT/dr = 0 upper |
591 |
|
|
C boundary condition. |
592 |
|
|
IF (implicitDiffusion) THEN |
593 |
jmc |
1.60 |
DO j=1-OLy,sNy+OLy |
594 |
|
|
DO i=1-OLx,sNx+OLx |
595 |
heimbach |
1.12 |
df(i,j) = 0. _d 0 |
596 |
adcroft |
1.1 |
ENDDO |
597 |
|
|
ENDDO |
598 |
|
|
ELSE |
599 |
jmc |
1.68 |
CALL GAD_DIFF_R(bi,bj,k,KappaR,TracerN,df,myThid) |
600 |
adcroft |
1.1 |
ENDIF |
601 |
|
|
|
602 |
jmc |
1.61 |
IF ( trUseDiffKr4 ) THEN |
603 |
jmc |
1.68 |
CALL GAD_BIHARM_R( bi,bj,k, diffKr4, TracerN, df, myThid ) |
604 |
jmc |
1.61 |
ENDIF |
605 |
|
|
|
606 |
adcroft |
1.1 |
#ifdef ALLOW_GMREDI |
607 |
|
|
C- GM/Redi flux in R |
608 |
jmc |
1.48 |
IF ( trUseGMRedi ) THEN |
609 |
adcroft |
1.1 |
C *note* should update GMREDI_RTRANSPORT to set df *aja* |
610 |
jmc |
1.68 |
CALL GMREDI_RTRANSPORT( |
611 |
jmc |
1.40 |
I iMin,iMax,jMin,jMax,bi,bj,k, |
612 |
jmc |
1.60 |
I TracerN,trIdentity, |
613 |
jmc |
1.40 |
U df, |
614 |
|
|
I myThid) |
615 |
adcroft |
1.1 |
ENDIF |
616 |
|
|
#endif |
617 |
|
|
|
618 |
jmc |
1.60 |
DO j=1-OLy,sNy+OLy |
619 |
|
|
DO i=1-OLx,sNx+OLx |
620 |
adcroft |
1.11 |
fVerT(i,j,kUp) = fVerT(i,j,kUp) + df(i,j)*maskUp(i,j) |
621 |
adcroft |
1.1 |
ENDDO |
622 |
|
|
ENDDO |
623 |
|
|
|
624 |
jmc |
1.32 |
#ifdef ALLOW_DIAGNOSTICS |
625 |
jmc |
1.41 |
C- Diagnostics of Tracer flux in R dir (mainly Diffusive terms), |
626 |
jmc |
1.32 |
C Explicit terms only & excluding advective terms: |
627 |
|
|
IF ( useDiagnostics .AND. |
628 |
jmc |
1.48 |
& (.NOT.implicitDiffusion .OR. trUseGMRedi) ) THEN |
629 |
jmc |
1.32 |
diagName = 'DFrE'//diagSufx |
630 |
jmc |
1.60 |
CALL DIAGNOSTICS_FILL( df, diagName, k,1, 2,bi,bj, myThid ) |
631 |
rpa |
1.66 |
#ifdef ALLOW_LAYERS |
632 |
rpa |
1.67 |
IF ( useLayers ) THEN |
633 |
|
|
CALL LAYERS_FILL_DFR( df, trIdentity, k, 1, 2,bi,bj, myThid ) |
634 |
|
|
ENDIF |
635 |
rpa |
1.66 |
#endif /* ALLOW_LAYERS */ |
636 |
jmc |
1.32 |
ENDIF |
637 |
|
|
#endif |
638 |
|
|
|
639 |
adcroft |
1.1 |
#ifdef ALLOW_KPP |
640 |
jmc |
1.29 |
C- Set non local KPP transport term (ghat): |
641 |
jmc |
1.48 |
IF ( trUseKPP .AND. k.GE.2 ) THEN |
642 |
jmc |
1.60 |
DO j=1-OLy,sNy+OLy |
643 |
|
|
DO i=1-OLx,sNx+OLx |
644 |
heimbach |
1.12 |
df(i,j) = 0. _d 0 |
645 |
adcroft |
1.1 |
ENDDO |
646 |
|
|
ENDDO |
647 |
jmc |
1.60 |
IF (trIdentity.EQ.GAD_TEMPERATURE) THEN |
648 |
adcroft |
1.1 |
CALL KPP_TRANSPORT_T( |
649 |
jmc |
1.47 |
I iMin,iMax,jMin,jMax,bi,bj,k,km1, |
650 |
|
|
O df, |
651 |
|
|
I myTime, myIter, myThid ) |
652 |
jmc |
1.60 |
ELSEIF (trIdentity.EQ.GAD_SALINITY) THEN |
653 |
adcroft |
1.1 |
CALL KPP_TRANSPORT_S( |
654 |
jmc |
1.47 |
I iMin,iMax,jMin,jMax,bi,bj,k,km1, |
655 |
|
|
O df, |
656 |
|
|
I myTime, myIter, myThid ) |
657 |
mlosch |
1.18 |
#ifdef ALLOW_PTRACERS |
658 |
jmc |
1.60 |
ELSEIF (trIdentity .GE. GAD_TR1) THEN |
659 |
mlosch |
1.18 |
CALL KPP_TRANSPORT_PTR( |
660 |
jmc |
1.47 |
I iMin,iMax,jMin,jMax,bi,bj,k,km1, |
661 |
jmc |
1.60 |
I trIdentity-GAD_TR1+1, |
662 |
jmc |
1.47 |
O df, |
663 |
|
|
I myTime, myIter, myThid ) |
664 |
mlosch |
1.18 |
#endif |
665 |
adcroft |
1.1 |
ELSE |
666 |
jmc |
1.54 |
WRITE(errorMessageUnit,*) |
667 |
jmc |
1.60 |
& 'tracer identity =', trIdentity, ' is not valid => STOP' |
668 |
jmc |
1.54 |
STOP 'ABNORMAL END: S/R GAD_CALC_RHS: invalid tracer identity' |
669 |
adcroft |
1.1 |
ENDIF |
670 |
jmc |
1.60 |
DO j=1-OLy,sNy+OLy |
671 |
|
|
DO i=1-OLx,sNx+OLx |
672 |
jmc |
1.43 |
fVerT(i,j,kUp) = fVerT(i,j,kUp) |
673 |
|
|
& + df(i,j)*maskUp(i,j)*rhoFacF(k) |
674 |
adcroft |
1.1 |
ENDDO |
675 |
|
|
ENDDO |
676 |
jmc |
1.54 |
#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 |
rpa |
1.66 |
#ifdef ALLOW_LAYERS |
685 |
rpa |
1.67 |
IF ( useLayers ) THEN |
686 |
|
|
CALL LAYERS_FILL_DFR( df, trIdentity, k, 1, 2,bi,bj, myThid ) |
687 |
|
|
ENDIF |
688 |
rpa |
1.66 |
#endif /* ALLOW_LAYERS */ |
689 |
jmc |
1.54 |
ENDIF |
690 |
|
|
#endif |
691 |
adcroft |
1.1 |
ENDIF |
692 |
jmc |
1.54 |
#endif /* ALLOW_KPP */ |
693 |
adcroft |
1.1 |
|
694 |
jahn |
1.51 |
#ifdef GAD_SMOLARKIEWICZ_HACK |
695 |
jahn |
1.52 |
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 |
jahn |
1.51 |
coj |
699 |
jahn |
1.52 |
coj Apply to all tracers except temperature |
700 |
jmc |
1.60 |
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 |
jahn |
1.52 |
coj Add outgoing fluxes |
705 |
jahn |
1.53 |
outFlux=deltaTLev(k)* |
706 |
jahn |
1.51 |
& _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 |
jmc |
1.65 |
trac = localT(i,j) |
714 |
jahn |
1.52 |
coj If they would reduce tracer by a fraction of more than |
715 |
|
|
coj SmolarkiewiczMaxFrac, scale them down |
716 |
jahn |
1.51 |
IF (outFlux.GT.0. _d 0 .AND. |
717 |
|
|
& outFlux.GT.SmolarkiewiczMaxFrac*trac) THEN |
718 |
jahn |
1.52 |
coj If tracer is already negative, scale flux to zero |
719 |
jahn |
1.51 |
fac = MAX(0. _d 0,SmolarkiewiczMaxFrac*trac/outFlux) |
720 |
jahn |
1.52 |
|
721 |
jahn |
1.51 |
IF (fZon(i+1,j).GT.0. _d 0) fZon(i+1,j)=fac*fZon(i+1,j) |
722 |
jmc |
1.54 |
IF (-fZon(i,j) .GT.0. _d 0) fZon(i,j) =fac*fZon(i,j) |
723 |
jahn |
1.51 |
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 |
jahn |
1.52 |
& 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 |
jahn |
1.53 |
coj thus it has an extra factor deltaTLev(k+1) |
733 |
|
|
gTrFac=deltaTLev(k+1) |
734 |
jahn |
1.52 |
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 |
jmc |
1.69 |
gTracer(i,j,k+1) = gTracer(i,j,k+1) |
746 |
jahn |
1.52 |
& +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 |
jmc |
1.69 |
gTracer(i,j,k+1) = gTracer(i,j,k+1) |
755 |
jahn |
1.52 |
& +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 |
jahn |
1.51 |
ENDIF |
761 |
jahn |
1.52 |
|
762 |
jahn |
1.51 |
ENDIF |
763 |
|
|
ENDDO |
764 |
|
|
ENDDO |
765 |
|
|
ENDIF |
766 |
|
|
#endif |
767 |
|
|
|
768 |
adcroft |
1.1 |
C-- Divergence of fluxes |
769 |
jmc |
1.43 |
C Anelastic: scale vertical fluxes by rhoFac and leave Horizontal fluxes unchanged |
770 |
mlosch |
1.59 |
C for Stevens OBC: keep only vertical diffusive contribution on boundaries |
771 |
jmc |
1.60 |
DO j=1-OLy,sNy+OLy-1 |
772 |
|
|
DO i=1-OLx,sNx+OLx-1 |
773 |
jmc |
1.69 |
gTracer(i,j,k) = gTracer(i,j,k) |
774 |
jmc |
1.43 |
& -_recip_hFacC(i,j,k,bi,bj)*recip_drF(k) |
775 |
|
|
& *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)*recip_rhoFacC(k) |
776 |
mlosch |
1.59 |
& *( (fZon(i+1,j)-fZon(i,j))*maskInC(i,j,bi,bj) |
777 |
|
|
& +(fMer(i,j+1)-fMer(i,j))*maskInC(i,j,bi,bj) |
778 |
jmc |
1.36 |
& +(fVerT(i,j,kDown)-fVerT(i,j,kUp))*rkSign |
779 |
jmc |
1.58 |
& -localT(i,j)*( (uTrans(i+1,j)-uTrans(i,j))*advFac |
780 |
|
|
& +(vTrans(i,j+1)-vTrans(i,j))*advFac |
781 |
jmc |
1.36 |
& +(rTransKp1(i,j)-rTrans(i,j))*rAdvFac |
782 |
mlosch |
1.59 |
& )*maskInC(i,j,bi,bj) |
783 |
adcroft |
1.1 |
& ) |
784 |
|
|
ENDDO |
785 |
|
|
ENDDO |
786 |
|
|
|
787 |
jmc |
1.27 |
#ifdef ALLOW_DEBUG |
788 |
jmc |
1.57 |
IF ( debugLevel .GE. debLevC |
789 |
jmc |
1.60 |
& .AND. trIdentity.EQ.GAD_TEMPERATURE |
790 |
jmc |
1.27 |
& .AND. k.EQ.2 .AND. myIter.EQ.1+nIter0 |
791 |
|
|
& .AND. nPx.EQ.1 .AND. nPy.EQ.1 |
792 |
|
|
& .AND. useCubedSphereExchange ) THEN |
793 |
|
|
CALL DEBUG_CS_CORNER_UV( ' fZon,fMer from GAD_CALC_RHS', |
794 |
|
|
& fZon,fMer, k, standardMessageUnit,bi,bj,myThid ) |
795 |
|
|
ENDIF |
796 |
|
|
#endif /* ALLOW_DEBUG */ |
797 |
jmc |
1.41 |
|
798 |
adcroft |
1.1 |
RETURN |
799 |
|
|
END |