/[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.69 - (show annotations) (download)
Fri Jun 24 22:18:47 2011 UTC (12 years, 10 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint63, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63a, checkpoint63b, checkpoint63c
Changes since 1.68: +3 -1 lines
avoid un-used variables

1 C $Header: /u/gcmpack/MITgcm/pkg/generic_advdiff/gad_advection.F,v 1.68 2011/06/21 16:13:52 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_AUTODIFF_TAMC
146 C msgBuf :: Informational/error message buffer
147 CHARACTER*(MAX_LEN_MBUF) msgBuf
148 #endif
149 #ifdef ALLOW_EXCH2
150 INTEGER myTile
151 #endif
152 #ifdef ALLOW_DIAGNOSTICS
153 CHARACTER*8 diagName
154 CHARACTER*4 diagSufx
155 LOGICAL doDiagAdvX, doDiagAdvY, doDiagAdvR
156 C- Functions:
157 CHARACTER*4 GAD_DIAG_SUFX
158 EXTERNAL GAD_DIAG_SUFX
159 LOGICAL DIAGNOSTICS_IS_ON
160 EXTERNAL DIAGNOSTICS_IS_ON
161 #endif
162 CEOP
163
164 #ifdef ALLOW_AUTODIFF_TAMC
165 act0 = tracerIdentity
166 max0 = maxpass
167 act1 = bi - myBxLo(myThid)
168 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
169 act2 = bj - myByLo(myThid)
170 max2 = myByHi(myThid) - myByLo(myThid) + 1
171 act3 = myThid - 1
172 max3 = nTx*nTy
173 act4 = ikey_dynamics - 1
174 igadkey = act0
175 & + act1*max0
176 & + act2*max0*max1
177 & + act3*max0*max1*max2
178 & + act4*max0*max1*max2*max3
179 IF (tracerIdentity.GT.maxpass) THEN
180 WRITE(msgBuf,'(A,2I3)')
181 & 'GAD_ADVECTION: maxpass < tracerIdentity ',
182 & maxpass, tracerIdentity
183 CALL PRINT_ERROR( msgBuf, myThid )
184 STOP 'ABNORMAL END: S/R GAD_ADVECTION'
185 ENDIF
186 #endif /* ALLOW_AUTODIFF_TAMC */
187
188 #ifdef ALLOW_DIAGNOSTICS
189 C-- Set diagnostics flags and suffix for the current tracer
190 doDiagAdvX = .FALSE.
191 doDiagAdvY = .FALSE.
192 doDiagAdvR = .FALSE.
193 IF ( useDiagnostics ) THEN
194 diagSufx = GAD_DIAG_SUFX( tracerIdentity, myThid )
195 diagName = 'ADVx'//diagSufx
196 doDiagAdvX = DIAGNOSTICS_IS_ON( diagName, myThid )
197 diagName = 'ADVy'//diagSufx
198 doDiagAdvY = DIAGNOSTICS_IS_ON( diagName, myThid )
199 diagName = 'ADVr'//diagSufx
200 doDiagAdvR = DIAGNOSTICS_IS_ON( diagName, myThid )
201 ENDIF
202 #endif
203
204 C-- Set up work arrays with valid (i.e. not NaN) values
205 C These inital values do not alter the numerical results. They
206 C just ensure that all memory references are to valid floating
207 C point numbers. This prevents spurious hardware signals due to
208 C uninitialised but inert locations.
209 DO j=1-OLy,sNy+OLy
210 DO i=1-OLx,sNx+OLx
211 xA(i,j) = 0. _d 0
212 yA(i,j) = 0. _d 0
213 uTrans(i,j) = 0. _d 0
214 vTrans(i,j) = 0. _d 0
215 rTrans(i,j) = 0. _d 0
216 fVerT(i,j,1) = 0. _d 0
217 fVerT(i,j,2) = 0. _d 0
218 rTransKp1(i,j)= 0. _d 0
219 #ifdef ALLOW_AUTODIFF_TAMC
220 localTij(i,j) = 0. _d 0
221 wfld(i,j) = 0. _d 0
222 #endif
223 ENDDO
224 ENDDO
225
226 C-- Set tile-specific parameters for horizontal fluxes
227 IF (useCubedSphereExchange) THEN
228 npass = 3
229 #ifdef ALLOW_AUTODIFF_TAMC
230 IF ( npass.GT.maxcube ) STOP 'maxcube needs to be = 3'
231 #endif
232 #ifdef ALLOW_EXCH2
233 myTile = W2_myTileList(bi,bj)
234 nCFace = exch2_myFace(myTile)
235 N_edge = exch2_isNedge(myTile).EQ.1
236 S_edge = exch2_isSedge(myTile).EQ.1
237 E_edge = exch2_isEedge(myTile).EQ.1
238 W_edge = exch2_isWedge(myTile).EQ.1
239 #else
240 nCFace = bi
241 N_edge = .TRUE.
242 S_edge = .TRUE.
243 E_edge = .TRUE.
244 W_edge = .TRUE.
245 #endif
246 ELSE
247 npass = 2
248 nCFace = 0
249 N_edge = .FALSE.
250 S_edge = .FALSE.
251 E_edge = .FALSE.
252 W_edge = .FALSE.
253 ENDIF
254
255 C-- Start of k loop for horizontal fluxes
256 DO k=1,Nr
257 #ifdef ALLOW_AUTODIFF_TAMC
258 kkey = (igadkey-1)*Nr + k
259 CADJ STORE tracer(:,:,k,bi,bj) =
260 CADJ & comlev1_bibj_k_gad, key=kkey, kind=isbyte
261 #endif /* ALLOW_AUTODIFF_TAMC */
262
263 C-- Get temporary terms used by tendency routines
264 CALL CALC_COMMON_FACTORS (
265 I uVel, vVel,
266 O uFld, vFld, uTrans, vTrans, xA, yA,
267 I k,bi,bj, myThid )
268
269 #ifdef ALLOW_GMREDI
270 C-- Residual transp = Bolus transp + Eulerian transp
271 IF (useGMRedi)
272 & CALL GMREDI_CALC_UVFLOW(
273 U uFld, vFld, uTrans, vTrans,
274 I k, bi, bj, myThid )
275 #endif /* ALLOW_GMREDI */
276
277 C-- Make local copy of tracer array and mask West & South
278 DO j=1-OLy,sNy+OLy
279 DO i=1-OLx,sNx+OLx
280 localTij(i,j)=tracer(i,j,k,bi,bj)
281 #ifdef ALLOW_OBCS
282 maskLocW(i,j) = _maskW(i,j,k,bi,bj)*maskInW(i,j,bi,bj)
283 maskLocS(i,j) = _maskS(i,j,k,bi,bj)*maskInS(i,j,bi,bj)
284 #else /* ALLOW_OBCS */
285 maskLocW(i,j) = _maskW(i,j,k,bi,bj)
286 maskLocS(i,j) = _maskS(i,j,k,bi,bj)
287 #endif /* ALLOW_OBCS */
288 ENDDO
289 ENDDO
290
291 cph-exch2#ifndef ALLOW_AUTODIFF_TAMC
292 IF (useCubedSphereExchange) THEN
293 withSigns = .FALSE.
294 CALL FILL_CS_CORNER_UV_RS(
295 & withSigns, maskLocW,maskLocS, bi,bj, myThid )
296 ENDIF
297 cph-exch2#endif
298
299 C-- Multiple passes for different directions on different tiles
300 C-- For cube need one pass for each of red, green and blue axes.
301 DO ipass=1,npass
302 #ifdef ALLOW_AUTODIFF_TAMC
303 passkey = ipass
304 & + (k-1) *maxpass
305 & + (igadkey-1)*maxpass*Nr
306 IF (npass .GT. maxpass) THEN
307 STOP 'GAD_ADVECTION: npass > maxcube. check tamc.h'
308 ENDIF
309 #endif /* ALLOW_AUTODIFF_TAMC */
310
311 interiorOnly = .FALSE.
312 overlapOnly = .FALSE.
313 IF (useCubedSphereExchange) THEN
314 C- CubedSphere : pass 3 times, with partial update of local tracer field
315 IF (ipass.EQ.1) THEN
316 overlapOnly = MOD(nCFace,3).EQ.0
317 interiorOnly = MOD(nCFace,3).NE.0
318 calc_fluxes_X = nCFace.EQ.6 .OR. nCFace.EQ.1 .OR. nCFace.EQ.2
319 calc_fluxes_Y = nCFace.EQ.3 .OR. nCFace.EQ.4 .OR. nCFace.EQ.5
320 ELSEIF (ipass.EQ.2) THEN
321 overlapOnly = MOD(nCFace,3).EQ.2
322 interiorOnly = MOD(nCFace,3).EQ.1
323 calc_fluxes_X = nCFace.EQ.2 .OR. nCFace.EQ.3 .OR. nCFace.EQ.4
324 calc_fluxes_Y = nCFace.EQ.5 .OR. nCFace.EQ.6 .OR. nCFace.EQ.1
325 ELSE
326 interiorOnly = .TRUE.
327 calc_fluxes_X = nCFace.EQ.5 .OR. nCFace.EQ.6
328 calc_fluxes_Y = nCFace.EQ.2 .OR. nCFace.EQ.3
329 ENDIF
330 ELSE
331 C- not CubedSphere
332 calc_fluxes_X = MOD(ipass,2).EQ.1
333 calc_fluxes_Y = .NOT.calc_fluxes_X
334 ENDIF
335
336 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
337 C-- X direction
338 C- Advective flux in X
339 DO j=1-Oly,sNy+Oly
340 DO i=1-Olx,sNx+Olx
341 af(i,j) = 0.
342 ENDDO
343 ENDDO
344 C
345 #ifdef ALLOW_AUTODIFF_TAMC
346 # ifndef DISABLE_MULTIDIM_ADVECTION
347 CADJ STORE localTij(:,:) =
348 CADJ & comlev1_bibj_k_gad_pass, key=passkey, kind=isbyte
349 CADJ STORE af(:,:) =
350 CADJ & comlev1_bibj_k_gad_pass, key=passkey, kind=isbyte
351 # endif
352 #endif /* ALLOW_AUTODIFF_TAMC */
353 C
354 IF (calc_fluxes_X) THEN
355
356 C- Do not compute fluxes if
357 C a) needed in overlap only
358 C and b) the overlap of myTile are not cube-face Edges
359 IF ( .NOT.overlapOnly .OR. N_edge .OR. S_edge ) THEN
360
361 C- Internal exchange for calculations in X
362 IF ( overlapOnly ) THEN
363 CALL FILL_CS_CORNER_TR_RL( 1, .FALSE.,
364 & localTij, bi,bj, myThid )
365 ENDIF
366
367 #ifdef ALLOW_AUTODIFF_TAMC
368 # ifndef DISABLE_MULTIDIM_ADVECTION
369 CADJ STORE localTij(:,:) =
370 CADJ & comlev1_bibj_k_gad_pass, key=passkey, kind=isbyte
371 # endif
372 #endif /* ALLOW_AUTODIFF_TAMC */
373
374 IF ( advectionScheme.EQ.ENUM_UPWIND_1RST
375 & .OR. advectionScheme.EQ.ENUM_DST2 ) THEN
376 CALL GAD_DST2U1_ADV_X( bi,bj,k, advectionScheme, .TRUE.,
377 I deltaTLev(k),uTrans,uFld,localTij,
378 O af, myThid )
379 ELSEIF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN
380 CALL GAD_FLUXLIMIT_ADV_X( bi,bj,k, .TRUE., deltaTLev(k),
381 I uTrans, uFld, maskLocW, localTij,
382 O af, myThid )
383 ELSEIF (advectionScheme.EQ.ENUM_DST3 ) THEN
384 CALL GAD_DST3_ADV_X( bi,bj,k, .TRUE., deltaTLev(k),
385 I uTrans, uFld, maskLocW, localTij,
386 O af, myThid )
387 ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
388 CALL GAD_DST3FL_ADV_X( bi,bj,k, .TRUE., deltaTLev(k),
389 I uTrans, uFld, maskLocW, localTij,
390 O af, myThid )
391 #ifndef ALLOW_AUTODIFF_TAMC
392 ELSEIF (advectionScheme.EQ.ENUM_OS7MP ) THEN
393 CALL GAD_OS7MP_ADV_X( bi,bj,k, .TRUE., deltaTLev(k),
394 I uTrans, uFld, maskLocW, localTij,
395 O af, myThid )
396 #endif
397 ELSE
398 STOP 'GAD_ADVECTION: adv. scheme incompatibale with multi-dim'
399 ENDIF
400
401 C- Internal exchange for next calculations in Y
402 IF ( overlapOnly .AND. ipass.EQ.1 ) THEN
403 CALL FILL_CS_CORNER_TR_RL( 2, .FALSE.,
404 & localTij, bi,bj, myThid )
405 ENDIF
406
407 C- Advective flux in X : done
408 ENDIF
409
410 C- Update the local tracer field where needed:
411 C use "maksInC" to prevent updating tracer field in OB regions
412
413 C update in overlap-Only
414 IF ( overlapOnly ) THEN
415 iMinUpd = 1-Olx+1
416 iMaxUpd = sNx+Olx-1
417 C- notes: these 2 lines below have no real effect (because recip_hFac=0
418 C in corner region) but safer to keep them.
419 IF ( W_edge ) iMinUpd = 1
420 IF ( E_edge ) iMaxUpd = sNx
421
422 IF ( S_edge ) THEN
423 DO j=1-Oly,0
424 DO i=iMinUpd,iMaxUpd
425 localTij(i,j) = localTij(i,j)
426 & -deltaTLev(k)*recip_rhoFacC(k)
427 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
428 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
429 & *( af(i+1,j)-af(i,j)
430 & -tracer(i,j,k,bi,bj)*(uTrans(i+1,j)-uTrans(i,j))
431 & )*maskInC(i,j,bi,bj)
432 ENDDO
433 ENDDO
434 ENDIF
435 IF ( N_edge ) THEN
436 DO j=sNy+1,sNy+Oly
437 DO i=iMinUpd,iMaxUpd
438 localTij(i,j) = localTij(i,j)
439 & -deltaTLev(k)*recip_rhoFacC(k)
440 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
441 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
442 & *( af(i+1,j)-af(i,j)
443 & -tracer(i,j,k,bi,bj)*(uTrans(i+1,j)-uTrans(i,j))
444 & )*maskInC(i,j,bi,bj)
445 ENDDO
446 ENDDO
447 ENDIF
448
449 ELSE
450 C do not only update the overlap
451 jMinUpd = 1-Oly
452 jMaxUpd = sNy+Oly
453 IF ( interiorOnly .AND. S_edge ) jMinUpd = 1
454 IF ( interiorOnly .AND. N_edge ) jMaxUpd = sNy
455 DO j=jMinUpd,jMaxUpd
456 DO i=1-Olx+1,sNx+Olx-1
457 localTij(i,j) = localTij(i,j)
458 & -deltaTLev(k)*recip_rhoFacC(k)
459 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
460 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
461 & *( af(i+1,j)-af(i,j)
462 & -tracer(i,j,k,bi,bj)*(uTrans(i+1,j)-uTrans(i,j))
463 & )*maskInC(i,j,bi,bj)
464 ENDDO
465 ENDDO
466 C- keep advective flux (for diagnostics)
467 DO j=1-Oly,sNy+Oly
468 DO i=1-Olx,sNx+Olx
469 afx(i,j) = af(i,j)
470 ENDDO
471 ENDDO
472
473 C- end if/else update overlap-Only
474 ENDIF
475
476 C-- End of X direction
477 ENDIF
478
479 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
480 C-- Y direction
481 cph-test
482 C- Advective flux in Y
483 DO j=1-Oly,sNy+Oly
484 DO i=1-Olx,sNx+Olx
485 af(i,j) = 0.
486 ENDDO
487 ENDDO
488 C
489 #ifdef ALLOW_AUTODIFF_TAMC
490 # ifndef DISABLE_MULTIDIM_ADVECTION
491 CADJ STORE localTij(:,:) =
492 CADJ & comlev1_bibj_k_gad_pass, key=passkey, kind=isbyte
493 CADJ STORE af(:,:) =
494 CADJ & comlev1_bibj_k_gad_pass, key=passkey, kind=isbyte
495 # endif
496 #endif /* ALLOW_AUTODIFF_TAMC */
497 C
498 IF (calc_fluxes_Y) THEN
499
500 C- Do not compute fluxes if
501 C a) needed in overlap only
502 C and b) the overlap of myTile are not cube-face edges
503 IF ( .NOT.overlapOnly .OR. E_edge .OR. W_edge ) THEN
504
505 C- Internal exchange for calculations in Y
506 IF ( overlapOnly ) THEN
507 CALL FILL_CS_CORNER_TR_RL( 2, .FALSE.,
508 & localTij, bi,bj, myThid )
509 ENDIF
510
511 C- Advective flux in Y
512 DO j=1-Oly,sNy+Oly
513 DO i=1-Olx,sNx+Olx
514 af(i,j) = 0.
515 ENDDO
516 ENDDO
517
518 #ifdef ALLOW_AUTODIFF_TAMC
519 #ifndef DISABLE_MULTIDIM_ADVECTION
520 CADJ STORE localTij(:,:) =
521 CADJ & comlev1_bibj_k_gad_pass, key=passkey, kind=isbyte
522 #endif
523 #endif /* ALLOW_AUTODIFF_TAMC */
524
525 IF ( advectionScheme.EQ.ENUM_UPWIND_1RST
526 & .OR. advectionScheme.EQ.ENUM_DST2 ) THEN
527 CALL GAD_DST2U1_ADV_Y( bi,bj,k, advectionScheme, .TRUE.,
528 I deltaTLev(k),vTrans,vFld,localTij,
529 O af, myThid )
530 ELSEIF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN
531 CALL GAD_FLUXLIMIT_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k),
532 I vTrans, vFld, maskLocS, localTij,
533 O af, myThid )
534 ELSEIF (advectionScheme.EQ.ENUM_DST3 ) THEN
535 CALL GAD_DST3_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k),
536 I vTrans, vFld, maskLocS, localTij,
537 O af, myThid )
538 ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
539 CALL GAD_DST3FL_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k),
540 I vTrans, vFld, maskLocS, localTij,
541 O af, myThid )
542 #ifndef ALLOW_AUTODIFF_TAMC
543 ELSEIF (advectionScheme.EQ.ENUM_OS7MP ) THEN
544 CALL GAD_OS7MP_ADV_Y( bi,bj,k, .TRUE., deltaTLev(k),
545 I vTrans, vFld, maskLocS, localTij,
546 O af, myThid )
547 #endif
548 ELSE
549 STOP 'GAD_ADVECTION: adv. scheme incompatibale with mutli-dim'
550 ENDIF
551
552 C- Internal exchange for next calculations in X
553 IF ( overlapOnly .AND. ipass.EQ.1 ) THEN
554 CALL FILL_CS_CORNER_TR_RL( 1, .FALSE.,
555 & localTij, bi,bj, myThid )
556 ENDIF
557
558 C- Advective flux in Y : done
559 ENDIF
560
561 C- Update the local tracer field where needed:
562 C use "maksInC" to prevent updating tracer field in OB regions
563
564 C update in overlap-Only
565 IF ( overlapOnly ) THEN
566 jMinUpd = 1-Oly+1
567 jMaxUpd = sNy+Oly-1
568 C- notes: these 2 lines below have no real effect (because recip_hFac=0
569 C in corner region) but safer to keep them.
570 IF ( S_edge ) jMinUpd = 1
571 IF ( N_edge ) jMaxUpd = sNy
572
573 IF ( W_edge ) THEN
574 DO j=jMinUpd,jMaxUpd
575 DO i=1-Olx,0
576 localTij(i,j) = localTij(i,j)
577 & -deltaTLev(k)*recip_rhoFacC(k)
578 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
579 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
580 & *( af(i,j+1)-af(i,j)
581 & -tracer(i,j,k,bi,bj)*(vTrans(i,j+1)-vTrans(i,j))
582 & )*maskInC(i,j,bi,bj)
583 ENDDO
584 ENDDO
585 ENDIF
586 IF ( E_edge ) THEN
587 DO j=jMinUpd,jMaxUpd
588 DO i=sNx+1,sNx+Olx
589 localTij(i,j) = localTij(i,j)
590 & -deltaTLev(k)*recip_rhoFacC(k)
591 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
592 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
593 & *( af(i,j+1)-af(i,j)
594 & -tracer(i,j,k,bi,bj)*(vTrans(i,j+1)-vTrans(i,j))
595 & )*maskInC(i,j,bi,bj)
596 ENDDO
597 ENDDO
598 ENDIF
599
600 ELSE
601 C do not only update the overlap
602 iMinUpd = 1-Olx
603 iMaxUpd = sNx+Olx
604 IF ( interiorOnly .AND. W_edge ) iMinUpd = 1
605 IF ( interiorOnly .AND. E_edge ) iMaxUpd = sNx
606 DO j=1-Oly+1,sNy+Oly-1
607 DO i=iMinUpd,iMaxUpd
608 localTij(i,j) = localTij(i,j)
609 & -deltaTLev(k)*recip_rhoFacC(k)
610 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
611 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
612 & *( af(i,j+1)-af(i,j)
613 & -tracer(i,j,k,bi,bj)*(vTrans(i,j+1)-vTrans(i,j))
614 & )*maskInC(i,j,bi,bj)
615 ENDDO
616 ENDDO
617 C- keep advective flux (for diagnostics)
618 DO j=1-Oly,sNy+Oly
619 DO i=1-Olx,sNx+Olx
620 afy(i,j) = af(i,j)
621 ENDDO
622 ENDDO
623
624 C end if/else update overlap-Only
625 ENDIF
626
627 C-- End of Y direction
628 ENDIF
629
630 C-- End of ipass loop
631 ENDDO
632
633 IF ( implicitAdvection ) THEN
634 C- explicit advection is done ; store tendency in gTracer:
635 DO j=1-Oly,sNy+Oly
636 DO i=1-Olx,sNx+Olx
637 gTracer(i,j,k,bi,bj)=
638 & (localTij(i,j)-tracer(i,j,k,bi,bj))/deltaTLev(k)
639 ENDDO
640 ENDDO
641 ELSE
642 C- horizontal advection done; store intermediate result in 3D array:
643 DO j=1-Oly,sNy+Oly
644 DO i=1-Olx,sNx+Olx
645 localTijk(i,j,k)=localTij(i,j)
646 ENDDO
647 ENDDO
648 ENDIF
649
650 #ifdef ALLOW_DIAGNOSTICS
651 IF ( doDiagAdvX ) THEN
652 diagName = 'ADVx'//diagSufx
653 CALL DIAGNOSTICS_FILL(afx,diagName, k,1, 2,bi,bj, myThid)
654 ENDIF
655 IF ( doDiagAdvY ) THEN
656 diagName = 'ADVy'//diagSufx
657 CALL DIAGNOSTICS_FILL(afy,diagName, k,1, 2,bi,bj, myThid)
658 ENDIF
659 #endif
660
661 #ifdef ALLOW_DEBUG
662 IF ( debugLevel .GE. debLevC
663 & .AND. tracerIdentity.EQ.GAD_TEMPERATURE
664 & .AND. k.LE.3 .AND. myIter.EQ.1+nIter0
665 & .AND. nPx.EQ.1 .AND. nPy.EQ.1
666 & .AND. useCubedSphereExchange ) THEN
667 CALL DEBUG_CS_CORNER_UV( ' afx,afy from GAD_ADVECTION',
668 & afx,afy, k, standardMessageUnit,bi,bj,myThid )
669 ENDIF
670 #endif /* ALLOW_DEBUG */
671
672 C-- End of K loop for horizontal fluxes
673 ENDDO
674
675 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
676
677 IF ( .NOT.implicitAdvection ) THEN
678 C-- Start of k loop for vertical flux
679 DO k=Nr,1,-1
680 #ifdef ALLOW_AUTODIFF_TAMC
681 kkey = (igadkey-1)*Nr + (Nr-k+1)
682 #endif /* ALLOW_AUTODIFF_TAMC */
683 C-- kUp Cycles through 1,2 to point to w-layer above
684 C-- kDown Cycles through 2,1 to point to w-layer below
685 kUp = 1+MOD(k+1,2)
686 kDown= 1+MOD(k,2)
687 c kp1=min(Nr,k+1)
688 kp1Msk=1.
689 if (k.EQ.Nr) kp1Msk=0.
690
691 #ifdef ALLOW_AUTODIFF_TAMC
692 CADJ STORE rtrans(:,:) =
693 CADJ & comlev1_bibj_k_gad, key=kkey, kind=isbyte
694 cphCADJ STORE wfld(:,:) =
695 cphCADJ & comlev1_bibj_k_gad, key=kkey, kind=isbyte
696 #endif
697
698 C-- Compute Vertical transport
699 #ifdef ALLOW_AIM
700 C- a hack to prevent Water-Vapor vert.transport into the stratospheric level Nr
701 IF ( k.EQ.1 .OR.
702 & (useAIM .AND. tracerIdentity.EQ.GAD_SALINITY .AND. k.EQ.Nr)
703 & ) THEN
704 #else
705 IF ( k.EQ.1 ) THEN
706 #endif
707
708 #ifdef ALLOW_AUTODIFF_TAMC
709 cphmultiCADJ STORE wfld(:,:) =
710 cphmultiCADJ & comlev1_bibj_k_gad, key=kkey, kind=isbyte
711 #endif /* ALLOW_AUTODIFF_TAMC */
712
713 C- Surface interface :
714 DO j=1-Oly,sNy+Oly
715 DO i=1-Olx,sNx+Olx
716 rTransKp1(i,j) = kp1Msk*rTrans(i,j)
717 wFld(i,j) = 0.
718 rTrans(i,j) = 0.
719 fVerT(i,j,kUp) = 0.
720 ENDDO
721 ENDDO
722
723 ELSE
724
725 #ifdef ALLOW_AUTODIFF_TAMC
726 cphmultiCADJ STORE wfld(:,:) =
727 cphmultiCADJ & comlev1_bibj_k_gad, key=kkey, kind=isbyte
728 #endif /* ALLOW_AUTODIFF_TAMC */
729
730 C- Interior interface :
731 DO j=1-Oly,sNy+Oly
732 DO i=1-Olx,sNx+Olx
733 rTransKp1(i,j) = kp1Msk*rTrans(i,j)
734 wFld(i,j) = wVel(i,j,k,bi,bj)
735 rTrans(i,j) = wVel(i,j,k,bi,bj)*rA(i,j,bi,bj)
736 & *deepFac2F(k)*rhoFacF(k)
737 & *maskC(i,j,k-1,bi,bj)
738 fVerT(i,j,kUp) = 0.
739 ENDDO
740 ENDDO
741
742 #ifdef ALLOW_GMREDI
743 C-- Residual transp = Bolus transp + Eulerian transp
744 IF (useGMRedi)
745 & CALL GMREDI_CALC_WFLOW(
746 U wFld, rTrans,
747 I k, bi, bj, myThid )
748 #endif /* ALLOW_GMREDI */
749
750 #ifdef ALLOW_AUTODIFF_TAMC
751 cphmultiCADJ STORE localTijk(:,:,k)
752 cphmultiCADJ & = comlev1_bibj_k_gad, key=kkey, kind=isbyte
753 cphmultiCADJ STORE rTrans(:,:)
754 cphmultiCADJ & = comlev1_bibj_k_gad, key=kkey, kind=isbyte
755 #endif /* ALLOW_AUTODIFF_TAMC */
756
757 C- Compute vertical advective flux in the interior:
758 IF ( vertAdvecScheme.EQ.ENUM_UPWIND_1RST
759 & .OR. vertAdvecScheme.EQ.ENUM_DST2 ) THEN
760 CALL GAD_DST2U1_ADV_R( bi,bj,k, advectionScheme,
761 I deltaTLev(k),rTrans,wFld,localTijk,
762 O fVerT(1-Olx,1-Oly,kUp), myThid )
763 ELSEIF( vertAdvecScheme.EQ.ENUM_FLUX_LIMIT) THEN
764 CALL GAD_FLUXLIMIT_ADV_R( bi,bj,k, deltaTLev(k),
765 I rTrans, wFld, localTijk,
766 O fVerT(1-Olx,1-Oly,kUp), myThid )
767 ELSEIF( vertAdvecScheme.EQ.ENUM_DST3 ) THEN
768 CALL GAD_DST3_ADV_R( bi,bj,k, deltaTLev(k),
769 I rTrans, wFld, localTijk,
770 O fVerT(1-Olx,1-Oly,kUp), myThid )
771 ELSEIF( vertAdvecScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
772 CALL GAD_DST3FL_ADV_R( bi,bj,k, deltaTLev(k),
773 I rTrans, wFld, localTijk,
774 O fVerT(1-Olx,1-Oly,kUp), myThid )
775 #ifndef ALLOW_AUTODIFF_TAMC
776 ELSEIF (vertAdvecScheme.EQ.ENUM_OS7MP ) THEN
777 CALL GAD_OS7MP_ADV_R( bi,bj,k, deltaTLev(k),
778 I rTrans, wFld, localTijk,
779 O fVerT(1-Olx,1-Oly,kUp), myThid )
780 #endif
781 ELSE
782 STOP 'GAD_ADVECTION: adv. scheme incompatibale with mutli-dim'
783 ENDIF
784
785 C- end Surface/Interior if bloc
786 ENDIF
787
788 #ifdef ALLOW_AUTODIFF_TAMC
789 cphmultiCADJ STORE rTrans(:,:)
790 cphmultiCADJ & = comlev1_bibj_k_gad, key=kkey, kind=isbyte
791 cphmultiCADJ STORE rTranskp1(:,:)
792 cphmultiCADJ & = comlev1_bibj_k_gad, key=kkey, kind=isbyte
793 cph --- following storing of fVerT is critical for correct
794 cph --- gradient with multiDimAdvection
795 cph --- Without it, kDown component is not properly recomputed
796 cph --- This is a TAF bug (and no warning available)
797 CADJ STORE fVerT(:,:,:)
798 CADJ & = comlev1_bibj_k_gad, key=kkey, kind=isbyte
799 #endif /* ALLOW_AUTODIFF_TAMC */
800
801 C-- Divergence of vertical fluxes
802 DO j=1-Oly,sNy+Oly
803 DO i=1-Olx,sNx+Olx
804 localTij(i,j) = localTijk(i,j,k)
805 & -deltaTLev(k)*recip_rhoFacC(k)
806 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
807 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
808 & *( fVerT(i,j,kDown)-fVerT(i,j,kUp)
809 & -tracer(i,j,k,bi,bj)*(rTransKp1(i,j)-rTrans(i,j))
810 & )*rkSign
811 gTracer(i,j,k,bi,bj)=
812 & (localTij(i,j)-tracer(i,j,k,bi,bj))/deltaTLev(k)
813 ENDDO
814 ENDDO
815
816 #ifdef ALLOW_DIAGNOSTICS
817 IF ( doDiagAdvR ) THEN
818 diagName = 'ADVr'//diagSufx
819 CALL DIAGNOSTICS_FILL( fVerT(1-Olx,1-Oly,kUp),
820 & diagName, k,1, 2,bi,bj, myThid)
821 ENDIF
822 #endif
823
824 C-- End of K loop for vertical flux
825 ENDDO
826 C-- end of if not.implicitAdvection block
827 ENDIF
828
829 RETURN
830 END

  ViewVC Help
Powered by ViewVC 1.1.22