/[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.66 - (show annotations) (download)
Mon Jun 6 15:44:00 2011 UTC (13 years ago) by jmc
Branch: MAIN
Changes since 1.65: +2 -2 lines
refine debugLevel criteria when printing messages

1 C $Header: /u/gcmpack/MITgcm/pkg/generic_advdiff/gad_advection.F,v 1.65 2011/05/20 21:08:56 jmc Exp $
2 C $Name: $
3
4 #include "GAD_OPTIONS.h"
5
6 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
7 CBOP
8 C !ROUTINE: GAD_ADVECTION
9
10 C !INTERFACE: ==========================================================
11 SUBROUTINE GAD_ADVECTION(
12 I implicitAdvection, advectionScheme, vertAdvecScheme,
13 I tracerIdentity, deltaTLev,
14 I uVel, vVel, wVel, tracer,
15 O gTracer,
16 I bi,bj, myTime,myIter,myThid)
17
18 C !DESCRIPTION:
19 C Calculates the tendency of a tracer due to advection.
20 C It uses the multi-dimensional method given in \ref{sect:multiDimAdvection}
21 C and can only be used for the non-linear advection schemes such as the
22 C direct-space-time method and flux-limiters.
23 C
24 C The algorithm is as follows:
25 C \begin{itemize}
26 C \item{$\theta^{(n+1/3)} = \theta^{(n)}
27 C - \Delta t \partial_x (u\theta^{(n)}) + \theta^{(n)} \partial_x u$}
28 C \item{$\theta^{(n+2/3)} = \theta^{(n+1/3)}
29 C - \Delta t \partial_y (v\theta^{(n+1/3)}) + \theta^{(n)} \partial_y v$}
30 C \item{$\theta^{(n+3/3)} = \theta^{(n+2/3)}
31 C - \Delta t \partial_r (w\theta^{(n+2/3)}) + \theta^{(n)} \partial_r w$}
32 C \item{$G_\theta = ( \theta^{(n+3/3)} - \theta^{(n)} )/\Delta t$}
33 C \end{itemize}
34 C
35 C The tendency (output) is over-written by this routine.
36
37 C !USES: ===============================================================
38 IMPLICIT NONE
39 #include "SIZE.h"
40 #include "EEPARAMS.h"
41 #include "PARAMS.h"
42 #include "GRID.h"
43 #include "GAD.h"
44 #ifdef ALLOW_AUTODIFF_TAMC
45 # include "tamc.h"
46 # include "tamc_keys.h"
47 # ifdef ALLOW_PTRACERS
48 # include "PTRACERS_SIZE.h"
49 # endif
50 #endif
51 #ifdef ALLOW_EXCH2
52 #include "W2_EXCH2_SIZE.h"
53 #include "W2_EXCH2_TOPOLOGY.h"
54 #endif /* ALLOW_EXCH2 */
55
56 C !INPUT PARAMETERS: ===================================================
57 C implicitAdvection :: implicit vertical advection (later on)
58 C advectionScheme :: advection scheme to use (Horizontal plane)
59 C vertAdvecScheme :: advection scheme to use (vertical direction)
60 C tracerIdentity :: tracer identifier (required only for OBCS)
61 C uVel :: velocity, zonal component
62 C vVel :: velocity, meridional component
63 C wVel :: velocity, vertical component
64 C tracer :: tracer field
65 C bi,bj :: tile indices
66 C myTime :: current time
67 C myIter :: iteration number
68 C myThid :: thread number
69 LOGICAL implicitAdvection
70 INTEGER advectionScheme, vertAdvecScheme
71 INTEGER tracerIdentity
72 _RL deltaTLev(Nr)
73 _RL uVel (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr,nSx,nSy)
74 _RL vVel (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr,nSx,nSy)
75 _RL wVel (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr,nSx,nSy)
76 _RL tracer(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr,nSx,nSy)
77 INTEGER bi,bj
78 _RL myTime
79 INTEGER myIter
80 INTEGER myThid
81
82 C !OUTPUT PARAMETERS: ==================================================
83 C gTracer :: tendency array
84 _RL gTracer(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nr,nSx,nSy)
85
86 C !LOCAL VARIABLES: ====================================================
87 C maskUp :: 2-D array for mask at W points
88 C maskLocW :: 2-D array for mask at West points
89 C maskLocS :: 2-D array for mask at South points
90 C [iMin,iMax]Upd :: loop range to update tracer field
91 C [jMin,jMax]Upd :: loop range to update tracer field
92 C i,j,k :: loop indices
93 C kUp :: index into 2 1/2D array, toggles between 1 and 2
94 C kDown :: index into 2 1/2D array, toggles between 2 and 1
95 C kp1 :: =k+1 for k<Nr, =Nr for k=Nr
96 C xA,yA :: areas of X and Y face of tracer cells
97 C uFld,vFld :: 2-D local copy of horizontal velocity, U,V components
98 C wFld :: 2-D local copy of vertical velocity
99 C uTrans,vTrans :: 2-D arrays of volume transports at U,V points
100 C rTrans :: 2-D arrays of volume transports at W points
101 C rTransKp1 :: vertical volume transport at interface k+1
102 C af :: 2-D array for horizontal advective flux
103 C afx :: 2-D array for horizontal advective flux, x direction
104 C afy :: 2-D array for horizontal advective flux, y direction
105 C fVerT :: 2 1/2D arrays for vertical advective flux
106 C localTij :: 2-D array, temporary local copy of tracer fld
107 C localTijk :: 3-D array, temporary local copy of tracer fld
108 C kp1Msk :: flag (0,1) for over-riding mask for W levels
109 C calc_fluxes_X :: logical to indicate to calculate fluxes in X dir
110 C calc_fluxes_Y :: logical to indicate to calculate fluxes in Y dir
111 C interiorOnly :: only update the interior of myTile, but not the edges
112 C overlapOnly :: only update the edges of myTile, but not the interior
113 C npass :: number of passes in multi-dimensional method
114 C ipass :: number of the current pass being made
115 C myTile :: variables used to determine which cube face
116 C nCFace :: owns a tile for cube grid runs using
117 C :: multi-dim advection.
118 C [N,S,E,W]_edge :: true if N,S,E,W edge of myTile is an Edge of the cube
119 c _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
120 _RS maskLocW(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
121 _RS maskLocS(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
122 INTEGER iMinUpd,iMaxUpd,jMinUpd,jMaxUpd
123 INTEGER i,j,k,kUp,kDown
124 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
125 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
126 _RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
127 _RL vFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
128 _RL wFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
129 _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
130 _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
131 _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
132 _RL rTransKp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
133 _RL af (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
134 _RL afx (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
135 _RL afy (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
136 _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
137 _RL localTij(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
138 _RL localTijk(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
139 _RL kp1Msk
140 LOGICAL calc_fluxes_X, calc_fluxes_Y, withSigns
141 LOGICAL interiorOnly, overlapOnly
142 INTEGER npass, ipass
143 INTEGER nCFace
144 LOGICAL N_edge, S_edge, E_edge, W_edge
145 #ifdef ALLOW_EXCH2
146 INTEGER myTile
147 #endif
148 #ifdef ALLOW_DIAGNOSTICS
149 CHARACTER*8 diagName
150 CHARACTER*4 diagSufx
151 LOGICAL doDiagAdvX, doDiagAdvY, doDiagAdvR
152 C- Functions:
153 CHARACTER*4 GAD_DIAG_SUFX
154 EXTERNAL GAD_DIAG_SUFX
155 LOGICAL DIAGNOSTICS_IS_ON
156 EXTERNAL DIAGNOSTICS_IS_ON
157 #endif
158 CEOP
159
160 #ifdef ALLOW_AUTODIFF_TAMC
161 act0 = tracerIdentity
162 max0 = maxpass
163 act1 = bi - myBxLo(myThid)
164 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
165 act2 = bj - myByLo(myThid)
166 max2 = myByHi(myThid) - myByLo(myThid) + 1
167 act3 = myThid - 1
168 max3 = nTx*nTy
169 act4 = ikey_dynamics - 1
170 igadkey = act0
171 & + act1*max0
172 & + act2*max0*max1
173 & + act3*max0*max1*max2
174 & + act4*max0*max1*max2*max3
175 IF (tracerIdentity.GT.maxpass) THEN
176 print *, 'ph-pass gad_advection ', maxpass, tracerIdentity
177 STOP 'maxpass seems smaller than tracerIdentity'
178 ENDIF
179 #endif /* ALLOW_AUTODIFF_TAMC */
180
181 #ifdef ALLOW_DIAGNOSTICS
182 C-- Set diagnostics flags and suffix for the current tracer
183 doDiagAdvX = .FALSE.
184 doDiagAdvY = .FALSE.
185 doDiagAdvR = .FALSE.
186 IF ( useDiagnostics ) THEN
187 diagSufx = GAD_DIAG_SUFX( tracerIdentity, myThid )
188 diagName = 'ADVx'//diagSufx
189 doDiagAdvX = DIAGNOSTICS_IS_ON( diagName, myThid )
190 diagName = 'ADVy'//diagSufx
191 doDiagAdvY = DIAGNOSTICS_IS_ON( diagName, myThid )
192 diagName = 'ADVr'//diagSufx
193 doDiagAdvR = DIAGNOSTICS_IS_ON( diagName, myThid )
194 ENDIF
195 #endif
196
197 C-- Set up work arrays with valid (i.e. not NaN) values
198 C These inital values do not alter the numerical results. They
199 C just ensure that all memory references are to valid floating
200 C point numbers. This prevents spurious hardware signals due to
201 C uninitialised but inert locations.
202 DO j=1-OLy,sNy+OLy
203 DO i=1-OLx,sNx+OLx
204 xA(i,j) = 0. _d 0
205 yA(i,j) = 0. _d 0
206 uTrans(i,j) = 0. _d 0
207 vTrans(i,j) = 0. _d 0
208 rTrans(i,j) = 0. _d 0
209 fVerT(i,j,1) = 0. _d 0
210 fVerT(i,j,2) = 0. _d 0
211 rTransKp1(i,j)= 0. _d 0
212 #ifdef ALLOW_AUTODIFF_TAMC
213 localTij(i,j) = 0. _d 0
214 wfld(i,j) = 0. _d 0
215 #endif
216 ENDDO
217 ENDDO
218
219 C-- Set tile-specific parameters for horizontal fluxes
220 IF (useCubedSphereExchange) THEN
221 npass = 3
222 #ifdef ALLOW_AUTODIFF_TAMC
223 IF ( npass.GT.maxcube ) STOP 'maxcube needs to be = 3'
224 #endif
225 #ifdef ALLOW_EXCH2
226 myTile = W2_myTileList(bi,bj)
227 nCFace = exch2_myFace(myTile)
228 N_edge = exch2_isNedge(myTile).EQ.1
229 S_edge = exch2_isSedge(myTile).EQ.1
230 E_edge = exch2_isEedge(myTile).EQ.1
231 W_edge = exch2_isWedge(myTile).EQ.1
232 #else
233 nCFace = bi
234 N_edge = .TRUE.
235 S_edge = .TRUE.
236 E_edge = .TRUE.
237 W_edge = .TRUE.
238 #endif
239 ELSE
240 npass = 2
241 nCFace = 0
242 N_edge = .FALSE.
243 S_edge = .FALSE.
244 E_edge = .FALSE.
245 W_edge = .FALSE.
246 ENDIF
247
248 C-- Start of k loop for horizontal fluxes
249 DO k=1,Nr
250 #ifdef ALLOW_AUTODIFF_TAMC
251 kkey = (igadkey-1)*Nr + k
252 CADJ STORE tracer(:,:,k,bi,bj) =
253 CADJ & comlev1_bibj_k_gad, key=kkey, kind=isbyte
254 #endif /* ALLOW_AUTODIFF_TAMC */
255
256 C-- Get temporary terms used by tendency routines
257 CALL CALC_COMMON_FACTORS (
258 I uVel, vVel,
259 O uFld, vFld, uTrans, vTrans, xA, yA,
260 I k,bi,bj, myThid )
261
262 #ifdef ALLOW_GMREDI
263 C-- Residual transp = Bolus transp + Eulerian transp
264 IF (useGMRedi)
265 & CALL GMREDI_CALC_UVFLOW(
266 U uFld, vFld, uTrans, vTrans,
267 I k, bi, bj, myThid )
268 #endif /* ALLOW_GMREDI */
269
270 C-- Make local copy of tracer array and mask West & South
271 DO j=1-OLy,sNy+OLy
272 DO i=1-OLx,sNx+OLx
273 localTij(i,j)=tracer(i,j,k,bi,bj)
274 #ifdef ALLOW_OBCS
275 maskLocW(i,j) = _maskW(i,j,k,bi,bj)*maskInW(i,j,bi,bj)
276 maskLocS(i,j) = _maskS(i,j,k,bi,bj)*maskInS(i,j,bi,bj)
277 #else /* ALLOW_OBCS */
278 maskLocW(i,j) = _maskW(i,j,k,bi,bj)
279 maskLocS(i,j) = _maskS(i,j,k,bi,bj)
280 #endif /* ALLOW_OBCS */
281 ENDDO
282 ENDDO
283
284 cph-exch2#ifndef ALLOW_AUTODIFF_TAMC
285 IF (useCubedSphereExchange) THEN
286 withSigns = .FALSE.
287 CALL FILL_CS_CORNER_UV_RS(
288 & withSigns, maskLocW,maskLocS, bi,bj, myThid )
289 ENDIF
290 cph-exch2#endif
291
292 C-- Multiple passes for different directions on different tiles
293 C-- For cube need one pass for each of red, green and blue axes.
294 DO ipass=1,npass
295 #ifdef ALLOW_AUTODIFF_TAMC
296 passkey = ipass
297 & + (k-1) *maxpass
298 & + (igadkey-1)*maxpass*Nr
299 IF (npass .GT. maxpass) THEN
300 STOP 'GAD_ADVECTION: npass > maxcube. check tamc.h'
301 ENDIF
302 #endif /* ALLOW_AUTODIFF_TAMC */
303
304 interiorOnly = .FALSE.
305 overlapOnly = .FALSE.
306 IF (useCubedSphereExchange) THEN
307 C- CubedSphere : pass 3 times, with partial update of local tracer field
308 IF (ipass.EQ.1) THEN
309 overlapOnly = MOD(nCFace,3).EQ.0
310 interiorOnly = MOD(nCFace,3).NE.0
311 calc_fluxes_X = nCFace.EQ.6 .OR. nCFace.EQ.1 .OR. nCFace.EQ.2
312 calc_fluxes_Y = nCFace.EQ.3 .OR. nCFace.EQ.4 .OR. nCFace.EQ.5
313 ELSEIF (ipass.EQ.2) THEN
314 overlapOnly = MOD(nCFace,3).EQ.2
315 interiorOnly = MOD(nCFace,3).EQ.1
316 calc_fluxes_X = nCFace.EQ.2 .OR. nCFace.EQ.3 .OR. nCFace.EQ.4
317 calc_fluxes_Y = nCFace.EQ.5 .OR. nCFace.EQ.6 .OR. nCFace.EQ.1
318 ELSE
319 interiorOnly = .TRUE.
320 calc_fluxes_X = nCFace.EQ.5 .OR. nCFace.EQ.6
321 calc_fluxes_Y = nCFace.EQ.2 .OR. nCFace.EQ.3
322 ENDIF
323 ELSE
324 C- not CubedSphere
325 calc_fluxes_X = MOD(ipass,2).EQ.1
326 calc_fluxes_Y = .NOT.calc_fluxes_X
327 ENDIF
328
329 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
330 C-- X direction
331 C- Advective flux in X
332 DO j=1-Oly,sNy+Oly
333 DO i=1-Olx,sNx+Olx
334 af(i,j) = 0.
335 ENDDO
336 ENDDO
337 C
338 #ifdef ALLOW_AUTODIFF_TAMC
339 # ifndef DISABLE_MULTIDIM_ADVECTION
340 CADJ STORE localTij(:,:) =
341 CADJ & comlev1_bibj_k_gad_pass, key=passkey, kind=isbyte
342 CADJ STORE af(:,:) =
343 CADJ & comlev1_bibj_k_gad_pass, key=passkey, kind=isbyte
344 # endif
345 #endif /* ALLOW_AUTODIFF_TAMC */
346 C
347 IF (calc_fluxes_X) THEN
348
349 C- Do not compute fluxes if
350 C a) needed in overlap only
351 C and b) the overlap of myTile are not cube-face Edges
352 IF ( .NOT.overlapOnly .OR. N_edge .OR. S_edge ) THEN
353
354 C- Internal exchange for calculations in X
355 IF ( overlapOnly ) THEN
356 CALL FILL_CS_CORNER_TR_RL( 1, .FALSE.,
357 & localTij, bi,bj, myThid )
358 ENDIF
359
360 #ifdef ALLOW_AUTODIFF_TAMC
361 # ifndef DISABLE_MULTIDIM_ADVECTION
362 CADJ STORE localTij(:,:) =
363 CADJ & comlev1_bibj_k_gad_pass, key=passkey, kind=isbyte
364 # endif
365 #endif /* ALLOW_AUTODIFF_TAMC */
366
367 IF ( advectionScheme.EQ.ENUM_UPWIND_1RST
368 & .OR. advectionScheme.EQ.ENUM_DST2 ) THEN
369 CALL GAD_DST2U1_ADV_X( bi,bj,k, advectionScheme, .TRUE.,
370 I deltaTLev(k),uTrans,uFld,localTij,
371 O af, myThid )
372 ELSEIF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN
373 CALL GAD_FLUXLIMIT_ADV_X( bi,bj,k, .TRUE., deltaTLev(k),
374 I uTrans, uFld, maskLocW, localTij,
375 O af, myThid )
376 ELSEIF (advectionScheme.EQ.ENUM_DST3 ) THEN
377 CALL GAD_DST3_ADV_X( bi,bj,k, .TRUE., deltaTLev(k),
378 I uTrans, uFld, maskLocW, localTij,
379 O af, myThid )
380 ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
381 CALL GAD_DST3FL_ADV_X( bi,bj,k, .TRUE., deltaTLev(k),
382 I uTrans, uFld, maskLocW, localTij,
383 O af, myThid )
384 #ifndef ALLOW_AUTODIFF_TAMC
385 ELSEIF (advectionScheme.EQ.ENUM_OS7MP ) THEN
386 CALL GAD_OS7MP_ADV_X( bi,bj,k, .TRUE., deltaTLev(k),
387 I uTrans, uFld, maskLocW, localTij,
388 O af, myThid )
389 #endif
390 ELSE
391 STOP 'GAD_ADVECTION: adv. scheme incompatibale with multi-dim'
392 ENDIF
393
394 C- Internal exchange for next calculations in Y
395 IF ( overlapOnly .AND. ipass.EQ.1 ) THEN
396 CALL FILL_CS_CORNER_TR_RL( 2, .FALSE.,
397 & localTij, bi,bj, myThid )
398 ENDIF
399
400 C- Advective flux in X : done
401 ENDIF
402
403 C- Update the local tracer field where needed:
404 C use "maksInC" to prevent updating tracer field in OB regions
405
406 C update in overlap-Only
407 IF ( overlapOnly ) THEN
408 iMinUpd = 1-Olx+1
409 iMaxUpd = sNx+Olx-1
410 C- notes: these 2 lines below have no real effect (because recip_hFac=0
411 C in corner region) but safer to keep them.
412 IF ( W_edge ) iMinUpd = 1
413 IF ( E_edge ) iMaxUpd = sNx
414
415 IF ( S_edge ) THEN
416 DO j=1-Oly,0
417 DO i=iMinUpd,iMaxUpd
418 localTij(i,j) = localTij(i,j)
419 & -deltaTLev(k)*recip_rhoFacC(k)
420 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
421 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
422 & *( af(i+1,j)-af(i,j)
423 & -tracer(i,j,k,bi,bj)*(uTrans(i+1,j)-uTrans(i,j))
424 & )*maskInC(i,j,bi,bj)
425 ENDDO
426 ENDDO
427 ENDIF
428 IF ( N_edge ) THEN
429 DO j=sNy+1,sNy+Oly
430 DO i=iMinUpd,iMaxUpd
431 localTij(i,j) = localTij(i,j)
432 & -deltaTLev(k)*recip_rhoFacC(k)
433 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
434 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
435 & *( af(i+1,j)-af(i,j)
436 & -tracer(i,j,k,bi,bj)*(uTrans(i+1,j)-uTrans(i,j))
437 & )*maskInC(i,j,bi,bj)
438 ENDDO
439 ENDDO
440 ENDIF
441
442 ELSE
443 C do not only update the overlap
444 jMinUpd = 1-Oly
445 jMaxUpd = sNy+Oly
446 IF ( interiorOnly .AND. S_edge ) jMinUpd = 1
447 IF ( interiorOnly .AND. N_edge ) jMaxUpd = sNy
448 DO j=jMinUpd,jMaxUpd
449 DO i=1-Olx+1,sNx+Olx-1
450 localTij(i,j) = localTij(i,j)
451 & -deltaTLev(k)*recip_rhoFacC(k)
452 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
453 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
454 & *( af(i+1,j)-af(i,j)
455 & -tracer(i,j,k,bi,bj)*(uTrans(i+1,j)-uTrans(i,j))
456 & )*maskInC(i,j,bi,bj)
457 ENDDO
458 ENDDO
459 C- keep advective flux (for diagnostics)
460 DO j=1-Oly,sNy+Oly
461 DO i=1-Olx,sNx+Olx
462 afx(i,j) = af(i,j)
463 ENDDO
464 ENDDO
465
466 C- end if/else update overlap-Only
467 ENDIF
468
469 C-- End of X direction
470 ENDIF
471
472 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
473 C-- Y direction
474 cph-test
475 C- Advective flux in Y
476 DO j=1-Oly,sNy+Oly
477 DO i=1-Olx,sNx+Olx
478 af(i,j) = 0.
479 ENDDO
480 ENDDO
481 C
482 #ifdef ALLOW_AUTODIFF_TAMC
483 # ifndef DISABLE_MULTIDIM_ADVECTION
484 CADJ STORE localTij(:,:) =
485 CADJ & comlev1_bibj_k_gad_pass, key=passkey, kind=isbyte
486 CADJ STORE af(:,:) =
487 CADJ & comlev1_bibj_k_gad_pass, key=passkey, kind=isbyte
488 # endif
489 #endif /* ALLOW_AUTODIFF_TAMC */
490 C
491 IF (calc_fluxes_Y) THEN
492
493 C- Do not compute fluxes if
494 C a) needed in overlap only
495 C and b) the overlap of myTile are not cube-face edges
496 IF ( .NOT.overlapOnly .OR. E_edge .OR. W_edge ) THEN
497
498 C- Internal exchange for calculations in Y
499 IF ( overlapOnly ) THEN
500 CALL FILL_CS_CORNER_TR_RL( 2, .FALSE.,
501 & localTij, bi,bj, myThid )
502 ENDIF
503
504 C- Advective flux in Y
505 DO j=1-Oly,sNy+Oly
506 DO i=1-Olx,sNx+Olx
507 af(i,j) = 0.
508 ENDDO
509 ENDDO
510
511 #ifdef ALLOW_AUTODIFF_TAMC
512 #ifndef DISABLE_MULTIDIM_ADVECTION
513 CADJ STORE localTij(:,:) =
514 CADJ & comlev1_bibj_k_gad_pass, key=passkey, kind=isbyte
515 #endif
516 #endif /* ALLOW_AUTODIFF_TAMC */
517
518 IF ( advectionScheme.EQ.ENUM_UPWIND_1RST
519 & .OR. advectionScheme.EQ.ENUM_DST2 ) THEN
520 CALL GAD_DST2U1_ADV_Y( bi,bj,k, advectionScheme, .TRUE.,
521 I deltaTLev(k),vTrans,vFld,localTij,
522 O af, myThid )
523 ELSEIF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN
524 CALL GAD_FLUXLIMIT_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k),
525 I vTrans, vFld, maskLocS, localTij,
526 O af, myThid )
527 ELSEIF (advectionScheme.EQ.ENUM_DST3 ) THEN
528 CALL GAD_DST3_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k),
529 I vTrans, vFld, maskLocS, localTij,
530 O af, myThid )
531 ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
532 CALL GAD_DST3FL_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k),
533 I vTrans, vFld, maskLocS, localTij,
534 O af, myThid )
535 #ifndef ALLOW_AUTODIFF_TAMC
536 ELSEIF (advectionScheme.EQ.ENUM_OS7MP ) THEN
537 CALL GAD_OS7MP_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k),
538 I vTrans, vFld, maskLocS, localTij,
539 O af, myThid )
540 #endif
541 ELSE
542 STOP 'GAD_ADVECTION: adv. scheme incompatibale with mutli-dim'
543 ENDIF
544
545 C- Internal exchange for next calculations in X
546 IF ( overlapOnly .AND. ipass.EQ.1 ) THEN
547 CALL FILL_CS_CORNER_TR_RL( 1, .FALSE.,
548 & localTij, bi,bj, myThid )
549 ENDIF
550
551 C- Advective flux in Y : done
552 ENDIF
553
554 C- Update the local tracer field where needed:
555 C use "maksInC" to prevent updating tracer field in OB regions
556
557 C update in overlap-Only
558 IF ( overlapOnly ) THEN
559 jMinUpd = 1-Oly+1
560 jMaxUpd = sNy+Oly-1
561 C- notes: these 2 lines below have no real effect (because recip_hFac=0
562 C in corner region) but safer to keep them.
563 IF ( S_edge ) jMinUpd = 1
564 IF ( N_edge ) jMaxUpd = sNy
565
566 IF ( W_edge ) THEN
567 DO j=jMinUpd,jMaxUpd
568 DO i=1-Olx,0
569 localTij(i,j) = localTij(i,j)
570 & -deltaTLev(k)*recip_rhoFacC(k)
571 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
572 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
573 & *( af(i,j+1)-af(i,j)
574 & -tracer(i,j,k,bi,bj)*(vTrans(i,j+1)-vTrans(i,j))
575 & )*maskInC(i,j,bi,bj)
576 ENDDO
577 ENDDO
578 ENDIF
579 IF ( E_edge ) THEN
580 DO j=jMinUpd,jMaxUpd
581 DO i=sNx+1,sNx+Olx
582 localTij(i,j) = localTij(i,j)
583 & -deltaTLev(k)*recip_rhoFacC(k)
584 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
585 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
586 & *( af(i,j+1)-af(i,j)
587 & -tracer(i,j,k,bi,bj)*(vTrans(i,j+1)-vTrans(i,j))
588 & )*maskInC(i,j,bi,bj)
589 ENDDO
590 ENDDO
591 ENDIF
592
593 ELSE
594 C do not only update the overlap
595 iMinUpd = 1-Olx
596 iMaxUpd = sNx+Olx
597 IF ( interiorOnly .AND. W_edge ) iMinUpd = 1
598 IF ( interiorOnly .AND. E_edge ) iMaxUpd = sNx
599 DO j=1-Oly+1,sNy+Oly-1
600 DO i=iMinUpd,iMaxUpd
601 localTij(i,j) = localTij(i,j)
602 & -deltaTLev(k)*recip_rhoFacC(k)
603 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
604 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
605 & *( af(i,j+1)-af(i,j)
606 & -tracer(i,j,k,bi,bj)*(vTrans(i,j+1)-vTrans(i,j))
607 & )*maskInC(i,j,bi,bj)
608 ENDDO
609 ENDDO
610 C- keep advective flux (for diagnostics)
611 DO j=1-Oly,sNy+Oly
612 DO i=1-Olx,sNx+Olx
613 afy(i,j) = af(i,j)
614 ENDDO
615 ENDDO
616
617 C end if/else update overlap-Only
618 ENDIF
619
620 C-- End of Y direction
621 ENDIF
622
623 C-- End of ipass loop
624 ENDDO
625
626 IF ( implicitAdvection ) THEN
627 C- explicit advection is done ; store tendency in gTracer:
628 DO j=1-Oly,sNy+Oly
629 DO i=1-Olx,sNx+Olx
630 gTracer(i,j,k,bi,bj)=
631 & (localTij(i,j)-tracer(i,j,k,bi,bj))/deltaTLev(k)
632 ENDDO
633 ENDDO
634 ELSE
635 C- horizontal advection done; store intermediate result in 3D array:
636 DO j=1-Oly,sNy+Oly
637 DO i=1-Olx,sNx+Olx
638 localTijk(i,j,k)=localTij(i,j)
639 ENDDO
640 ENDDO
641 ENDIF
642
643 #ifdef ALLOW_DIAGNOSTICS
644 IF ( doDiagAdvX ) THEN
645 diagName = 'ADVx'//diagSufx
646 CALL DIAGNOSTICS_FILL(afx,diagName, k,1, 2,bi,bj, myThid)
647 ENDIF
648 IF ( doDiagAdvY ) THEN
649 diagName = 'ADVy'//diagSufx
650 CALL DIAGNOSTICS_FILL(afy,diagName, k,1, 2,bi,bj, myThid)
651 ENDIF
652 #endif
653
654 #ifdef ALLOW_DEBUG
655 IF ( debugLevel .GE. debLevC
656 & .AND. tracerIdentity.EQ.GAD_TEMPERATURE
657 & .AND. k.LE.3 .AND. myIter.EQ.1+nIter0
658 & .AND. nPx.EQ.1 .AND. nPy.EQ.1
659 & .AND. useCubedSphereExchange ) THEN
660 CALL DEBUG_CS_CORNER_UV( ' afx,afy from GAD_ADVECTION',
661 & afx,afy, k, standardMessageUnit,bi,bj,myThid )
662 ENDIF
663 #endif /* ALLOW_DEBUG */
664
665 C-- End of K loop for horizontal fluxes
666 ENDDO
667
668 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
669
670 IF ( .NOT.implicitAdvection ) THEN
671 C-- Start of k loop for vertical flux
672 DO k=Nr,1,-1
673 #ifdef ALLOW_AUTODIFF_TAMC
674 kkey = (igadkey-1)*Nr + (Nr-k+1)
675 #endif /* ALLOW_AUTODIFF_TAMC */
676 C-- kUp Cycles through 1,2 to point to w-layer above
677 C-- kDown Cycles through 2,1 to point to w-layer below
678 kUp = 1+MOD(k+1,2)
679 kDown= 1+MOD(k,2)
680 c kp1=min(Nr,k+1)
681 kp1Msk=1.
682 if (k.EQ.Nr) kp1Msk=0.
683
684 #ifdef ALLOW_AUTODIFF_TAMC
685 CADJ STORE rtrans(:,:) =
686 CADJ & comlev1_bibj_k_gad, key=kkey, kind=isbyte
687 cphCADJ STORE wfld(:,:) =
688 cphCADJ & comlev1_bibj_k_gad, key=kkey, kind=isbyte
689 #endif
690
691 C-- Compute Vertical transport
692 #ifdef ALLOW_AIM
693 C- a hack to prevent Water-Vapor vert.transport into the stratospheric level Nr
694 IF ( k.EQ.1 .OR.
695 & (useAIM .AND. tracerIdentity.EQ.GAD_SALINITY .AND. k.EQ.Nr)
696 & ) THEN
697 #else
698 IF ( k.EQ.1 ) THEN
699 #endif
700
701 #ifdef ALLOW_AUTODIFF_TAMC
702 cphmultiCADJ STORE wfld(:,:) =
703 cphmultiCADJ & comlev1_bibj_k_gad, key=kkey, kind=isbyte
704 #endif /* ALLOW_AUTODIFF_TAMC */
705
706 C- Surface interface :
707 DO j=1-Oly,sNy+Oly
708 DO i=1-Olx,sNx+Olx
709 rTransKp1(i,j) = kp1Msk*rTrans(i,j)
710 wFld(i,j) = 0.
711 rTrans(i,j) = 0.
712 fVerT(i,j,kUp) = 0.
713 ENDDO
714 ENDDO
715
716 ELSE
717
718 #ifdef ALLOW_AUTODIFF_TAMC
719 cphmultiCADJ STORE wfld(:,:) =
720 cphmultiCADJ & comlev1_bibj_k_gad, key=kkey, kind=isbyte
721 #endif /* ALLOW_AUTODIFF_TAMC */
722
723 C- Interior interface :
724 DO j=1-Oly,sNy+Oly
725 DO i=1-Olx,sNx+Olx
726 rTransKp1(i,j) = kp1Msk*rTrans(i,j)
727 wFld(i,j) = wVel(i,j,k,bi,bj)
728 rTrans(i,j) = wVel(i,j,k,bi,bj)*rA(i,j,bi,bj)
729 & *deepFac2F(k)*rhoFacF(k)
730 & *maskC(i,j,k-1,bi,bj)
731 fVerT(i,j,kUp) = 0.
732 ENDDO
733 ENDDO
734
735 #ifdef ALLOW_GMREDI
736 C-- Residual transp = Bolus transp + Eulerian transp
737 IF (useGMRedi)
738 & CALL GMREDI_CALC_WFLOW(
739 U wFld, rTrans,
740 I k, bi, bj, myThid )
741 #endif /* ALLOW_GMREDI */
742
743 #ifdef ALLOW_AUTODIFF_TAMC
744 cphmultiCADJ STORE localTijk(:,:,k)
745 cphmultiCADJ & = comlev1_bibj_k_gad, key=kkey, kind=isbyte
746 cphmultiCADJ STORE rTrans(:,:)
747 cphmultiCADJ & = comlev1_bibj_k_gad, key=kkey, kind=isbyte
748 #endif /* ALLOW_AUTODIFF_TAMC */
749
750 C- Compute vertical advective flux in the interior:
751 IF ( vertAdvecScheme.EQ.ENUM_UPWIND_1RST
752 & .OR. vertAdvecScheme.EQ.ENUM_DST2 ) THEN
753 CALL GAD_DST2U1_ADV_R( bi,bj,k, advectionScheme,
754 I deltaTLev(k),rTrans,wFld,localTijk,
755 O fVerT(1-Olx,1-Oly,kUp), myThid )
756 ELSEIF( vertAdvecScheme.EQ.ENUM_FLUX_LIMIT) THEN
757 CALL GAD_FLUXLIMIT_ADV_R( bi,bj,k, deltaTLev(k),
758 I rTrans, wFld, localTijk,
759 O fVerT(1-Olx,1-Oly,kUp), myThid )
760 ELSEIF( vertAdvecScheme.EQ.ENUM_DST3 ) THEN
761 CALL GAD_DST3_ADV_R( bi,bj,k, deltaTLev(k),
762 I rTrans, wFld, localTijk,
763 O fVerT(1-Olx,1-Oly,kUp), myThid )
764 ELSEIF( vertAdvecScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
765 CALL GAD_DST3FL_ADV_R( bi,bj,k, deltaTLev(k),
766 I rTrans, wFld, localTijk,
767 O fVerT(1-Olx,1-Oly,kUp), myThid )
768 #ifndef ALLOW_AUTODIFF_TAMC
769 ELSEIF (vertAdvecScheme.EQ.ENUM_OS7MP ) THEN
770 CALL GAD_OS7MP_ADV_R( bi,bj,k, deltaTLev(k),
771 I rTrans, wFld, localTijk,
772 O fVerT(1-Olx,1-Oly,kUp), myThid )
773 #endif
774 ELSE
775 STOP 'GAD_ADVECTION: adv. scheme incompatibale with mutli-dim'
776 ENDIF
777
778 C- end Surface/Interior if bloc
779 ENDIF
780
781 #ifdef ALLOW_AUTODIFF_TAMC
782 cphmultiCADJ STORE rTrans(:,:)
783 cphmultiCADJ & = comlev1_bibj_k_gad, key=kkey, kind=isbyte
784 cphmultiCADJ STORE rTranskp1(:,:)
785 cphmultiCADJ & = comlev1_bibj_k_gad, key=kkey, kind=isbyte
786 cph --- following storing of fVerT is critical for correct
787 cph --- gradient with multiDimAdvection
788 cph --- Without it, kDown component is not properly recomputed
789 cph --- This is a TAF bug (and no warning available)
790 CADJ STORE fVerT(:,:,:)
791 CADJ & = comlev1_bibj_k_gad, key=kkey, kind=isbyte
792 #endif /* ALLOW_AUTODIFF_TAMC */
793
794 C-- Divergence of vertical fluxes
795 DO j=1-Oly,sNy+Oly
796 DO i=1-Olx,sNx+Olx
797 localTij(i,j) = localTijk(i,j,k)
798 & -deltaTLev(k)*recip_rhoFacC(k)
799 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
800 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
801 & *( fVerT(i,j,kDown)-fVerT(i,j,kUp)
802 & -tracer(i,j,k,bi,bj)*(rTransKp1(i,j)-rTrans(i,j))
803 & )*rkSign
804 gTracer(i,j,k,bi,bj)=
805 & (localTij(i,j)-tracer(i,j,k,bi,bj))/deltaTLev(k)
806 ENDDO
807 ENDDO
808
809 #ifdef ALLOW_DIAGNOSTICS
810 IF ( doDiagAdvR ) THEN
811 diagName = 'ADVr'//diagSufx
812 CALL DIAGNOSTICS_FILL( fVerT(1-Olx,1-Oly,kUp),
813 & diagName, k,1, 2,bi,bj, myThid)
814 ENDIF
815 #endif
816
817 C-- End of K loop for vertical flux
818 ENDDO
819 C-- end of if not.implicitAdvection block
820 ENDIF
821
822 RETURN
823 END

  ViewVC Help
Powered by ViewVC 1.1.22