/[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.35 - (hide annotations) (download)
Sun Apr 3 16:01:21 2005 UTC (19 years, 2 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57g_pre, checkpoint57g_post, checkpoint57i_post, checkpoint57f_post, checkpoint57h_pre, checkpoint57h_post, checkpoint57h_done, checkpoint57f_pre
Changes since 1.34: +33 -7 lines
We fool the adjoint:
When using DST3FL(=33) we revert to unlimited DST3(=30) in reverse mode
if we set inAdExact=.FALSE.

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

  ViewVC Help
Powered by ViewVC 1.1.22