/[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.39 - (hide annotations) (download)
Sun Feb 26 01:56:27 2006 UTC (18 years, 3 months ago) by jmc
Branch: MAIN
Changes since 1.38: +30 -4 lines
just a try to fix the adjoint Pb (duplicated input argument)

1 jmc 1.39 C $Header: /u/gcmpack/MITgcm/pkg/generic_advdiff/gad_calc_rhs.F,v 1.38 2005/11/06 22:14:02 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 adcroft 1.1 SUBROUTINE GAD_CALC_RHS(
11     I bi,bj,iMin,iMax,jMin,jMax,k,kM1,kUp,kDown,
12 jmc 1.23 I xA,yA,uTrans,vTrans,rTrans,rTransKp1,maskUp,
13     I uVel, vVel, wVel,
14 jmc 1.38 I diffKh, diffK4, KappaR, Tracer, TracAB,
15 jmc 1.26 I tracerIdentity, advectionScheme, vertAdvecScheme,
16 jmc 1.23 I calcAdvection, implicitAdvection,
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     C kup :: index into 2 1/2D array, toggles between 1|2
58     C kdown :: index into 2 1/2D array, toggles between 2|1
59     C kp1 :: =k+1 for k<Nr, =Nr for k=Nr
60     C xA,yA :: areas of X and Y face of tracer cells
61     C uTrans,vTrans :: 2-D arrays of volume transports at U,V points
62     C rTrans :: 2-D arrays of volume transports at W points
63     C rTransKp1 :: 2-D array of volume trans at W pts, interf k+1
64     C maskUp :: 2-D array for mask at W points
65     C uVel,vVel,wVel :: 3 components of the velcity field (3-D array)
66     C diffKh :: horizontal diffusion coefficient
67     C diffK4 :: bi-harmonic diffusion coefficient
68 jmc 1.30 C KappaR :: 2-D array for vertical diffusion coefficient, interf k
69 jmc 1.39 C Tracer :: tracer field @ time-step n (Note: only used
70     C if applying AB on tracer field rather than on tendency Gtr)
71     C TracAB :: current tracer field (@ time-step n if applying AB on Gtr
72     C or extrapolated fwd in time to n+1/2 if applying AB on Tr)
73 edhill 1.24 C tracerIdentity :: tracer identifier (required for KPP,GM)
74 jmc 1.26 C advectionScheme :: advection scheme to use (Horizontal plane)
75     C vertAdvecScheme :: advection scheme to use (Vertical direction)
76 edhill 1.24 C calcAdvection :: =False if Advec computed with multiDim scheme
77     C implicitAdvection:: =True if vertical Advec computed implicitly
78 jmc 1.27 C myTime :: current time
79     C myIter :: iteration number
80 edhill 1.24 C myThid :: thread number
81 adcroft 1.11 INTEGER bi,bj,iMin,iMax,jMin,jMax
82 adcroft 1.1 INTEGER k,kUp,kDown,kM1
83     _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
84     _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
85     _RL uTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
86     _RL vTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
87     _RL rTrans(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
88 jmc 1.23 _RL rTransKp1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
89 adcroft 1.1 _RS maskUp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
90 jmc 1.23 _RL uVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
91     _RL vVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
92     _RL wVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
93 adcroft 1.1 _RL diffKh, diffK4
94 jmc 1.30 _RL KappaR(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95 adcroft 1.1 _RL Tracer(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
96 jmc 1.38 _RL TracAB(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
97 adcroft 1.1 INTEGER tracerIdentity
98 jmc 1.26 INTEGER advectionScheme, vertAdvecScheme
99 jmc 1.14 LOGICAL calcAdvection
100 jmc 1.23 LOGICAL implicitAdvection
101 jmc 1.27 _RL myTime
102     INTEGER myIter, myThid
103 adcroft 1.11
104     C !OUTPUT PARAMETERS: ==================================================
105 jmc 1.38 C gTracer :: tendency array
106 edhill 1.24 C fVerT :: 2 1/2D arrays for vertical advective flux
107 adcroft 1.11 _RL gTracer(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
108 adcroft 1.1 _RL fVerT (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
109    
110 adcroft 1.11 C !LOCAL VARIABLES: ====================================================
111 edhill 1.24 C i,j :: loop indices
112     C df4 :: used for storing del^2 T for bi-harmonic term
113     C fZon :: zonal flux
114 jmc 1.32 C fMer :: meridional flux
115 edhill 1.24 C af :: advective flux
116     C df :: diffusive flux
117     C localT :: local copy of tracer field
118 jmc 1.38 C locABT :: local copy of (AB-extrapolated) tracer field
119 jmc 1.32 #ifdef ALLOW_DIAGNOSTICS
120     CHARACTER*8 diagName
121     CHARACTER*4 GAD_DIAG_SUFX, diagSufx
122     EXTERNAL GAD_DIAG_SUFX
123     #endif
124 adcroft 1.1 INTEGER i,j
125     _RL df4 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
126     _RL fZon (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
127     _RL fMer (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
128     _RL af (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
129     _RL df (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
130     _RL localT(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
131 jmc 1.38 _RL locABT(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
132 jmc 1.23 _RL advFac, rAdvFac
133 adcroft 1.11 CEOP
134 adcroft 1.1
135     #ifdef ALLOW_AUTODIFF_TAMC
136     C-- only the kUp part of fverT is set in this subroutine
137     C-- the kDown is still required
138     fVerT(1,1,kDown) = fVerT(1,1,kDown)
139     #endif
140 heimbach 1.13
141 jmc 1.32 #ifdef ALLOW_DIAGNOSTICS
142     C-- Set diagnostic suffix for the current tracer
143     IF ( useDiagnostics ) THEN
144     diagSufx = GAD_DIAG_SUFX( tracerIdentity, myThid )
145     ENDIF
146     #endif
147    
148 jmc 1.23 advFac = 0. _d 0
149     IF (calcAdvection) advFac = 1. _d 0
150 jmc 1.36 rAdvFac = rkSign*advFac
151 jmc 1.23 IF (implicitAdvection) rAdvFac = 0. _d 0
152    
153 adcroft 1.1 DO j=1-OLy,sNy+OLy
154     DO i=1-OLx,sNx+OLx
155 heimbach 1.12 fZon(i,j) = 0. _d 0
156     fMer(i,j) = 0. _d 0
157     fVerT(i,j,kUp) = 0. _d 0
158 heimbach 1.13 df(i,j) = 0. _d 0
159     df4(i,j) = 0. _d 0
160 adcroft 1.1 ENDDO
161     ENDDO
162    
163     C-- Make local copy of tracer array
164 jmc 1.39 #ifdef ALLOW_ADAMSBASHFORTH_3
165 adcroft 1.1 DO j=1-OLy,sNy+OLy
166     DO i=1-OLx,sNx+OLx
167 jmc 1.37 localT(i,j)=Tracer(i,j,k,bi,bj)
168 jmc 1.38 locABT(i,j)=TracAB(i,j,k,bi,bj)
169 adcroft 1.1 ENDDO
170     ENDDO
171 jmc 1.39 #else /* ALLOW_ADAMSBASHFORTH_3 */
172     DO j=1-OLy,sNy+OLy
173     DO i=1-OLx,sNx+OLx
174     localT(i,j)=TracAB(i,j,k,bi,bj)
175     locABT(i,j)=TracAB(i,j,k,bi,bj)
176     ENDDO
177     ENDDO
178     #endif /* ALLOW_ADAMSBASHFORTH_3 */
179 adcroft 1.1
180 adcroft 1.8 C-- Unless we have already calculated the advection terms we initialize
181     C the tendency to zero.
182 jmc 1.20 C <== now done earlier at the beginning of thermodynamics.
183     c IF (calcAdvection) THEN
184     c DO j=1-Oly,sNy+Oly
185     c DO i=1-Olx,sNx+Olx
186     c gTracer(i,j,k,bi,bj)=0. _d 0
187     c ENDDO
188     c ENDDO
189     c ENDIF
190 adcroft 1.1
191     C-- Pre-calculate del^2 T if bi-harmonic coefficient is non-zero
192     IF (diffK4 .NE. 0.) THEN
193     CALL GAD_GRAD_X(bi,bj,k,xA,localT,fZon,myThid)
194     CALL GAD_GRAD_Y(bi,bj,k,yA,localT,fMer,myThid)
195     CALL GAD_DEL2(bi,bj,k,fZon,fMer,df4,myThid)
196     ENDIF
197    
198     C-- Initialize net flux in X direction
199     DO j=1-Oly,sNy+Oly
200     DO i=1-Olx,sNx+Olx
201 heimbach 1.12 fZon(i,j) = 0. _d 0
202 adcroft 1.1 ENDDO
203     ENDDO
204    
205     C- Advective flux in X
206 jmc 1.14 IF (calcAdvection) THEN
207 jmc 1.32 IF (advectionScheme.EQ.ENUM_CENTERED_2ND) THEN
208 jmc 1.38 CALL GAD_C2_ADV_X(bi,bj,k,uTrans,locABT,af,myThid)
209 jmc 1.37 ELSEIF ( advectionScheme.EQ.ENUM_UPWIND_1RST
210     & .OR. advectionScheme.EQ.ENUM_DST2 ) THEN
211     CALL GAD_DST2U1_ADV_X( bi,bj,k, advectionScheme,
212 jmc 1.38 I dTtracerLev(k), uTrans, uVel, locABT,
213 jmc 1.37 O af, myThid )
214 jmc 1.32 ELSEIF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN
215     CALL GAD_FLUXLIMIT_ADV_X( bi,bj,k, dTtracerLev(k),
216 jmc 1.38 I uTrans, uVel, maskW(1-Olx,1-Oly,k,bi,bj), locABT,
217 jmc 1.32 O af, myThid )
218     ELSEIF (advectionScheme.EQ.ENUM_UPWIND_3RD ) THEN
219 jmc 1.38 CALL GAD_U3_ADV_X(bi,bj,k,uTrans,locABT,af,myThid)
220 jmc 1.32 ELSEIF (advectionScheme.EQ.ENUM_CENTERED_4TH) THEN
221 jmc 1.38 CALL GAD_C4_ADV_X(bi,bj,k,uTrans,locABT,af,myThid)
222 jmc 1.32 ELSEIF (advectionScheme.EQ.ENUM_DST3 ) THEN
223     CALL GAD_DST3_ADV_X( bi,bj,k, dTtracerLev(k),
224 jmc 1.38 I uTrans, uVel, maskW(1-Olx,1-Oly,k,bi,bj), locABT,
225 jmc 1.32 O af, myThid )
226     ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
227 heimbach 1.35 IF ( inAdMode ) THEN
228     cph This block is to trick the adjoint:
229     cph IF inAdExact=.FALSE., we want to use DST3
230     cph with limiters in forward, but without limiters in reverse.
231     CALL GAD_DST3_ADV_X( bi,bj,k, dTtracerLev(k),
232 jmc 1.38 I uTrans, uVel, maskW(1-Olx,1-Oly,k,bi,bj), locABT,
233 heimbach 1.35 O af, myThid )
234     ELSE
235 jmc 1.32 CALL GAD_DST3FL_ADV_X( bi,bj,k, dTtracerLev(k),
236 jmc 1.38 I uTrans, uVel, maskW(1-Olx,1-Oly,k,bi,bj), locABT,
237 heimbach 1.35 O af, myThid )
238     ENDIF
239 jmc 1.32 ELSE
240     STOP 'GAD_CALC_RHS: Bad advectionScheme (X)'
241     ENDIF
242     DO j=1-Oly,sNy+Oly
243     DO i=1-Olx,sNx+Olx
244     fZon(i,j) = fZon(i,j) + af(i,j)
245     ENDDO
246     ENDDO
247     #ifdef ALLOW_DIAGNOSTICS
248     IF ( useDiagnostics ) THEN
249     diagName = 'ADVx'//diagSufx
250 jmc 1.33 CALL DIAGNOSTICS_FILL(af,diagName, k,1, 2,bi,bj, myThid)
251 jmc 1.32 ENDIF
252     #endif
253 adcroft 1.8 ENDIF
254 adcroft 1.1
255     C- Diffusive flux in X
256     IF (diffKh.NE.0.) THEN
257     CALL GAD_DIFF_X(bi,bj,k,xA,diffKh,localT,df,myThid)
258     ELSE
259 adcroft 1.5 DO j=1-Oly,sNy+Oly
260     DO i=1-Olx,sNx+Olx
261 heimbach 1.12 df(i,j) = 0. _d 0
262 adcroft 1.1 ENDDO
263     ENDDO
264     ENDIF
265    
266 jmc 1.32 C- Add bi-harmonic diffusive flux in X
267     IF (diffK4 .NE. 0.) THEN
268     CALL GAD_BIHARM_X(bi,bj,k,xA,df4,diffK4,df,myThid)
269     ENDIF
270    
271 adcroft 1.1 #ifdef ALLOW_GMREDI
272     C- GM/Redi flux in X
273     IF (useGMRedi) THEN
274 jmc 1.38 C *note* should update GMREDI_XTRANSPORT to set df *aja*
275 adcroft 1.1 CALL GMREDI_XTRANSPORT(
276     I iMin,iMax,jMin,jMax,bi,bj,K,
277 jmc 1.39 #ifdef ALLOW_ADAMSBASHFORTH_3
278 heimbach 1.15 I xA,Tracer,tracerIdentity,
279 jmc 1.39 #else
280     I xA,TracAB,tracerIdentity,
281     #endif
282 adcroft 1.1 U df,
283     I myThid)
284     ENDIF
285     #endif
286 adcroft 1.5 DO j=1-Oly,sNy+Oly
287     DO i=1-Olx,sNx+Olx
288 adcroft 1.1 fZon(i,j) = fZon(i,j) + df(i,j)
289     ENDDO
290     ENDDO
291    
292 jmc 1.32 #ifdef ALLOW_DIAGNOSTICS
293     C- Diagnostics of Tracer flux in X dir (mainly Diffusive term),
294     C excluding advective terms:
295     IF ( useDiagnostics .AND.
296     & (diffKh.NE.0. .OR. diffK4 .NE.0. .OR. useGMRedi) ) THEN
297     diagName = 'DIFx'//diagSufx
298 jmc 1.33 CALL DIAGNOSTICS_FILL(df,diagName, k,1, 2,bi,bj, myThid)
299 adcroft 1.1 ENDIF
300 jmc 1.32 #endif
301 adcroft 1.1
302     C-- Initialize net flux in Y direction
303     DO j=1-Oly,sNy+Oly
304     DO i=1-Olx,sNx+Olx
305 heimbach 1.12 fMer(i,j) = 0. _d 0
306 adcroft 1.1 ENDDO
307     ENDDO
308    
309     C- Advective flux in Y
310 jmc 1.14 IF (calcAdvection) THEN
311 jmc 1.32 IF (advectionScheme.EQ.ENUM_CENTERED_2ND) THEN
312 jmc 1.38 CALL GAD_C2_ADV_Y(bi,bj,k,vTrans,locABT,af,myThid)
313 jmc 1.37 ELSEIF ( advectionScheme.EQ.ENUM_UPWIND_1RST
314     & .OR. advectionScheme.EQ.ENUM_DST2 ) THEN
315     CALL GAD_DST2U1_ADV_Y( bi,bj,k, advectionScheme,
316 jmc 1.38 I dTtracerLev(k), vTrans, vVel, locABT,
317 jmc 1.37 O af, myThid )
318 jmc 1.32 ELSEIF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN
319     CALL GAD_FLUXLIMIT_ADV_Y( bi,bj,k, dTtracerLev(k),
320 jmc 1.38 I vTrans, vVel, maskS(1-Olx,1-Oly,k,bi,bj), locABT,
321 jmc 1.32 O af, myThid )
322     ELSEIF (advectionScheme.EQ.ENUM_UPWIND_3RD ) THEN
323 jmc 1.38 CALL GAD_U3_ADV_Y(bi,bj,k,vTrans,locABT,af,myThid)
324 jmc 1.32 ELSEIF (advectionScheme.EQ.ENUM_CENTERED_4TH) THEN
325 jmc 1.38 CALL GAD_C4_ADV_Y(bi,bj,k,vTrans,locABT,af,myThid)
326 jmc 1.32 ELSEIF (advectionScheme.EQ.ENUM_DST3 ) THEN
327     CALL GAD_DST3_ADV_Y( bi,bj,k, dTtracerLev(k),
328 jmc 1.38 I vTrans, vVel, maskS(1-Olx,1-Oly,k,bi,bj), locABT,
329 jmc 1.32 O af, myThid )
330     ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
331 heimbach 1.35 IF ( inAdMode ) THEN
332     cph This block is to trick the adjoint:
333     cph IF inAdExact=.FALSE., we want to use DST3
334     cph with limiters in forward, but without limiters in reverse.
335     CALL GAD_DST3_ADV_Y( bi,bj,k, dTtracerLev(k),
336 jmc 1.38 I vTrans, vVel, maskS(1-Olx,1-Oly,k,bi,bj), locABT,
337 heimbach 1.35 O af, myThid )
338     ELSE
339 jmc 1.32 CALL GAD_DST3FL_ADV_Y( bi,bj,k, dTtracerLev(k),
340 jmc 1.38 I vTrans, vVel, maskS(1-Olx,1-Oly,k,bi,bj), locABT,
341 heimbach 1.35 O af, myThid )
342     ENDIF
343 jmc 1.32 ELSE
344     STOP 'GAD_CALC_RHS: Bad advectionScheme (Y)'
345     ENDIF
346     DO j=1-Oly,sNy+Oly
347     DO i=1-Olx,sNx+Olx
348     fMer(i,j) = fMer(i,j) + af(i,j)
349     ENDDO
350     ENDDO
351     #ifdef ALLOW_DIAGNOSTICS
352     IF ( useDiagnostics ) THEN
353     diagName = 'ADVy'//diagSufx
354 jmc 1.33 CALL DIAGNOSTICS_FILL(af,diagName, k,1, 2,bi,bj, myThid)
355 jmc 1.32 ENDIF
356     #endif
357 adcroft 1.8 ENDIF
358 adcroft 1.1
359     C- Diffusive flux in Y
360     IF (diffKh.NE.0.) THEN
361     CALL GAD_DIFF_Y(bi,bj,k,yA,diffKh,localT,df,myThid)
362     ELSE
363     DO j=1-Oly,sNy+Oly
364     DO i=1-Olx,sNx+Olx
365 heimbach 1.12 df(i,j) = 0. _d 0
366 adcroft 1.1 ENDDO
367     ENDDO
368     ENDIF
369    
370 jmc 1.32 C- Add bi-harmonic flux in Y
371     IF (diffK4 .NE. 0.) THEN
372     CALL GAD_BIHARM_Y(bi,bj,k,yA,df4,diffK4,df,myThid)
373     ENDIF
374    
375 adcroft 1.1 #ifdef ALLOW_GMREDI
376     C- GM/Redi flux in Y
377     IF (useGMRedi) THEN
378 jmc 1.38 C *note* should update GMREDI_YTRANSPORT to set df *aja*
379 adcroft 1.1 CALL GMREDI_YTRANSPORT(
380     I iMin,iMax,jMin,jMax,bi,bj,K,
381 jmc 1.39 #ifdef ALLOW_ADAMSBASHFORTH_3
382 heimbach 1.15 I yA,Tracer,tracerIdentity,
383 jmc 1.39 #else
384     I yA,TracAB,tracerIdentity,
385     #endif
386 adcroft 1.1 U df,
387     I myThid)
388     ENDIF
389     #endif
390     DO j=1-Oly,sNy+Oly
391     DO i=1-Olx,sNx+Olx
392     fMer(i,j) = fMer(i,j) + df(i,j)
393     ENDDO
394     ENDDO
395    
396 jmc 1.32 #ifdef ALLOW_DIAGNOSTICS
397     C- Diagnostics of Tracer flux in Y dir (mainly Diffusive terms),
398     C excluding advective terms:
399     IF ( useDiagnostics .AND.
400     & (diffKh.NE.0. .OR. diffK4 .NE.0. .OR. useGMRedi) ) THEN
401     diagName = 'DIFy'//diagSufx
402 jmc 1.33 CALL DIAGNOSTICS_FILL(df,diagName, k,1, 2,bi,bj, myThid)
403 adcroft 1.1 ENDIF
404 jmc 1.32 #endif
405 adcroft 1.1
406 jmc 1.16 C-- Compute vertical flux fVerT(kUp) at interface k (between k-1 & k):
407 adcroft 1.1 C- Advective flux in R
408 jmc 1.25 #ifdef ALLOW_AIM
409     C- a hack to prevent Water-Vapor vert.transport into the stratospheric level Nr
410     IF (calcAdvection .AND. .NOT.implicitAdvection .AND. K.GE.2 .AND.
411     & (.NOT.useAIM .OR.tracerIdentity.NE.GAD_SALINITY .OR.K.LT.Nr)
412     & ) THEN
413     #else
414 jmc 1.23 IF (calcAdvection .AND. .NOT.implicitAdvection .AND. K.GE.2) THEN
415 jmc 1.25 #endif
416 jmc 1.2 C- Compute vertical advective flux in the interior:
417 jmc 1.32 IF (vertAdvecScheme.EQ.ENUM_CENTERED_2ND) THEN
418 jmc 1.38 CALL GAD_C2_ADV_R(bi,bj,k,rTrans,TracAB,af,myThid)
419 jmc 1.37 ELSEIF ( vertAdvecScheme.EQ.ENUM_UPWIND_1RST
420     & .OR. vertAdvecScheme.EQ.ENUM_DST2 ) THEN
421     CALL GAD_DST2U1_ADV_R( bi,bj,k, vertAdvecScheme,
422 jmc 1.38 I dTtracerLev(k),rTrans,wVel,TracAB(1-Olx,1-Oly,1,bi,bj),
423 jmc 1.37 O af, myThid )
424 jmc 1.32 ELSEIF (vertAdvecScheme.EQ.ENUM_FLUX_LIMIT) THEN
425 jmc 1.37 CALL GAD_FLUXLIMIT_ADV_R( bi,bj,k,
426 jmc 1.38 I dTtracerLev(k),rTrans,wVel,TracAB(1-Olx,1-Oly,1,bi,bj),
427 jmc 1.37 O af, myThid )
428 jmc 1.32 ELSEIF (vertAdvecScheme.EQ.ENUM_UPWIND_3RD ) THEN
429 jmc 1.38 CALL GAD_U3_ADV_R(bi,bj,k,rTrans,TracAB,af,myThid)
430 jmc 1.32 ELSEIF (vertAdvecScheme.EQ.ENUM_CENTERED_4TH) THEN
431 jmc 1.38 CALL GAD_C4_ADV_R(bi,bj,k,rTrans,TracAB,af,myThid)
432 jmc 1.32 ELSEIF (vertAdvecScheme.EQ.ENUM_DST3 ) THEN
433 jmc 1.37 CALL GAD_DST3_ADV_R( bi,bj,k,
434 jmc 1.38 I dTtracerLev(k),rTrans,wVel,TracAB(1-Olx,1-Oly,1,bi,bj),
435 jmc 1.37 O af, myThid )
436 jmc 1.32 ELSEIF (vertAdvecScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
437 heimbach 1.35 cph This block is to trick the adjoint:
438     cph IF inAdExact=.FALSE., we want to use DST3
439     cph with limiters in forward, but without limiters in reverse.
440     IF ( inAdMode ) THEN
441 jmc 1.37 CALL GAD_DST3_ADV_R( bi,bj,k,
442 jmc 1.38 I dTtracerLev(k),rTrans,wVel,TracAB(1-Olx,1-Oly,1,bi,bj),
443 jmc 1.37 O af, myThid )
444 heimbach 1.35 ELSE
445 jmc 1.37 CALL GAD_DST3FL_ADV_R( bi,bj,k,
446 jmc 1.38 I dTtracerLev(k),rTrans,wVel,TracAB(1-Olx,1-Oly,1,bi,bj),
447 jmc 1.37 O af, myThid )
448 heimbach 1.35 ENDIF
449 jmc 1.32 ELSE
450     STOP 'GAD_CALC_RHS: Bad vertAdvecScheme (R)'
451     ENDIF
452 jmc 1.23 C- add the advective flux to fVerT
453 jmc 1.32 DO j=1-Oly,sNy+Oly
454     DO i=1-Olx,sNx+Olx
455     fVerT(i,j,kUp) = fVerT(i,j,kUp) + af(i,j)
456     ENDDO
457 jmc 1.2 ENDDO
458 jmc 1.32 #ifdef ALLOW_DIAGNOSTICS
459     IF ( useDiagnostics ) THEN
460     diagName = 'ADVr'//diagSufx
461 jmc 1.33 CALL DIAGNOSTICS_FILL(af,diagName, k,1, 2,bi,bj, myThid)
462 jmc 1.34 C- note: needs to explicitly increment the counter since DIAGNOSTICS_FILL
463     C does it only if k=1 (never the case here)
464     IF ( k.EQ.2 ) CALL DIAGNOSTICS_COUNT(diagName,bi,bj,myThid)
465 jmc 1.32 ENDIF
466     #endif
467 adcroft 1.8 ENDIF
468 adcroft 1.1
469     C- Diffusive flux in R
470     C Note: For K=1 then KM1=1 and this gives a dT/dr = 0 upper
471     C boundary condition.
472     IF (implicitDiffusion) THEN
473 adcroft 1.5 DO j=1-Oly,sNy+Oly
474     DO i=1-Olx,sNx+Olx
475 heimbach 1.12 df(i,j) = 0. _d 0
476 adcroft 1.1 ENDDO
477     ENDDO
478     ELSE
479 jmc 1.39 #ifdef ALLOW_ADAMSBASHFORTH_3
480 jmc 1.37 CALL GAD_DIFF_R(bi,bj,k,KappaR,Tracer,df,myThid)
481 jmc 1.39 #else
482     CALL GAD_DIFF_R(bi,bj,k,KappaR,TracAB,df,myThid)
483     #endif
484 adcroft 1.1 ENDIF
485    
486     #ifdef ALLOW_GMREDI
487     C- GM/Redi flux in R
488     IF (useGMRedi) THEN
489     C *note* should update GMREDI_RTRANSPORT to set df *aja*
490     CALL GMREDI_RTRANSPORT(
491     I iMin,iMax,jMin,jMax,bi,bj,K,
492 jmc 1.39 #ifdef ALLOW_ADAMSBASHFORTH_3
493 heimbach 1.15 I Tracer,tracerIdentity,
494 jmc 1.39 #else
495     I TracAB,tracerIdentity,
496     #endif
497 adcroft 1.1 U df,
498     I myThid)
499     ENDIF
500     #endif
501    
502 adcroft 1.5 DO j=1-Oly,sNy+Oly
503     DO i=1-Olx,sNx+Olx
504 adcroft 1.11 fVerT(i,j,kUp) = fVerT(i,j,kUp) + df(i,j)*maskUp(i,j)
505 adcroft 1.1 ENDDO
506     ENDDO
507    
508 jmc 1.32 #ifdef ALLOW_DIAGNOSTICS
509     C- Diagnostics of Tracer flux in R dir (mainly Diffusive terms),
510     C Explicit terms only & excluding advective terms:
511     IF ( useDiagnostics .AND.
512     & (.NOT.implicitDiffusion .OR. useGMRedi) ) THEN
513     diagName = 'DFrE'//diagSufx
514 jmc 1.33 CALL DIAGNOSTICS_FILL(df,diagName, k,1, 2,bi,bj, myThid)
515 jmc 1.32 ENDIF
516     #endif
517    
518 adcroft 1.1 #ifdef ALLOW_KPP
519 jmc 1.29 C- Set non local KPP transport term (ghat):
520     IF ( useKPP .AND. k.GE.2 ) THEN
521 adcroft 1.5 DO j=1-Oly,sNy+Oly
522     DO i=1-Olx,sNx+Olx
523 heimbach 1.12 df(i,j) = 0. _d 0
524 adcroft 1.1 ENDDO
525     ENDDO
526     IF (tracerIdentity.EQ.GAD_TEMPERATURE) THEN
527     CALL KPP_TRANSPORT_T(
528     I iMin,iMax,jMin,jMax,bi,bj,k,km1,
529 jmc 1.29 O df )
530 adcroft 1.1 ELSEIF (tracerIdentity.EQ.GAD_SALINITY) THEN
531     CALL KPP_TRANSPORT_S(
532     I iMin,iMax,jMin,jMax,bi,bj,k,km1,
533 jmc 1.29 O df )
534 mlosch 1.18 #ifdef ALLOW_PTRACERS
535 dimitri 1.22 ELSEIF (tracerIdentity .GE. GAD_TR1) THEN
536 mlosch 1.18 CALL KPP_TRANSPORT_PTR(
537     I iMin,iMax,jMin,jMax,bi,bj,k,km1,
538 jmc 1.29 I tracerIdentity-GAD_TR1+1,
539     O df )
540 mlosch 1.18 #endif
541 adcroft 1.1 ELSE
542 mlosch 1.18 PRINT*,'invalid tracer indentity: ', tracerIdentity
543 adcroft 1.1 STOP 'GAD_CALC_RHS: Ooops'
544     ENDIF
545 adcroft 1.5 DO j=1-Oly,sNy+Oly
546     DO i=1-Olx,sNx+Olx
547 adcroft 1.11 fVerT(i,j,kUp) = fVerT(i,j,kUp) + df(i,j)*maskUp(i,j)
548 adcroft 1.1 ENDDO
549     ENDDO
550     ENDIF
551     #endif
552    
553     C-- Divergence of fluxes
554 adcroft 1.10 DO j=1-Oly,sNy+Oly-1
555     DO i=1-Olx,sNx+Olx-1
556 adcroft 1.8 gTracer(i,j,k,bi,bj)=gTracer(i,j,k,bi,bj)
557 jmc 1.23 & -_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)*recip_rA(i,j,bi,bj)
558     & *( (fZon(i+1,j)-fZon(i,j))
559     & +(fMer(i,j+1)-fMer(i,j))
560 jmc 1.36 & +(fVerT(i,j,kDown)-fVerT(i,j,kUp))*rkSign
561 jmc 1.23 & -localT(i,j)*( (uTrans(i+1,j)-uTrans(i,j))
562     & +(vTrans(i,j+1)-vTrans(i,j))
563 jmc 1.36 & +(rTransKp1(i,j)-rTrans(i,j))*rAdvFac
564 jmc 1.23 & )*advFac
565 adcroft 1.1 & )
566     ENDDO
567     ENDDO
568    
569 jmc 1.27 #ifdef ALLOW_DEBUG
570     IF ( debugLevel .GE. debLevB
571 jmc 1.28 & .AND. tracerIdentity.EQ.GAD_TEMPERATURE
572 jmc 1.27 & .AND. k.EQ.2 .AND. myIter.EQ.1+nIter0
573     & .AND. nPx.EQ.1 .AND. nPy.EQ.1
574     & .AND. useCubedSphereExchange ) THEN
575     CALL DEBUG_CS_CORNER_UV( ' fZon,fMer from GAD_CALC_RHS',
576     & fZon,fMer, k, standardMessageUnit,bi,bj,myThid )
577     ENDIF
578     #endif /* ALLOW_DEBUG */
579    
580 adcroft 1.1 RETURN
581     END

  ViewVC Help
Powered by ViewVC 1.1.22