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

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

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


Revision 1.77 - (show annotations) (download)
Sun Mar 13 01:44:02 2016 UTC (8 years, 1 month ago) by jmc
Branch: MAIN
CVS Tags: checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65v, checkpoint65w, checkpoint65u, checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, HEAD
Changes since 1.76: +106 -5 lines
- from Darren: add PPM and PQM advection schemes (number 40-42 and 50-52)
  with 2 types of limiter (see: Engwirda & Kelley, submit. to JCP);
  Note (from Darren): unlimited PPM/PQM scheme (40 & 50) are just for
  testing and not for actual use.

1 C $Header: /u/gcmpack/MITgcm/pkg/generic_advdiff/gad_advection.F,v 1.76 2015/06/03 13:39:22 rpa Exp $
2 C $Name: $
3
4 #include "GAD_OPTIONS.h"
5 #ifdef ALLOW_AUTODIFF
6 # include "AUTODIFF_OPTIONS.h"
7 #endif
8
9 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
10 CBOP
11 C !ROUTINE: GAD_ADVECTION
12
13 C !INTERFACE: ==========================================================
14 SUBROUTINE GAD_ADVECTION(
15 I implicitAdvection, advectionScheme, vertAdvecScheme,
16 I trIdentity, deltaTLev,
17 I uFld, vFld, wFld, tracer,
18 O gTracer,
19 I bi,bj, myTime,myIter,myThid)
20
21 C !DESCRIPTION:
22 C Calculates the tendency of a tracer due to advection.
23 C It uses the multi-dimensional method given in \ref{sect:multiDimAdvection}
24 C and can only be used for the non-linear advection schemes such as the
25 C direct-space-time method and flux-limiters.
26 C
27 C The algorithm is as follows:
28 C \begin{itemize}
29 C \item{$\theta^{(n+1/3)} = \theta^{(n)}
30 C - \Delta t \partial_x (u\theta^{(n)}) + \theta^{(n)} \partial_x u$}
31 C \item{$\theta^{(n+2/3)} = \theta^{(n+1/3)}
32 C - \Delta t \partial_y (v\theta^{(n+1/3)}) + \theta^{(n)} \partial_y v$}
33 C \item{$\theta^{(n+3/3)} = \theta^{(n+2/3)}
34 C - \Delta t \partial_r (w\theta^{(n+2/3)}) + \theta^{(n)} \partial_r w$}
35 C \item{$G_\theta = ( \theta^{(n+3/3)} - \theta^{(n)} )/\Delta t$}
36 C \end{itemize}
37 C
38 C The tendency (output) is over-written by this routine.
39
40 C !USES: ===============================================================
41 IMPLICIT NONE
42 #include "SIZE.h"
43 #include "EEPARAMS.h"
44 #include "PARAMS.h"
45 #include "GRID.h"
46 #include "GAD.h"
47 #ifdef ALLOW_AUTODIFF
48 # include "tamc.h"
49 # include "tamc_keys.h"
50 # ifdef ALLOW_PTRACERS
51 # include "PTRACERS_SIZE.h"
52 # endif
53 #endif
54 #ifdef ALLOW_EXCH2
55 #include "W2_EXCH2_SIZE.h"
56 #include "W2_EXCH2_TOPOLOGY.h"
57 #endif /* ALLOW_EXCH2 */
58
59 C !INPUT PARAMETERS: ===================================================
60 C implicitAdvection :: implicit vertical advection (later on)
61 C advectionScheme :: advection scheme to use (Horizontal plane)
62 C vertAdvecScheme :: advection scheme to use (vertical direction)
63 C trIdentity :: tracer identifier
64 C uFld :: Advection velocity field, zonal component
65 C vFld :: Advection velocity field, meridional component
66 C wFld :: Advection velocity field, vertical component
67 C tracer :: tracer field
68 C bi,bj :: tile indices
69 C myTime :: current time
70 C myIter :: iteration number
71 C myThid :: thread number
72 LOGICAL implicitAdvection
73 INTEGER advectionScheme, vertAdvecScheme
74 INTEGER trIdentity
75 _RL deltaTLev(Nr)
76 _RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
77 _RL vFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
78 _RL wFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
79 _RL tracer(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
80 INTEGER bi,bj
81 _RL myTime
82 INTEGER myIter
83 INTEGER myThid
84
85 C !OUTPUT PARAMETERS: ==================================================
86 C gTracer :: tendency array
87 _RL gTracer(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
88
89 C !FUNCTIONS: ==========================================================
90 #ifdef ALLOW_DIAGNOSTICS
91 CHARACTER*4 GAD_DIAG_SUFX
92 EXTERNAL GAD_DIAG_SUFX
93 LOGICAL DIAGNOSTICS_IS_ON
94 EXTERNAL DIAGNOSTICS_IS_ON
95 #endif
96
97 C !LOCAL VARIABLES: ====================================================
98 C maskUp :: 2-D array for mask at W points
99 C maskLocW :: 2-D array for mask at West points
100 C maskLocS :: 2-D array for mask at South points
101 C [iMin,iMax]Upd :: loop range to update tracer field
102 C [jMin,jMax]Upd :: loop range to update tracer field
103 C i,j,k :: loop indices
104 C kUp :: index into 2 1/2D array, toggles between 1 and 2
105 C kDown :: index into 2 1/2D array, toggles between 2 and 1
106 C xA,yA :: areas of X and Y face of tracer cells
107 C uTrans,vTrans :: 2-D arrays of volume transports at U,V points
108 C rTrans :: 2-D arrays of volume transports at W points
109 C rTransKp :: vertical volume transport at interface k+1
110 C rTran3d :: 3-D array of volume transport at W points
111 C afr :: 3-D array of vertical advective flux
112 C af :: 2-D array for horizontal advective flux
113 C afx :: 2-D array for horizontal advective flux, x direction
114 C afy :: 2-D array for horizontal advective flux, y direction
115 C fVerT :: 2 1/2D arrays for vertical advective flux
116 C localTij :: 2-D array, temporary local copy of tracer field
117 C localT3d :: 3-D array, temporary local copy of tracer field
118 C kp1Msk :: flag (0,1) for over-riding mask for W levels
119 C calc_fluxes_X :: logical to indicate to calculate fluxes in X dir
120 C calc_fluxes_Y :: logical to indicate to calculate fluxes in Y dir
121 C interiorOnly :: only update the interior of myTile, but not the edges
122 C overlapOnly :: only update the edges of myTile, but not the interior
123 C npass :: number of passes in multi-dimensional method
124 C ipass :: number of the current pass being made
125 C myTile :: variables used to determine which cube face
126 C nCFace :: owns a tile for cube grid runs using
127 C :: multi-dim advection.
128 C [N,S,E,W]_edge :: true if N,S,E,W edge of myTile is an Edge of the cube
129 c _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
130 _RS maskLocW(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
131 _RS maskLocS(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
132 INTEGER iMinUpd,iMaxUpd,jMinUpd,jMaxUpd
133 INTEGER i,j,k,kUp,kDown
134 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
135 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
136 _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
137 _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
138 _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
139 _RL rTransKp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
140 _RL rTran3d (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
141 _RL afr (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
142 _RL af (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
143 _RL afx (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
144 _RL afy (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
145 _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
146 _RL localTij(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
147 _RL localT3d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
148 #ifdef GAD_MULTIDIM_COMPRESSIBLE
149 _RL tmpTrac
150 _RL localVol(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
151 _RL locVol3d(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
152 #endif
153 _RL kp1Msk
154 LOGICAL calc_fluxes_X, calc_fluxes_Y, withSigns
155 LOGICAL interiorOnly, overlapOnly
156 INTEGER npass, ipass
157 INTEGER nCFace
158 LOGICAL N_edge, S_edge, E_edge, W_edge
159 #ifdef ALLOW_AUTODIFF_TAMC
160 C msgBuf :: Informational/error message buffer
161 CHARACTER*(MAX_LEN_MBUF) msgBuf
162 #endif
163 #ifdef ALLOW_EXCH2
164 INTEGER myTile
165 #endif
166 #ifdef ALLOW_DIAGNOSTICS
167 CHARACTER*8 diagName
168 CHARACTER*4 diagSufx
169 LOGICAL doDiagAdvX, doDiagAdvY, doDiagAdvR
170 #endif /* ALLOW_DIAGNOSTICS */
171 CEOP
172
173 #ifdef ALLOW_AUTODIFF_TAMC
174 act0 = trIdentity
175 max0 = maxpass
176 act1 = bi - myBxLo(myThid)
177 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
178 act2 = bj - myByLo(myThid)
179 max2 = myByHi(myThid) - myByLo(myThid) + 1
180 act3 = myThid - 1
181 max3 = nTx*nTy
182 act4 = ikey_dynamics - 1
183 igadkey = act0
184 & + act1*max0
185 & + act2*max0*max1
186 & + act3*max0*max1*max2
187 & + act4*max0*max1*max2*max3
188 IF (trIdentity.GT.maxpass) THEN
189 WRITE(msgBuf,'(A,2I3)')
190 & 'GAD_ADVECTION: maxpass < trIdentity ',
191 & maxpass, trIdentity
192 CALL PRINT_ERROR( msgBuf, myThid )
193 STOP 'ABNORMAL END: S/R GAD_ADVECTION'
194 ENDIF
195 #endif /* ALLOW_AUTODIFF_TAMC */
196
197 #ifdef ALLOW_DIAGNOSTICS
198 C-- Set diagnostics flags and suffix for the current tracer
199 doDiagAdvX = .FALSE.
200 doDiagAdvY = .FALSE.
201 doDiagAdvR = .FALSE.
202 IF ( useDiagnostics ) THEN
203 diagSufx = GAD_DIAG_SUFX( trIdentity, myThid )
204 diagName = 'ADVx'//diagSufx
205 doDiagAdvX = DIAGNOSTICS_IS_ON( diagName, myThid )
206 diagName = 'ADVy'//diagSufx
207 doDiagAdvY = DIAGNOSTICS_IS_ON( diagName, myThid )
208 diagName = 'ADVr'//diagSufx
209 doDiagAdvR = DIAGNOSTICS_IS_ON( diagName, myThid )
210 ENDIF
211 #endif /* ALLOW_DIAGNOSTICS */
212
213 C-- Set up work arrays with valid (i.e. not NaN) values
214 C These inital values do not alter the numerical results. They
215 C just ensure that all memory references are to valid floating
216 C point numbers. This prevents spurious hardware signals due to
217 C uninitialised but inert locations.
218 DO j=1-OLy,sNy+OLy
219 DO i=1-OLx,sNx+OLx
220 C- xA,yA,vFld,uTrans,vTrans are set over the full domain
221 C => no need for extra initialisation
222 c xA(i,j) = 0. _d 0
223 c yA(i,j) = 0. _d 0
224 c uTrans(i,j) = 0. _d 0
225 c vTrans(i,j) = 0. _d 0
226 C- rTransKp is set over the full domain: no need for extra initialisation
227 c rTransKp(i,j)= 0. _d 0
228 C- rTrans and fVerT need to be initialised to zero:
229 rTrans(i,j) = 0. _d 0
230 fVerT(i,j,1) = 0. _d 0
231 fVerT(i,j,2) = 0. _d 0
232 #ifdef ALLOW_AUTODIFF
233 # ifdef GAD_MULTIDIM_COMPRESSIBLE
234 localVol(i,j) = 0. _d 0
235 # endif
236 localTij(i,j) = 0. _d 0
237 #endif /* ALLOW_AUTODIFF */
238 ENDDO
239 ENDDO
240
241 C-- Set tile-specific parameters for horizontal fluxes
242 IF (useCubedSphereExchange) THEN
243 npass = 3
244 #ifdef ALLOW_AUTODIFF_TAMC
245 IF ( npass.GT.maxcube ) STOP 'maxcube needs to be = 3'
246 #endif
247 #ifdef ALLOW_EXCH2
248 myTile = W2_myTileList(bi,bj)
249 nCFace = exch2_myFace(myTile)
250 N_edge = exch2_isNedge(myTile).EQ.1
251 S_edge = exch2_isSedge(myTile).EQ.1
252 E_edge = exch2_isEedge(myTile).EQ.1
253 W_edge = exch2_isWedge(myTile).EQ.1
254 #else
255 nCFace = bi
256 N_edge = .TRUE.
257 S_edge = .TRUE.
258 E_edge = .TRUE.
259 W_edge = .TRUE.
260 #endif
261 ELSE
262 npass = 2
263 nCFace = 0
264 N_edge = .FALSE.
265 S_edge = .FALSE.
266 E_edge = .FALSE.
267 W_edge = .FALSE.
268 ENDIF
269
270 C-- Start of k loop for horizontal fluxes
271 DO k=1,Nr
272 #ifdef ALLOW_AUTODIFF_TAMC
273 kkey = (igadkey-1)*Nr + k
274 CADJ STORE tracer(:,:,k,bi,bj) =
275 CADJ & comlev1_bibj_k_gad, key=kkey, kind=isbyte
276 #endif /* ALLOW_AUTODIFF_TAMC */
277
278 C-- Get temporary terms used by tendency routines
279 DO j=1-OLy,sNy+OLy
280 DO i=1-OLx,sNx+OLx
281 xA(i,j) = _dyG(i,j,bi,bj)*deepFacC(k)
282 & *drF(k)*_hFacW(i,j,k,bi,bj)
283 yA(i,j) = _dxG(i,j,bi,bj)*deepFacC(k)
284 & *drF(k)*_hFacS(i,j,k,bi,bj)
285 ENDDO
286 ENDDO
287 C-- Calculate "volume transports" through tracer cell faces.
288 C anelastic: scaled by rhoFacC (~ mass transport)
289 DO j=1-OLy,sNy+OLy
290 DO i=1-OLx,sNx+OLx
291 uTrans(i,j) = uFld(i,j,k)*xA(i,j)*rhoFacC(k)
292 vTrans(i,j) = vFld(i,j,k)*yA(i,j)*rhoFacC(k)
293 ENDDO
294 ENDDO
295
296 C-- Make local copy of tracer array and mask West & South
297 DO j=1-OLy,sNy+OLy
298 DO i=1-OLx,sNx+OLx
299 localTij(i,j) = tracer(i,j,k,bi,bj)
300 #ifdef GAD_MULTIDIM_COMPRESSIBLE
301 localVol(i,j) = rA(i,j,bi,bj)*deepFac2C(k)
302 & *rhoFacC(k)*drF(k)*hFacC(i,j,k,bi,bj)
303 & + ( oneRS - maskC(i,j,k,bi,bj) )
304 #endif
305 #ifdef ALLOW_OBCS
306 maskLocW(i,j) = _maskW(i,j,k,bi,bj)*maskInW(i,j,bi,bj)
307 maskLocS(i,j) = _maskS(i,j,k,bi,bj)*maskInS(i,j,bi,bj)
308 #else /* ALLOW_OBCS */
309 maskLocW(i,j) = _maskW(i,j,k,bi,bj)
310 maskLocS(i,j) = _maskS(i,j,k,bi,bj)
311 #endif /* ALLOW_OBCS */
312 ENDDO
313 ENDDO
314
315 IF (useCubedSphereExchange) THEN
316 withSigns = .FALSE.
317 CALL FILL_CS_CORNER_UV_RS(
318 & withSigns, maskLocW,maskLocS, bi,bj, myThid )
319 ENDIF
320
321 C-- Multiple passes for different directions on different tiles
322 C-- For cube need one pass for each of red, green and blue axes.
323 DO ipass=1,npass
324 #ifdef ALLOW_AUTODIFF_TAMC
325 passkey = ipass
326 & + (k-1) *maxpass
327 & + (igadkey-1)*maxpass*Nr
328 IF (npass .GT. maxpass) THEN
329 STOP 'GAD_ADVECTION: npass > maxcube. check tamc.h'
330 ENDIF
331 #endif /* ALLOW_AUTODIFF_TAMC */
332
333 interiorOnly = .FALSE.
334 overlapOnly = .FALSE.
335 IF (useCubedSphereExchange) THEN
336 C- CubedSphere : pass 3 times, with partial update of local tracer field
337 IF (ipass.EQ.1) THEN
338 overlapOnly = MOD(nCFace,3).EQ.0
339 interiorOnly = MOD(nCFace,3).NE.0
340 calc_fluxes_X = nCFace.EQ.6 .OR. nCFace.EQ.1 .OR. nCFace.EQ.2
341 calc_fluxes_Y = nCFace.EQ.3 .OR. nCFace.EQ.4 .OR. nCFace.EQ.5
342 ELSEIF (ipass.EQ.2) THEN
343 overlapOnly = MOD(nCFace,3).EQ.2
344 interiorOnly = MOD(nCFace,3).EQ.1
345 calc_fluxes_X = nCFace.EQ.2 .OR. nCFace.EQ.3 .OR. nCFace.EQ.4
346 calc_fluxes_Y = nCFace.EQ.5 .OR. nCFace.EQ.6 .OR. nCFace.EQ.1
347 ELSE
348 interiorOnly = .TRUE.
349 calc_fluxes_X = nCFace.EQ.5 .OR. nCFace.EQ.6
350 calc_fluxes_Y = nCFace.EQ.2 .OR. nCFace.EQ.3
351 ENDIF
352 ELSE
353 C- not CubedSphere
354 calc_fluxes_X = MOD(ipass,2).EQ.1
355 calc_fluxes_Y = .NOT.calc_fluxes_X
356 ENDIF
357
358 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
359 C-- X direction
360
361 #ifdef ALLOW_AUTODIFF
362 C- Always reset advective flux in X
363 DO j=1-OLy,sNy+OLy
364 DO i=1-OLx,sNx+OLx
365 af(i,j) = 0.
366 ENDDO
367 ENDDO
368 # ifndef DISABLE_MULTIDIM_ADVECTION
369 CADJ STORE localTij(:,:) =
370 CADJ & comlev1_bibj_k_gad_pass, key=passkey, kind=isbyte
371 CADJ STORE af(:,:) =
372 CADJ & comlev1_bibj_k_gad_pass, key=passkey, kind=isbyte
373 # endif
374 #endif /* ALLOW_AUTODIFF */
375
376 IF (calc_fluxes_X) THEN
377
378 C- Do not compute fluxes if
379 C a) needed in overlap only
380 C and b) the overlap of myTile are not cube-face Edges
381 IF ( .NOT.overlapOnly .OR. N_edge .OR. S_edge ) THEN
382
383 C- Internal exchange for calculations in X
384 IF ( overlapOnly ) THEN
385 CALL FILL_CS_CORNER_TR_RL( 1, .FALSE.,
386 & localTij, bi,bj, myThid )
387 ENDIF
388
389 C- Advective flux in X
390 #ifndef ALLOW_AUTODIFF
391 DO j=1-OLy,sNy+OLy
392 DO i=1-OLx,sNx+OLx
393 af(i,j) = 0.
394 ENDDO
395 ENDDO
396 #else /* ALLOW_AUTODIFF */
397 # ifndef DISABLE_MULTIDIM_ADVECTION
398 CADJ STORE localTij(:,:) =
399 CADJ & comlev1_bibj_k_gad_pass, key=passkey, kind=isbyte
400 # endif
401 #endif /* ALLOW_AUTODIFF */
402
403 IF ( advectionScheme.EQ.ENUM_UPWIND_1RST
404 & .OR. advectionScheme.EQ.ENUM_DST2 ) THEN
405 CALL GAD_DST2U1_ADV_X( bi,bj,k, advectionScheme, .TRUE.,
406 I deltaTLev(k),uTrans,uFld(1-OLx,1-OLy,k), localTij,
407 O af, myThid )
408 ELSEIF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN
409 CALL GAD_FLUXLIMIT_ADV_X( bi,bj,k, .TRUE., deltaTLev(k),
410 I uTrans, uFld(1-OLx,1-OLy,k), maskLocW, localTij,
411 O af, myThid )
412 ELSEIF (advectionScheme.EQ.ENUM_DST3 ) THEN
413 CALL GAD_DST3_ADV_X( bi,bj,k, .TRUE., deltaTLev(k),
414 I uTrans, uFld(1-OLx,1-OLy,k), maskLocW, localTij,
415 O af, myThid )
416 ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
417 CALL GAD_DST3FL_ADV_X( bi,bj,k, .TRUE., deltaTLev(k),
418 I uTrans, uFld(1-OLx,1-OLy,k), maskLocW, localTij,
419 O af, myThid )
420 #ifndef ALLOW_AUTODIFF
421 ELSEIF (advectionScheme.EQ.ENUM_OS7MP ) THEN
422 CALL GAD_OS7MP_ADV_X( bi,bj,k, .TRUE., deltaTLev(k),
423 I uTrans, uFld(1-OLx,1-OLy,k), maskLocW, localTij,
424 O af, myThid )
425 ELSEIF (advectionScheme.EQ.ENUM_PPM_NULL_LIMIT .OR.
426 & advectionScheme.EQ.ENUM_PPM_MONO_LIMIT .OR.
427 & advectionScheme.EQ.ENUM_PPM_WENO_LIMIT) THEN
428 CALL GAD_PPM_ADV_X( advectionScheme, bi, bj, k , .TRUE.,
429 I deltaTLev(k), uFld(1-OLx,1-OLy,k), uTrans, localTij,
430 O af, myThid )
431 ELSEIF (advectionScheme.EQ.ENUM_PQM_NULL_LIMIT .OR.
432 & advectionScheme.EQ.ENUM_PQM_MONO_LIMIT .OR.
433 & advectionScheme.EQ.ENUM_PQM_WENO_LIMIT) THEN
434 CALL GAD_PQM_ADV_X( advectionScheme, bi, bj, k , .TRUE.,
435 I deltaTLev(k), uFld(1-OLx,1-OLy,k), uTrans, localTij,
436 O af, myThid )
437 #endif /* ndef ALLOW_AUTODIFF */
438 ELSE
439 STOP 'GAD_ADVECTION: adv. scheme incompatibale with multi-dim'
440 ENDIF
441
442 #ifdef ALLOW_OBCS
443 IF ( useOBCS ) THEN
444 C- replace advective flux with 1st order upwind scheme estimate
445 CALL OBCS_U1_ADV_TRACER( .TRUE., trIdentity, bi, bj, k,
446 I maskLocW, uTrans, localTij,
447 U af, myThid )
448 ENDIF
449 #endif /* ALLOW_OBCS */
450
451 C- Internal exchange for next calculations in Y
452 IF ( overlapOnly .AND. ipass.EQ.1 ) THEN
453 CALL FILL_CS_CORNER_TR_RL( 2, .FALSE.,
454 & localTij, bi,bj, myThid )
455 ENDIF
456
457 C- Advective flux in X : done
458 ENDIF
459
460 C- Update the local tracer field where needed:
461 C use "maksInC" to prevent updating tracer field in OB regions
462 #ifdef ALLOW_AUTODIFF_TAMC
463 # ifdef GAD_MULTIDIM_COMPRESSIBLE
464 CADJ STORE localVol(:,:) =
465 CADJ & comlev1_bibj_k_gad_pass, key=passkey, kind=isbyte
466 CADJ STORE localTij(:,:) =
467 CADJ & comlev1_bibj_k_gad_pass, key=passkey, kind=isbyte
468 # endif
469 #endif /* ALLOW_AUTODIFF_TAMC */
470
471 C update in overlap-Only
472 IF ( overlapOnly ) THEN
473 iMinUpd = 1-OLx+1
474 iMaxUpd = sNx+OLx-1
475 C- notes: these 2 lines below have no real effect (because recip_hFac=0
476 C in corner region) but safer to keep them.
477 IF ( W_edge ) iMinUpd = 1
478 IF ( E_edge ) iMaxUpd = sNx
479
480 IF ( S_edge ) THEN
481 DO j=1-OLy,0
482 DO i=iMinUpd,iMaxUpd
483 #ifdef GAD_MULTIDIM_COMPRESSIBLE
484 tmpTrac = localTij(i,j)*localVol(i,j)
485 & -deltaTLev(k)*( af(i+1,j) - af(i,j) )
486 & *maskInC(i,j,bi,bj)
487 localVol(i,j) = localVol(i,j)
488 & -deltaTLev(k)*( uTrans(i+1,j) - uTrans(i,j) )
489 & *maskInC(i,j,bi,bj)
490 localTij(i,j) = tmpTrac/localVol(i,j)
491 #else /* GAD_MULTIDIM_COMPRESSIBLE */
492 localTij(i,j) = localTij(i,j)
493 & -deltaTLev(k)*recip_rhoFacC(k)
494 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
495 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
496 & *( af(i+1,j)-af(i,j)
497 & -tracer(i,j,k,bi,bj)*(uTrans(i+1,j)-uTrans(i,j))
498 & )*maskInC(i,j,bi,bj)
499 #endif /* GAD_MULTIDIM_COMPRESSIBLE */
500 ENDDO
501 ENDDO
502 ENDIF
503 IF ( N_edge ) THEN
504 DO j=sNy+1,sNy+OLy
505 DO i=iMinUpd,iMaxUpd
506 #ifdef GAD_MULTIDIM_COMPRESSIBLE
507 tmpTrac = localTij(i,j)*localVol(i,j)
508 & -deltaTLev(k)*( af(i+1,j) - af(i,j) )
509 & *maskInC(i,j,bi,bj)
510 localVol(i,j) = localVol(i,j)
511 & -deltaTLev(k)*( uTrans(i+1,j) - uTrans(i,j) )
512 & *maskInC(i,j,bi,bj)
513 localTij(i,j) = tmpTrac/localVol(i,j)
514 #else /* GAD_MULTIDIM_COMPRESSIBLE */
515 localTij(i,j) = localTij(i,j)
516 & -deltaTLev(k)*recip_rhoFacC(k)
517 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
518 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
519 & *( af(i+1,j)-af(i,j)
520 & -tracer(i,j,k,bi,bj)*(uTrans(i+1,j)-uTrans(i,j))
521 & )*maskInC(i,j,bi,bj)
522 #endif /* GAD_MULTIDIM_COMPRESSIBLE */
523 ENDDO
524 ENDDO
525 ENDIF
526
527 ELSE
528 C do not only update the overlap
529 jMinUpd = 1-OLy
530 jMaxUpd = sNy+OLy
531 IF ( interiorOnly .AND. S_edge ) jMinUpd = 1
532 IF ( interiorOnly .AND. N_edge ) jMaxUpd = sNy
533 DO j=jMinUpd,jMaxUpd
534 DO i=1-OLx+1,sNx+OLx-1
535 #ifdef GAD_MULTIDIM_COMPRESSIBLE
536 tmpTrac = localTij(i,j)*localVol(i,j)
537 & -deltaTLev(k)*( af(i+1,j) - af(i,j) )
538 & *maskInC(i,j,bi,bj)
539 localVol(i,j) = localVol(i,j)
540 & -deltaTLev(k)*( uTrans(i+1,j) - uTrans(i,j) )
541 & *maskInC(i,j,bi,bj)
542 localTij(i,j) = tmpTrac/localVol(i,j)
543 #else /* GAD_MULTIDIM_COMPRESSIBLE */
544 localTij(i,j) = localTij(i,j)
545 & -deltaTLev(k)*recip_rhoFacC(k)
546 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
547 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
548 & *( af(i+1,j)-af(i,j)
549 & -tracer(i,j,k,bi,bj)*(uTrans(i+1,j)-uTrans(i,j))
550 & )*maskInC(i,j,bi,bj)
551 #endif /* GAD_MULTIDIM_COMPRESSIBLE */
552 ENDDO
553 ENDDO
554 C- keep advective flux (for diagnostics)
555 DO j=1-OLy,sNy+OLy
556 DO i=1-OLx,sNx+OLx
557 afx(i,j) = af(i,j)
558 ENDDO
559 ENDDO
560
561 C- end if/else update overlap-Only
562 ENDIF
563
564 C-- End of X direction
565 ENDIF
566
567 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
568 C-- Y direction
569
570 #ifdef ALLOW_AUTODIFF
571 C- Always reset advective flux in Y
572 DO j=1-OLy,sNy+OLy
573 DO i=1-OLx,sNx+OLx
574 af(i,j) = 0.
575 ENDDO
576 ENDDO
577 # ifndef DISABLE_MULTIDIM_ADVECTION
578 CADJ STORE localTij(:,:) =
579 CADJ & comlev1_bibj_k_gad_pass, key=passkey, kind=isbyte
580 CADJ STORE af(:,:) =
581 CADJ & comlev1_bibj_k_gad_pass, key=passkey, kind=isbyte
582 # endif
583 #endif /* ALLOW_AUTODIFF */
584
585 IF (calc_fluxes_Y) THEN
586
587 C- Do not compute fluxes if
588 C a) needed in overlap only
589 C and b) the overlap of myTile are not cube-face edges
590 IF ( .NOT.overlapOnly .OR. E_edge .OR. W_edge ) THEN
591
592 C- Internal exchange for calculations in Y
593 IF ( overlapOnly ) THEN
594 CALL FILL_CS_CORNER_TR_RL( 2, .FALSE.,
595 & localTij, bi,bj, myThid )
596 ENDIF
597
598 C- Advective flux in Y
599 #ifndef ALLOW_AUTODIFF
600 DO j=1-OLy,sNy+OLy
601 DO i=1-OLx,sNx+OLx
602 af(i,j) = 0.
603 ENDDO
604 ENDDO
605 #else /* ALLOW_AUTODIFF */
606 #ifndef DISABLE_MULTIDIM_ADVECTION
607 CADJ STORE localTij(:,:) =
608 CADJ & comlev1_bibj_k_gad_pass, key=passkey, kind=isbyte
609 #endif
610 #endif /* ALLOW_AUTODIFF */
611
612 IF ( advectionScheme.EQ.ENUM_UPWIND_1RST
613 & .OR. advectionScheme.EQ.ENUM_DST2 ) THEN
614 CALL GAD_DST2U1_ADV_Y( bi,bj,k, advectionScheme, .TRUE.,
615 I deltaTLev(k),vTrans,vFld(1-OLx,1-OLy,k), localTij,
616 O af, myThid )
617 ELSEIF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN
618 CALL GAD_FLUXLIMIT_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k),
619 I vTrans, vFld(1-OLx,1-OLy,k), maskLocS, localTij,
620 O af, myThid )
621 ELSEIF (advectionScheme.EQ.ENUM_DST3 ) THEN
622 CALL GAD_DST3_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k),
623 I vTrans, vFld(1-OLx,1-OLy,k), maskLocS, localTij,
624 O af, myThid )
625 ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
626 CALL GAD_DST3FL_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k),
627 I vTrans, vFld(1-OLx,1-OLy,k), maskLocS, localTij,
628 O af, myThid )
629 #ifndef ALLOW_AUTODIFF
630 ELSEIF (advectionScheme.EQ.ENUM_OS7MP ) THEN
631 CALL GAD_OS7MP_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k),
632 I vTrans, vFld(1-OLx,1-OLy,k), maskLocS, localTij,
633 O af, myThid )
634 ELSEIF (advectionScheme.EQ.ENUM_PPM_NULL_LIMIT .OR.
635 & advectionScheme.EQ.ENUM_PPM_MONO_LIMIT .OR.
636 & advectionScheme.EQ.ENUM_PPM_WENO_LIMIT) THEN
637 CALL GAD_PPM_ADV_Y(advectionScheme, bi, bj, k , .TRUE.,
638 I deltaTLev(k), vFld(1-OLX,1-OLy,k), vTrans, localTij,
639 O af, myThid )
640 ELSEIF (advectionScheme.EQ.ENUM_PQM_NULL_LIMIT .OR.
641 & advectionScheme.EQ.ENUM_PQM_MONO_LIMIT .OR.
642 & advectionScheme.EQ.ENUM_PQM_WENO_LIMIT) THEN
643 CALL GAD_PQM_ADV_Y(advectionScheme, bi, bj, k , .TRUE.,
644 I deltaTLev(k), vFld(1-OLX,1-OLy,k), vTrans, localTij,
645 O af, myThid )
646 #endif /* ndef ALLOW_AUTODIFF */
647 ELSE
648 STOP 'GAD_ADVECTION: adv. scheme incompatibale with mutli-dim'
649 ENDIF
650
651 #ifdef ALLOW_OBCS
652 IF ( useOBCS ) THEN
653 C- replace advective flux with 1st order upwind scheme estimate
654 CALL OBCS_U1_ADV_TRACER( .FALSE., trIdentity, bi, bj, k,
655 I maskLocS, vTrans, localTij,
656 U af, myThid )
657 ENDIF
658 #endif /* ALLOW_OBCS */
659
660 C- Internal exchange for next calculations in X
661 IF ( overlapOnly .AND. ipass.EQ.1 ) THEN
662 CALL FILL_CS_CORNER_TR_RL( 1, .FALSE.,
663 & localTij, bi,bj, myThid )
664 ENDIF
665
666 C- Advective flux in Y : done
667 ENDIF
668
669 C- Update the local tracer field where needed:
670 C use "maksInC" to prevent updating tracer field in OB regions
671 #ifdef ALLOW_AUTODIFF_TAMC
672 # ifdef GAD_MULTIDIM_COMPRESSIBLE
673 CADJ STORE localVol(:,:) =
674 CADJ & comlev1_bibj_k_gad_pass, key=passkey, kind=isbyte
675 CADJ STORE localTij(:,:) =
676 CADJ & comlev1_bibj_k_gad_pass, key=passkey, kind=isbyte
677 # endif
678 #endif /* ALLOW_AUTODIFF_TAMC */
679
680 C update in overlap-Only
681 IF ( overlapOnly ) THEN
682 jMinUpd = 1-OLy+1
683 jMaxUpd = sNy+OLy-1
684 C- notes: these 2 lines below have no real effect (because recip_hFac=0
685 C in corner region) but safer to keep them.
686 IF ( S_edge ) jMinUpd = 1
687 IF ( N_edge ) jMaxUpd = sNy
688
689 IF ( W_edge ) THEN
690 DO j=jMinUpd,jMaxUpd
691 DO i=1-OLx,0
692 #ifdef GAD_MULTIDIM_COMPRESSIBLE
693 tmpTrac = localTij(i,j)*localVol(i,j)
694 & -deltaTLev(k)*( af(i,j+1) - af(i,j) )
695 & *maskInC(i,j,bi,bj)
696 localVol(i,j) = localVol(i,j)
697 & -deltaTLev(k)*( vTrans(i,j+1) - vTrans(i,j) )
698 & *maskInC(i,j,bi,bj)
699 localTij(i,j) = tmpTrac/localVol(i,j)
700 #else /* GAD_MULTIDIM_COMPRESSIBLE */
701 localTij(i,j) = localTij(i,j)
702 & -deltaTLev(k)*recip_rhoFacC(k)
703 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
704 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
705 & *( af(i,j+1)-af(i,j)
706 & -tracer(i,j,k,bi,bj)*(vTrans(i,j+1)-vTrans(i,j))
707 & )*maskInC(i,j,bi,bj)
708 #endif /* GAD_MULTIDIM_COMPRESSIBLE */
709 ENDDO
710 ENDDO
711 ENDIF
712 IF ( E_edge ) THEN
713 DO j=jMinUpd,jMaxUpd
714 DO i=sNx+1,sNx+OLx
715 #ifdef GAD_MULTIDIM_COMPRESSIBLE
716 tmpTrac = localTij(i,j)*localVol(i,j)
717 & -deltaTLev(k)*( af(i,j+1) - af(i,j) )
718 & *maskInC(i,j,bi,bj)
719 localVol(i,j) = localVol(i,j)
720 & -deltaTLev(k)*( vTrans(i,j+1) - vTrans(i,j) )
721 & *maskInC(i,j,bi,bj)
722 localTij(i,j) = tmpTrac/localVol(i,j)
723 #else /* GAD_MULTIDIM_COMPRESSIBLE */
724 localTij(i,j) = localTij(i,j)
725 & -deltaTLev(k)*recip_rhoFacC(k)
726 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
727 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
728 & *( af(i,j+1)-af(i,j)
729 & -tracer(i,j,k,bi,bj)*(vTrans(i,j+1)-vTrans(i,j))
730 & )*maskInC(i,j,bi,bj)
731 #endif /* GAD_MULTIDIM_COMPRESSIBLE */
732 ENDDO
733 ENDDO
734 ENDIF
735
736 ELSE
737 C do not only update the overlap
738 iMinUpd = 1-OLx
739 iMaxUpd = sNx+OLx
740 IF ( interiorOnly .AND. W_edge ) iMinUpd = 1
741 IF ( interiorOnly .AND. E_edge ) iMaxUpd = sNx
742 DO j=1-OLy+1,sNy+OLy-1
743 DO i=iMinUpd,iMaxUpd
744 #ifdef GAD_MULTIDIM_COMPRESSIBLE
745 tmpTrac = localTij(i,j)*localVol(i,j)
746 & -deltaTLev(k)*( af(i,j+1) - af(i,j) )
747 & *maskInC(i,j,bi,bj)
748 localVol(i,j) = localVol(i,j)
749 & -deltaTLev(k)*( vTrans(i,j+1) - vTrans(i,j) )
750 & *maskInC(i,j,bi,bj)
751 localTij(i,j) = tmpTrac/localVol(i,j)
752 #else /* GAD_MULTIDIM_COMPRESSIBLE */
753 localTij(i,j) = localTij(i,j)
754 & -deltaTLev(k)*recip_rhoFacC(k)
755 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
756 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
757 & *( af(i,j+1)-af(i,j)
758 & -tracer(i,j,k,bi,bj)*(vTrans(i,j+1)-vTrans(i,j))
759 & )*maskInC(i,j,bi,bj)
760 #endif /* GAD_MULTIDIM_COMPRESSIBLE */
761 ENDDO
762 ENDDO
763 C- keep advective flux (for diagnostics)
764 DO j=1-OLy,sNy+OLy
765 DO i=1-OLx,sNx+OLx
766 afy(i,j) = af(i,j)
767 ENDDO
768 ENDDO
769
770 C end if/else update overlap-Only
771 ENDIF
772
773 C-- End of Y direction
774 ENDIF
775
776 C-- End of ipass loop
777 ENDDO
778
779 IF ( implicitAdvection ) THEN
780 C- explicit advection is done ; store tendency in gTracer:
781 #ifdef GAD_MULTIDIM_COMPRESSIBLE
782 STOP 'GAD_ADVECTION: missing code for implicitAdvection'
783 #endif /* GAD_MULTIDIM_COMPRESSIBLE */
784 DO j=1-OLy,sNy+OLy
785 DO i=1-OLx,sNx+OLx
786 gTracer(i,j,k) =
787 & ( localTij(i,j) - tracer(i,j,k,bi,bj) )/deltaTLev(k)
788 ENDDO
789 ENDDO
790 ELSE
791 C- horizontal advection done; store intermediate result in 3D array:
792 DO j=1-OLy,sNy+OLy
793 DO i=1-OLx,sNx+OLx
794 #ifdef GAD_MULTIDIM_COMPRESSIBLE
795 locVol3d(i,j,k) = localVol(i,j)
796 #endif /* GAD_MULTIDIM_COMPRESSIBLE */
797 localT3d(i,j,k) = localTij(i,j)
798 ENDDO
799 ENDDO
800 ENDIF
801
802 #ifdef ALLOW_DIAGNOSTICS
803 IF ( doDiagAdvX ) THEN
804 diagName = 'ADVx'//diagSufx
805 CALL DIAGNOSTICS_FILL( afx, diagName, k,1, 2,bi,bj, myThid )
806 ENDIF
807 IF ( doDiagAdvY ) THEN
808 diagName = 'ADVy'//diagSufx
809 CALL DIAGNOSTICS_FILL( afy, diagName, k,1, 2,bi,bj, myThid )
810 ENDIF
811 #ifdef ALLOW_LAYERS
812 IF ( useLayers ) THEN
813 CALL LAYERS_FILL(afx,trIdentity,'AFX',k,1,2,bi,bj,myThid)
814 CALL LAYERS_FILL(afy,trIdentity,'AFY',k,1,2,bi,bj,myThid)
815 ENDIF
816 #endif /* ALLOW_LAYERS */
817 #endif /* ALLOW_DIAGNOSTICS */
818
819 #ifdef ALLOW_DEBUG
820 IF ( debugLevel .GE. debLevC
821 & .AND. trIdentity.EQ.GAD_TEMPERATURE
822 & .AND. k.LE.3 .AND. myIter.EQ.1+nIter0
823 & .AND. nPx.EQ.1 .AND. nPy.EQ.1
824 & .AND. useCubedSphereExchange ) THEN
825 CALL DEBUG_CS_CORNER_UV( ' afx,afy from GAD_ADVECTION',
826 & afx,afy, k, standardMessageUnit,bi,bj,myThid )
827 ENDIF
828 #endif /* ALLOW_DEBUG */
829
830 C-- End of K loop for horizontal fluxes
831 ENDDO
832
833 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
834
835 IF ( .NOT.implicitAdvection ) THEN
836
837 #ifndef ALLOW_AUTODIFF
838 IF (vertAdvecScheme.EQ.ENUM_PPM_NULL_LIMIT .OR.
839 & vertAdvecScheme.EQ.ENUM_PPM_MONO_LIMIT .OR.
840 & vertAdvecScheme.EQ.ENUM_PPM_WENO_LIMIT) THEN
841 C-- PPM-style vertical advection
842 DO k=1,Nr
843 IF (k.EQ.1) THEN
844 C-- vertical transport: surface
845 DO j=1-OLy,sNy+OLy
846 DO i=1-OLx,sNx+OLx
847 rTran3d(i,j,k) = 0. _d 0
848 ENDDO
849 ENDDO
850 ELSE
851 C-- vertical transport: interior
852 DO j=1-OLy,sNy+OLy
853 DO i=1-OLx,sNx+OLx
854 rTran3d(i,j,k) = wFld(i,j,k)*rA(i,j,bi,bj)
855 & *deepFac2F(k)*rhoFacF(k)
856 & *maskC(i,j,k-1,bi,bj)
857 ENDDO
858 ENDDO
859 ENDIF
860 ENDDO
861 C-- calc. PPM vertical flux data
862 CALL GAD_PPM_ADV_R( vertAdvecScheme, bi, bj,
863 I deltaTLev, wFld, rTran3d, localT3d,
864 O afr, myThid )
865 ENDIF
866
867 IF (vertAdvecScheme.EQ.ENUM_PQM_NULL_LIMIT .OR.
868 & vertAdvecScheme.EQ.ENUM_PQM_MONO_LIMIT .OR.
869 & vertAdvecScheme.EQ.ENUM_PQM_WENO_LIMIT) THEN
870 C-- PQM-style vertical advection
871 DO k=1,Nr
872 IF (k.EQ.1) THEN
873 C-- vertical transport: surface
874 DO j=1-OLy,sNy+OLy
875 DO i=1-OLx,sNx+OLx
876 rTran3d(i,j,k) = 0. _d 0
877 ENDDO
878 ENDDO
879 ELSE
880 C-- vertical transport: interior
881 DO j=1-OLy,sNy+OLy
882 DO i=1-OLx,sNx+OLx
883 rTran3d(i,j,k) = wFld(i,j,k)*rA(i,j,bi,bj)
884 & *deepFac2F(k)*rhoFacF(k)
885 & *maskC(i,j,k-1,bi,bj)
886 ENDDO
887 ENDDO
888 ENDIF
889 ENDDO
890 C-- calc. PQM vertical flux data
891 CALL GAD_PQM_ADV_R( vertAdvecScheme, bi, bj,
892 I deltaTLev, wFld, rTran3d, localT3d,
893 O afr, myThid )
894 ENDIF
895 #endif /* ndef ALLOW_AUTODIFF */
896
897 C-- Start of k loop for vertical flux
898 DO k=Nr,1,-1
899 #ifdef ALLOW_AUTODIFF_TAMC
900 kkey = (igadkey-1)*Nr + (Nr-k+1)
901 #endif /* ALLOW_AUTODIFF_TAMC */
902 C-- kUp Cycles through 1,2 to point to w-layer above
903 C-- kDown Cycles through 2,1 to point to w-layer below
904 kUp = 1+MOD(k+1,2)
905 kDown= 1+MOD(k,2)
906 kp1Msk=1.
907 IF (k.EQ.Nr) kp1Msk=0.
908
909 #ifdef ALLOW_AUTODIFF_TAMC
910 CADJ STORE rtrans(:,:) =
911 CADJ & comlev1_bibj_k_gad, key=kkey, kind=isbyte
912 #endif
913
914 C-- Compute Vertical transport
915 #ifdef ALLOW_AIM
916 C- a hack to prevent Water-Vapor vert.transport into the stratospheric level Nr
917 IF ( k.EQ.1 .OR.
918 & (useAIM .AND. trIdentity.EQ.GAD_SALINITY .AND. k.EQ.Nr)
919 & ) THEN
920 #else
921 IF ( k.EQ.1 ) THEN
922 #endif
923
924 C- Surface interface :
925 DO j=1-OLy,sNy+OLy
926 DO i=1-OLx,sNx+OLx
927 rTransKp(i,j) = kp1Msk*rTrans(i,j)
928 rTrans(i,j) = 0.
929 fVerT(i,j,kUp) = 0.
930 ENDDO
931 ENDDO
932
933 ELSE
934
935 C- Interior interface :
936 DO j=1-OLy,sNy+OLy
937 DO i=1-OLx,sNx+OLx
938 rTransKp(i,j) = kp1Msk*rTrans(i,j)
939 rTrans(i,j) = wFld(i,j,k)*rA(i,j,bi,bj)
940 & *deepFac2F(k)*rhoFacF(k)
941 & *maskC(i,j,k-1,bi,bj)
942 fVerT(i,j,kUp) = 0.
943 ENDDO
944 ENDDO
945
946 #ifdef ALLOW_AUTODIFF_TAMC
947 cphmultiCADJ STORE localT3d(:,:,k)
948 cphmultiCADJ & = comlev1_bibj_k_gad, key=kkey, kind=isbyte
949 cphmultiCADJ STORE rTrans(:,:)
950 cphmultiCADJ & = comlev1_bibj_k_gad, key=kkey, kind=isbyte
951 #endif /* ALLOW_AUTODIFF_TAMC */
952
953 C- Compute vertical advective flux in the interior:
954 IF ( vertAdvecScheme.EQ.ENUM_UPWIND_1RST
955 & .OR. vertAdvecScheme.EQ.ENUM_DST2 ) THEN
956 CALL GAD_DST2U1_ADV_R( bi,bj,k, advectionScheme,
957 I deltaTLev(k),rTrans,wFld(1-OLx,1-OLy,k),localT3d,
958 O fVerT(1-OLx,1-OLy,kUp), myThid )
959 ELSEIF( vertAdvecScheme.EQ.ENUM_FLUX_LIMIT) THEN
960 CALL GAD_FLUXLIMIT_ADV_R( bi,bj,k, deltaTLev(k),
961 I rTrans, wFld(1-OLx,1-OLy,k), localT3d,
962 O fVerT(1-OLx,1-OLy,kUp), myThid )
963 ELSEIF( vertAdvecScheme.EQ.ENUM_DST3 ) THEN
964 CALL GAD_DST3_ADV_R( bi,bj,k, deltaTLev(k),
965 I rTrans, wFld(1-OLx,1-OLy,k), localT3d,
966 O fVerT(1-OLx,1-OLy,kUp), myThid )
967 ELSEIF( vertAdvecScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
968 CALL GAD_DST3FL_ADV_R( bi,bj,k, deltaTLev(k),
969 I rTrans, wFld(1-OLx,1-OLy,k), localT3d,
970 O fVerT(1-OLx,1-OLy,kUp), myThid )
971 #ifndef ALLOW_AUTODIFF
972 ELSEIF (vertAdvecScheme.EQ.ENUM_OS7MP ) THEN
973 CALL GAD_OS7MP_ADV_R( bi,bj,k, deltaTLev(k),
974 I rTrans, wFld(1-OLx,1-OLy,k), localT3d,
975 O fVerT(1-OLx,1-OLy,kUp), myThid )
976 ELSEIF (vertAdvecScheme.EQ.ENUM_PPM_NULL_LIMIT .OR.
977 & vertAdvecScheme.EQ.ENUM_PPM_MONO_LIMIT .OR.
978 & vertAdvecScheme.EQ.ENUM_PPM_WENO_LIMIT .OR.
979 & vertAdvecScheme.EQ.ENUM_PQM_NULL_LIMIT .OR.
980 & vertAdvecScheme.EQ.ENUM_PQM_MONO_LIMIT .OR.
981 & vertAdvecScheme.EQ.ENUM_PQM_WENO_LIMIT) THEN
982 C- copy level from 3d flux data
983 DO j = 1-OLy,sNy+OLy
984 DO i = 1-OLx,sNx+OLx
985 fVerT(i,j,kUp) = afr(i,j,k)
986 ENDDO
987 ENDDO
988 #endif /* ndef ALLOW_AUTODIFF */
989 ELSE
990 STOP 'GAD_ADVECTION: adv. scheme incompatibale with mutli-dim'
991 ENDIF
992
993 C- end Surface/Interior if bloc
994 ENDIF
995
996 #ifdef ALLOW_AUTODIFF_TAMC
997 cphmultiCADJ STORE rTrans(:,:)
998 cphmultiCADJ & = comlev1_bibj_k_gad, key=kkey, kind=isbyte
999 cphmultiCADJ STORE rTranskp(:,:)
1000 cphmultiCADJ & = comlev1_bibj_k_gad, key=kkey, kind=isbyte
1001 cph --- following storing of fVerT is critical for correct
1002 cph --- gradient with multiDimAdvection
1003 cph --- Without it, kDown component is not properly recomputed
1004 cph --- This is a TAF bug (and no warning available)
1005 CADJ STORE fVerT(:,:,:)
1006 CADJ & = comlev1_bibj_k_gad, key=kkey, kind=isbyte
1007 #endif /* ALLOW_AUTODIFF_TAMC */
1008
1009 C-- Divergence of vertical fluxes
1010 #ifdef GAD_MULTIDIM_COMPRESSIBLE
1011 DO j=1-OLy,sNy+OLy
1012 DO i=1-OLx,sNx+OLx
1013 tmpTrac = localT3d(i,j,k)*locVol3d(i,j,k)
1014 & -deltaTLev(k)*( fVerT(i,j,kDown)-fVerT(i,j,kUp) )
1015 & *rkSign*maskInC(i,j,bi,bj)
1016 localVol(i,j) = locVol3d(i,j,k)
1017 & -deltaTLev(k)*( rTransKp(i,j) - rTrans(i,j) )
1018 & *rkSign*maskInC(i,j,bi,bj)
1019 C- localTij only needed for Variance Bugget: can be move there
1020 localTij(i,j) = tmpTrac/localVol(i,j)
1021 C-- without rescaling of tendencies:
1022 c gTracer(i,j,k) =
1023 c & ( localTij(i,j) - tracer(i,j,k,bi,bj) )/deltaTLev(k)
1024 C-- Non-Lin Free-Surf: consistent with rescaling of tendencies
1025 C (in FREESURF_RESCALE_G) and RealFreshFlux/addMass.
1026 C Also valid for linear Free-Surf (r & r* coords) w/wout RealFreshFlux
1027 C and consistent with linFSConserveTr and "surfExpan_" monitor.
1028 gTracer(i,j,k) =
1029 & ( tmpTrac - tracer(i,j,k,bi,bj)*localVol(i,j) )
1030 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
1031 & *recip_drF(k)*_recip_hFacC(i,j,k,bi,bj)
1032 & *recip_rhoFacC(k)
1033 & /deltaTLev(k)
1034 ENDDO
1035 ENDDO
1036 #else /* GAD_MULTIDIM_COMPRESSIBLE */
1037 DO j=1-OLy,sNy+OLy
1038 DO i=1-OLx,sNx+OLx
1039 localTij(i,j) = localT3d(i,j,k)
1040 & -deltaTLev(k)*recip_rhoFacC(k)
1041 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
1042 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
1043 & *( fVerT(i,j,kDown)-fVerT(i,j,kUp)
1044 & -tracer(i,j,k,bi,bj)*(rTransKp(i,j)-rTrans(i,j))
1045 & )*rkSign*maskInC(i,j,bi,bj)
1046 gTracer(i,j,k) =
1047 & ( localTij(i,j) - tracer(i,j,k,bi,bj) )/deltaTLev(k)
1048 ENDDO
1049 ENDDO
1050 #endif /* GAD_MULTIDIM_COMPRESSIBLE */
1051
1052 #ifdef ALLOW_DIAGNOSTICS
1053 IF ( doDiagAdvR ) THEN
1054 diagName = 'ADVr'//diagSufx
1055 CALL DIAGNOSTICS_FILL( fVerT(1-OLx,1-OLy,kUp),
1056 & diagName, k,1, 2,bi,bj, myThid )
1057 ENDIF
1058 #ifdef ALLOW_LAYERS
1059 IF ( useLayers ) THEN
1060 CALL LAYERS_FILL( fVerT(1-OLx,1-OLy,kUp), trIdentity,
1061 & 'AFR', k, 1, 2,bi,bj, myThid)
1062 ENDIF
1063 #endif /* ALLOW_LAYERS */
1064 #endif /* ALLOW_DIAGNOSTICS */
1065
1066 C-- End of K loop for vertical flux
1067 ENDDO
1068 C-- end of if not.implicitAdvection block
1069 ENDIF
1070
1071 RETURN
1072 END

  ViewVC Help
Powered by ViewVC 1.1.22