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

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

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


Revision 1.52 - (show annotations) (download)
Wed Apr 23 18:32:20 2008 UTC (16 years, 1 month ago) by jahn
Branch: MAIN
CVS Tags: checkpoint60, checkpoint61, checkpoint59q, checkpoint59r, checkpoint61f, checkpoint61n, checkpoint61q, checkpoint61e, checkpoint61g, checkpoint61d, checkpoint61b, checkpoint61c, checkpoint61a, checkpoint61o, checkpoint61l, checkpoint61m, checkpoint61j, checkpoint61k, checkpoint61h, checkpoint61i, checkpoint61r, checkpoint61p
Changes since 1.51: +47 -22 lines
fix faulty Smolarkiewicz (positive Redi) hack

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

  ViewVC Help
Powered by ViewVC 1.1.22