/[MITgcm]/MITgcm_contrib/jscott/code_rafmod/gad_calc_rhs_raf.F
ViewVC logotype

Annotation of /MITgcm_contrib/jscott/code_rafmod/gad_calc_rhs_raf.F

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


Revision 1.2 - (hide annotations) (download)
Thu Sep 3 20:40:01 2009 UTC (15 years, 11 months ago) by jscott
Branch: MAIN
CVS Tags: HEAD
Changes since 1.1: +117 -35 lines
update code for crude ML horiz mixing scheme

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

  ViewVC Help
Powered by ViewVC 1.1.22