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