/[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.62 - (show annotations) (download)
Wed Jul 4 20:21:45 2012 UTC (11 years, 10 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64, checkpoint63p, checkpoint63q, checkpoint63r, checkpoint63s, checkpoint64q, checkpoint64p, checkpoint64r, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f
Changes since 1.61: +33 -30 lines
move flag "inAdMode" from PARAMS.h to AUTODIFF_PARAMS.h

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

  ViewVC Help
Powered by ViewVC 1.1.22