/[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.49 - (show annotations) (download)
Sun Sep 23 17:12:16 2007 UTC (16 years, 8 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint59h
Changes since 1.48: +6 -25 lines
add trUseGMRedi & trUseKPP to the argument list

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

  ViewVC Help
Powered by ViewVC 1.1.22