/[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.50 - (show annotations) (download)
Sat Oct 20 14:56:43 2007 UTC (16 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint59p, checkpoint59m, checkpoint59l, checkpoint59o, checkpoint59n, checkpoint59i, checkpoint59k, checkpoint59j
Changes since 1.49: +1 -7 lines
Adding 7th Order One Step method (OS7MP) to adjoint.
Gradient check OK, but long-term stability to be seen.

1 C $Header: /u/gcmpack/MITgcm/pkg/generic_advdiff/gad_calc_rhs.F,v 1.49 2007/09/23 17:12:16 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 ELSEIF (advectionScheme.EQ.ENUM_OS7MP ) THEN
246 CALL GAD_OS7MP_ADV_X( bi,bj,k, .TRUE., dTtracerLev(k),
247 I uTrans, uFld, maskW(1-Olx,1-Oly,k,bi,bj), locABT,
248 O af, myThid )
249 ELSE
250 STOP 'GAD_CALC_RHS: Bad advectionScheme (X)'
251 ENDIF
252 DO j=1-Oly,sNy+Oly
253 DO i=1-Olx,sNx+Olx
254 fZon(i,j) = fZon(i,j) + af(i,j)
255 ENDDO
256 ENDDO
257 #ifdef ALLOW_DIAGNOSTICS
258 IF ( useDiagnostics ) THEN
259 diagName = 'ADVx'//diagSufx
260 CALL DIAGNOSTICS_FILL(af,diagName, k,1, 2,bi,bj, myThid)
261 ENDIF
262 #endif
263 ENDIF
264
265 C- Diffusive flux in X
266 IF (diffKh.NE.0.) THEN
267 CALL GAD_DIFF_X(bi,bj,k,xA,diffKh,localT,df,myThid)
268 ELSE
269 DO j=1-Oly,sNy+Oly
270 DO i=1-Olx,sNx+Olx
271 df(i,j) = 0. _d 0
272 ENDDO
273 ENDDO
274 ENDIF
275
276 C- Add bi-harmonic diffusive flux in X
277 IF (diffK4 .NE. 0.) THEN
278 CALL GAD_BIHARM_X(bi,bj,k,xA,df4,diffK4,df,myThid)
279 ENDIF
280
281 #ifdef ALLOW_GMREDI
282 C- GM/Redi flux in X
283 IF ( trUseGMRedi ) THEN
284 C *note* should update GMREDI_XTRANSPORT to set df *aja*
285 IF ( applyAB_onTracer ) THEN
286 CALL GMREDI_XTRANSPORT(
287 I iMin,iMax,jMin,jMax,bi,bj,k,
288 I xA,TracerN,tracerIdentity,
289 U df,
290 I myThid)
291 ELSE
292 CALL GMREDI_XTRANSPORT(
293 I iMin,iMax,jMin,jMax,bi,bj,k,
294 I xA,TracAB, tracerIdentity,
295 U df,
296 I myThid)
297 ENDIF
298 ENDIF
299 #endif
300 C anelastic: advect.fluxes are scaled by rhoFac but hor.diff. flx are not
301 DO j=1-Oly,sNy+Oly
302 DO i=1-Olx,sNx+Olx
303 fZon(i,j) = fZon(i,j) + df(i,j)*rhoFacC(k)
304 ENDDO
305 ENDDO
306
307 #ifdef ALLOW_DIAGNOSTICS
308 C- Diagnostics of Tracer flux in X dir (mainly Diffusive term),
309 C excluding advective terms:
310 IF ( useDiagnostics .AND.
311 & (diffKh.NE.0. .OR. diffK4 .NE.0. .OR. trUseGMRedi) ) THEN
312 diagName = 'DFxE'//diagSufx
313 CALL DIAGNOSTICS_FILL(df,diagName, k,1, 2,bi,bj, myThid)
314 ENDIF
315 #endif
316
317 C-- Initialize net flux in Y direction
318 DO j=1-Oly,sNy+Oly
319 DO i=1-Olx,sNx+Olx
320 fMer(i,j) = 0. _d 0
321 ENDDO
322 ENDDO
323
324 C- Advective flux in Y
325 IF (calcAdvection) THEN
326 IF (advectionScheme.EQ.ENUM_CENTERED_2ND) THEN
327 CALL GAD_C2_ADV_Y(bi,bj,k,vTrans,locABT,af,myThid)
328 ELSEIF ( advectionScheme.EQ.ENUM_UPWIND_1RST
329 & .OR. advectionScheme.EQ.ENUM_DST2 ) THEN
330 CALL GAD_DST2U1_ADV_Y( bi,bj,k, advectionScheme, .TRUE.,
331 I dTtracerLev(k), vTrans, vFld, locABT,
332 O af, myThid )
333 ELSEIF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN
334 CALL GAD_FLUXLIMIT_ADV_Y( bi,bj,k, .TRUE., dTtracerLev(k),
335 I vTrans, vFld, maskS(1-Olx,1-Oly,k,bi,bj), locABT,
336 O af, myThid )
337 ELSEIF (advectionScheme.EQ.ENUM_UPWIND_3RD ) THEN
338 CALL GAD_U3_ADV_Y(bi,bj,k,vTrans,locABT,af,myThid)
339 ELSEIF (advectionScheme.EQ.ENUM_CENTERED_4TH) THEN
340 CALL GAD_C4_ADV_Y(bi,bj,k,vTrans,locABT,af,myThid)
341 ELSEIF (advectionScheme.EQ.ENUM_DST3 ) THEN
342 CALL GAD_DST3_ADV_Y( bi,bj,k, .TRUE., dTtracerLev(k),
343 I vTrans, vFld, maskS(1-Olx,1-Oly,k,bi,bj), locABT,
344 O af, myThid )
345 ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
346 IF ( inAdMode ) THEN
347 cph This block is to trick the adjoint:
348 cph IF inAdExact=.FALSE., we want to use DST3
349 cph with limiters in forward, but without limiters in reverse.
350 CALL GAD_DST3_ADV_Y( bi,bj,k, .TRUE., dTtracerLev(k),
351 I vTrans, vFld, maskS(1-Olx,1-Oly,k,bi,bj), locABT,
352 O af, myThid )
353 ELSE
354 CALL GAD_DST3FL_ADV_Y( bi,bj,k, .TRUE., dTtracerLev(k),
355 I vTrans, vFld, maskS(1-Olx,1-Oly,k,bi,bj), locABT,
356 O af, myThid )
357 ENDIF
358 ELSEIF (advectionScheme.EQ.ENUM_OS7MP ) THEN
359 CALL GAD_OS7MP_ADV_Y( bi,bj,k, .TRUE., dTtracerLev(k),
360 I vTrans, vFld, maskS(1-Olx,1-Oly,k,bi,bj), locABT,
361 O af, myThid )
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 ( trUseGMRedi ) 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. trUseGMRedi) ) 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 ELSEIF (vertAdvecScheme.EQ.ENUM_OS7MP ) THEN
474 CALL GAD_OS7MP_ADV_R( bi,bj,k,
475 I dTtracerLev(k),rTrans,wFld,TracAB(1-Olx,1-Oly,1,bi,bj),
476 O af, myThid )
477 ELSE
478 STOP 'GAD_CALC_RHS: Bad vertAdvecScheme (R)'
479 ENDIF
480 C- add the advective flux to fVerT
481 DO j=1-Oly,sNy+Oly
482 DO i=1-Olx,sNx+Olx
483 fVerT(i,j,kUp) = fVerT(i,j,kUp) + af(i,j)
484 ENDDO
485 ENDDO
486 #ifdef ALLOW_DIAGNOSTICS
487 IF ( useDiagnostics ) THEN
488 diagName = 'ADVr'//diagSufx
489 CALL DIAGNOSTICS_FILL(af,diagName, k,1, 2,bi,bj, myThid)
490 C- note: needs to explicitly increment the counter since DIAGNOSTICS_FILL
491 C does it only if k=1 (never the case here)
492 IF ( k.EQ.2 ) CALL DIAGNOSTICS_COUNT(diagName,bi,bj,myThid)
493 ENDIF
494 #endif
495 ENDIF
496
497 C- Diffusive flux in R
498 C Note: For K=1 then KM1=1 and this gives a dT/dr = 0 upper
499 C boundary condition.
500 IF (implicitDiffusion) THEN
501 DO j=1-Oly,sNy+Oly
502 DO i=1-Olx,sNx+Olx
503 df(i,j) = 0. _d 0
504 ENDDO
505 ENDDO
506 ELSE
507 IF ( applyAB_onTracer ) THEN
508 CALL GAD_DIFF_R(bi,bj,k,KappaR,TracerN,df,myThid)
509 ELSE
510 CALL GAD_DIFF_R(bi,bj,k,KappaR,TracAB, df,myThid)
511 ENDIF
512 ENDIF
513
514 #ifdef ALLOW_GMREDI
515 C- GM/Redi flux in R
516 IF ( trUseGMRedi ) THEN
517 C *note* should update GMREDI_RTRANSPORT to set df *aja*
518 IF ( applyAB_onTracer ) THEN
519 CALL GMREDI_RTRANSPORT(
520 I iMin,iMax,jMin,jMax,bi,bj,k,
521 I TracerN,tracerIdentity,
522 U df,
523 I myThid)
524 ELSE
525 CALL GMREDI_RTRANSPORT(
526 I iMin,iMax,jMin,jMax,bi,bj,k,
527 I TracAB, tracerIdentity,
528 U df,
529 I myThid)
530 ENDIF
531 ENDIF
532 #endif
533
534 DO j=1-Oly,sNy+Oly
535 DO i=1-Olx,sNx+Olx
536 fVerT(i,j,kUp) = fVerT(i,j,kUp) + df(i,j)*maskUp(i,j)
537 ENDDO
538 ENDDO
539
540 #ifdef ALLOW_DIAGNOSTICS
541 C- Diagnostics of Tracer flux in R dir (mainly Diffusive terms),
542 C Explicit terms only & excluding advective terms:
543 IF ( useDiagnostics .AND.
544 & (.NOT.implicitDiffusion .OR. trUseGMRedi) ) THEN
545 diagName = 'DFrE'//diagSufx
546 CALL DIAGNOSTICS_FILL(df,diagName, k,1, 2,bi,bj, myThid)
547 ENDIF
548 #endif
549
550 #ifdef ALLOW_KPP
551 C- Set non local KPP transport term (ghat):
552 IF ( trUseKPP .AND. k.GE.2 ) THEN
553 DO j=1-Oly,sNy+Oly
554 DO i=1-Olx,sNx+Olx
555 df(i,j) = 0. _d 0
556 ENDDO
557 ENDDO
558 IF (tracerIdentity.EQ.GAD_TEMPERATURE) THEN
559 CALL KPP_TRANSPORT_T(
560 I iMin,iMax,jMin,jMax,bi,bj,k,km1,
561 O df,
562 I myTime, myIter, myThid )
563 ELSEIF (tracerIdentity.EQ.GAD_SALINITY) THEN
564 CALL KPP_TRANSPORT_S(
565 I iMin,iMax,jMin,jMax,bi,bj,k,km1,
566 O df,
567 I myTime, myIter, myThid )
568 #ifdef ALLOW_PTRACERS
569 ELSEIF (tracerIdentity .GE. GAD_TR1) THEN
570 CALL KPP_TRANSPORT_PTR(
571 I iMin,iMax,jMin,jMax,bi,bj,k,km1,
572 I tracerIdentity-GAD_TR1+1,
573 O df,
574 I myTime, myIter, myThid )
575 #endif
576 ELSE
577 PRINT*,'invalid tracer indentity: ', tracerIdentity
578 STOP 'GAD_CALC_RHS: Ooops'
579 ENDIF
580 DO j=1-Oly,sNy+Oly
581 DO i=1-Olx,sNx+Olx
582 fVerT(i,j,kUp) = fVerT(i,j,kUp)
583 & + df(i,j)*maskUp(i,j)*rhoFacF(k)
584 ENDDO
585 ENDDO
586 ENDIF
587 #endif
588
589 C-- Divergence of fluxes
590 C Anelastic: scale vertical fluxes by rhoFac and leave Horizontal fluxes unchanged
591 DO j=1-Oly,sNy+Oly-1
592 DO i=1-Olx,sNx+Olx-1
593 gTracer(i,j,k,bi,bj)=gTracer(i,j,k,bi,bj)
594 & -_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)
595 & *recip_rA(i,j,bi,bj)*recip_deepFac2C(k)*recip_rhoFacC(k)
596 & *( (fZon(i+1,j)-fZon(i,j))
597 & +(fMer(i,j+1)-fMer(i,j))
598 & +(fVerT(i,j,kDown)-fVerT(i,j,kUp))*rkSign
599 & -localT(i,j)*( (uTrans(i+1,j)-uTrans(i,j))
600 & +(vTrans(i,j+1)-vTrans(i,j))
601 & +(rTransKp1(i,j)-rTrans(i,j))*rAdvFac
602 & )*advFac
603 & )
604 ENDDO
605 ENDDO
606
607 #ifdef ALLOW_DEBUG
608 IF ( debugLevel .GE. debLevB
609 & .AND. tracerIdentity.EQ.GAD_TEMPERATURE
610 & .AND. k.EQ.2 .AND. myIter.EQ.1+nIter0
611 & .AND. nPx.EQ.1 .AND. nPy.EQ.1
612 & .AND. useCubedSphereExchange ) THEN
613 CALL DEBUG_CS_CORNER_UV( ' fZon,fMer from GAD_CALC_RHS',
614 & fZon,fMer, k, standardMessageUnit,bi,bj,myThid )
615 ENDIF
616 #endif /* ALLOW_DEBUG */
617
618 RETURN
619 END

  ViewVC Help
Powered by ViewVC 1.1.22