/[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.70 - (hide annotations) (download)
Mon Aug 18 12:22:46 2014 UTC (9 years, 9 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint65c
Changes since 1.69: +26 -29 lines
- change Tracer argument (drop bi,bj indices) in S/R GAD_CALC_RHS,
  GAD_C2_ADV_R, GAD_U3_ADV_R, GAD_C4_ADV_R, GAD_DIFF_R, GAD_BIHARM_R
  + also in GMREDI_X/Y/RTRANSPORT ; and update corresponding calls in
  S/R temp/salt/ptracers_integrate.F

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

  ViewVC Help
Powered by ViewVC 1.1.22