/[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.67 - (show annotations) (download)
Wed Jun 4 22:18:55 2014 UTC (9 years, 11 months ago) by rpa
Branch: MAIN
CVS Tags: checkpoint65, checkpoint64y, checkpoint64z
Changes since 1.66: +13 -5 lines
fixed bugs in layers

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

  ViewVC Help
Powered by ViewVC 1.1.22