/[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.34 - (hide annotations) (download)
Thu Feb 17 00:03:07 2005 UTC (19 years, 3 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57e_post, eckpoint57e_pre
Changes since 1.33: +4 -1 lines
fix diagnostics of vertical advective flux.

1 jmc 1.34 C $Header: /u/gcmpack/MITgcm/pkg/generic_advdiff/gad_calc_rhs.F,v 1.33 2004/12/20 19:10:13 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     CALL GAD_DST3FL_ADV_X( bi,bj,k, dTtracerLev(k),
207     I uTrans, uVel, maskW(1-Olx,1-Oly,k,bi,bj), localT,
208     O af, myThid )
209     ELSE
210     STOP 'GAD_CALC_RHS: Bad advectionScheme (X)'
211     ENDIF
212     DO j=1-Oly,sNy+Oly
213     DO i=1-Olx,sNx+Olx
214     fZon(i,j) = fZon(i,j) + af(i,j)
215     ENDDO
216     ENDDO
217     #ifdef ALLOW_DIAGNOSTICS
218     IF ( useDiagnostics ) THEN
219     diagName = 'ADVx'//diagSufx
220 jmc 1.33 CALL DIAGNOSTICS_FILL(af,diagName, k,1, 2,bi,bj, myThid)
221 jmc 1.32 ENDIF
222     #endif
223 adcroft 1.8 ENDIF
224 adcroft 1.1
225     C- Diffusive flux in X
226     IF (diffKh.NE.0.) THEN
227     CALL GAD_DIFF_X(bi,bj,k,xA,diffKh,localT,df,myThid)
228     ELSE
229 adcroft 1.5 DO j=1-Oly,sNy+Oly
230     DO i=1-Olx,sNx+Olx
231 heimbach 1.12 df(i,j) = 0. _d 0
232 adcroft 1.1 ENDDO
233     ENDDO
234     ENDIF
235    
236 jmc 1.32 C- Add bi-harmonic diffusive flux in X
237     IF (diffK4 .NE. 0.) THEN
238     CALL GAD_BIHARM_X(bi,bj,k,xA,df4,diffK4,df,myThid)
239     ENDIF
240    
241 adcroft 1.1 #ifdef ALLOW_GMREDI
242     C- GM/Redi flux in X
243     IF (useGMRedi) THEN
244     C *note* should update GMREDI_XTRANSPORT to use localT and set df *aja*
245     CALL GMREDI_XTRANSPORT(
246     I iMin,iMax,jMin,jMax,bi,bj,K,
247 heimbach 1.15 I xA,Tracer,tracerIdentity,
248 adcroft 1.1 U df,
249     I myThid)
250     ENDIF
251     #endif
252 adcroft 1.5 DO j=1-Oly,sNy+Oly
253     DO i=1-Olx,sNx+Olx
254 adcroft 1.1 fZon(i,j) = fZon(i,j) + df(i,j)
255     ENDDO
256     ENDDO
257    
258 jmc 1.32 #ifdef ALLOW_DIAGNOSTICS
259     C- Diagnostics of Tracer flux in X dir (mainly Diffusive term),
260     C excluding advective terms:
261     IF ( useDiagnostics .AND.
262     & (diffKh.NE.0. .OR. diffK4 .NE.0. .OR. useGMRedi) ) THEN
263     diagName = 'DIFx'//diagSufx
264 jmc 1.33 CALL DIAGNOSTICS_FILL(df,diagName, k,1, 2,bi,bj, myThid)
265 adcroft 1.1 ENDIF
266 jmc 1.32 #endif
267 adcroft 1.1
268     C-- Initialize net flux in Y direction
269     DO j=1-Oly,sNy+Oly
270     DO i=1-Olx,sNx+Olx
271 heimbach 1.12 fMer(i,j) = 0. _d 0
272 adcroft 1.1 ENDDO
273     ENDDO
274    
275     C- Advective flux in Y
276 jmc 1.14 IF (calcAdvection) THEN
277 jmc 1.32 IF (advectionScheme.EQ.ENUM_CENTERED_2ND) THEN
278     CALL GAD_C2_ADV_Y(bi,bj,k,vTrans,localT,af,myThid)
279     ELSEIF (advectionScheme.EQ.ENUM_FLUX_LIMIT) THEN
280     CALL GAD_FLUXLIMIT_ADV_Y( bi,bj,k, dTtracerLev(k),
281     I vTrans, vVel, maskS(1-Olx,1-Oly,k,bi,bj), localT,
282     O af, myThid )
283     ELSEIF (advectionScheme.EQ.ENUM_UPWIND_3RD ) THEN
284     CALL GAD_U3_ADV_Y(bi,bj,k,vTrans,localT,af,myThid)
285     ELSEIF (advectionScheme.EQ.ENUM_CENTERED_4TH) THEN
286     CALL GAD_C4_ADV_Y(bi,bj,k,vTrans,localT,af,myThid)
287     ELSEIF (advectionScheme.EQ.ENUM_DST3 ) THEN
288     CALL GAD_DST3_ADV_Y( bi,bj,k, dTtracerLev(k),
289     I vTrans, vVel, maskS(1-Olx,1-Oly,k,bi,bj), localT,
290     O af, myThid )
291     ELSEIF (advectionScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
292     CALL GAD_DST3FL_ADV_Y( bi,bj,k, dTtracerLev(k),
293     I vTrans, vVel, maskS(1-Olx,1-Oly,k,bi,bj), localT,
294     O af, myThid )
295     ELSE
296     STOP 'GAD_CALC_RHS: Bad advectionScheme (Y)'
297     ENDIF
298     DO j=1-Oly,sNy+Oly
299     DO i=1-Olx,sNx+Olx
300     fMer(i,j) = fMer(i,j) + af(i,j)
301     ENDDO
302     ENDDO
303     #ifdef ALLOW_DIAGNOSTICS
304     IF ( useDiagnostics ) THEN
305     diagName = 'ADVy'//diagSufx
306 jmc 1.33 CALL DIAGNOSTICS_FILL(af,diagName, k,1, 2,bi,bj, myThid)
307 jmc 1.32 ENDIF
308     #endif
309 adcroft 1.8 ENDIF
310 adcroft 1.1
311     C- Diffusive flux in Y
312     IF (diffKh.NE.0.) THEN
313     CALL GAD_DIFF_Y(bi,bj,k,yA,diffKh,localT,df,myThid)
314     ELSE
315     DO j=1-Oly,sNy+Oly
316     DO i=1-Olx,sNx+Olx
317 heimbach 1.12 df(i,j) = 0. _d 0
318 adcroft 1.1 ENDDO
319     ENDDO
320     ENDIF
321    
322 jmc 1.32 C- Add bi-harmonic flux in Y
323     IF (diffK4 .NE. 0.) THEN
324     CALL GAD_BIHARM_Y(bi,bj,k,yA,df4,diffK4,df,myThid)
325     ENDIF
326    
327 adcroft 1.1 #ifdef ALLOW_GMREDI
328     C- GM/Redi flux in Y
329     IF (useGMRedi) THEN
330 heimbach 1.7 C *note* should update GMREDI_YTRANSPORT to use localT and set df *aja*
331 adcroft 1.1 CALL GMREDI_YTRANSPORT(
332     I iMin,iMax,jMin,jMax,bi,bj,K,
333 heimbach 1.15 I yA,Tracer,tracerIdentity,
334 adcroft 1.1 U df,
335     I myThid)
336     ENDIF
337     #endif
338     DO j=1-Oly,sNy+Oly
339     DO i=1-Olx,sNx+Olx
340     fMer(i,j) = fMer(i,j) + df(i,j)
341     ENDDO
342     ENDDO
343    
344 jmc 1.32 #ifdef ALLOW_DIAGNOSTICS
345     C- Diagnostics of Tracer flux in Y dir (mainly Diffusive terms),
346     C excluding advective terms:
347     IF ( useDiagnostics .AND.
348     & (diffKh.NE.0. .OR. diffK4 .NE.0. .OR. useGMRedi) ) THEN
349     diagName = 'DIFy'//diagSufx
350 jmc 1.33 CALL DIAGNOSTICS_FILL(df,diagName, k,1, 2,bi,bj, myThid)
351 adcroft 1.1 ENDIF
352 jmc 1.32 #endif
353 adcroft 1.1
354 jmc 1.16 C-- Compute vertical flux fVerT(kUp) at interface k (between k-1 & k):
355 adcroft 1.1 C- Advective flux in R
356 jmc 1.25 #ifdef ALLOW_AIM
357     C- a hack to prevent Water-Vapor vert.transport into the stratospheric level Nr
358     IF (calcAdvection .AND. .NOT.implicitAdvection .AND. K.GE.2 .AND.
359     & (.NOT.useAIM .OR.tracerIdentity.NE.GAD_SALINITY .OR.K.LT.Nr)
360     & ) THEN
361     #else
362 jmc 1.23 IF (calcAdvection .AND. .NOT.implicitAdvection .AND. K.GE.2) THEN
363 jmc 1.25 #endif
364 jmc 1.2 C- Compute vertical advective flux in the interior:
365 jmc 1.32 IF (vertAdvecScheme.EQ.ENUM_CENTERED_2ND) THEN
366     CALL GAD_C2_ADV_R(bi,bj,k,rTrans,tracer,af,myThid)
367     ELSEIF (vertAdvecScheme.EQ.ENUM_FLUX_LIMIT) THEN
368     CALL GAD_FLUXLIMIT_ADV_R(
369     & bi,bj,k,dTtracerLev(k),rTrans,wVel,tracer,af,myThid)
370     ELSEIF (vertAdvecScheme.EQ.ENUM_UPWIND_3RD ) THEN
371     CALL GAD_U3_ADV_R(bi,bj,k,rTrans,tracer,af,myThid)
372     ELSEIF (vertAdvecScheme.EQ.ENUM_CENTERED_4TH) THEN
373     CALL GAD_C4_ADV_R(bi,bj,k,rTrans,tracer,af,myThid)
374     ELSEIF (vertAdvecScheme.EQ.ENUM_DST3 ) THEN
375     CALL GAD_DST3_ADV_R(
376     & bi,bj,k,dTtracerLev(k),rTrans,wVel,tracer,af,myThid)
377     ELSEIF (vertAdvecScheme.EQ.ENUM_DST3_FLUX_LIMIT ) THEN
378     CALL GAD_DST3FL_ADV_R(
379     & bi,bj,k,dTtracerLev(k),rTrans,wVel,tracer,af,myThid)
380     ELSE
381     STOP 'GAD_CALC_RHS: Bad vertAdvecScheme (R)'
382     ENDIF
383 jmc 1.23 C- add the advective flux to fVerT
384 jmc 1.32 DO j=1-Oly,sNy+Oly
385     DO i=1-Olx,sNx+Olx
386     fVerT(i,j,kUp) = fVerT(i,j,kUp) + af(i,j)
387     ENDDO
388 jmc 1.2 ENDDO
389 jmc 1.32 #ifdef ALLOW_DIAGNOSTICS
390     IF ( useDiagnostics ) THEN
391     diagName = 'ADVr'//diagSufx
392 jmc 1.33 CALL DIAGNOSTICS_FILL(af,diagName, k,1, 2,bi,bj, myThid)
393 jmc 1.34 C- note: needs to explicitly increment the counter since DIAGNOSTICS_FILL
394     C does it only if k=1 (never the case here)
395     IF ( k.EQ.2 ) CALL DIAGNOSTICS_COUNT(diagName,bi,bj,myThid)
396 jmc 1.32 ENDIF
397     #endif
398 adcroft 1.8 ENDIF
399 adcroft 1.1
400     C- Diffusive flux in R
401     C Note: For K=1 then KM1=1 and this gives a dT/dr = 0 upper
402     C boundary condition.
403     IF (implicitDiffusion) THEN
404 adcroft 1.5 DO j=1-Oly,sNy+Oly
405     DO i=1-Olx,sNx+Olx
406 heimbach 1.12 df(i,j) = 0. _d 0
407 adcroft 1.1 ENDDO
408     ENDDO
409     ELSE
410 jmc 1.30 CALL GAD_DIFF_R(bi,bj,k,KappaR,tracer,df,myThid)
411 adcroft 1.1 ENDIF
412    
413     #ifdef ALLOW_GMREDI
414     C- GM/Redi flux in R
415     IF (useGMRedi) THEN
416     C *note* should update GMREDI_RTRANSPORT to set df *aja*
417     CALL GMREDI_RTRANSPORT(
418     I iMin,iMax,jMin,jMax,bi,bj,K,
419 heimbach 1.15 I Tracer,tracerIdentity,
420 adcroft 1.1 U df,
421     I myThid)
422     ENDIF
423     #endif
424    
425 adcroft 1.5 DO j=1-Oly,sNy+Oly
426     DO i=1-Olx,sNx+Olx
427 adcroft 1.11 fVerT(i,j,kUp) = fVerT(i,j,kUp) + df(i,j)*maskUp(i,j)
428 adcroft 1.1 ENDDO
429     ENDDO
430    
431 jmc 1.32 #ifdef ALLOW_DIAGNOSTICS
432     C- Diagnostics of Tracer flux in R dir (mainly Diffusive terms),
433     C Explicit terms only & excluding advective terms:
434     IF ( useDiagnostics .AND.
435     & (.NOT.implicitDiffusion .OR. useGMRedi) ) THEN
436     diagName = 'DFrE'//diagSufx
437 jmc 1.33 CALL DIAGNOSTICS_FILL(df,diagName, k,1, 2,bi,bj, myThid)
438 jmc 1.32 ENDIF
439     #endif
440    
441 adcroft 1.1 #ifdef ALLOW_KPP
442 jmc 1.29 C- Set non local KPP transport term (ghat):
443     IF ( useKPP .AND. k.GE.2 ) THEN
444 adcroft 1.5 DO j=1-Oly,sNy+Oly
445     DO i=1-Olx,sNx+Olx
446 heimbach 1.12 df(i,j) = 0. _d 0
447 adcroft 1.1 ENDDO
448     ENDDO
449     IF (tracerIdentity.EQ.GAD_TEMPERATURE) THEN
450     CALL KPP_TRANSPORT_T(
451     I iMin,iMax,jMin,jMax,bi,bj,k,km1,
452 jmc 1.29 O df )
453 adcroft 1.1 ELSEIF (tracerIdentity.EQ.GAD_SALINITY) THEN
454     CALL KPP_TRANSPORT_S(
455     I iMin,iMax,jMin,jMax,bi,bj,k,km1,
456 jmc 1.29 O df )
457 mlosch 1.18 #ifdef ALLOW_PTRACERS
458 dimitri 1.22 ELSEIF (tracerIdentity .GE. GAD_TR1) THEN
459 mlosch 1.18 CALL KPP_TRANSPORT_PTR(
460     I iMin,iMax,jMin,jMax,bi,bj,k,km1,
461 jmc 1.29 I tracerIdentity-GAD_TR1+1,
462     O df )
463 mlosch 1.18 #endif
464 adcroft 1.1 ELSE
465 mlosch 1.18 PRINT*,'invalid tracer indentity: ', tracerIdentity
466 adcroft 1.1 STOP 'GAD_CALC_RHS: Ooops'
467     ENDIF
468 adcroft 1.5 DO j=1-Oly,sNy+Oly
469     DO i=1-Olx,sNx+Olx
470 adcroft 1.11 fVerT(i,j,kUp) = fVerT(i,j,kUp) + df(i,j)*maskUp(i,j)
471 adcroft 1.1 ENDDO
472     ENDDO
473     ENDIF
474     #endif
475    
476     C-- Divergence of fluxes
477 adcroft 1.10 DO j=1-Oly,sNy+Oly-1
478     DO i=1-Olx,sNx+Olx-1
479 adcroft 1.8 gTracer(i,j,k,bi,bj)=gTracer(i,j,k,bi,bj)
480 jmc 1.23 & -_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)*recip_rA(i,j,bi,bj)
481     & *( (fZon(i+1,j)-fZon(i,j))
482     & +(fMer(i,j+1)-fMer(i,j))
483     & +(fVerT(i,j,kUp)-fVerT(i,j,kDown))*rkFac
484     & -localT(i,j)*( (uTrans(i+1,j)-uTrans(i,j))
485     & +(vTrans(i,j+1)-vTrans(i,j))
486     & +(rTrans(i,j)-rTransKp1(i,j))*rAdvFac
487     & )*advFac
488 adcroft 1.1 & )
489     ENDDO
490     ENDDO
491    
492 jmc 1.27 #ifdef ALLOW_DEBUG
493     IF ( debugLevel .GE. debLevB
494 jmc 1.28 & .AND. tracerIdentity.EQ.GAD_TEMPERATURE
495 jmc 1.27 & .AND. k.EQ.2 .AND. myIter.EQ.1+nIter0
496     & .AND. nPx.EQ.1 .AND. nPy.EQ.1
497     & .AND. useCubedSphereExchange ) THEN
498     CALL DEBUG_CS_CORNER_UV( ' fZon,fMer from GAD_CALC_RHS',
499     & fZon,fMer, k, standardMessageUnit,bi,bj,myThid )
500     ENDIF
501     #endif /* ALLOW_DEBUG */
502    
503 adcroft 1.1 RETURN
504     END

  ViewVC Help
Powered by ViewVC 1.1.22