/[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.60 - (show annotations) (download)
Fri Mar 9 20:13:44 2012 UTC (12 years, 3 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63k
Changes since 1.59: +107 -84 lines
- allow to switch to upwind 1rst order advection scheme for the advective
  flux computaion at the open-boundary.

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

  ViewVC Help
Powered by ViewVC 1.1.22