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

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

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


Revision 1.41 - (hide annotations) (download)
Sun Jun 18 23:27:44 2006 UTC (18 years, 3 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint58m_post, checkpoint58o_post, checkpoint58p_post, checkpoint58q_post, checkpoint58r_post, checkpoint58n_post, checkpoint58k_post, checkpoint58l_post
Changes since 1.40: +41 -40 lines
make a local copy of velocity to pass (like u,v,r_Trans) to tracer advection S/R

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

  ViewVC Help
Powered by ViewVC 1.1.22