/[MITgcm]/MITgcm/pkg/generic_advdiff/gad_calc_rhs.F
ViewVC logotype

Contents of /MITgcm/pkg/generic_advdiff/gad_calc_rhs.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.47 - (show annotations) (download)
Thu May 3 21:34:39 2007 UTC (17 years, 1 month ago) by jmc
Branch: MAIN
CVS Tags: checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59c, checkpoint59b
Changes since 1.46: +11 -8 lines
add aguments to KPP_TRANSPORT calls.

1 C $Header: /u/gcmpack/MITgcm/pkg/generic_advdiff/gad_calc_rhs.F,v 1.46 2007/04/04 01:39:06 jmc Exp $
2 C $Name: $
3
4 #include "GAD_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: GAD_CALC_RHS
8
9 C !INTERFACE: ==========================================================
10 SUBROUTINE GAD_CALC_RHS(
11 I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
12 I xA, yA, maskUp, uFld, vFld, wFld,
13 I uTrans, vTrans, rTrans, rTransKp1,
14 I diffKh, diffK4, KappaR, TracerN, TracAB,
15 I tracerIdentity, advectionScheme, vertAdvecScheme,
16 I calcAdvection, implicitAdvection, applyAB_onTracer,
17 U fVerT, gTracer,
18 I myTime, myIter, myThid )
19
20 C !DESCRIPTION:
21 C Calculates the tendency of a tracer due to advection and diffusion.
22 C It calculates the fluxes in each direction indepentently and then
23 C sets the tendency to the divergence of these fluxes. The advective
24 C fluxes are only calculated here when using the linear advection schemes
25 C otherwise only the diffusive and parameterized fluxes are calculated.
26 C
27 C Contributions to the flux are calculated and added:
28 C \begin{equation*}
29 C {\bf F} = {\bf F}_{adv} + {\bf F}_{diff} +{\bf F}_{GM} + {\bf F}_{KPP}
30 C \end{equation*}
31 C
32 C The tendency is the divergence of the fluxes:
33 C \begin{equation*}
34 C G_\theta = G_\theta + \nabla \cdot {\bf F}
35 C \end{equation*}
36 C
37 C The tendency is assumed to contain data on entry.
38
39 C !USES: ===============================================================
40 IMPLICIT NONE
41 #include "SIZE.h"
42 #include "EEPARAMS.h"
43 #include "PARAMS.h"
44 #include "GRID.h"
45 #include "SURFACE.h"
46 #include "GAD.h"
47
48 #ifdef ALLOW_AUTODIFF_TAMC
49 #include "tamc.h"
50 #include "tamc_keys.h"
51 #endif /* ALLOW_AUTODIFF_TAMC */
52
53 C !INPUT PARAMETERS: ===================================================
54 C bi,bj :: tile indices
55 C iMin,iMax :: loop range for called routines
56 C jMin,jMax :: loop range for called routines
57 C k :: vertical index
58 C kM1 :: =k-1 for k>1, =1 for k=1
59 C kUp :: index into 2 1/2D array, toggles between 1|2
60 C kDown :: index into 2 1/2D array, toggles between 2|1
61 C xA,yA :: areas of X and Y face of tracer cells
62 C maskUp :: 2-D array for mask at W points
63 C uFld,vFld,wFld :: Local copy of velocity field (3 components)
64 C uTrans,vTrans :: 2-D arrays of volume transports at U,V points
65 C rTrans :: 2-D arrays of volume transports at W points
66 C rTransKp1 :: 2-D array of volume trans at W pts, interf k+1
67 C diffKh :: horizontal diffusion coefficient
68 C diffK4 :: bi-harmonic diffusion coefficient
69 C KappaR :: 2-D array for vertical diffusion coefficient, interf k
70 C TracerN :: tracer field @ time-step n (Note: only used
71 C if applying AB on tracer field rather than on tendency gTr)
72 C TracAB :: current tracer field (@ time-step n if applying AB on gTr
73 C or extrapolated fwd in time to n+1/2 if applying AB on Tr)
74 C tracerIdentity :: tracer identifier (required for KPP,GM)
75 C advectionScheme :: advection scheme to use (Horizontal plane)
76 C vertAdvecScheme :: advection scheme to use (Vertical direction)
77 C calcAdvection :: =False if Advec computed with multiDim scheme
78 C implicitAdvection:: =True if vertical Advec computed implicitly
79 C applyAB_onTracer :: apply Adams-Bashforth on Tracer (rather than on gTr)
80 C myTime :: current time
81 C myIter :: iteration number
82 C myThid :: thread number
83 INTEGER bi,bj,iMin,iMax,jMin,jMax
84 INTEGER k,kUp,kDown,kM1
85 _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86 _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
87 _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
88 _RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
89 _RL vFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
90 _RL wFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
91 _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
92 _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
93 _RL rTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94 _RL rTransKp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95 _RL diffKh, diffK4
96 _RL KappaR(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
97 _RL TracerN(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
98 _RL TracAB (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
99 INTEGER tracerIdentity
100 INTEGER advectionScheme, vertAdvecScheme
101 LOGICAL calcAdvection
102 LOGICAL implicitAdvection, applyAB_onTracer
103 _RL myTime
104 INTEGER myIter, myThid
105
106 C !OUTPUT PARAMETERS: ==================================================
107 C gTracer :: tendency array
108 C fVerT :: 2 1/2D arrays for vertical advective flux
109 _RL gTracer(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
110 _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
111
112 C !LOCAL VARIABLES: ====================================================
113 C i,j :: loop indices
114 C df4 :: used for storing del^2 T for bi-harmonic term
115 C fZon :: zonal flux
116 C fMer :: meridional flux
117 C af :: advective flux
118 C df :: diffusive flux
119 C localT :: local copy of tracer field
120 C locABT :: local copy of (AB-extrapolated) tracer field
121 #ifdef ALLOW_DIAGNOSTICS
122 CHARACTER*8 diagName
123 CHARACTER*4 GAD_DIAG_SUFX, diagSufx
124 EXTERNAL GAD_DIAG_SUFX
125 #endif
126 INTEGER i,j
127 _RL df4 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
128 _RL fZon (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
129 _RL fMer (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
130 _RL af (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
131 _RL df (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
132 _RL localT(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
133 _RL locABT(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
134 _RL advFac, rAdvFac
135 CEOP
136
137 #ifdef ALLOW_AUTODIFF_TAMC
138 C-- only the kUp part of fverT is set in this subroutine
139 C-- the kDown is still required
140 fVerT(1,1,kDown) = fVerT(1,1,kDown)
141 #endif
142
143 #ifdef ALLOW_DIAGNOSTICS
144 C-- Set diagnostic suffix for the current tracer
145 IF ( useDiagnostics ) THEN
146 diagSufx = GAD_DIAG_SUFX( tracerIdentity, myThid )
147 ENDIF
148 #endif
149
150 advFac = 0. _d 0
151 IF (calcAdvection) advFac = 1. _d 0
152 rAdvFac = rkSign*advFac
153 IF (implicitAdvection) rAdvFac = 0. _d 0
154
155 DO j=1-OLy,sNy+OLy
156 DO i=1-OLx,sNx+OLx
157 fZon(i,j) = 0. _d 0
158 fMer(i,j) = 0. _d 0
159 fVerT(i,j,kUp) = 0. _d 0
160 df(i,j) = 0. _d 0
161 df4(i,j) = 0. _d 0
162 ENDDO
163 ENDDO
164
165 C-- Make local copy of tracer array
166 IF ( applyAB_onTracer ) THEN
167 DO j=1-OLy,sNy+OLy
168 DO i=1-OLx,sNx+OLx
169 localT(i,j)=TracerN(i,j,k,bi,bj)
170 locABT(i,j)= TracAB(i,j,k,bi,bj)
171 ENDDO
172 ENDDO
173 ELSE
174 DO j=1-OLy,sNy+OLy
175 DO i=1-OLx,sNx+OLx
176 localT(i,j)= TracAB(i,j,k,bi,bj)
177 locABT(i,j)= TracAB(i,j,k,bi,bj)
178 ENDDO
179 ENDDO
180 ENDIF
181
182 C-- Unless we have already calculated the advection terms we initialize
183 C the tendency to zero.
184 C <== now done earlier at the beginning of thermodynamics.
185 c IF (calcAdvection) THEN
186 c DO j=1-Oly,sNy+Oly
187 c DO i=1-Olx,sNx+Olx
188 c gTracer(i,j,k,bi,bj)=0. _d 0
189 c ENDDO
190 c ENDDO
191 c ENDIF
192
193 C-- Pre-calculate del^2 T if bi-harmonic coefficient is non-zero
194 IF (diffK4 .NE. 0.) THEN
195 CALL GAD_GRAD_X(bi,bj,k,xA,localT,fZon,myThid)
196 CALL GAD_GRAD_Y(bi,bj,k,yA,localT,fMer,myThid)
197 CALL GAD_DEL2(bi,bj,k,fZon,fMer,df4,myThid)
198 ENDIF
199
200 C-- Initialize net flux in X direction
201 DO j=1-Oly,sNy+Oly
202 DO i=1-Olx,sNx+Olx
203 fZon(i,j) = 0. _d 0
204 ENDDO
205 ENDDO
206
207 C- Advective flux in X
208 IF (calcAdvection) THEN
209 IF (advectionScheme.EQ.ENUM_CENTERED_2ND) THEN
210 CALL GAD_C2_ADV_X(bi,bj,k,uTrans,locABT,af,myThid)
211 ELSEIF ( advectionScheme.EQ.ENUM_UPWIND_1RST
212 & .OR. advectionScheme.EQ.ENUM_DST2 ) THEN
213 CALL GAD_DST2U1_ADV_X( bi,bj,k, advectionScheme, .TRUE.,
214 I dTtracerLev(k), uTrans, uFld, locABT,
215 O af, myThid )
216 ELSEIF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN
217 CALL GAD_FLUXLIMIT_ADV_X( bi,bj,k, .TRUE., dTtracerLev(k),
218 I uTrans, uFld, maskW(1-Olx,1-Oly,k,bi,bj), locABT,
219 O af, myThid )
220 ELSEIF (advectionScheme.EQ.ENUM_UPWIND_3RD ) THEN
221 CALL GAD_U3_ADV_X(bi,bj,k,uTrans,locABT,af,myThid)
222 ELSEIF (advectionScheme.EQ.ENUM_CENTERED_4TH) THEN
223 CALL GAD_C4_ADV_X(bi,bj,k,uTrans,locABT,af,myThid)
224 ELSEIF (advectionScheme.EQ.ENUM_DST3 ) THEN
225 CALL GAD_DST3_ADV_X( bi,bj,k, .TRUE., dTtracerLev(k),
226 I uTrans, uFld, maskW(1-Olx,1-Oly,k,bi,bj), locABT,
227 O af, myThid )
228 ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
229 IF ( inAdMode ) THEN
230 cph This block is to trick the adjoint:
231 cph IF inAdExact=.FALSE., we want to use DST3
232 cph with limiters in forward, but without limiters in reverse.
233 CALL GAD_DST3_ADV_X( bi,bj,k, .TRUE., dTtracerLev(k),
234 I uTrans, uFld, maskW(1-Olx,1-Oly,k,bi,bj), locABT,
235 O af, myThid )
236 ELSE
237 CALL GAD_DST3FL_ADV_X( bi,bj,k, .TRUE., dTtracerLev(k),
238 I uTrans, uFld, maskW(1-Olx,1-Oly,k,bi,bj), locABT,
239 O af, myThid )
240 ENDIF
241 #ifndef ALLOW_AUTODIFF_TAMC
242 ELSEIF (advectionScheme.EQ.ENUM_OS7MP ) THEN
243 CALL GAD_OS7MP_ADV_X( bi,bj,k, .TRUE., dTtracerLev(k),
244 I uTrans, uFld, maskW(1-Olx,1-Oly,k,bi,bj), locABT,
245 O af, myThid )
246 #endif
247 ELSE
248 STOP 'GAD_CALC_RHS: Bad advectionScheme (X)'
249 ENDIF
250 DO j=1-Oly,sNy+Oly
251 DO i=1-Olx,sNx+Olx
252 fZon(i,j) = fZon(i,j) + af(i,j)
253 ENDDO
254 ENDDO
255 #ifdef ALLOW_DIAGNOSTICS
256 IF ( useDiagnostics ) THEN
257 diagName = 'ADVx'//diagSufx
258 CALL DIAGNOSTICS_FILL(af,diagName, k,1, 2,bi,bj, myThid)
259 ENDIF
260 #endif
261 ENDIF
262
263 C- Diffusive flux in X
264 IF (diffKh.NE.0.) THEN
265 CALL GAD_DIFF_X(bi,bj,k,xA,diffKh,localT,df,myThid)
266 ELSE
267 DO j=1-Oly,sNy+Oly
268 DO i=1-Olx,sNx+Olx
269 df(i,j) = 0. _d 0
270 ENDDO
271 ENDDO
272 ENDIF
273
274 C- Add bi-harmonic diffusive flux in X
275 IF (diffK4 .NE. 0.) THEN
276 CALL GAD_BIHARM_X(bi,bj,k,xA,df4,diffK4,df,myThid)
277 ENDIF
278
279 #ifdef ALLOW_GMREDI
280 C- GM/Redi flux in X
281 IF (useGMRedi) THEN
282 C *note* should update GMREDI_XTRANSPORT to set df *aja*
283 IF ( applyAB_onTracer ) THEN
284 CALL GMREDI_XTRANSPORT(
285 I iMin,iMax,jMin,jMax,bi,bj,k,
286 I xA,TracerN,tracerIdentity,
287 U df,
288 I myThid)
289 ELSE
290 CALL GMREDI_XTRANSPORT(
291 I iMin,iMax,jMin,jMax,bi,bj,k,
292 I xA,TracAB, tracerIdentity,
293 U df,
294 I myThid)
295 ENDIF
296 ENDIF
297 #endif
298 C anelastic: advect.fluxes are scaled by rhoFac but hor.diff. flx are not
299 DO j=1-Oly,sNy+Oly
300 DO i=1-Olx,sNx+Olx
301 fZon(i,j) = fZon(i,j) + df(i,j)*rhoFacC(k)
302 ENDDO
303 ENDDO
304
305 #ifdef ALLOW_DIAGNOSTICS
306 C- Diagnostics of Tracer flux in X dir (mainly Diffusive term),
307 C excluding advective terms:
308 IF ( useDiagnostics .AND.
309 & (diffKh.NE.0. .OR. diffK4 .NE.0. .OR. useGMRedi) ) THEN
310 diagName = 'DFxE'//diagSufx
311 CALL DIAGNOSTICS_FILL(df,diagName, k,1, 2,bi,bj, myThid)
312 ENDIF
313 #endif
314
315 C-- Initialize net flux in Y direction
316 DO j=1-Oly,sNy+Oly
317 DO i=1-Olx,sNx+Olx
318 fMer(i,j) = 0. _d 0
319 ENDDO
320 ENDDO
321
322 C- Advective flux in Y
323 IF (calcAdvection) THEN
324 IF (advectionScheme.EQ.ENUM_CENTERED_2ND) THEN
325 CALL GAD_C2_ADV_Y(bi,bj,k,vTrans,locABT,af,myThid)
326 ELSEIF ( advectionScheme.EQ.ENUM_UPWIND_1RST
327 & .OR. advectionScheme.EQ.ENUM_DST2 ) THEN
328 CALL GAD_DST2U1_ADV_Y( bi,bj,k, advectionScheme, .TRUE.,
329 I dTtracerLev(k), vTrans, vFld, locABT,
330 O af, myThid )
331 ELSEIF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN
332 CALL GAD_FLUXLIMIT_ADV_Y( bi,bj,k, .TRUE., dTtracerLev(k),
333 I vTrans, vFld, maskS(1-Olx,1-Oly,k,bi,bj), locABT,
334 O af, myThid )
335 ELSEIF (advectionScheme.EQ.ENUM_UPWIND_3RD ) THEN
336 CALL GAD_U3_ADV_Y(bi,bj,k,vTrans,locABT,af,myThid)
337 ELSEIF (advectionScheme.EQ.ENUM_CENTERED_4TH) THEN
338 CALL GAD_C4_ADV_Y(bi,bj,k,vTrans,locABT,af,myThid)
339 ELSEIF (advectionScheme.EQ.ENUM_DST3 ) THEN
340 CALL GAD_DST3_ADV_Y( bi,bj,k, .TRUE., dTtracerLev(k),
341 I vTrans, vFld, maskS(1-Olx,1-Oly,k,bi,bj), locABT,
342 O af, myThid )
343 ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
344 IF ( inAdMode ) THEN
345 cph This block is to trick the adjoint:
346 cph IF inAdExact=.FALSE., we want to use DST3
347 cph with limiters in forward, but without limiters in reverse.
348 CALL GAD_DST3_ADV_Y( bi,bj,k, .TRUE., dTtracerLev(k),
349 I vTrans, vFld, maskS(1-Olx,1-Oly,k,bi,bj), locABT,
350 O af, myThid )
351 ELSE
352 CALL GAD_DST3FL_ADV_Y( bi,bj,k, .TRUE., dTtracerLev(k),
353 I vTrans, vFld, maskS(1-Olx,1-Oly,k,bi,bj), locABT,
354 O af, myThid )
355 ENDIF
356 #ifndef ALLOW_AUTODIFF_TAMC
357 ELSEIF (advectionScheme.EQ.ENUM_OS7MP ) THEN
358 CALL GAD_OS7MP_ADV_Y( bi,bj,k, .TRUE., dTtracerLev(k),
359 I vTrans, vFld, maskS(1-Olx,1-Oly,k,bi,bj), locABT,
360 O af, myThid )
361 #endif
362 ELSE
363 STOP 'GAD_CALC_RHS: Bad advectionScheme (Y)'
364 ENDIF
365 DO j=1-Oly,sNy+Oly
366 DO i=1-Olx,sNx+Olx
367 fMer(i,j) = fMer(i,j) + af(i,j)
368 ENDDO
369 ENDDO
370 #ifdef ALLOW_DIAGNOSTICS
371 IF ( useDiagnostics ) THEN
372 diagName = 'ADVy'//diagSufx
373 CALL DIAGNOSTICS_FILL(af,diagName, k,1, 2,bi,bj, myThid)
374 ENDIF
375 #endif
376 ENDIF
377
378 C- Diffusive flux in Y
379 IF (diffKh.NE.0.) THEN
380 CALL GAD_DIFF_Y(bi,bj,k,yA,diffKh,localT,df,myThid)
381 ELSE
382 DO j=1-Oly,sNy+Oly
383 DO i=1-Olx,sNx+Olx
384 df(i,j) = 0. _d 0
385 ENDDO
386 ENDDO
387 ENDIF
388
389 C- Add bi-harmonic flux in Y
390 IF (diffK4 .NE. 0.) THEN
391 CALL GAD_BIHARM_Y(bi,bj,k,yA,df4,diffK4,df,myThid)
392 ENDIF
393
394 #ifdef ALLOW_GMREDI
395 C- GM/Redi flux in Y
396 IF (useGMRedi) THEN
397 C *note* should update GMREDI_YTRANSPORT to set df *aja*
398 IF ( applyAB_onTracer ) THEN
399 CALL GMREDI_YTRANSPORT(
400 I iMin,iMax,jMin,jMax,bi,bj,k,
401 I yA,TracerN,tracerIdentity,
402 U df,
403 I myThid)
404 ELSE
405 CALL GMREDI_YTRANSPORT(
406 I iMin,iMax,jMin,jMax,bi,bj,k,
407 I yA,TracAB, tracerIdentity,
408 U df,
409 I myThid)
410 ENDIF
411 ENDIF
412 #endif
413 C anelastic: advect.fluxes are scaled by rhoFac but hor.diff. flx are not
414 DO j=1-Oly,sNy+Oly
415 DO i=1-Olx,sNx+Olx
416 fMer(i,j) = fMer(i,j) + df(i,j)*rhoFacC(k)
417 ENDDO
418 ENDDO
419
420 #ifdef ALLOW_DIAGNOSTICS
421 C- Diagnostics of Tracer flux in Y dir (mainly Diffusive terms),
422 C excluding advective terms:
423 IF ( useDiagnostics .AND.
424 & (diffKh.NE.0. .OR. diffK4 .NE.0. .OR. useGMRedi) ) THEN
425 diagName = 'DFyE'//diagSufx
426 CALL DIAGNOSTICS_FILL(df,diagName, k,1, 2,bi,bj, myThid)
427 ENDIF
428 #endif
429
430 C-- Compute vertical flux fVerT(kUp) at interface k (between k-1 & k):
431 C- Advective flux in R
432 #ifdef ALLOW_AIM
433 C- a hack to prevent Water-Vapor vert.transport into the stratospheric level Nr
434 IF (calcAdvection .AND. .NOT.implicitAdvection .AND. k.GE.2 .AND.
435 & (.NOT.useAIM .OR.tracerIdentity.NE.GAD_SALINITY .OR.k.LT.Nr)
436 & ) THEN
437 #else
438 IF (calcAdvection .AND. .NOT.implicitAdvection .AND. k.GE.2) THEN
439 #endif
440 C- Compute vertical advective flux in the interior:
441 IF (vertAdvecScheme.EQ.ENUM_CENTERED_2ND) THEN
442 CALL GAD_C2_ADV_R(bi,bj,k,rTrans,TracAB,af,myThid)
443 ELSEIF ( vertAdvecScheme.EQ.ENUM_UPWIND_1RST
444 & .OR. vertAdvecScheme.EQ.ENUM_DST2 ) THEN
445 CALL GAD_DST2U1_ADV_R( bi,bj,k, vertAdvecScheme,
446 I dTtracerLev(k),rTrans,wFld,TracAB(1-Olx,1-Oly,1,bi,bj),
447 O af, myThid )
448 ELSEIF (vertAdvecScheme.EQ.ENUM_FLUX_LIMIT) THEN
449 CALL GAD_FLUXLIMIT_ADV_R( bi,bj,k,
450 I dTtracerLev(k),rTrans,wFld,TracAB(1-Olx,1-Oly,1,bi,bj),
451 O af, myThid )
452 ELSEIF (vertAdvecScheme.EQ.ENUM_UPWIND_3RD ) THEN
453 CALL GAD_U3_ADV_R(bi,bj,k,rTrans,TracAB,af,myThid)
454 ELSEIF (vertAdvecScheme.EQ.ENUM_CENTERED_4TH) THEN
455 CALL GAD_C4_ADV_R(bi,bj,k,rTrans,TracAB,af,myThid)
456 ELSEIF (vertAdvecScheme.EQ.ENUM_DST3 ) THEN
457 CALL GAD_DST3_ADV_R( bi,bj,k,
458 I dTtracerLev(k),rTrans,wFld,TracAB(1-Olx,1-Oly,1,bi,bj),
459 O af, myThid )
460 ELSEIF (vertAdvecScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
461 cph This block is to trick the adjoint:
462 cph IF inAdExact=.FALSE., we want to use DST3
463 cph with limiters in forward, but without limiters in reverse.
464 IF ( inAdMode ) THEN
465 CALL GAD_DST3_ADV_R( bi,bj,k,
466 I dTtracerLev(k),rTrans,wFld,TracAB(1-Olx,1-Oly,1,bi,bj),
467 O af, myThid )
468 ELSE
469 CALL GAD_DST3FL_ADV_R( bi,bj,k,
470 I dTtracerLev(k),rTrans,wFld,TracAB(1-Olx,1-Oly,1,bi,bj),
471 O af, myThid )
472 ENDIF
473 #ifndef ALLOW_AUTODIFF_TAMC
474 ELSEIF (vertAdvecScheme.EQ.ENUM_OS7MP ) THEN
475 CALL GAD_OS7MP_ADV_R( bi,bj,k,
476 I dTtracerLev(k),rTrans,wFld,TracAB(1-Olx,1-Oly,1,bi,bj),
477 O af, myThid )
478 #endif
479 ELSE
480 STOP 'GAD_CALC_RHS: Bad vertAdvecScheme (R)'
481 ENDIF
482 C- add the advective flux to fVerT
483 DO j=1-Oly,sNy+Oly
484 DO i=1-Olx,sNx+Olx
485 fVerT(i,j,kUp) = fVerT(i,j,kUp) + af(i,j)
486 ENDDO
487 ENDDO
488 #ifdef ALLOW_DIAGNOSTICS
489 IF ( useDiagnostics ) THEN
490 diagName = 'ADVr'//diagSufx
491 CALL DIAGNOSTICS_FILL(af,diagName, k,1, 2,bi,bj, myThid)
492 C- note: needs to explicitly increment the counter since DIAGNOSTICS_FILL
493 C does it only if k=1 (never the case here)
494 IF ( k.EQ.2 ) CALL DIAGNOSTICS_COUNT(diagName,bi,bj,myThid)
495 ENDIF
496 #endif
497 ENDIF
498
499 C- Diffusive flux in R
500 C Note: For K=1 then KM1=1 and this gives a dT/dr = 0 upper
501 C boundary condition.
502 IF (implicitDiffusion) THEN
503 DO j=1-Oly,sNy+Oly
504 DO i=1-Olx,sNx+Olx
505 df(i,j) = 0. _d 0
506 ENDDO
507 ENDDO
508 ELSE
509 IF ( applyAB_onTracer ) THEN
510 CALL GAD_DIFF_R(bi,bj,k,KappaR,TracerN,df,myThid)
511 ELSE
512 CALL GAD_DIFF_R(bi,bj,k,KappaR,TracAB, df,myThid)
513 ENDIF
514 ENDIF
515
516 #ifdef ALLOW_GMREDI
517 C- GM/Redi flux in R
518 IF (useGMRedi) THEN
519 C *note* should update GMREDI_RTRANSPORT to set df *aja*
520 IF ( applyAB_onTracer ) THEN
521 CALL GMREDI_RTRANSPORT(
522 I iMin,iMax,jMin,jMax,bi,bj,k,
523 I TracerN,tracerIdentity,
524 U df,
525 I myThid)
526 ELSE
527 CALL GMREDI_RTRANSPORT(
528 I iMin,iMax,jMin,jMax,bi,bj,k,
529 I TracAB, tracerIdentity,
530 U df,
531 I myThid)
532 ENDIF
533 ENDIF
534 #endif
535
536 DO j=1-Oly,sNy+Oly
537 DO i=1-Olx,sNx+Olx
538 fVerT(i,j,kUp) = fVerT(i,j,kUp) + df(i,j)*maskUp(i,j)
539 ENDDO
540 ENDDO
541
542 #ifdef ALLOW_DIAGNOSTICS
543 C- Diagnostics of Tracer flux in R dir (mainly Diffusive terms),
544 C Explicit terms only & excluding advective terms:
545 IF ( useDiagnostics .AND.
546 & (.NOT.implicitDiffusion .OR. useGMRedi) ) THEN
547 diagName = 'DFrE'//diagSufx
548 CALL DIAGNOSTICS_FILL(df,diagName, k,1, 2,bi,bj, myThid)
549 ENDIF
550 #endif
551
552 #ifdef ALLOW_KPP
553 C- Set non local KPP transport term (ghat):
554 IF ( useKPP .AND. k.GE.2 ) THEN
555 DO j=1-Oly,sNy+Oly
556 DO i=1-Olx,sNx+Olx
557 df(i,j) = 0. _d 0
558 ENDDO
559 ENDDO
560 IF (tracerIdentity.EQ.GAD_TEMPERATURE) THEN
561 CALL KPP_TRANSPORT_T(
562 I iMin,iMax,jMin,jMax,bi,bj,k,km1,
563 O df,
564 I myTime, myIter, myThid )
565 ELSEIF (tracerIdentity.EQ.GAD_SALINITY) THEN
566 CALL KPP_TRANSPORT_S(
567 I iMin,iMax,jMin,jMax,bi,bj,k,km1,
568 O df,
569 I myTime, myIter, myThid )
570 #ifdef ALLOW_PTRACERS
571 ELSEIF (tracerIdentity .GE. GAD_TR1) THEN
572 CALL KPP_TRANSPORT_PTR(
573 I iMin,iMax,jMin,jMax,bi,bj,k,km1,
574 I tracerIdentity-GAD_TR1+1,
575 O df,
576 I myTime, myIter, myThid )
577 #endif
578 ELSE
579 PRINT*,'invalid tracer indentity: ', tracerIdentity
580 STOP 'GAD_CALC_RHS: Ooops'
581 ENDIF
582 DO j=1-Oly,sNy+Oly
583 DO i=1-Olx,sNx+Olx
584 fVerT(i,j,kUp) = fVerT(i,j,kUp)
585 & + df(i,j)*maskUp(i,j)*rhoFacF(k)
586 ENDDO
587 ENDDO
588 ENDIF
589 #endif
590
591 C-- Divergence of fluxes
592 C Anelastic: scale vertical fluxes by rhoFac and leave Horizontal fluxes unchanged
593 DO j=1-Oly,sNy+Oly-1
594 DO i=1-Olx,sNx+Olx-1
595 gTracer(i,j,k,bi,bj)=gTracer(i,j,k,bi,bj)
596 & -_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
597 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)*recip_rhoFacC(k)
598 & *( (fZon(i+1,j)-fZon(i,j))
599 & +(fMer(i,j+1)-fMer(i,j))
600 & +(fVerT(i,j,kDown)-fVerT(i,j,kUp))*rkSign
601 & -localT(i,j)*( (uTrans(i+1,j)-uTrans(i,j))
602 & +(vTrans(i,j+1)-vTrans(i,j))
603 & +(rTransKp1(i,j)-rTrans(i,j))*rAdvFac
604 & )*advFac
605 & )
606 ENDDO
607 ENDDO
608
609 #ifdef ALLOW_DEBUG
610 IF ( debugLevel .GE. debLevB
611 & .AND. tracerIdentity.EQ.GAD_TEMPERATURE
612 & .AND. k.EQ.2 .AND. myIter.EQ.1+nIter0
613 & .AND. nPx.EQ.1 .AND. nPy.EQ.1
614 & .AND. useCubedSphereExchange ) THEN
615 CALL DEBUG_CS_CORNER_UV( ' fZon,fMer from GAD_CALC_RHS',
616 & fZon,fMer, k, standardMessageUnit,bi,bj,myThid )
617 ENDIF
618 #endif /* ALLOW_DEBUG */
619
620 RETURN
621 END

  ViewVC Help
Powered by ViewVC 1.1.22