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

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

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


Revision 1.48 - (show annotations) (download)
Tue Sep 18 15:27:56 2007 UTC (16 years, 8 months ago) by jmc
Branch: MAIN
Changes since 1.47: +32 -9 lines
poor handeling of PTRACERS_useGMRedi & PTRACERS_useKPP (with no equivalent
flag for T & S) is prone to serious bugs, and rarely (never ?) tested.
Temporary fix (but trUseGMRedi & trUseKPP should be passed as argument).

1 C $Header: /u/gcmpack/MITgcm/pkg/generic_advdiff/gad_calc_rhs.F,v 1.47 2007/05/03 21:34:39 jmc Exp $
2 C $Name: $
3
4 #include "GAD_OPTIONS.h"
5 #ifdef ALLOW_PTRACERS
6 # include "PTRACERS_OPTIONS.h"
7 #endif
8
9 CBOP
10 C !ROUTINE: GAD_CALC_RHS
11
12 C !INTERFACE: ==========================================================
13 SUBROUTINE GAD_CALC_RHS(
14 I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
15 I xA, yA, maskUp, uFld, vFld, wFld,
16 I uTrans, vTrans, rTrans, rTransKp1,
17 I diffKh, diffK4, KappaR, TracerN, TracAB,
18 I tracerIdentity, advectionScheme, vertAdvecScheme,
19 I calcAdvection, implicitAdvection, applyAB_onTracer,
20 U fVerT, gTracer,
21 I myTime, myIter, myThid )
22
23 C !DESCRIPTION:
24 C Calculates the tendency of a tracer due to advection and diffusion.
25 C It calculates the fluxes in each direction indepentently and then
26 C sets the tendency to the divergence of these fluxes. The advective
27 C fluxes are only calculated here when using the linear advection schemes
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 C The tendency is the divergence of the fluxes:
36 C \begin{equation*}
37 C G_\theta = G_\theta + \nabla \cdot {\bf F}
38 C \end{equation*}
39 C
40 C The tendency is assumed to contain data on entry.
41
42 C !USES: ===============================================================
43 IMPLICIT NONE
44 #include "SIZE.h"
45 #include "EEPARAMS.h"
46 #include "PARAMS.h"
47 #include "GRID.h"
48 #include "SURFACE.h"
49 #include "GAD.h"
50 #ifdef ALLOW_PTRACERS
51 # include "PTRACERS_SIZE.h"
52 # include "PTRACERS.h"
53 #endif
54
55 #ifdef ALLOW_AUTODIFF_TAMC
56 #include "tamc.h"
57 #include "tamc_keys.h"
58 #endif /* ALLOW_AUTODIFF_TAMC */
59
60 C !INPUT PARAMETERS: ===================================================
61 C bi,bj :: tile indices
62 C iMin,iMax :: loop range for called routines
63 C jMin,jMax :: loop range for called routines
64 C k :: vertical index
65 C kM1 :: =k-1 for k>1, =1 for k=1
66 C kUp :: index into 2 1/2D array, toggles between 1|2
67 C kDown :: index into 2 1/2D array, toggles between 2|1
68 C xA,yA :: areas of X and Y face of tracer cells
69 C maskUp :: 2-D array for mask at W points
70 C uFld,vFld,wFld :: Local copy of velocity field (3 components)
71 C uTrans,vTrans :: 2-D arrays of volume transports at U,V points
72 C rTrans :: 2-D arrays of volume transports at W points
73 C rTransKp1 :: 2-D array of volume trans at W pts, interf k+1
74 C diffKh :: horizontal diffusion coefficient
75 C diffK4 :: bi-harmonic diffusion coefficient
76 C KappaR :: 2-D array for vertical diffusion coefficient, interf k
77 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 C or extrapolated fwd in time to n+1/2 if applying AB on Tr)
81 C tracerIdentity :: tracer identifier (required for KPP,GM)
82 C advectionScheme :: advection scheme to use (Horizontal plane)
83 C vertAdvecScheme :: advection scheme to use (Vertical direction)
84 C calcAdvection :: =False if Advec computed with multiDim scheme
85 C--- needs to pass those 2 flags as argument
86 C trUseGMRedi :: true if this tracer uses GM-Redi
87 C trUseKPP :: true if this tracer uses KPP
88 C---
89 C implicitAdvection:: =True if vertical Advec computed implicitly
90 C applyAB_onTracer :: apply Adams-Bashforth on Tracer (rather than on gTr)
91 C myTime :: current time
92 C myIter :: iteration number
93 C myThid :: thread number
94 INTEGER bi,bj,iMin,iMax,jMin,jMax
95 INTEGER k,kUp,kDown,kM1
96 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
98 _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
99 _RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
100 _RL vFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
101 _RL wFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
102 _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
103 _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
104 _RL rTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
105 _RL rTransKp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
106 _RL diffKh, diffK4
107 _RL KappaR(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
108 _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 INTEGER tracerIdentity
111 INTEGER advectionScheme, vertAdvecScheme
112 LOGICAL calcAdvection, trUseGMRedi, trUseKPP
113 LOGICAL implicitAdvection, applyAB_onTracer
114 _RL myTime
115 INTEGER myIter, myThid
116
117 C !OUTPUT PARAMETERS: ==================================================
118 C gTracer :: tendency array
119 C fVerT :: 2 1/2D arrays for vertical advective flux
120 _RL gTracer(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
121 _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
122
123 C !LOCAL VARIABLES: ====================================================
124 C i,j :: loop indices
125 C df4 :: used for storing del^2 T for bi-harmonic term
126 C fZon :: zonal flux
127 C fMer :: meridional flux
128 C af :: advective flux
129 C df :: diffusive flux
130 C localT :: local copy of tracer field
131 C locABT :: local copy of (AB-extrapolated) tracer field
132 #ifdef ALLOW_DIAGNOSTICS
133 CHARACTER*8 diagName
134 CHARACTER*4 GAD_DIAG_SUFX, diagSufx
135 EXTERNAL GAD_DIAG_SUFX
136 #endif
137 INTEGER i,j
138 _RL df4 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
139 _RL fZon (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
140 _RL fMer (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 _RL locABT(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
145 _RL advFac, rAdvFac
146 CEOP
147
148 #ifdef ALLOW_AUTODIFF_TAMC
149 C-- only the kUp part of fverT is set in this subroutine
150 C-- the kDown is still required
151 fVerT(1,1,kDown) = fVerT(1,1,kDown)
152 #endif
153
154 trUseGMRedi = useGMRedi
155 trUseKPP = useKPP
156 #ifdef ALLOW_PTRACERS
157 IF ( usePTRACERS ) THEN
158 i = MAX( 1, tracerIdentity - GAD_TR1 + 1 )
159 trUseGMRedi = useGMRedi .AND. ( tracerIdentity.LT.GAD_TR1
160 & .OR. PTRACERS_useGMRedi(i) )
161 trUseKPP = useKPP .AND. ( tracerIdentity.LT.GAD_TR1
162 & .OR. PTRACERS_useKPP(i) )
163 ENDIF
164 #endif
165
166 #ifdef ALLOW_DIAGNOSTICS
167 C-- Set diagnostic suffix for the current tracer
168 IF ( useDiagnostics ) THEN
169 diagSufx = GAD_DIAG_SUFX( tracerIdentity, myThid )
170 ENDIF
171 #endif
172
173 advFac = 0. _d 0
174 IF (calcAdvection) advFac = 1. _d 0
175 rAdvFac = rkSign*advFac
176 IF (implicitAdvection) rAdvFac = 0. _d 0
177
178 DO j=1-OLy,sNy+OLy
179 DO i=1-OLx,sNx+OLx
180 fZon(i,j) = 0. _d 0
181 fMer(i,j) = 0. _d 0
182 fVerT(i,j,kUp) = 0. _d 0
183 df(i,j) = 0. _d 0
184 df4(i,j) = 0. _d 0
185 ENDDO
186 ENDDO
187
188 C-- Make local copy of tracer array
189 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
205 C-- Unless we have already calculated the advection terms we initialize
206 C the tendency to zero.
207 C <== now done earlier at the beginning of thermodynamics.
208 c IF (calcAdvection) THEN
209 c DO j=1-Oly,sNy+Oly
210 c DO i=1-Olx,sNx+Olx
211 c gTracer(i,j,k,bi,bj)=0. _d 0
212 c ENDDO
213 c ENDDO
214 c ENDIF
215
216 C-- Pre-calculate del^2 T if bi-harmonic coefficient is non-zero
217 IF (diffK4 .NE. 0.) THEN
218 CALL GAD_GRAD_X(bi,bj,k,xA,localT,fZon,myThid)
219 CALL GAD_GRAD_Y(bi,bj,k,yA,localT,fMer,myThid)
220 CALL GAD_DEL2(bi,bj,k,fZon,fMer,df4,myThid)
221 ENDIF
222
223 C-- Initialize net flux in X direction
224 DO j=1-Oly,sNy+Oly
225 DO i=1-Olx,sNx+Olx
226 fZon(i,j) = 0. _d 0
227 ENDDO
228 ENDDO
229
230 C- Advective flux in X
231 IF (calcAdvection) THEN
232 IF (advectionScheme.EQ.ENUM_CENTERED_2ND) THEN
233 CALL GAD_C2_ADV_X(bi,bj,k,uTrans,locABT,af,myThid)
234 ELSEIF ( advectionScheme.EQ.ENUM_UPWIND_1RST
235 & .OR. advectionScheme.EQ.ENUM_DST2 ) THEN
236 CALL GAD_DST2U1_ADV_X( bi,bj,k, advectionScheme, .TRUE.,
237 I dTtracerLev(k), uTrans, uFld, locABT,
238 O af, myThid )
239 ELSEIF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN
240 CALL GAD_FLUXLIMIT_ADV_X( bi,bj,k, .TRUE., dTtracerLev(k),
241 I uTrans, uFld, maskW(1-Olx,1-Oly,k,bi,bj), locABT,
242 O af, myThid )
243 ELSEIF (advectionScheme.EQ.ENUM_UPWIND_3RD ) THEN
244 CALL GAD_U3_ADV_X(bi,bj,k,uTrans,locABT,af,myThid)
245 ELSEIF (advectionScheme.EQ.ENUM_CENTERED_4TH) THEN
246 CALL GAD_C4_ADV_X(bi,bj,k,uTrans,locABT,af,myThid)
247 ELSEIF (advectionScheme.EQ.ENUM_DST3 ) THEN
248 CALL GAD_DST3_ADV_X( bi,bj,k, .TRUE., dTtracerLev(k),
249 I uTrans, uFld, maskW(1-Olx,1-Oly,k,bi,bj), locABT,
250 O af, myThid )
251 ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
252 IF ( inAdMode ) THEN
253 cph This block is to trick the adjoint:
254 cph IF inAdExact=.FALSE., we want to use DST3
255 cph with limiters in forward, but without limiters in reverse.
256 CALL GAD_DST3_ADV_X( bi,bj,k, .TRUE., dTtracerLev(k),
257 I uTrans, uFld, maskW(1-Olx,1-Oly,k,bi,bj), locABT,
258 O af, myThid )
259 ELSE
260 CALL GAD_DST3FL_ADV_X( bi,bj,k, .TRUE., dTtracerLev(k),
261 I uTrans, uFld, maskW(1-Olx,1-Oly,k,bi,bj), locABT,
262 O af, myThid )
263 ENDIF
264 #ifndef ALLOW_AUTODIFF_TAMC
265 ELSEIF (advectionScheme.EQ.ENUM_OS7MP ) THEN
266 CALL GAD_OS7MP_ADV_X( bi,bj,k, .TRUE., dTtracerLev(k),
267 I uTrans, uFld, maskW(1-Olx,1-Oly,k,bi,bj), locABT,
268 O af, myThid )
269 #endif
270 ELSE
271 STOP 'GAD_CALC_RHS: Bad advectionScheme (X)'
272 ENDIF
273 DO j=1-Oly,sNy+Oly
274 DO i=1-Olx,sNx+Olx
275 fZon(i,j) = fZon(i,j) + af(i,j)
276 ENDDO
277 ENDDO
278 #ifdef ALLOW_DIAGNOSTICS
279 IF ( useDiagnostics ) THEN
280 diagName = 'ADVx'//diagSufx
281 CALL DIAGNOSTICS_FILL(af,diagName, k,1, 2,bi,bj, myThid)
282 ENDIF
283 #endif
284 ENDIF
285
286 C- Diffusive flux in X
287 IF (diffKh.NE.0.) THEN
288 CALL GAD_DIFF_X(bi,bj,k,xA,diffKh,localT,df,myThid)
289 ELSE
290 DO j=1-Oly,sNy+Oly
291 DO i=1-Olx,sNx+Olx
292 df(i,j) = 0. _d 0
293 ENDDO
294 ENDDO
295 ENDIF
296
297 C- Add bi-harmonic diffusive flux in X
298 IF (diffK4 .NE. 0.) THEN
299 CALL GAD_BIHARM_X(bi,bj,k,xA,df4,diffK4,df,myThid)
300 ENDIF
301
302 #ifdef ALLOW_GMREDI
303 C- GM/Redi flux in X
304 IF ( trUseGMRedi ) THEN
305 C *note* should update GMREDI_XTRANSPORT to set df *aja*
306 IF ( applyAB_onTracer ) THEN
307 CALL GMREDI_XTRANSPORT(
308 I iMin,iMax,jMin,jMax,bi,bj,k,
309 I xA,TracerN,tracerIdentity,
310 U df,
311 I myThid)
312 ELSE
313 CALL GMREDI_XTRANSPORT(
314 I iMin,iMax,jMin,jMax,bi,bj,k,
315 I xA,TracAB, tracerIdentity,
316 U df,
317 I myThid)
318 ENDIF
319 ENDIF
320 #endif
321 C anelastic: advect.fluxes are scaled by rhoFac but hor.diff. flx are not
322 DO j=1-Oly,sNy+Oly
323 DO i=1-Olx,sNx+Olx
324 fZon(i,j) = fZon(i,j) + df(i,j)*rhoFacC(k)
325 ENDDO
326 ENDDO
327
328 #ifdef ALLOW_DIAGNOSTICS
329 C- Diagnostics of Tracer flux in X dir (mainly Diffusive term),
330 C excluding advective terms:
331 IF ( useDiagnostics .AND.
332 & (diffKh.NE.0. .OR. diffK4 .NE.0. .OR. trUseGMRedi) ) THEN
333 diagName = 'DFxE'//diagSufx
334 CALL DIAGNOSTICS_FILL(df,diagName, k,1, 2,bi,bj, myThid)
335 ENDIF
336 #endif
337
338 C-- Initialize net flux in Y direction
339 DO j=1-Oly,sNy+Oly
340 DO i=1-Olx,sNx+Olx
341 fMer(i,j) = 0. _d 0
342 ENDDO
343 ENDDO
344
345 C- Advective flux in Y
346 IF (calcAdvection) THEN
347 IF (advectionScheme.EQ.ENUM_CENTERED_2ND) THEN
348 CALL GAD_C2_ADV_Y(bi,bj,k,vTrans,locABT,af,myThid)
349 ELSEIF ( advectionScheme.EQ.ENUM_UPWIND_1RST
350 & .OR. advectionScheme.EQ.ENUM_DST2 ) THEN
351 CALL GAD_DST2U1_ADV_Y( bi,bj,k, advectionScheme, .TRUE.,
352 I dTtracerLev(k), vTrans, vFld, locABT,
353 O af, myThid )
354 ELSEIF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN
355 CALL GAD_FLUXLIMIT_ADV_Y( bi,bj,k, .TRUE., dTtracerLev(k),
356 I vTrans, vFld, maskS(1-Olx,1-Oly,k,bi,bj), locABT,
357 O af, myThid )
358 ELSEIF (advectionScheme.EQ.ENUM_UPWIND_3RD ) THEN
359 CALL GAD_U3_ADV_Y(bi,bj,k,vTrans,locABT,af,myThid)
360 ELSEIF (advectionScheme.EQ.ENUM_CENTERED_4TH) THEN
361 CALL GAD_C4_ADV_Y(bi,bj,k,vTrans,locABT,af,myThid)
362 ELSEIF (advectionScheme.EQ.ENUM_DST3 ) THEN
363 CALL GAD_DST3_ADV_Y( bi,bj,k, .TRUE., dTtracerLev(k),
364 I vTrans, vFld, maskS(1-Olx,1-Oly,k,bi,bj), locABT,
365 O af, myThid )
366 ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
367 IF ( inAdMode ) THEN
368 cph This block is to trick the adjoint:
369 cph IF inAdExact=.FALSE., we want to use DST3
370 cph with limiters in forward, but without limiters in reverse.
371 CALL GAD_DST3_ADV_Y( bi,bj,k, .TRUE., dTtracerLev(k),
372 I vTrans, vFld, maskS(1-Olx,1-Oly,k,bi,bj), locABT,
373 O af, myThid )
374 ELSE
375 CALL GAD_DST3FL_ADV_Y( bi,bj,k, .TRUE., dTtracerLev(k),
376 I vTrans, vFld, maskS(1-Olx,1-Oly,k,bi,bj), locABT,
377 O af, myThid )
378 ENDIF
379 #ifndef ALLOW_AUTODIFF_TAMC
380 ELSEIF (advectionScheme.EQ.ENUM_OS7MP ) THEN
381 CALL GAD_OS7MP_ADV_Y( bi,bj,k, .TRUE., dTtracerLev(k),
382 I vTrans, vFld, maskS(1-Olx,1-Oly,k,bi,bj), locABT,
383 O af, myThid )
384 #endif
385 ELSE
386 STOP 'GAD_CALC_RHS: Bad advectionScheme (Y)'
387 ENDIF
388 DO j=1-Oly,sNy+Oly
389 DO i=1-Olx,sNx+Olx
390 fMer(i,j) = fMer(i,j) + af(i,j)
391 ENDDO
392 ENDDO
393 #ifdef ALLOW_DIAGNOSTICS
394 IF ( useDiagnostics ) THEN
395 diagName = 'ADVy'//diagSufx
396 CALL DIAGNOSTICS_FILL(af,diagName, k,1, 2,bi,bj, myThid)
397 ENDIF
398 #endif
399 ENDIF
400
401 C- Diffusive flux in Y
402 IF (diffKh.NE.0.) THEN
403 CALL GAD_DIFF_Y(bi,bj,k,yA,diffKh,localT,df,myThid)
404 ELSE
405 DO j=1-Oly,sNy+Oly
406 DO i=1-Olx,sNx+Olx
407 df(i,j) = 0. _d 0
408 ENDDO
409 ENDDO
410 ENDIF
411
412 C- Add bi-harmonic flux in Y
413 IF (diffK4 .NE. 0.) THEN
414 CALL GAD_BIHARM_Y(bi,bj,k,yA,df4,diffK4,df,myThid)
415 ENDIF
416
417 #ifdef ALLOW_GMREDI
418 C- GM/Redi flux in Y
419 IF ( trUseGMRedi ) THEN
420 C *note* should update GMREDI_YTRANSPORT to set df *aja*
421 IF ( applyAB_onTracer ) THEN
422 CALL GMREDI_YTRANSPORT(
423 I iMin,iMax,jMin,jMax,bi,bj,k,
424 I yA,TracerN,tracerIdentity,
425 U df,
426 I myThid)
427 ELSE
428 CALL GMREDI_YTRANSPORT(
429 I iMin,iMax,jMin,jMax,bi,bj,k,
430 I yA,TracAB, tracerIdentity,
431 U df,
432 I myThid)
433 ENDIF
434 ENDIF
435 #endif
436 C anelastic: advect.fluxes are scaled by rhoFac but hor.diff. flx are not
437 DO j=1-Oly,sNy+Oly
438 DO i=1-Olx,sNx+Olx
439 fMer(i,j) = fMer(i,j) + df(i,j)*rhoFacC(k)
440 ENDDO
441 ENDDO
442
443 #ifdef ALLOW_DIAGNOSTICS
444 C- Diagnostics of Tracer flux in Y dir (mainly Diffusive terms),
445 C excluding advective terms:
446 IF ( useDiagnostics .AND.
447 & (diffKh.NE.0. .OR. diffK4 .NE.0. .OR. trUseGMRedi) ) THEN
448 diagName = 'DFyE'//diagSufx
449 CALL DIAGNOSTICS_FILL(df,diagName, k,1, 2,bi,bj, myThid)
450 ENDIF
451 #endif
452
453 C-- Compute vertical flux fVerT(kUp) at interface k (between k-1 & k):
454 C- Advective flux in R
455 #ifdef ALLOW_AIM
456 C- a hack to prevent Water-Vapor vert.transport into the stratospheric level Nr
457 IF (calcAdvection .AND. .NOT.implicitAdvection .AND. k.GE.2 .AND.
458 & (.NOT.useAIM .OR.tracerIdentity.NE.GAD_SALINITY .OR.k.LT.Nr)
459 & ) THEN
460 #else
461 IF (calcAdvection .AND. .NOT.implicitAdvection .AND. k.GE.2) THEN
462 #endif
463 C- Compute vertical advective flux in the interior:
464 IF (vertAdvecScheme.EQ.ENUM_CENTERED_2ND) THEN
465 CALL GAD_C2_ADV_R(bi,bj,k,rTrans,TracAB,af,myThid)
466 ELSEIF ( vertAdvecScheme.EQ.ENUM_UPWIND_1RST
467 & .OR. vertAdvecScheme.EQ.ENUM_DST2 ) THEN
468 CALL GAD_DST2U1_ADV_R( bi,bj,k, vertAdvecScheme,
469 I dTtracerLev(k),rTrans,wFld,TracAB(1-Olx,1-Oly,1,bi,bj),
470 O af, myThid )
471 ELSEIF (vertAdvecScheme.EQ.ENUM_FLUX_LIMIT) THEN
472 CALL GAD_FLUXLIMIT_ADV_R( bi,bj,k,
473 I dTtracerLev(k),rTrans,wFld,TracAB(1-Olx,1-Oly,1,bi,bj),
474 O af, myThid )
475 ELSEIF (vertAdvecScheme.EQ.ENUM_UPWIND_3RD ) THEN
476 CALL GAD_U3_ADV_R(bi,bj,k,rTrans,TracAB,af,myThid)
477 ELSEIF (vertAdvecScheme.EQ.ENUM_CENTERED_4TH) THEN
478 CALL GAD_C4_ADV_R(bi,bj,k,rTrans,TracAB,af,myThid)
479 ELSEIF (vertAdvecScheme.EQ.ENUM_DST3 ) THEN
480 CALL GAD_DST3_ADV_R( bi,bj,k,
481 I dTtracerLev(k),rTrans,wFld,TracAB(1-Olx,1-Oly,1,bi,bj),
482 O af, myThid )
483 ELSEIF (vertAdvecScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
484 cph This block is to trick the adjoint:
485 cph IF inAdExact=.FALSE., we want to use DST3
486 cph with limiters in forward, but without limiters in reverse.
487 IF ( inAdMode ) THEN
488 CALL GAD_DST3_ADV_R( bi,bj,k,
489 I dTtracerLev(k),rTrans,wFld,TracAB(1-Olx,1-Oly,1,bi,bj),
490 O af, myThid )
491 ELSE
492 CALL GAD_DST3FL_ADV_R( bi,bj,k,
493 I dTtracerLev(k),rTrans,wFld,TracAB(1-Olx,1-Oly,1,bi,bj),
494 O af, myThid )
495 ENDIF
496 #ifndef ALLOW_AUTODIFF_TAMC
497 ELSEIF (vertAdvecScheme.EQ.ENUM_OS7MP ) THEN
498 CALL GAD_OS7MP_ADV_R( bi,bj,k,
499 I dTtracerLev(k),rTrans,wFld,TracAB(1-Olx,1-Oly,1,bi,bj),
500 O af, myThid )
501 #endif
502 ELSE
503 STOP 'GAD_CALC_RHS: Bad vertAdvecScheme (R)'
504 ENDIF
505 C- add the advective flux to fVerT
506 DO j=1-Oly,sNy+Oly
507 DO i=1-Olx,sNx+Olx
508 fVerT(i,j,kUp) = fVerT(i,j,kUp) + af(i,j)
509 ENDDO
510 ENDDO
511 #ifdef ALLOW_DIAGNOSTICS
512 IF ( useDiagnostics ) THEN
513 diagName = 'ADVr'//diagSufx
514 CALL DIAGNOSTICS_FILL(af,diagName, k,1, 2,bi,bj, myThid)
515 C- note: needs to explicitly increment the counter since DIAGNOSTICS_FILL
516 C does it only if k=1 (never the case here)
517 IF ( k.EQ.2 ) CALL DIAGNOSTICS_COUNT(diagName,bi,bj,myThid)
518 ENDIF
519 #endif
520 ENDIF
521
522 C- Diffusive flux in R
523 C Note: For K=1 then KM1=1 and this gives a dT/dr = 0 upper
524 C boundary condition.
525 IF (implicitDiffusion) THEN
526 DO j=1-Oly,sNy+Oly
527 DO i=1-Olx,sNx+Olx
528 df(i,j) = 0. _d 0
529 ENDDO
530 ENDDO
531 ELSE
532 IF ( applyAB_onTracer ) THEN
533 CALL GAD_DIFF_R(bi,bj,k,KappaR,TracerN,df,myThid)
534 ELSE
535 CALL GAD_DIFF_R(bi,bj,k,KappaR,TracAB, df,myThid)
536 ENDIF
537 ENDIF
538
539 #ifdef ALLOW_GMREDI
540 C- GM/Redi flux in R
541 IF ( trUseGMRedi ) THEN
542 C *note* should update GMREDI_RTRANSPORT to set df *aja*
543 IF ( applyAB_onTracer ) THEN
544 CALL GMREDI_RTRANSPORT(
545 I iMin,iMax,jMin,jMax,bi,bj,k,
546 I TracerN,tracerIdentity,
547 U df,
548 I myThid)
549 ELSE
550 CALL GMREDI_RTRANSPORT(
551 I iMin,iMax,jMin,jMax,bi,bj,k,
552 I TracAB, tracerIdentity,
553 U df,
554 I myThid)
555 ENDIF
556 ENDIF
557 #endif
558
559 DO j=1-Oly,sNy+Oly
560 DO i=1-Olx,sNx+Olx
561 fVerT(i,j,kUp) = fVerT(i,j,kUp) + df(i,j)*maskUp(i,j)
562 ENDDO
563 ENDDO
564
565 #ifdef ALLOW_DIAGNOSTICS
566 C- Diagnostics of Tracer flux in R dir (mainly Diffusive terms),
567 C Explicit terms only & excluding advective terms:
568 IF ( useDiagnostics .AND.
569 & (.NOT.implicitDiffusion .OR. trUseGMRedi) ) THEN
570 diagName = 'DFrE'//diagSufx
571 CALL DIAGNOSTICS_FILL(df,diagName, k,1, 2,bi,bj, myThid)
572 ENDIF
573 #endif
574
575 #ifdef ALLOW_KPP
576 C- Set non local KPP transport term (ghat):
577 IF ( trUseKPP .AND. k.GE.2 ) THEN
578 DO j=1-Oly,sNy+Oly
579 DO i=1-Olx,sNx+Olx
580 df(i,j) = 0. _d 0
581 ENDDO
582 ENDDO
583 IF (tracerIdentity.EQ.GAD_TEMPERATURE) THEN
584 CALL KPP_TRANSPORT_T(
585 I iMin,iMax,jMin,jMax,bi,bj,k,km1,
586 O df,
587 I myTime, myIter, myThid )
588 ELSEIF (tracerIdentity.EQ.GAD_SALINITY) THEN
589 CALL KPP_TRANSPORT_S(
590 I iMin,iMax,jMin,jMax,bi,bj,k,km1,
591 O df,
592 I myTime, myIter, myThid )
593 #ifdef ALLOW_PTRACERS
594 ELSEIF (tracerIdentity .GE. GAD_TR1) THEN
595 CALL KPP_TRANSPORT_PTR(
596 I iMin,iMax,jMin,jMax,bi,bj,k,km1,
597 I tracerIdentity-GAD_TR1+1,
598 O df,
599 I myTime, myIter, myThid )
600 #endif
601 ELSE
602 PRINT*,'invalid tracer indentity: ', tracerIdentity
603 STOP 'GAD_CALC_RHS: Ooops'
604 ENDIF
605 DO j=1-Oly,sNy+Oly
606 DO i=1-Olx,sNx+Olx
607 fVerT(i,j,kUp) = fVerT(i,j,kUp)
608 & + df(i,j)*maskUp(i,j)*rhoFacF(k)
609 ENDDO
610 ENDDO
611 ENDIF
612 #endif
613
614 C-- Divergence of fluxes
615 C Anelastic: scale vertical fluxes by rhoFac and leave Horizontal fluxes unchanged
616 DO j=1-Oly,sNy+Oly-1
617 DO i=1-Olx,sNx+Olx-1
618 gTracer(i,j,k,bi,bj)=gTracer(i,j,k,bi,bj)
619 & -_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
620 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)*recip_rhoFacC(k)
621 & *( (fZon(i+1,j)-fZon(i,j))
622 & +(fMer(i,j+1)-fMer(i,j))
623 & +(fVerT(i,j,kDown)-fVerT(i,j,kUp))*rkSign
624 & -localT(i,j)*( (uTrans(i+1,j)-uTrans(i,j))
625 & +(vTrans(i,j+1)-vTrans(i,j))
626 & +(rTransKp1(i,j)-rTrans(i,j))*rAdvFac
627 & )*advFac
628 & )
629 ENDDO
630 ENDDO
631
632 #ifdef ALLOW_DEBUG
633 IF ( debugLevel .GE. debLevB
634 & .AND. tracerIdentity.EQ.GAD_TEMPERATURE
635 & .AND. k.EQ.2 .AND. myIter.EQ.1+nIter0
636 & .AND. nPx.EQ.1 .AND. nPy.EQ.1
637 & .AND. useCubedSphereExchange ) THEN
638 CALL DEBUG_CS_CORNER_UV( ' fZon,fMer from GAD_CALC_RHS',
639 & fZon,fMer, k, standardMessageUnit,bi,bj,myThid )
640 ENDIF
641 #endif /* ALLOW_DEBUG */
642
643 RETURN
644 END

  ViewVC Help
Powered by ViewVC 1.1.22