/[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.66 - (show annotations) (download)
Wed Jun 4 14:48:32 2014 UTC (10 years ago) by rpa
Branch: MAIN
Changes since 1.65: +13 -1 lines
Adding new capabilities to layers to enable calculation of diapycnal flux.

1 C $Header: /u/gcmpack/MITgcm/pkg/generic_advdiff/gad_calc_rhs.F,v 1.65 2014/05/23 19:59:47 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 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 CALL LAYERS_FILL_DFX( df, trIdentity, k, 1, 2,bi,bj, myThid )
343 #endif /* ALLOW_LAYERS */
344 ENDIF
345 #endif
346
347 C-- Initialize net flux in Y direction
348 DO j=1-OLy,sNy+OLy
349 DO i=1-OLx,sNx+OLx
350 fMer(i,j) = 0. _d 0
351 ENDDO
352 ENDDO
353
354 C- Advective flux in Y
355 IF (calcAdvection) THEN
356 IF (advectionScheme.EQ.ENUM_CENTERED_2ND) THEN
357 CALL GAD_C2_ADV_Y( bi,bj,k, vTrans, locABT, af, myThid )
358 ELSEIF ( advectionScheme.EQ.ENUM_UPWIND_1RST
359 & .OR. advectionScheme.EQ.ENUM_DST2 ) THEN
360 CALL GAD_DST2U1_ADV_Y( bi,bj,k, advectionScheme, .TRUE.,
361 I deltaTLev(k), vTrans, vFld, locABT,
362 O af, myThid )
363 ELSE
364 DO j=1-OLy,sNy+OLy
365 DO i=1-OLx,sNx+OLx
366 #ifdef ALLOW_OBCS
367 maskLocS(i,j) = _maskS(i,j,k,bi,bj)*maskInS(i,j,bi,bj)
368 #else /* ALLOW_OBCS */
369 maskLocS(i,j) = _maskS(i,j,k,bi,bj)
370 #endif /* ALLOW_OBCS */
371 ENDDO
372 ENDDO
373 IF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN
374 CALL GAD_FLUXLIMIT_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k),
375 I vTrans, vFld, maskLocS, locABT,
376 O af, myThid )
377 ELSEIF (advectionScheme.EQ.ENUM_UPWIND_3RD ) THEN
378 CALL GAD_U3_ADV_Y( bi,bj,k, vTrans, maskLocS, locABT,
379 O af, myThid )
380 ELSEIF (advectionScheme.EQ.ENUM_CENTERED_4TH) THEN
381 CALL GAD_C4_ADV_Y( bi,bj,k, vTrans, maskLocS, locABT,
382 O af, myThid )
383 #ifdef ALLOW_AUTODIFF
384 ELSEIF( advectionScheme.EQ.ENUM_DST3 .OR.
385 & (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT .AND. inAdMode)
386 & ) THEN
387 cph This block is to trick the adjoint:
388 cph If inAdExact=.FALSE., we want to use DST3
389 cph with limiters in forward, but without limiters in reverse.
390 #else /* ALLOW_AUTODIFF */
391 ELSEIF( advectionScheme.EQ.ENUM_DST3 ) THEN
392 #endif /* ALLOW_AUTODIFF */
393 CALL GAD_DST3_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k),
394 I vTrans, vFld, maskLocS, locABT,
395 O af, myThid )
396 ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
397 CALL GAD_DST3FL_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k),
398 I vTrans, vFld, maskLocS, locABT,
399 O af, myThid )
400 #ifndef ALLOW_AUTODIFF
401 ELSEIF (advectionScheme.EQ.ENUM_OS7MP ) THEN
402 CALL GAD_OS7MP_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k),
403 I vTrans, vFld, maskLocS, locABT,
404 O af, myThid )
405 #endif
406 ELSE
407 STOP 'GAD_CALC_RHS: Bad advectionScheme (Y)'
408 ENDIF
409 ENDIF
410 #ifdef ALLOW_OBCS
411 IF ( useOBCS ) THEN
412 C- replace advective flux with 1st order upwind scheme estimate
413 CALL OBCS_U1_ADV_TRACER( .FALSE., trIdentity, bi, bj, k,
414 I maskS(1-OLx,1-OLy,k,bi,bj),
415 I vTrans, locABT,
416 U af, myThid )
417 ENDIF
418 #endif /* ALLOW_OBCS */
419 DO j=1-OLy,sNy+OLy
420 DO i=1-OLx,sNx+OLx
421 fMer(i,j) = fMer(i,j) + af(i,j)
422 ENDDO
423 ENDDO
424 #ifdef ALLOW_DIAGNOSTICS
425 IF ( useDiagnostics ) THEN
426 diagName = 'ADVy'//diagSufx
427 CALL DIAGNOSTICS_FILL( af, diagName, k,1, 2,bi,bj, myThid )
428 ENDIF
429 #endif
430 ENDIF
431
432 C- Diffusive flux in Y
433 IF (diffKh.NE.0.) THEN
434 CALL GAD_DIFF_Y(bi,bj,k,yA,diffKh,localT,df,myThid)
435 ELSE
436 DO j=1-OLy,sNy+OLy
437 DO i=1-OLx,sNx+OLx
438 df(i,j) = 0. _d 0
439 ENDDO
440 ENDDO
441 ENDIF
442
443 C- Add bi-harmonic flux in Y
444 IF (diffK4 .NE. 0.) THEN
445 CALL GAD_BIHARM_Y(bi,bj,k,yA,df4,diffK4,df,myThid)
446 ENDIF
447
448 #ifdef ALLOW_GMREDI
449 C- GM/Redi flux in Y
450 IF ( trUseGMRedi ) THEN
451 C *note* should update GMREDI_YTRANSPORT to set df *aja*
452 IF ( applyAB_onTracer ) THEN
453 CALL GMREDI_YTRANSPORT(
454 I iMin,iMax,jMin,jMax,bi,bj,k,
455 I yA,TracerN,trIdentity,
456 U df,
457 I myThid)
458 ELSE
459 CALL GMREDI_YTRANSPORT(
460 I iMin,iMax,jMin,jMax,bi,bj,k,
461 I yA,TracAB, trIdentity,
462 U df,
463 I myThid)
464 ENDIF
465 ENDIF
466 #endif
467 C anelastic: advect.fluxes are scaled by rhoFac but hor.diff. flx are not
468 DO j=1-OLy,sNy+OLy
469 DO i=1-OLx,sNx+OLx
470 fMer(i,j) = fMer(i,j) + df(i,j)*rhoFacC(k)
471 ENDDO
472 ENDDO
473
474 #ifdef ALLOW_DIAGNOSTICS
475 C- Diagnostics of Tracer flux in Y dir (mainly Diffusive terms),
476 C excluding advective terms:
477 IF ( useDiagnostics .AND.
478 & (diffKh.NE.0. .OR. diffK4 .NE.0. .OR. trUseGMRedi) ) THEN
479 diagName = 'DFyE'//diagSufx
480 CALL DIAGNOSTICS_FILL( df, diagName, k,1, 2,bi,bj, myThid )
481 #ifdef ALLOW_LAYERS
482 CALL LAYERS_FILL_DFY( df, trIdentity,k, 1, 2,bi,bj, myThid )
483 #endif /* ALLOW_LAYERS */
484 ENDIF
485 #endif
486
487 C-- Compute vertical flux fVerT(kUp) at interface k (between k-1 & k):
488 C- Advective flux in R
489 #ifdef ALLOW_AIM
490 C- a hack to prevent Water-Vapor vert.transport into the stratospheric level Nr
491 IF (calcAdvection .AND. .NOT.implicitAdvection .AND. k.GE.2 .AND.
492 & (.NOT.useAIM .OR. trIdentity.NE.GAD_SALINITY .OR. k.LT.Nr)
493 & ) THEN
494 #else
495 IF (calcAdvection .AND. .NOT.implicitAdvection .AND. k.GE.2) THEN
496 #endif
497 C- Compute vertical advective flux in the interior:
498 IF (vertAdvecScheme.EQ.ENUM_CENTERED_2ND) THEN
499 CALL GAD_C2_ADV_R( bi,bj,k, rTrans, TracAB, af, myThid )
500 ELSEIF ( vertAdvecScheme.EQ.ENUM_UPWIND_1RST
501 & .OR. vertAdvecScheme.EQ.ENUM_DST2 ) THEN
502 CALL GAD_DST2U1_ADV_R( bi,bj,k,vertAdvecScheme,deltaTLev(k),
503 I rTrans, wFld, TracAB(1-OLx,1-OLy,1,bi,bj),
504 O af, myThid )
505 ELSEIF (vertAdvecScheme.EQ.ENUM_FLUX_LIMIT) THEN
506 CALL GAD_FLUXLIMIT_ADV_R( bi,bj,k, deltaTLev(k),
507 I rTrans, wFld, TracAB(1-OLx,1-OLy,1,bi,bj),
508 O af, myThid )
509 ELSEIF (vertAdvecScheme.EQ.ENUM_UPWIND_3RD ) THEN
510 CALL GAD_U3_ADV_R( bi,bj,k, rTrans, TracAB, af, myThid )
511 ELSEIF (vertAdvecScheme.EQ.ENUM_CENTERED_4TH) THEN
512 CALL GAD_C4_ADV_R( bi,bj,k, rTrans, TracAB, af, myThid )
513 #ifdef ALLOW_AUTODIFF
514 ELSEIF( vertAdvecScheme.EQ.ENUM_DST3 .OR.
515 & (vertAdvecScheme.EQ.ENUM_DST3_FLUX_LIMIT .AND. inAdMode)
516 & ) THEN
517 cph This block is to trick the adjoint:
518 cph If inAdExact=.FALSE., we want to use DST3
519 cph with limiters in forward, but without limiters in reverse.
520 #else /* ALLOW_AUTODIFF */
521 ELSEIF( vertAdvecScheme.EQ.ENUM_DST3 ) THEN
522 #endif /* ALLOW_AUTODIFF */
523 CALL GAD_DST3_ADV_R( bi,bj,k, deltaTLev(k),
524 I rTrans, wFld, TracAB(1-OLx,1-OLy,1,bi,bj),
525 O af, myThid )
526 ELSEIF (vertAdvecScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
527 CALL GAD_DST3FL_ADV_R( bi,bj,k, deltaTLev(k),
528 I rTrans, wFld, TracAB(1-OLx,1-OLy,1,bi,bj),
529 O af, myThid )
530 #ifndef ALLOW_AUTODIFF
531 ELSEIF (vertAdvecScheme.EQ.ENUM_OS7MP ) THEN
532 CALL GAD_OS7MP_ADV_R( bi,bj,k, deltaTLev(k),
533 I rTrans, wFld, TracAB(1-OLx,1-OLy,1,bi,bj),
534 O af, myThid )
535 #endif
536 ELSE
537 STOP 'GAD_CALC_RHS: Bad vertAdvecScheme (R)'
538 ENDIF
539 C- add the advective flux to fVerT
540 DO j=1-OLy,sNy+OLy
541 DO i=1-OLx,sNx+OLx
542 fVerT(i,j,kUp) = fVerT(i,j,kUp) + af(i,j)*maskInC(i,j,bi,bj)
543 ENDDO
544 ENDDO
545 #ifdef ALLOW_DIAGNOSTICS
546 IF ( useDiagnostics ) THEN
547 diagName = 'ADVr'//diagSufx
548 CALL DIAGNOSTICS_FILL( af, diagName, k,1, 2,bi,bj, myThid )
549 C- note: needs to explicitly increment the counter since DIAGNOSTICS_FILL
550 C does it only if k=1 (never the case here)
551 IF ( k.EQ.2 ) CALL DIAGNOSTICS_COUNT(diagName,bi,bj,myThid)
552 ENDIF
553 #endif
554 ENDIF
555
556 C- Diffusive flux in R
557 C Note: For K=1 then KM1=1 and this gives a dT/dr = 0 upper
558 C boundary condition.
559 IF (implicitDiffusion) THEN
560 DO j=1-OLy,sNy+OLy
561 DO i=1-OLx,sNx+OLx
562 df(i,j) = 0. _d 0
563 ENDDO
564 ENDDO
565 ELSE
566 IF ( applyAB_onTracer ) THEN
567 CALL GAD_DIFF_R(bi,bj,k,KappaR,TracerN,df,myThid)
568 ELSE
569 CALL GAD_DIFF_R(bi,bj,k,KappaR,TracAB, df,myThid)
570 ENDIF
571 ENDIF
572
573 IF ( trUseDiffKr4 ) THEN
574 IF ( applyAB_onTracer ) THEN
575 CALL GAD_BIHARM_R( bi,bj,k, diffKr4, TracerN, df, myThid )
576 ELSE
577 CALL GAD_BIHARM_R( bi,bj,k, diffKr4, TracAB, df, myThid )
578 ENDIF
579 ENDIF
580
581 #ifdef ALLOW_GMREDI
582 C- GM/Redi flux in R
583 IF ( trUseGMRedi ) THEN
584 C *note* should update GMREDI_RTRANSPORT to set df *aja*
585 IF ( applyAB_onTracer ) THEN
586 CALL GMREDI_RTRANSPORT(
587 I iMin,iMax,jMin,jMax,bi,bj,k,
588 I TracerN,trIdentity,
589 U df,
590 I myThid)
591 ELSE
592 CALL GMREDI_RTRANSPORT(
593 I iMin,iMax,jMin,jMax,bi,bj,k,
594 I TracAB, trIdentity,
595 U df,
596 I myThid)
597 ENDIF
598 ENDIF
599 #endif
600
601 DO j=1-OLy,sNy+OLy
602 DO i=1-OLx,sNx+OLx
603 fVerT(i,j,kUp) = fVerT(i,j,kUp) + df(i,j)*maskUp(i,j)
604 ENDDO
605 ENDDO
606
607 #ifdef ALLOW_DIAGNOSTICS
608 C- Diagnostics of Tracer flux in R dir (mainly Diffusive terms),
609 C Explicit terms only & excluding advective terms:
610 IF ( useDiagnostics .AND.
611 & (.NOT.implicitDiffusion .OR. trUseGMRedi) ) THEN
612 diagName = 'DFrE'//diagSufx
613 CALL DIAGNOSTICS_FILL( df, diagName, k,1, 2,bi,bj, myThid )
614 #ifdef ALLOW_LAYERS
615 CALL LAYERS_FILL_DFR( df, trIdentity, k, 1, 2,bi,bj, myThid )
616 #endif /* ALLOW_LAYERS */
617 ENDIF
618 #endif
619
620 #ifdef ALLOW_KPP
621 C- Set non local KPP transport term (ghat):
622 IF ( trUseKPP .AND. k.GE.2 ) THEN
623 DO j=1-OLy,sNy+OLy
624 DO i=1-OLx,sNx+OLx
625 df(i,j) = 0. _d 0
626 ENDDO
627 ENDDO
628 IF (trIdentity.EQ.GAD_TEMPERATURE) THEN
629 CALL KPP_TRANSPORT_T(
630 I iMin,iMax,jMin,jMax,bi,bj,k,km1,
631 O df,
632 I myTime, myIter, myThid )
633 ELSEIF (trIdentity.EQ.GAD_SALINITY) THEN
634 CALL KPP_TRANSPORT_S(
635 I iMin,iMax,jMin,jMax,bi,bj,k,km1,
636 O df,
637 I myTime, myIter, myThid )
638 #ifdef ALLOW_PTRACERS
639 ELSEIF (trIdentity .GE. GAD_TR1) THEN
640 CALL KPP_TRANSPORT_PTR(
641 I iMin,iMax,jMin,jMax,bi,bj,k,km1,
642 I trIdentity-GAD_TR1+1,
643 O df,
644 I myTime, myIter, myThid )
645 #endif
646 ELSE
647 WRITE(errorMessageUnit,*)
648 & 'tracer identity =', trIdentity, ' is not valid => STOP'
649 STOP 'ABNORMAL END: S/R GAD_CALC_RHS: invalid tracer identity'
650 ENDIF
651 DO j=1-OLy,sNy+OLy
652 DO i=1-OLx,sNx+OLx
653 fVerT(i,j,kUp) = fVerT(i,j,kUp)
654 & + df(i,j)*maskUp(i,j)*rhoFacF(k)
655 ENDDO
656 ENDDO
657 #ifdef ALLOW_DIAGNOSTICS
658 C- Diagnostics of Non-Local Tracer (vertical) flux
659 IF ( useDiagnostics ) THEN
660 diagName = 'KPPg'//diagSufx
661 CALL DIAGNOSTICS_FILL( df, diagName, k,1, 2,bi,bj, myThid )
662 C- note: needs to explicitly increment the counter since DIAGNOSTICS_FILL
663 C does it only if k=1 (never the case here)
664 IF ( k.EQ.2 ) CALL DIAGNOSTICS_COUNT(diagName,bi,bj,myThid)
665 #ifdef ALLOW_LAYERS
666 CALL LAYERS_FILL_DFR( df, trIdentity, k, 1, 2,bi,bj, myThid )
667 #endif /* ALLOW_LAYERS */
668 ENDIF
669 #endif
670 ENDIF
671 #endif /* ALLOW_KPP */
672
673 #ifdef GAD_SMOLARKIEWICZ_HACK
674 coj Hack to make redi (and everything else in this s/r) positive
675 coj (see Smolarkiewicz MWR 1989 and Bott MWR 1989).
676 coj Only works if 'down' is k+1 and k loop in thermodynamics is k=Nr,1,-1
677 coj
678 coj Apply to all tracers except temperature
679 IF ( trIdentity.NE.GAD_TEMPERATURE .AND.
680 & trIdentity.NE.GAD_SALINITY ) THEN
681 DO j=1-OLy,sNy+OLy-1
682 DO i=1-OLx,sNx+OLx-1
683 coj Add outgoing fluxes
684 outFlux=deltaTLev(k)*
685 & _recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
686 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)*recip_rhoFacC(k)
687 & *( MAX(0. _d 0,fZon(i+1,j)) + MAX(0. _d 0,-fZon(i,j))
688 & +MAX(0. _d 0,fMer(i,j+1)) + MAX(0. _d 0,-fMer(i,j))
689 & +MAX(0. _d 0,fVerT(i,j,kDown)*rkSign)
690 & +MAX(0. _d 0,-fVerT(i,j,kUp)*rkSign)
691 & )
692 trac = localT(i,j)
693 coj If they would reduce tracer by a fraction of more than
694 coj SmolarkiewiczMaxFrac, scale them down
695 IF (outFlux.GT.0. _d 0 .AND.
696 & outFlux.GT.SmolarkiewiczMaxFrac*trac) THEN
697 coj If tracer is already negative, scale flux to zero
698 fac = MAX(0. _d 0,SmolarkiewiczMaxFrac*trac/outFlux)
699
700 IF (fZon(i+1,j).GT.0. _d 0) fZon(i+1,j)=fac*fZon(i+1,j)
701 IF (-fZon(i,j) .GT.0. _d 0) fZon(i,j) =fac*fZon(i,j)
702 IF (fMer(i,j+1).GT.0. _d 0) fMer(i,j+1)=fac*fMer(i,j+1)
703 IF (-fMer(i,j) .GT.0. _d 0) fMer(i,j) =fac*fMer(i,j)
704 IF (-fVerT(i,j,kUp)*rkSign .GT.0. _d 0)
705 & fVerT(i,j,kUp)=fac*fVerT(i,j,kUp)
706
707 IF (k.LT.Nr .AND. fVerT(i,j,kDown)*rkSign.GT.0. _d 0) THEN
708 coj Down flux is special: it has already been applied in lower layer,
709 coj so we have to readjust this.
710 coj Note: for k+1, gTracer is now the updated tracer, not the tendency!
711 coj thus it has an extra factor deltaTLev(k+1)
712 gTrFac=deltaTLev(k+1)
713 coj Other factors that have been applied to gTracer since the last call:
714 #ifdef NONLIN_FRSURF
715 IF (nonlinFreeSurf.GT.0) THEN
716 IF (select_rStar.GT.0) THEN
717 #ifndef DISABLE_RSTAR_CODE
718 gTrFac = gTrFac/rStarExpC(i,j,bi,bj)
719 #endif /* DISABLE_RSTAR_CODE */
720 ENDIF
721 ENDIF
722 #endif /* NONLIN_FRSURF */
723 coj Now: undo down flux, ...
724 gTracer(i,j,k+1,bi,bj)=gTracer(i,j,k+1,bi,bj)
725 & +gTrFac
726 & *_recip_hFacC(i,j,k+1,bi,bj)*recip_drF(k+1)
727 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k+1)
728 & *recip_rhoFacC(k+1)
729 & *( -fVerT(i,j,kDown)*rkSign )
730 coj ... scale ...
731 fVerT(i,j,kDown)=fac*fVerT(i,j,kDown)
732 coj ... and reapply
733 gTracer(i,j,k+1,bi,bj)=gTracer(i,j,k+1,bi,bj)
734 & +gTrFac
735 & *_recip_hFacC(i,j,k+1,bi,bj)*recip_drF(k+1)
736 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k+1)
737 & *recip_rhoFacC(k+1)
738 & *( fVerT(i,j,kDown)*rkSign )
739 ENDIF
740
741 ENDIF
742 ENDDO
743 ENDDO
744 ENDIF
745 #endif
746
747 C-- Divergence of fluxes
748 C Anelastic: scale vertical fluxes by rhoFac and leave Horizontal fluxes unchanged
749 C for Stevens OBC: keep only vertical diffusive contribution on boundaries
750 DO j=1-OLy,sNy+OLy-1
751 DO i=1-OLx,sNx+OLx-1
752 gTracer(i,j,k,bi,bj)=gTracer(i,j,k,bi,bj)
753 & -_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
754 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)*recip_rhoFacC(k)
755 & *( (fZon(i+1,j)-fZon(i,j))*maskInC(i,j,bi,bj)
756 & +(fMer(i,j+1)-fMer(i,j))*maskInC(i,j,bi,bj)
757 & +(fVerT(i,j,kDown)-fVerT(i,j,kUp))*rkSign
758 & -localT(i,j)*( (uTrans(i+1,j)-uTrans(i,j))*advFac
759 & +(vTrans(i,j+1)-vTrans(i,j))*advFac
760 & +(rTransKp1(i,j)-rTrans(i,j))*rAdvFac
761 & )*maskInC(i,j,bi,bj)
762 & )
763 ENDDO
764 ENDDO
765
766 #ifdef ALLOW_DEBUG
767 IF ( debugLevel .GE. debLevC
768 & .AND. trIdentity.EQ.GAD_TEMPERATURE
769 & .AND. k.EQ.2 .AND. myIter.EQ.1+nIter0
770 & .AND. nPx.EQ.1 .AND. nPy.EQ.1
771 & .AND. useCubedSphereExchange ) THEN
772 CALL DEBUG_CS_CORNER_UV( ' fZon,fMer from GAD_CALC_RHS',
773 & fZon,fMer, k, standardMessageUnit,bi,bj,myThid )
774 ENDIF
775 #endif /* ALLOW_DEBUG */
776
777 RETURN
778 END

  ViewVC Help
Powered by ViewVC 1.1.22