/[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.58 - (show annotations) (download)
Wed Oct 22 00:28:47 2008 UTC (15 years, 7 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint61f, checkpoint61e, checkpoint61g, checkpoint61h, checkpoint61i
Changes since 1.57: +5 -5 lines
changes in FILL_CS_CORNER_TR_RL argument list.

1 C $Header: /u/gcmpack/MITgcm/pkg/generic_advdiff/gad_advection.F,v 1.57 2008/02/08 20:03:36 jmc 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 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)
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, byte=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 maskLocW(i,j)=_maskW(i,j,k,bi,bj)
275 maskLocS(i,j)=_maskS(i,j,k,bi,bj)
276 ENDDO
277 ENDDO
278
279 cph-exch2#ifndef ALLOW_AUTODIFF_TAMC
280 IF (useCubedSphereExchange) THEN
281 withSigns = .FALSE.
282 CALL FILL_CS_CORNER_UV_RS(
283 & withSigns, maskLocW,maskLocS, bi,bj, myThid )
284 ENDIF
285 cph-exch2#endif
286
287 C-- Multiple passes for different directions on different tiles
288 C-- For cube need one pass for each of red, green and blue axes.
289 DO ipass=1,npass
290 #ifdef ALLOW_AUTODIFF_TAMC
291 passkey = ipass
292 & + (k-1) *maxpass
293 & + (igadkey-1)*maxpass*Nr
294 IF (npass .GT. maxpass) THEN
295 STOP 'GAD_ADVECTION: npass > maxcube. check tamc.h'
296 ENDIF
297 #endif /* ALLOW_AUTODIFF_TAMC */
298
299 interiorOnly = .FALSE.
300 overlapOnly = .FALSE.
301 IF (useCubedSphereExchange) THEN
302 #ifdef MULTIDIM_OLD_VERSION
303 C- CubedSphere : pass 3 times, with full update of local tracer field
304 IF (ipass.EQ.1) THEN
305 calc_fluxes_X = nCFace.EQ.1 .OR. nCFace.EQ.2
306 calc_fluxes_Y = nCFace.EQ.4 .OR. nCFace.EQ.5
307 ELSEIF (ipass.EQ.2) THEN
308 calc_fluxes_X = nCFace.EQ.3 .OR. nCFace.EQ.4
309 calc_fluxes_Y = nCFace.EQ.6 .OR. nCFace.EQ.1
310 #else /* MULTIDIM_OLD_VERSION */
311 C- CubedSphere : pass 3 times, with partial update of local tracer field
312 IF (ipass.EQ.1) THEN
313 overlapOnly = MOD(nCFace,3).EQ.0
314 interiorOnly = MOD(nCFace,3).NE.0
315 calc_fluxes_X = nCFace.EQ.6 .OR. nCFace.EQ.1 .OR. nCFace.EQ.2
316 calc_fluxes_Y = nCFace.EQ.3 .OR. nCFace.EQ.4 .OR. nCFace.EQ.5
317 ELSEIF (ipass.EQ.2) THEN
318 overlapOnly = MOD(nCFace,3).EQ.2
319 interiorOnly = MOD(nCFace,3).EQ.1
320 calc_fluxes_X = nCFace.EQ.2 .OR. nCFace.EQ.3 .OR. nCFace.EQ.4
321 calc_fluxes_Y = nCFace.EQ.5 .OR. nCFace.EQ.6 .OR. nCFace.EQ.1
322 #endif /* MULTIDIM_OLD_VERSION */
323 ELSE
324 interiorOnly = .TRUE.
325 calc_fluxes_X = nCFace.EQ.5 .OR. nCFace.EQ.6
326 calc_fluxes_Y = nCFace.EQ.2 .OR. nCFace.EQ.3
327 ENDIF
328 ELSE
329 C- not CubedSphere
330 calc_fluxes_X = MOD(ipass,2).EQ.1
331 calc_fluxes_Y = .NOT.calc_fluxes_X
332 ENDIF
333
334 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
335 C-- X direction
336 C- Advective flux in X
337 DO j=1-Oly,sNy+Oly
338 DO i=1-Olx,sNx+Olx
339 af(i,j) = 0.
340 ENDDO
341 ENDDO
342 C
343 #ifdef ALLOW_AUTODIFF_TAMC
344 # ifndef DISABLE_MULTIDIM_ADVECTION
345 CADJ STORE localTij(:,:) =
346 CADJ & comlev1_bibj_k_gad_pass, key=passkey, byte=isbyte
347 CADJ STORE af(:,:) =
348 CADJ & comlev1_bibj_k_gad_pass, key=passkey, byte=isbyte
349 # endif
350 #endif /* ALLOW_AUTODIFF_TAMC */
351 C
352 IF (calc_fluxes_X) THEN
353
354 C- Do not compute fluxes if
355 C a) needed in overlap only
356 C and b) the overlap of myTile are not cube-face Edges
357 IF ( .NOT.overlapOnly .OR. N_edge .OR. S_edge ) THEN
358
359 C- Internal exchange for calculations in X
360 #ifdef MULTIDIM_OLD_VERSION
361 IF ( useCubedSphereExchange ) THEN
362 #else
363 IF ( overlapOnly ) THEN
364 #endif
365 CALL FILL_CS_CORNER_TR_RL( 1, .FALSE.,
366 & localTij, bi,bj, myThid )
367 ENDIF
368
369 #ifdef ALLOW_AUTODIFF_TAMC
370 # ifndef DISABLE_MULTIDIM_ADVECTION
371 CADJ STORE localTij(:,:) =
372 CADJ & comlev1_bibj_k_gad_pass, key=passkey, byte=isbyte
373 # endif
374 #endif /* ALLOW_AUTODIFF_TAMC */
375
376 IF ( advectionScheme.EQ.ENUM_UPWIND_1RST
377 & .OR. advectionScheme.EQ.ENUM_DST2 ) THEN
378 CALL GAD_DST2U1_ADV_X( bi,bj,k, advectionScheme, .TRUE.,
379 I dTtracerLev(k),uTrans,uFld,localTij,
380 O af, myThid )
381 ELSEIF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN
382 CALL GAD_FLUXLIMIT_ADV_X( bi,bj,k, .TRUE., dTtracerLev(k),
383 I uTrans, uFld, maskLocW, localTij,
384 O af, myThid )
385 ELSEIF (advectionScheme.EQ.ENUM_DST3 ) THEN
386 CALL GAD_DST3_ADV_X( bi,bj,k, .TRUE., dTtracerLev(k),
387 I uTrans, uFld, maskLocW, localTij,
388 O af, myThid )
389 ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
390 CALL GAD_DST3FL_ADV_X( bi,bj,k, .TRUE., dTtracerLev(k),
391 I uTrans, uFld, maskLocW, localTij,
392 O af, myThid )
393 #ifndef ALLOW_AUTODIFF_TAMC
394 ELSEIF (advectionScheme.EQ.ENUM_OS7MP ) THEN
395 CALL GAD_OS7MP_ADV_X( bi,bj,k, .TRUE., dTtracerLev(k),
396 I uTrans, uFld, maskLocW, localTij,
397 O af, myThid )
398 #endif
399 ELSE
400 STOP 'GAD_ADVECTION: adv. scheme incompatibale with multi-dim'
401 ENDIF
402
403 C- Internal exchange for next calculations in Y
404 IF ( overlapOnly .AND. ipass.EQ.1 ) THEN
405 CALL FILL_CS_CORNER_TR_RL( 2, .FALSE.,
406 & localTij, bi,bj, myThid )
407 ENDIF
408
409 C- Advective flux in X : done
410 ENDIF
411
412 C- Update the local tracer field where needed:
413
414 C update in overlap-Only
415 IF ( overlapOnly ) THEN
416 iMinUpd = 1-Olx+1
417 iMaxUpd = sNx+Olx-1
418 C- notes: these 2 lines below have no real effect (because recip_hFac=0
419 C in corner region) but safer to keep them.
420 IF ( W_edge ) iMinUpd = 1
421 IF ( E_edge ) iMaxUpd = sNx
422
423 IF ( S_edge ) THEN
424 DO j=1-Oly,0
425 DO i=iMinUpd,iMaxUpd
426 localTij(i,j) = localTij(i,j)
427 & -dTtracerLev(k)*recip_rhoFacC(k)
428 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
429 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
430 & *( af(i+1,j)-af(i,j)
431 & -tracer(i,j,k,bi,bj)*(uTrans(i+1,j)-uTrans(i,j))
432 & )
433 ENDDO
434 ENDDO
435 ENDIF
436 IF ( N_edge ) THEN
437 DO j=sNy+1,sNy+Oly
438 DO i=iMinUpd,iMaxUpd
439 localTij(i,j) = localTij(i,j)
440 & -dTtracerLev(k)*recip_rhoFacC(k)
441 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
442 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
443 & *( af(i+1,j)-af(i,j)
444 & -tracer(i,j,k,bi,bj)*(uTrans(i+1,j)-uTrans(i,j))
445 & )
446 ENDDO
447 ENDDO
448 ENDIF
449
450 ELSE
451 C do not only update the overlap
452 jMinUpd = 1-Oly
453 jMaxUpd = sNy+Oly
454 IF ( interiorOnly .AND. S_edge ) jMinUpd = 1
455 IF ( interiorOnly .AND. N_edge ) jMaxUpd = sNy
456 DO j=jMinUpd,jMaxUpd
457 DO i=1-Olx+1,sNx+Olx-1
458 localTij(i,j) = localTij(i,j)
459 & -dTtracerLev(k)*recip_rhoFacC(k)
460 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
461 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
462 & *( af(i+1,j)-af(i,j)
463 & -tracer(i,j,k,bi,bj)*(uTrans(i+1,j)-uTrans(i,j))
464 & )
465 ENDDO
466 ENDDO
467 C- keep advective flux (for diagnostics)
468 DO j=1-Oly,sNy+Oly
469 DO i=1-Olx,sNx+Olx
470 afx(i,j) = af(i,j)
471 ENDDO
472 ENDDO
473
474 #ifdef ALLOW_OBCS
475 C- Apply open boundary conditions
476 IF ( useOBCS ) THEN
477 IF (tracerIdentity.EQ.GAD_TEMPERATURE) THEN
478 CALL OBCS_APPLY_TLOC( bi, bj, k, localTij, myThid )
479 ELSEIF (tracerIdentity.EQ.GAD_SALINITY) THEN
480 CALL OBCS_APPLY_SLOC( bi, bj, k, localTij, myThid )
481 #ifdef ALLOW_PTRACERS
482 ELSEIF (tracerIdentity.GE.GAD_TR1) THEN
483 CALL OBCS_APPLY_PTRACER( bi, bj, k,
484 & tracerIdentity-GAD_TR1+1, localTij, myThid )
485 #endif /* ALLOW_PTRACERS */
486 ENDIF
487 ENDIF
488 #endif /* ALLOW_OBCS */
489
490 C- end if/else update overlap-Only
491 ENDIF
492
493 C-- End of X direction
494 ENDIF
495
496 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
497 C-- Y direction
498 cph-test
499 C- Advective flux in Y
500 DO j=1-Oly,sNy+Oly
501 DO i=1-Olx,sNx+Olx
502 af(i,j) = 0.
503 ENDDO
504 ENDDO
505 C
506 #ifdef ALLOW_AUTODIFF_TAMC
507 # ifndef DISABLE_MULTIDIM_ADVECTION
508 CADJ STORE localTij(:,:) =
509 CADJ & comlev1_bibj_k_gad_pass, key=passkey, byte=isbyte
510 CADJ STORE af(:,:) =
511 CADJ & comlev1_bibj_k_gad_pass, key=passkey, byte=isbyte
512 # endif
513 #endif /* ALLOW_AUTODIFF_TAMC */
514 C
515 IF (calc_fluxes_Y) THEN
516
517 C- Do not compute fluxes if
518 C a) needed in overlap only
519 C and b) the overlap of myTile are not cube-face edges
520 IF ( .NOT.overlapOnly .OR. E_edge .OR. W_edge ) THEN
521
522 C- Internal exchange for calculations in Y
523 #ifdef MULTIDIM_OLD_VERSION
524 IF ( useCubedSphereExchange ) THEN
525 #else
526 IF ( overlapOnly ) THEN
527 #endif
528 CALL FILL_CS_CORNER_TR_RL( 2, .FALSE.,
529 & localTij, bi,bj, myThid )
530 ENDIF
531
532 C- Advective flux in Y
533 DO j=1-Oly,sNy+Oly
534 DO i=1-Olx,sNx+Olx
535 af(i,j) = 0.
536 ENDDO
537 ENDDO
538
539 #ifdef ALLOW_AUTODIFF_TAMC
540 #ifndef DISABLE_MULTIDIM_ADVECTION
541 CADJ STORE localTij(:,:) =
542 CADJ & comlev1_bibj_k_gad_pass, key=passkey, byte=isbyte
543 #endif
544 #endif /* ALLOW_AUTODIFF_TAMC */
545
546 IF ( advectionScheme.EQ.ENUM_UPWIND_1RST
547 & .OR. advectionScheme.EQ.ENUM_DST2 ) THEN
548 CALL GAD_DST2U1_ADV_Y( bi,bj,k, advectionScheme, .TRUE.,
549 I dTtracerLev(k),vTrans,vFld,localTij,
550 O af, myThid )
551 ELSEIF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN
552 CALL GAD_FLUXLIMIT_ADV_Y( bi,bj,k, .TRUE., dTtracerLev(k),
553 I vTrans, vFld, maskLocS, localTij,
554 O af, myThid )
555 ELSEIF (advectionScheme.EQ.ENUM_DST3 ) THEN
556 CALL GAD_DST3_ADV_Y( bi,bj,k, .TRUE., dTtracerLev(k),
557 I vTrans, vFld, maskLocS, localTij,
558 O af, myThid )
559 ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
560 CALL GAD_DST3FL_ADV_Y( bi,bj,k, .TRUE., dTtracerLev(k),
561 I vTrans, vFld, maskLocS, localTij,
562 O af, myThid )
563 #ifndef ALLOW_AUTODIFF_TAMC
564 ELSEIF (advectionScheme.EQ.ENUM_OS7MP ) THEN
565 CALL GAD_OS7MP_ADV_Y( bi,bj,k, .TRUE., dTtracerLev(k),
566 I vTrans, vFld, maskLocS, localTij,
567 O af, myThid )
568 #endif
569 ELSE
570 STOP 'GAD_ADVECTION: adv. scheme incompatibale with mutli-dim'
571 ENDIF
572
573 C- Internal exchange for next calculations in X
574 IF ( overlapOnly .AND. ipass.EQ.1 ) THEN
575 CALL FILL_CS_CORNER_TR_RL( 1, .FALSE.,
576 & localTij, bi,bj, myThid )
577 ENDIF
578
579 C- Advective flux in Y : done
580 ENDIF
581
582 C- Update the local tracer field where needed:
583
584 C update in overlap-Only
585 IF ( overlapOnly ) THEN
586 jMinUpd = 1-Oly+1
587 jMaxUpd = sNy+Oly-1
588 C- notes: these 2 lines below have no real effect (because recip_hFac=0
589 C in corner region) but safer to keep them.
590 IF ( S_edge ) jMinUpd = 1
591 IF ( N_edge ) jMaxUpd = sNy
592
593 IF ( W_edge ) THEN
594 DO j=jMinUpd,jMaxUpd
595 DO i=1-Olx,0
596 localTij(i,j) = localTij(i,j)
597 & -dTtracerLev(k)*recip_rhoFacC(k)
598 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
599 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
600 & *( af(i,j+1)-af(i,j)
601 & -tracer(i,j,k,bi,bj)*(vTrans(i,j+1)-vTrans(i,j))
602 & )
603 ENDDO
604 ENDDO
605 ENDIF
606 IF ( E_edge ) THEN
607 DO j=jMinUpd,jMaxUpd
608 DO i=sNx+1,sNx+Olx
609 localTij(i,j) = localTij(i,j)
610 & -dTtracerLev(k)*recip_rhoFacC(k)
611 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
612 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
613 & *( af(i,j+1)-af(i,j)
614 & -tracer(i,j,k,bi,bj)*(vTrans(i,j+1)-vTrans(i,j))
615 & )
616 ENDDO
617 ENDDO
618 ENDIF
619
620 ELSE
621 C do not only update the overlap
622 iMinUpd = 1-Olx
623 iMaxUpd = sNx+Olx
624 IF ( interiorOnly .AND. W_edge ) iMinUpd = 1
625 IF ( interiorOnly .AND. E_edge ) iMaxUpd = sNx
626 DO j=1-Oly+1,sNy+Oly-1
627 DO i=iMinUpd,iMaxUpd
628 localTij(i,j) = localTij(i,j)
629 & -dTtracerLev(k)*recip_rhoFacC(k)
630 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
631 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
632 & *( af(i,j+1)-af(i,j)
633 & -tracer(i,j,k,bi,bj)*(vTrans(i,j+1)-vTrans(i,j))
634 & )
635 ENDDO
636 ENDDO
637 C- keep advective flux (for diagnostics)
638 DO j=1-Oly,sNy+Oly
639 DO i=1-Olx,sNx+Olx
640 afy(i,j) = af(i,j)
641 ENDDO
642 ENDDO
643
644 #ifdef ALLOW_OBCS
645 C- Apply open boundary conditions
646 IF (useOBCS) THEN
647 IF (tracerIdentity.EQ.GAD_TEMPERATURE) THEN
648 CALL OBCS_APPLY_TLOC( bi, bj, k, localTij, myThid )
649 ELSEIF (tracerIdentity.EQ.GAD_SALINITY) THEN
650 CALL OBCS_APPLY_SLOC( bi, bj, k, localTij, myThid )
651 #ifdef ALLOW_PTRACERS
652 ELSEIF (tracerIdentity.GE.GAD_TR1) THEN
653 CALL OBCS_APPLY_PTRACER( bi, bj, k,
654 & tracerIdentity-GAD_TR1+1, localTij, myThid )
655 #endif /* ALLOW_PTRACERS */
656 ENDIF
657 ENDIF
658 #endif /* ALLOW_OBCS */
659
660 C end if/else update overlap-Only
661 ENDIF
662
663 C-- End of Y direction
664 ENDIF
665
666 C-- End of ipass loop
667 ENDDO
668
669 IF ( implicitAdvection ) THEN
670 C- explicit advection is done ; store tendency in gTracer:
671 DO j=1-Oly,sNy+Oly
672 DO i=1-Olx,sNx+Olx
673 gTracer(i,j,k,bi,bj)=
674 & (localTij(i,j)-tracer(i,j,k,bi,bj))/dTtracerLev(k)
675 ENDDO
676 ENDDO
677 ELSE
678 C- horizontal advection done; store intermediate result in 3D array:
679 DO j=1-Oly,sNy+Oly
680 DO i=1-Olx,sNx+Olx
681 localTijk(i,j,k)=localTij(i,j)
682 ENDDO
683 ENDDO
684 ENDIF
685
686 #ifdef ALLOW_DIAGNOSTICS
687 IF ( doDiagAdvX ) THEN
688 diagName = 'ADVx'//diagSufx
689 CALL DIAGNOSTICS_FILL(afx,diagName, k,1, 2,bi,bj, myThid)
690 ENDIF
691 IF ( doDiagAdvY ) THEN
692 diagName = 'ADVy'//diagSufx
693 CALL DIAGNOSTICS_FILL(afy,diagName, k,1, 2,bi,bj, myThid)
694 ENDIF
695 #endif
696
697 #ifdef ALLOW_DEBUG
698 IF ( debugLevel .GE. debLevB
699 & .AND. tracerIdentity.EQ.GAD_TEMPERATURE
700 & .AND. k.LE.3 .AND. myIter.EQ.1+nIter0
701 & .AND. nPx.EQ.1 .AND. nPy.EQ.1
702 & .AND. useCubedSphereExchange ) THEN
703 CALL DEBUG_CS_CORNER_UV( ' afx,afy from GAD_ADVECTION',
704 & afx,afy, k, standardMessageUnit,bi,bj,myThid )
705 ENDIF
706 #endif /* ALLOW_DEBUG */
707
708 C-- End of K loop for horizontal fluxes
709 ENDDO
710
711 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
712
713 IF ( .NOT.implicitAdvection ) THEN
714 C-- Start of k loop for vertical flux
715 DO k=Nr,1,-1
716 #ifdef ALLOW_AUTODIFF_TAMC
717 kkey = (igadkey-1)*Nr + (Nr-k+1)
718 #endif /* ALLOW_AUTODIFF_TAMC */
719 C-- kUp Cycles through 1,2 to point to w-layer above
720 C-- kDown Cycles through 2,1 to point to w-layer below
721 kUp = 1+MOD(k+1,2)
722 kDown= 1+MOD(k,2)
723 c kp1=min(Nr,k+1)
724 kp1Msk=1.
725 if (k.EQ.Nr) kp1Msk=0.
726
727 #ifdef ALLOW_AUTODIFF_TAMC
728 CADJ STORE rtrans(:,:) =
729 CADJ & comlev1_bibj_k_gad, key=kkey, byte=isbyte
730 cphCADJ STORE wfld(:,:) =
731 cphCADJ & comlev1_bibj_k_gad, key=kkey, byte=isbyte
732 #endif
733
734 C-- Compute Vertical transport
735 #ifdef ALLOW_AIM
736 C- a hack to prevent Water-Vapor vert.transport into the stratospheric level Nr
737 IF ( k.EQ.1 .OR.
738 & (useAIM .AND. tracerIdentity.EQ.GAD_SALINITY .AND. k.EQ.Nr)
739 & ) THEN
740 #else
741 IF ( k.EQ.1 ) THEN
742 #endif
743
744 #ifdef ALLOW_AUTODIFF_TAMC
745 cphmultiCADJ STORE wfld(:,:) =
746 cphmultiCADJ & comlev1_bibj_k_gad, key=kkey, byte=isbyte
747 #endif /* ALLOW_AUTODIFF_TAMC */
748
749 C- Surface interface :
750 DO j=1-Oly,sNy+Oly
751 DO i=1-Olx,sNx+Olx
752 rTransKp1(i,j) = kp1Msk*rTrans(i,j)
753 wFld(i,j) = 0.
754 rTrans(i,j) = 0.
755 fVerT(i,j,kUp) = 0.
756 ENDDO
757 ENDDO
758
759 ELSE
760
761 #ifdef ALLOW_AUTODIFF_TAMC
762 cphmultiCADJ STORE wfld(:,:) =
763 cphmultiCADJ & comlev1_bibj_k_gad, key=kkey, byte=isbyte
764 #endif /* ALLOW_AUTODIFF_TAMC */
765
766 C- Interior interface :
767 DO j=1-Oly,sNy+Oly
768 DO i=1-Olx,sNx+Olx
769 rTransKp1(i,j) = kp1Msk*rTrans(i,j)
770 wFld(i,j) = wVel(i,j,k,bi,bj)
771 rTrans(i,j) = wVel(i,j,k,bi,bj)*rA(i,j,bi,bj)
772 & *deepFac2F(k)*rhoFacF(k)
773 & *maskC(i,j,k-1,bi,bj)
774 fVerT(i,j,kUp) = 0.
775 ENDDO
776 ENDDO
777
778 #ifdef ALLOW_GMREDI
779 C-- Residual transp = Bolus transp + Eulerian transp
780 IF (useGMRedi)
781 & CALL GMREDI_CALC_WFLOW(
782 U wFld, rTrans,
783 I k, bi, bj, myThid )
784 #endif /* ALLOW_GMREDI */
785
786 #ifdef ALLOW_AUTODIFF_TAMC
787 cphmultiCADJ STORE localTijk(:,:,k)
788 cphmultiCADJ & = comlev1_bibj_k_gad, key=kkey, byte=isbyte
789 cphmultiCADJ STORE rTrans(:,:)
790 cphmultiCADJ & = comlev1_bibj_k_gad, key=kkey, byte=isbyte
791 #endif /* ALLOW_AUTODIFF_TAMC */
792
793 C- Compute vertical advective flux in the interior:
794 IF ( vertAdvecScheme.EQ.ENUM_UPWIND_1RST
795 & .OR. vertAdvecScheme.EQ.ENUM_DST2 ) THEN
796 CALL GAD_DST2U1_ADV_R( bi,bj,k, advectionScheme,
797 I dTtracerLev(k),rTrans,wFld,localTijk,
798 O fVerT(1-Olx,1-Oly,kUp), myThid )
799 ELSEIF( vertAdvecScheme.EQ.ENUM_FLUX_LIMIT) THEN
800 CALL GAD_FLUXLIMIT_ADV_R( bi,bj,k, dTtracerLev(k),
801 I rTrans, wFld, localTijk,
802 O fVerT(1-Olx,1-Oly,kUp), myThid )
803 ELSEIF( vertAdvecScheme.EQ.ENUM_DST3 ) THEN
804 CALL GAD_DST3_ADV_R( bi,bj,k, dTtracerLev(k),
805 I rTrans, wFld, localTijk,
806 O fVerT(1-Olx,1-Oly,kUp), myThid )
807 ELSEIF( vertAdvecScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
808 CALL GAD_DST3FL_ADV_R( bi,bj,k, dTtracerLev(k),
809 I rTrans, wFld, localTijk,
810 O fVerT(1-Olx,1-Oly,kUp), myThid )
811 #ifndef ALLOW_AUTODIFF_TAMC
812 ELSEIF (vertAdvecScheme.EQ.ENUM_OS7MP ) THEN
813 CALL GAD_OS7MP_ADV_R( bi,bj,k, dTtracerLev(k),
814 I rTrans, wFld, localTijk,
815 O fVerT(1-Olx,1-Oly,kUp), myThid )
816 #endif
817 ELSE
818 STOP 'GAD_ADVECTION: adv. scheme incompatibale with mutli-dim'
819 ENDIF
820
821 C- end Surface/Interior if bloc
822 ENDIF
823
824 #ifdef ALLOW_AUTODIFF_TAMC
825 cphmultiCADJ STORE rTrans(:,:)
826 cphmultiCADJ & = comlev1_bibj_k_gad, key=kkey, byte=isbyte
827 cphmultiCADJ STORE rTranskp1(:,:)
828 cphmultiCADJ & = comlev1_bibj_k_gad, key=kkey, byte=isbyte
829 cph --- following storing of fVerT is critical for correct
830 cph --- gradient with multiDimAdvection
831 cph --- Without it, kDown component is not properly recomputed
832 cph --- This is a TAF bug (and no warning available)
833 CADJ STORE fVerT(:,:,:)
834 CADJ & = comlev1_bibj_k_gad, key=kkey, byte=isbyte
835 #endif /* ALLOW_AUTODIFF_TAMC */
836
837 C-- Divergence of vertical fluxes
838 DO j=1-Oly,sNy+Oly
839 DO i=1-Olx,sNx+Olx
840 localTij(i,j) = localTijk(i,j,k)
841 & -dTtracerLev(k)*recip_rhoFacC(k)
842 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
843 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
844 & *( fVerT(i,j,kDown)-fVerT(i,j,kUp)
845 & -tracer(i,j,k,bi,bj)*(rTransKp1(i,j)-rTrans(i,j))
846 & )*rkSign
847 gTracer(i,j,k,bi,bj)=
848 & (localTij(i,j)-tracer(i,j,k,bi,bj))/dTtracerLev(k)
849 ENDDO
850 ENDDO
851
852 #ifdef ALLOW_DIAGNOSTICS
853 IF ( doDiagAdvR ) THEN
854 diagName = 'ADVr'//diagSufx
855 CALL DIAGNOSTICS_FILL( fVerT(1-Olx,1-Oly,kUp),
856 & diagName, k,1, 2,bi,bj, myThid)
857 ENDIF
858 #endif
859
860 C-- End of K loop for vertical flux
861 ENDDO
862 C-- end of if not.implicitAdvection block
863 ENDIF
864
865 RETURN
866 END

  ViewVC Help
Powered by ViewVC 1.1.22