/[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.53 - (show annotations) (download)
Fri Oct 19 19:18:01 2007 UTC (16 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint59m, checkpoint59l, checkpoint59i, checkpoint59k, checkpoint59j
Changes since 1.52: +7 -1 lines
Final fix for multiDimAdvection adjoint (fVerT(:,:.kDown))

1 C $Header: /u/gcmpack/MITgcm/pkg/generic_advdiff/gad_advection.F,v 1.52 2007/10/19 15:52:08 heimbach 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, :: loop range for called routines
91 C jMin,jMax :: loop range for called routines
92 C [iMin,iMax]Upd :: loop range to update tracer field
93 C [jMin,jMax]Upd :: loop range to update tracer field
94 C i,j,k :: loop indices
95 C kUp :: index into 2 1/2D array, toggles between 1 and 2
96 C kDown :: index into 2 1/2D array, toggles between 2 and 1
97 C kp1 :: =k+1 for k<Nr, =Nr for k=Nr
98 C xA,yA :: areas of X and Y face of tracer cells
99 C uFld,vFld :: 2-D local copy of horizontal velocity, U,V components
100 C wFld :: 2-D local copy of vertical velocity
101 C uTrans,vTrans :: 2-D arrays of volume transports at U,V points
102 C rTrans :: 2-D arrays of volume transports at W points
103 C rTransKp1 :: vertical volume transport at interface k+1
104 C af :: 2-D array for horizontal advective flux
105 C afx :: 2-D array for horizontal advective flux, x direction
106 C afy :: 2-D array for horizontal advective flux, y direction
107 C fVerT :: 2 1/2D arrays for vertical advective flux
108 C localTij :: 2-D array, temporary local copy of tracer fld
109 C localTijk :: 3-D array, temporary local copy of tracer fld
110 C kp1Msk :: flag (0,1) for over-riding mask for W levels
111 C calc_fluxes_X :: logical to indicate to calculate fluxes in X dir
112 C calc_fluxes_Y :: logical to indicate to calculate fluxes in Y dir
113 C interiorOnly :: only update the interior of myTile, but not the edges
114 C overlapOnly :: only update the edges of myTile, but not the interior
115 C nipass :: number of passes in multi-dimensional method
116 C ipass :: number of the current pass being made
117 C myTile :: variables used to determine which cube face
118 C nCFace :: owns a tile for cube grid runs using
119 C :: multi-dim advection.
120 C [N,S,E,W]_edge :: true if N,S,E,W edge of myTile is an Edge of the cube
121 c _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
122 _RS maskLocW(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
123 _RS maskLocS(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
124 INTEGER iMin,iMax,jMin,jMax
125 INTEGER iMinUpd,iMaxUpd,jMinUpd,jMaxUpd
126 INTEGER i,j,k,kUp,kDown
127 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
128 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
129 _RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
130 _RL vFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
131 _RL wFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
132 _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
133 _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
134 _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
135 _RL rTransKp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
136 _RL af (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
137 _RL afx (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
138 _RL afy (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
139 _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
140 _RL localTij(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
141 _RL localTijk(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
142 _RL kp1Msk
143 LOGICAL calc_fluxes_X, calc_fluxes_Y, withSigns
144 LOGICAL interiorOnly, overlapOnly
145 INTEGER nipass,ipass
146 INTEGER nCFace
147 LOGICAL N_edge, S_edge, E_edge, W_edge
148 #ifdef ALLOW_EXCH2
149 INTEGER myTile
150 #endif
151 #ifdef ALLOW_DIAGNOSTICS
152 CHARACTER*8 diagName
153 CHARACTER*4 GAD_DIAG_SUFX, diagSufx
154 EXTERNAL GAD_DIAG_SUFX
155 #endif
156 CEOP
157
158 #ifdef ALLOW_AUTODIFF_TAMC
159 act0 = tracerIdentity
160 max0 = maxpass
161 act1 = bi - myBxLo(myThid)
162 max1 = myBxHi(myThid) - myBxLo(myThid) + 1
163 act2 = bj - myByLo(myThid)
164 max2 = myByHi(myThid) - myByLo(myThid) + 1
165 act3 = myThid - 1
166 max3 = nTx*nTy
167 act4 = ikey_dynamics - 1
168 igadkey = act0
169 & + act1*max0
170 & + act2*max0*max1
171 & + act3*max0*max1*max2
172 & + act4*max0*max1*max2*max3
173 if (tracerIdentity.GT.maxpass) then
174 print *, 'ph-pass gad_advection ', maxpass, tracerIdentity
175 STOP 'maxpass seems smaller than tracerIdentity'
176 endif
177 #endif /* ALLOW_AUTODIFF_TAMC */
178
179 #ifdef ALLOW_DIAGNOSTICS
180 C-- Set diagnostic suffix for the current tracer
181 IF ( useDiagnostics ) THEN
182 diagSufx = GAD_DIAG_SUFX( tracerIdentity, myThid )
183 ENDIF
184 #endif
185
186 C-- Set up work arrays with valid (i.e. not NaN) values
187 C These inital values do not alter the numerical results. They
188 C just ensure that all memory references are to valid floating
189 C point numbers. This prevents spurious hardware signals due to
190 C uninitialised but inert locations.
191 DO j=1-OLy,sNy+OLy
192 DO i=1-OLx,sNx+OLx
193 xA(i,j) = 0. _d 0
194 yA(i,j) = 0. _d 0
195 uTrans(i,j) = 0. _d 0
196 vTrans(i,j) = 0. _d 0
197 rTrans(i,j) = 0. _d 0
198 fVerT(i,j,1) = 0. _d 0
199 fVerT(i,j,2) = 0. _d 0
200 rTransKp1(i,j)= 0. _d 0
201 #ifdef ALLOW_AUTODIFF_TAMC
202 localTij(i,j) = 0. _d 0
203 wfld(i,j) = 0. _d 0
204 #endif
205 ENDDO
206 ENDDO
207
208 C-- Set tile-specific parameters for horizontal fluxes
209 IF (useCubedSphereExchange) THEN
210 nipass=3
211 #ifdef ALLOW_AUTODIFF_TAMC
212 IF ( nipass.GT.maxcube ) STOP 'maxcube needs to be = 3'
213 #endif
214 #ifdef ALLOW_EXCH2
215 myTile = W2_myTileList(bi)
216 nCFace = exch2_myFace(myTile)
217 N_edge = exch2_isNedge(myTile).EQ.1
218 S_edge = exch2_isSedge(myTile).EQ.1
219 E_edge = exch2_isEedge(myTile).EQ.1
220 W_edge = exch2_isWedge(myTile).EQ.1
221 #else
222 nCFace = bi
223 N_edge = .TRUE.
224 S_edge = .TRUE.
225 E_edge = .TRUE.
226 W_edge = .TRUE.
227 #endif
228 ELSE
229 nipass=2
230 nCFace = bi
231 N_edge = .FALSE.
232 S_edge = .FALSE.
233 E_edge = .FALSE.
234 W_edge = .FALSE.
235 ENDIF
236
237 iMin = 1-OLx
238 iMax = sNx+OLx
239 jMin = 1-OLy
240 jMax = sNy+OLy
241
242 C-- Start of k loop for horizontal fluxes
243 DO k=1,Nr
244 #ifdef ALLOW_AUTODIFF_TAMC
245 kkey = (igadkey-1)*Nr + k
246 CADJ STORE tracer(:,:,k,bi,bj) =
247 CADJ & comlev1_bibj_k_gad, key=kkey, byte=isbyte
248 #endif /* ALLOW_AUTODIFF_TAMC */
249
250 C-- Get temporary terms used by tendency routines
251 CALL CALC_COMMON_FACTORS (
252 I uVel, vVel,
253 O uFld, vFld, uTrans, vTrans, xA, yA,
254 I k,bi,bj, myThid )
255
256 #ifdef ALLOW_GMREDI
257 C-- Residual transp = Bolus transp + Eulerian transp
258 IF (useGMRedi)
259 & CALL GMREDI_CALC_UVFLOW(
260 U uFld, vFld, uTrans, vTrans,
261 I k, bi, bj, myThid )
262 #endif /* ALLOW_GMREDI */
263
264 C-- Make local copy of tracer array and mask West & South
265 DO j=1-OLy,sNy+OLy
266 DO i=1-OLx,sNx+OLx
267 localTij(i,j)=tracer(i,j,k,bi,bj)
268 maskLocW(i,j)=maskW(i,j,k,bi,bj)
269 maskLocS(i,j)=maskS(i,j,k,bi,bj)
270 ENDDO
271 ENDDO
272
273 cph-exch2#ifndef ALLOW_AUTODIFF_TAMC
274 IF (useCubedSphereExchange) THEN
275 withSigns = .FALSE.
276 CALL FILL_CS_CORNER_UV_RS(
277 & withSigns, maskLocW,maskLocS, bi,bj, myThid )
278 ENDIF
279 cph-exch2#endif
280
281 C-- Multiple passes for different directions on different tiles
282 C-- For cube need one pass for each of red, green and blue axes.
283 DO ipass=1,nipass
284 #ifdef ALLOW_AUTODIFF_TAMC
285 passkey = ipass
286 & + (k-1) *maxpass
287 & + (igadkey-1)*maxpass*Nr
288 IF (nipass .GT. maxpass) THEN
289 STOP 'GAD_ADVECTION: nipass > maxcube. check tamc.h'
290 ENDIF
291 #endif /* ALLOW_AUTODIFF_TAMC */
292
293 interiorOnly = .FALSE.
294 overlapOnly = .FALSE.
295 IF (useCubedSphereExchange) THEN
296 #ifdef MULTIDIM_OLD_VERSION
297 C- CubedSphere : pass 3 times, with full update of local tracer field
298 IF (ipass.EQ.1) THEN
299 calc_fluxes_X = nCFace.EQ.1 .OR. nCFace.EQ.2
300 calc_fluxes_Y = nCFace.EQ.4 .OR. nCFace.EQ.5
301 ELSEIF (ipass.EQ.2) THEN
302 calc_fluxes_X = nCFace.EQ.3 .OR. nCFace.EQ.4
303 calc_fluxes_Y = nCFace.EQ.6 .OR. nCFace.EQ.1
304 #else /* MULTIDIM_OLD_VERSION */
305 C- CubedSphere : pass 3 times, with partial update of local tracer field
306 IF (ipass.EQ.1) THEN
307 overlapOnly = MOD(nCFace,3).EQ.0
308 interiorOnly = MOD(nCFace,3).NE.0
309 calc_fluxes_X = nCFace.EQ.6 .OR. nCFace.EQ.1 .OR. nCFace.EQ.2
310 calc_fluxes_Y = nCFace.EQ.3 .OR. nCFace.EQ.4 .OR. nCFace.EQ.5
311 ELSEIF (ipass.EQ.2) THEN
312 overlapOnly = MOD(nCFace,3).EQ.2
313 calc_fluxes_X = nCFace.EQ.2 .OR. nCFace.EQ.3 .OR. nCFace.EQ.4
314 calc_fluxes_Y = nCFace.EQ.5 .OR. nCFace.EQ.6 .OR. nCFace.EQ.1
315 #endif /* MULTIDIM_OLD_VERSION */
316 ELSE
317 calc_fluxes_X = nCFace.EQ.5 .OR. nCFace.EQ.6
318 calc_fluxes_Y = nCFace.EQ.2 .OR. nCFace.EQ.3
319 ENDIF
320 ELSE
321 C- not CubedSphere
322 calc_fluxes_X = MOD(ipass,2).EQ.1
323 calc_fluxes_Y = .NOT.calc_fluxes_X
324 ENDIF
325
326 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
327 C-- X direction
328 C- Advective flux in X
329 DO j=1-Oly,sNy+Oly
330 DO i=1-Olx,sNx+Olx
331 af(i,j) = 0.
332 ENDDO
333 ENDDO
334 C
335 #ifdef ALLOW_AUTODIFF_TAMC
336 # ifndef DISABLE_MULTIDIM_ADVECTION
337 CADJ STORE localTij(:,:) =
338 CADJ & comlev1_bibj_k_gad_pass, key=passkey, byte=isbyte
339 CADJ STORE af(:,:) =
340 CADJ & comlev1_bibj_k_gad_pass, key=passkey, byte=isbyte
341 # endif
342 #endif /* ALLOW_AUTODIFF_TAMC */
343 C
344 IF (calc_fluxes_X) THEN
345
346 C- Do not compute fluxes if
347 C a) needed in overlap only
348 C and b) the overlap of myTile are not cube-face Edges
349 IF ( .NOT.overlapOnly .OR. N_edge .OR. S_edge ) THEN
350
351 cph-exch2#ifndef ALLOW_AUTODIFF_TAMC
352 C- Internal exchange for calculations in X
353 #ifdef MULTIDIM_OLD_VERSION
354 IF ( useCubedSphereExchange ) THEN
355 #else
356 IF ( useCubedSphereExchange .AND.
357 & ( overlapOnly .OR. ipass.EQ.1 ) ) THEN
358 #endif
359 CALL FILL_CS_CORNER_TR_RL( .TRUE., .FALSE.,
360 & localTij, bi,bj, myThid )
361 ENDIF
362 cph-exch2#endif
363
364 #ifdef ALLOW_AUTODIFF_TAMC
365 # ifndef DISABLE_MULTIDIM_ADVECTION
366 CADJ STORE localTij(:,:) =
367 CADJ & comlev1_bibj_k_gad_pass, key=passkey, byte=isbyte
368 # endif
369 #endif /* ALLOW_AUTODIFF_TAMC */
370
371 IF ( advectionScheme.EQ.ENUM_UPWIND_1RST
372 & .OR. advectionScheme.EQ.ENUM_DST2 ) THEN
373 CALL GAD_DST2U1_ADV_X( bi,bj,k, advectionScheme, .TRUE.,
374 I dTtracerLev(k),uTrans,uFld,localTij,
375 O af, myThid )
376 ELSEIF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN
377 CALL GAD_FLUXLIMIT_ADV_X( bi,bj,k, .TRUE., dTtracerLev(k),
378 I uTrans, uFld, maskLocW, localTij,
379 O af, myThid )
380 ELSEIF (advectionScheme.EQ.ENUM_DST3 ) THEN
381 CALL GAD_DST3_ADV_X( bi,bj,k, .TRUE., dTtracerLev(k),
382 I uTrans, uFld, maskLocW, localTij,
383 O af, myThid )
384 ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
385 CALL GAD_DST3FL_ADV_X( bi,bj,k, .TRUE., dTtracerLev(k),
386 I uTrans, uFld, maskLocW, localTij,
387 O af, myThid )
388 #ifndef ALLOW_AUTODIFF_TAMC
389 ELSEIF (advectionScheme.EQ.ENUM_OS7MP ) THEN
390 CALL GAD_OS7MP_ADV_X( bi,bj,k, .TRUE., dTtracerLev(k),
391 I uTrans, uFld, maskLocW, localTij,
392 O af, myThid )
393 #endif
394 ELSE
395 STOP 'GAD_ADVECTION: adv. scheme incompatibale with multi-dim'
396 ENDIF
397
398 C- Advective flux in X : done
399 ENDIF
400
401 cph-exch2#ifndef ALLOW_AUTODIFF_TAMC
402 C- Internal exchange for next calculations in Y
403 IF ( overlapOnly .AND. ipass.EQ.1 ) THEN
404 CALL FILL_CS_CORNER_TR_RL(.FALSE., .FALSE.,
405 & localTij, bi,bj, myThid )
406 ENDIF
407 cph-exch2#endif
408
409 C- Update the local tracer field where needed:
410
411 C update in overlap-Only
412 IF ( overlapOnly ) THEN
413 iMinUpd = 1-Olx+1
414 iMaxUpd = sNx+Olx-1
415 C- notes: these 2 lines below have no real effect (because recip_hFac=0
416 C in corner region) but safer to keep them.
417 IF ( W_edge ) iMinUpd = 1
418 IF ( E_edge ) iMaxUpd = sNx
419
420 IF ( S_edge ) THEN
421 DO j=1-Oly,0
422 DO i=iMinUpd,iMaxUpd
423 localTij(i,j) = localTij(i,j)
424 & -dTtracerLev(k)*recip_rhoFacC(k)
425 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
426 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
427 & *( af(i+1,j)-af(i,j)
428 & -tracer(i,j,k,bi,bj)*(uTrans(i+1,j)-uTrans(i,j))
429 & )
430 ENDDO
431 ENDDO
432 ENDIF
433 IF ( N_edge ) THEN
434 DO j=sNy+1,sNy+Oly
435 DO i=iMinUpd,iMaxUpd
436 localTij(i,j) = localTij(i,j)
437 & -dTtracerLev(k)*recip_rhoFacC(k)
438 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
439 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
440 & *( af(i+1,j)-af(i,j)
441 & -tracer(i,j,k,bi,bj)*(uTrans(i+1,j)-uTrans(i,j))
442 & )
443 ENDDO
444 ENDDO
445 ENDIF
446
447 ELSE
448 C do not only update the overlap
449 jMinUpd = 1-Oly
450 jMaxUpd = sNy+Oly
451 IF ( interiorOnly .AND. S_edge ) jMinUpd = 1
452 IF ( interiorOnly .AND. N_edge ) jMaxUpd = sNy
453 DO j=jMinUpd,jMaxUpd
454 DO i=1-Olx+1,sNx+Olx-1
455 localTij(i,j) = localTij(i,j)
456 & -dTtracerLev(k)*recip_rhoFacC(k)
457 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
458 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
459 & *( af(i+1,j)-af(i,j)
460 & -tracer(i,j,k,bi,bj)*(uTrans(i+1,j)-uTrans(i,j))
461 & )
462 ENDDO
463 ENDDO
464 C- keep advective flux (for diagnostics)
465 DO j=1-Oly,sNy+Oly
466 DO i=1-Olx,sNx+Olx
467 afx(i,j) = af(i,j)
468 ENDDO
469 ENDDO
470
471 #ifdef ALLOW_OBCS
472 C- Apply open boundary conditions
473 IF ( useOBCS ) THEN
474 IF (tracerIdentity.EQ.GAD_TEMPERATURE) THEN
475 CALL OBCS_APPLY_TLOC( bi, bj, k, localTij, myThid )
476 ELSEIF (tracerIdentity.EQ.GAD_SALINITY) THEN
477 CALL OBCS_APPLY_SLOC( bi, bj, k, localTij, myThid )
478 #ifdef ALLOW_PTRACERS
479 ELSEIF (tracerIdentity.GE.GAD_TR1) THEN
480 CALL OBCS_APPLY_PTRACER( bi, bj, k,
481 & tracerIdentity-GAD_TR1+1, localTij, myThid )
482 #endif /* ALLOW_PTRACERS */
483 ENDIF
484 ENDIF
485 #endif /* ALLOW_OBCS */
486
487 C- end if/else update overlap-Only
488 ENDIF
489
490 C-- End of X direction
491 ENDIF
492
493 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
494 C-- Y direction
495 cph-test
496 C- Advective flux in Y
497 DO j=1-Oly,sNy+Oly
498 DO i=1-Olx,sNx+Olx
499 af(i,j) = 0.
500 ENDDO
501 ENDDO
502 C
503 #ifdef ALLOW_AUTODIFF_TAMC
504 # ifndef DISABLE_MULTIDIM_ADVECTION
505 CADJ STORE localTij(:,:) =
506 CADJ & comlev1_bibj_k_gad_pass, key=passkey, byte=isbyte
507 CADJ STORE af(:,:) =
508 CADJ & comlev1_bibj_k_gad_pass, key=passkey, byte=isbyte
509 # endif
510 #endif /* ALLOW_AUTODIFF_TAMC */
511 C
512 IF (calc_fluxes_Y) THEN
513
514 C- Do not compute fluxes if
515 C a) needed in overlap only
516 C and b) the overlap of myTile are not cube-face edges
517 IF ( .NOT.overlapOnly .OR. E_edge .OR. W_edge ) THEN
518
519 cph-exch2#ifndef ALLOW_AUTODIFF_TAMC
520 C- Internal exchange for calculations in Y
521 #ifdef MULTIDIM_OLD_VERSION
522 IF ( useCubedSphereExchange ) THEN
523 #else
524 IF ( useCubedSphereExchange .AND.
525 & ( overlapOnly .OR. ipass.EQ.1 ) ) THEN
526 #endif
527 CALL FILL_CS_CORNER_TR_RL(.FALSE., .FALSE.,
528 & localTij, bi,bj, myThid )
529 ENDIF
530 cph-exch2#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- Advective flux in Y : done
574 ENDIF
575
576 cph-exch2#ifndef ALLOW_AUTODIFF_TAMC
577 C- Internal exchange for next calculations in X
578 IF ( overlapOnly .AND. ipass.EQ.1 ) THEN
579 CALL FILL_CS_CORNER_TR_RL( .TRUE., .FALSE.,
580 & localTij, bi,bj, myThid )
581 ENDIF
582 cph-exch2#endif
583
584 C- Update the local tracer field where needed:
585
586 C update in overlap-Only
587 IF ( overlapOnly ) THEN
588 jMinUpd = 1-Oly+1
589 jMaxUpd = sNy+Oly-1
590 C- notes: these 2 lines below have no real effect (because recip_hFac=0
591 C in corner region) but safer to keep them.
592 IF ( S_edge ) jMinUpd = 1
593 IF ( N_edge ) jMaxUpd = sNy
594
595 IF ( W_edge ) THEN
596 DO j=jMinUpd,jMaxUpd
597 DO i=1-Olx,0
598 localTij(i,j) = localTij(i,j)
599 & -dTtracerLev(k)*recip_rhoFacC(k)
600 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
601 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
602 & *( af(i,j+1)-af(i,j)
603 & -tracer(i,j,k,bi,bj)*(vTrans(i,j+1)-vTrans(i,j))
604 & )
605 ENDDO
606 ENDDO
607 ENDIF
608 IF ( E_edge ) THEN
609 DO j=jMinUpd,jMaxUpd
610 DO i=sNx+1,sNx+Olx
611 localTij(i,j) = localTij(i,j)
612 & -dTtracerLev(k)*recip_rhoFacC(k)
613 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
614 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
615 & *( af(i,j+1)-af(i,j)
616 & -tracer(i,j,k,bi,bj)*(vTrans(i,j+1)-vTrans(i,j))
617 & )
618 ENDDO
619 ENDDO
620 ENDIF
621
622 ELSE
623 C do not only update the overlap
624 iMinUpd = 1-Olx
625 iMaxUpd = sNx+Olx
626 IF ( interiorOnly .AND. W_edge ) iMinUpd = 1
627 IF ( interiorOnly .AND. E_edge ) iMaxUpd = sNx
628 DO j=1-Oly+1,sNy+Oly-1
629 DO i=iMinUpd,iMaxUpd
630 localTij(i,j) = localTij(i,j)
631 & -dTtracerLev(k)*recip_rhoFacC(k)
632 & *_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
633 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)
634 & *( af(i,j+1)-af(i,j)
635 & -tracer(i,j,k,bi,bj)*(vTrans(i,j+1)-vTrans(i,j))
636 & )
637 ENDDO
638 ENDDO
639 C- keep advective flux (for diagnostics)
640 DO j=1-Oly,sNy+Oly
641 DO i=1-Olx,sNx+Olx
642 afy(i,j) = af(i,j)
643 ENDDO
644 ENDDO
645
646 #ifdef ALLOW_OBCS
647 C- Apply open boundary conditions
648 IF (useOBCS) THEN
649 IF (tracerIdentity.EQ.GAD_TEMPERATURE) THEN
650 CALL OBCS_APPLY_TLOC( bi, bj, k, localTij, myThid )
651 ELSEIF (tracerIdentity.EQ.GAD_SALINITY) THEN
652 CALL OBCS_APPLY_SLOC( bi, bj, k, localTij, myThid )
653 #ifdef ALLOW_PTRACERS
654 ELSEIF (tracerIdentity.GE.GAD_TR1) THEN
655 CALL OBCS_APPLY_PTRACER( bi, bj, k,
656 & tracerIdentity-GAD_TR1+1, localTij, myThid )
657 #endif /* ALLOW_PTRACERS */
658 ENDIF
659 ENDIF
660 #endif /* ALLOW_OBCS */
661
662 C end if/else update overlap-Only
663 ENDIF
664
665 C-- End of Y direction
666 ENDIF
667
668 C-- End of ipass loop
669 ENDDO
670
671 IF ( implicitAdvection ) THEN
672 C- explicit advection is done ; store tendency in gTracer:
673 DO j=1-Oly,sNy+Oly
674 DO i=1-Olx,sNx+Olx
675 gTracer(i,j,k,bi,bj)=
676 & (localTij(i,j)-tracer(i,j,k,bi,bj))/dTtracerLev(k)
677 ENDDO
678 ENDDO
679 ELSE
680 C- horizontal advection done; store intermediate result in 3D array:
681 DO j=1-Oly,sNy+Oly
682 DO i=1-Olx,sNx+Olx
683 localTijk(i,j,k)=localTij(i,j)
684 ENDDO
685 ENDDO
686 ENDIF
687
688 #ifdef ALLOW_DIAGNOSTICS
689 IF ( useDiagnostics ) THEN
690 diagName = 'ADVx'//diagSufx
691 CALL DIAGNOSTICS_FILL(afx,diagName, k,1, 2,bi,bj, myThid)
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 ( useDiagnostics ) 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