/[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.55 - (show annotations) (download)
Thu Feb 7 23:59:02 2008 UTC (16 years, 3 months ago) by jmc
Branch: MAIN
Changes since 1.54: +48 -64 lines
cleaner version: avoid updating if unnecessary.

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

  ViewVC Help
Powered by ViewVC 1.1.22