/[MITgcm]/MITgcm/model/src/salt_integrate.F
ViewVC logotype

Annotation of /MITgcm/model/src/salt_integrate.F

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


Revision 1.18 - (hide annotations) (download)
Sun Oct 9 18:13:09 2016 UTC (7 years, 8 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, HEAD
Changes since 1.17: +4 -2 lines
- with INCLUDE_IMPLVERTADV_CODE defined, also call MOM_U,V_IMPLICIT_R &
  GAD_IMPLICIT_R  (instead of IMPLDIFF) when just implicitViscosity and
  implicitDiffusion (respectively) are used (even without momImplVertAdv
  or temp,salt,PTRACERS_ImplVertAdv).

1 jmc 1.18 C $Header: /u/gcmpack/MITgcm/model/src/salt_integrate.F,v 1.17 2014/09/05 21:07:14 jmc Exp $
2 jmc 1.7 C $Name: $
3 jmc 1.1
4     #include "PACKAGES_CONFIG.h"
5     #include "CPP_OPTIONS.h"
6 jmc 1.7 #ifdef ALLOW_AUTODIFF
7     # include "AUTODIFF_OPTIONS.h"
8     #endif
9 jmc 1.4 #ifdef ALLOW_GENERIC_ADVDIFF
10     # include "GAD_OPTIONS.h"
11     #endif
12 jmc 1.1
13     CBOP
14     C !ROUTINE: SALT_INTEGRATE
15     C !INTERFACE:
16     SUBROUTINE SALT_INTEGRATE(
17 jmc 1.5 I bi, bj, recip_hFac,
18 jmc 1.4 I uFld, vFld, wFld,
19     U KappaRk,
20 jmc 1.1 I myTime, myIter, myThid )
21     C !DESCRIPTION: \bv
22     C *==========================================================*
23     C | SUBROUTINE SALT_INTEGRATE
24 jmc 1.12 C | o Calculate tendency for salinity and integrates
25     C | forward in time. The salinity array is updated here
26     C | while adjustments (filters, conv.adjustment) are applied
27     C | later, in S/R TRACERS_CORRECTION_STEP.
28 jmc 1.1 C *==========================================================*
29 jmc 1.8 C | A procedure called APPLY_FORCING_S is called from
30 jmc 1.1 C | here. These procedures can be used to add per problem
31     C | E-P flux source terms.
32     C | Note: Although it is slightly counter-intuitive the
33     C | EXTERNAL_FORCING routine is not the place to put
34     C | file I/O. Instead files that are required to
35     C | calculate the external source terms are generally
36     C | read during the model main loop. This makes the
37     C | logistics of multi-processing simpler and also
38     C | makes the adjoint generation simpler. It also
39     C | allows for I/O to overlap computation where that
40     C | is supported by hardware.
41     C | Aside from the problem specific term the code here
42     C | forms the tendency terms due to advection and mixing
43     C | The baseline implementation here uses a centered
44     C | difference form for the advection term and a tensorial
45     C | divergence of a flux form for the diffusive term. The
46 jmc 1.4 C | diffusive term is formulated so that isopycnal mixing
47     C | and GM-style subgrid-scale terms can be incorporated by
48 jmc 1.1 C | simply setting the diffusion tensor terms appropriately.
49     C *==========================================================*
50     C \ev
51    
52     C !USES:
53     IMPLICIT NONE
54     C == GLobal variables ==
55     #include "SIZE.h"
56     #include "EEPARAMS.h"
57     #include "PARAMS.h"
58 jmc 1.5 #include "GRID.h"
59     #include "DYNVARS.h"
60 jmc 1.1 #include "RESTART.h"
61     #ifdef ALLOW_GENERIC_ADVDIFF
62 jmc 1.4 # include "GAD.h"
63     # include "GAD_SOM_VARS.h"
64 jmc 1.1 #endif
65 jmc 1.7 #ifdef ALLOW_AUTODIFF
66 jmc 1.1 # include "tamc.h"
67     # include "tamc_keys.h"
68     #endif
69    
70     C !INPUT/OUTPUT PARAMETERS:
71     C == Routine arguments ==
72 jmc 1.5 C bi, bj, :: tile indices
73     C recip_hFac :: reciprocal of cell open-depth factor (@ next iter)
74     C uFld,vFld :: Local copy of horizontal velocity field
75     C wFld :: Local copy of vertical velocity field
76     C KappaRk :: Vertical diffusion for Salinity
77     C myTime :: current time
78     C myIter :: current iteration number
79     C myThid :: my Thread Id. number
80 jmc 1.4 INTEGER bi, bj
81 jmc 1.5 _RS recip_hFac(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
82     _RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
83     _RL vFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
84     _RL wFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
85     _RL KappaRk (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
86 jmc 1.1 _RL myTime
87     INTEGER myIter
88     INTEGER myThid
89     CEOP
90    
91     #ifdef ALLOW_GENERIC_ADVDIFF
92 jmc 1.10 #ifdef ALLOW_DIAGNOSTICS
93     C !FUNCTIONS:
94     LOGICAL DIAGNOSTICS_IS_ON
95     EXTERNAL DIAGNOSTICS_IS_ON
96     #endif /* ALLOW_DIAGNOSTICS */
97    
98 jmc 1.1 C !LOCAL VARIABLES:
99 jmc 1.5 C iMin, iMax :: 1rst index loop range
100     C jMin, jMax :: 2nd index loop range
101     C k :: vertical index
102     C kM1 :: =k-1 for k>1, =1 for k=1
103     C kUp :: index into 2 1/2D array, toggles between 1|2
104     C kDown :: index into 2 1/2D array, toggles between 2|1
105     C xA :: Tracer cell face area normal to X
106     C yA :: Tracer cell face area normal to X
107     C maskUp :: Land/water mask for Wvel points (interface k)
108     C uTrans :: Zonal volume transport through cell face
109     C vTrans :: Meridional volume transport through cell face
110     C rTrans :: Vertical volume transport at interface k
111     C rTransKp :: Vertical volume transport at inteface k+1
112     C fZon :: Flux of salt (S) in the zonal direction
113     C fMer :: Flux of salt (S) in the meridional direction
114     C fVer :: Flux of salt (S) in the vertical direction
115     C at the upper(U) and lower(D) faces of a cell.
116 jmc 1.14 C gS_loc :: Salinity tendency (local to this S/R)
117 jmc 1.10 C gsForc :: Salinity forcing tendency
118 jmc 1.8 C gs_AB :: Adams-Bashforth salinity tendency increment
119 jmc 1.4 INTEGER iMin, iMax, jMin, jMax
120 jmc 1.1 INTEGER i, j, k
121     INTEGER kUp, kDown, kM1
122     _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
123     _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
124     _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
125     _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
126     _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
127     _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
128     _RL rTransKp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
129 jmc 1.5 _RL fZon (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
130     _RL fMer (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
131     _RL fVer (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
132 jmc 1.14 _RL gS_loc (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
133 jmc 1.10 _RL gsForc (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
134 jmc 1.1 _RL gs_AB (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
135 jmc 1.10 #ifdef ALLOW_DIAGNOSTICS
136     LOGICAL diagForcing, diagAB_tend
137     #endif
138 jmc 1.1 LOGICAL calcAdvection
139     INTEGER iterNb
140     #ifdef ALLOW_ADAMSBASHFORTH_3
141 jmc 1.11 INTEGER m2
142 jmc 1.1 #endif
143 jmc 1.4 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
144    
145     iterNb = myIter
146     IF (staggerTimeStep) iterNb = myIter - 1
147 jmc 1.1
148 jmc 1.17 C- Loop ranges for daughter routines
149     c iMin = 1
150     c iMax = sNx
151     c jMin = 1
152     c jMax = sNy
153     C Regarding model dynamics, only needs to get correct tracer tendency
154     C (gS_loc) in tile interior (1:sNx,1:sNy);
155     C However, for some diagnostics, we may want to get valid tendency
156     C extended over 1 point in tile halo region (0:sNx+1,0:sNy=1).
157     iMin = 0
158     iMax = sNx+1
159     jMin = 0
160     jMax = sNy+1
161    
162 jmc 1.10 #ifdef ALLOW_DIAGNOSTICS
163     diagForcing = .FALSE.
164     diagAB_tend = .FALSE.
165     IF ( useDiagnostics .AND. saltForcing )
166     & diagForcing = DIAGNOSTICS_IS_ON( 'gS_Forc ', myThid )
167     IF ( useDiagnostics .AND. AdamsBashforthGs )
168     & diagAB_tend = DIAGNOSTICS_IS_ON( 'AB_gS ', myThid )
169     #endif
170    
171 jmc 1.1 #ifdef ALLOW_AUTODIFF_TAMC
172     act1 = bi - myBxLo(myThid)
173     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
174     act2 = bj - myByLo(myThid)
175     max2 = myByHi(myThid) - myByLo(myThid) + 1
176     act3 = myThid - 1
177     max3 = nTx*nTy
178     act4 = ikey_dynamics - 1
179     itdkey = (act1 + 1) + act2*max1
180     & + act3*max1*max2
181     & + act4*max1*max2*max3
182     #endif /* ALLOW_AUTODIFF_TAMC */
183    
184 jmc 1.16 C- Apply AB on S :
185     IF ( AdamsBashforth_S ) THEN
186     C compute S^n+1/2 (stored in gsNm) extrapolating S forward in time
187     #ifdef ALLOW_ADAMSBASHFORTH_3
188     c m1 = 1 + MOD(iterNb+1,2)
189     c m2 = 1 + MOD( iterNb ,2)
190     CALL ADAMS_BASHFORTH3(
191     I bi, bj, 0, Nr,
192     I salt(1-OLx,1-OLy,1,bi,bj),
193     U gsNm, gs_AB,
194     I saltStartAB, iterNb, myThid )
195     #else /* ALLOW_ADAMSBASHFORTH_3 */
196     CALL ADAMS_BASHFORTH2(
197     I bi, bj, 0, Nr,
198     I salt(1-OLx,1-OLy,1,bi,bj),
199     U gsNm1(1-OLx,1-OLy,1,bi,bj), gs_AB,
200     I saltStartAB, iterNb, myThid )
201     #endif /* ALLOW_ADAMSBASHFORTH_3 */
202     ENDIF
203    
204 jmc 1.4 C- Tracer tendency needs to be set to zero (moved here from gad_calc_rhs):
205     DO k=1,Nr
206     DO j=1-OLy,sNy+OLy
207     DO i=1-OLx,sNx+OLx
208 jmc 1.14 gS_loc(i,j,k) = 0. _d 0
209 jmc 1.4 ENDDO
210     ENDDO
211     ENDDO
212 jmc 1.1 DO j=1-OLy,sNy+OLy
213     DO i=1-OLx,sNx+OLx
214 jmc 1.5 fVer(i,j,1) = 0. _d 0
215     fVer(i,j,2) = 0. _d 0
216 jmc 1.1 ENDDO
217     ENDDO
218 jmc 1.4 #ifdef ALLOW_AUTODIFF
219     DO k=1,Nr
220     DO j=1-OLy,sNy+OLy
221     DO i=1-OLx,sNx+OLx
222     kappaRk(i,j,k) = 0. _d 0
223     ENDDO
224     ENDDO
225     ENDDO
226 jmc 1.16 CADJ STORE wFld(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
227     CADJ STORE salt(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
228     # ifdef ALLOW_ADAMSBASHFORTH_3
229     CADJ STORE gsNm(:,:,:,bi,bj,1) = comlev1_bibj, key=itdkey, byte=isbyte
230     CADJ STORE gsNm(:,:,:,bi,bj,2) = comlev1_bibj, key=itdkey, byte=isbyte
231     # else
232     CADJ STORE gsNm1(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
233     # endif
234 jmc 1.4 #endif /* ALLOW_AUTODIFF */
235 jmc 1.1
236 jmc 1.4 #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
237 jmc 1.5 CALL CALC_3D_DIFFUSIVITY(
238 jmc 1.4 I bi, bj, iMin, iMax, jMin, jMax,
239     I GAD_SALINITY, useGMredi, useKPP,
240     O kappaRk,
241     I myThid )
242     #endif /* INCLUDE_CALC_DIFFUSIVITY_CALL */
243    
244     #ifndef DISABLE_MULTIDIM_ADVECTION
245     C-- Some advection schemes are better calculated using a multi-dimensional
246     C method in the absence of any other terms and, if used, is done here.
247     C
248     C The CPP flag DISABLE_MULTIDIM_ADVECTION is currently unset in GAD_OPTIONS.h
249     C The default is to use multi-dimensinal advection for non-linear advection
250     C schemes. However, for the sake of efficiency of the adjoint it is necessary
251     C to be able to exclude this scheme to avoid excessive storage and
252     C recomputation. It *is* differentiable, if you need it.
253     C Edit GAD_OPTIONS.h and #define DISABLE_MULTIDIM_ADVECTION to
254     C disable this section of code.
255     #ifdef GAD_ALLOW_TS_SOM_ADV
256     # ifdef ALLOW_AUTODIFF_TAMC
257     CADJ STORE som_S = comlev1_bibj, key=itdkey, byte=isbyte
258     # endif
259     IF ( saltSOM_Advection ) THEN
260     # ifdef ALLOW_DEBUG
261     IF (debugMode) CALL DEBUG_CALL('GAD_SOM_ADVECT',myThid)
262     # endif
263     CALL GAD_SOM_ADVECT(
264 jmc 1.8 I saltImplVertAdv,
265     I saltAdvScheme, saltVertAdvScheme, GAD_SALINITY,
266     I dTtracerLev, uFld, vFld, wFld, salt,
267 jmc 1.4 U som_S,
268 jmc 1.14 O gS_loc,
269 jmc 1.4 I bi, bj, myTime, myIter, myThid )
270     ELSEIF (saltMultiDimAdvec) THEN
271     #else /* GAD_ALLOW_TS_SOM_ADV */
272     IF (saltMultiDimAdvec) THEN
273     #endif /* GAD_ALLOW_TS_SOM_ADV */
274     # ifdef ALLOW_DEBUG
275     IF (debugMode) CALL DEBUG_CALL('GAD_ADVECTION',myThid)
276     # endif
277     CALL GAD_ADVECTION(
278 jmc 1.8 I saltImplVertAdv,
279     I saltAdvScheme, saltVertAdvScheme, GAD_SALINITY,
280     I dTtracerLev, uFld, vFld, wFld, salt,
281 jmc 1.14 O gS_loc,
282 jmc 1.4 I bi, bj, myTime, myIter, myThid )
283     ENDIF
284     #endif /* DISABLE_MULTIDIM_ADVECTION */
285    
286     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
287    
288     C- Start vertical index (k) loop (Nr:1)
289     calcAdvection = saltAdvection .AND. .NOT.saltMultiDimAdvec
290 jmc 1.1 DO k=Nr,1,-1
291     #ifdef ALLOW_AUTODIFF_TAMC
292     kkey = (itdkey-1)*Nr + k
293     #endif
294     kM1 = MAX(1,k-1)
295     kUp = 1+MOD(k+1,2)
296     kDown= 1+MOD(k,2)
297    
298     #ifdef ALLOW_AUTODIFF_TAMC
299 jmc 1.5 CADJ STORE fVer(:,:,:) = comlev1_bibj_k, key=kkey,
300 jmc 1.1 CADJ & byte=isbyte, kind = isbyte
301 jmc 1.14 CADJ STORE gS_loc(:,:,k) = comlev1_bibj_k, key=kkey,
302 jmc 1.1 CADJ & byte=isbyte, kind = isbyte
303     # ifdef ALLOW_ADAMSBASHFORTH_3
304     CADJ STORE gsNm(:,:,k,bi,bj,1) = comlev1_bibj_k, key=kkey,
305     CADJ & byte=isbyte, kind = isbyte
306     CADJ STORE gsNm(:,:,k,bi,bj,2) = comlev1_bibj_k, key=kkey,
307     CADJ & kind = isbyte
308     # else
309     CADJ STORE gsNm1(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey,
310     CADJ & byte=isbyte, kind = isbyte
311     # endif
312     #endif /* ALLOW_AUTODIFF_TAMC */
313     CALL CALC_ADV_FLOW(
314     I uFld, vFld, wFld,
315     U rTrans,
316     O uTrans, vTrans, rTransKp,
317     O maskUp, xA, yA,
318     I k, bi, bj, myThid )
319    
320 jmc 1.10 C-- Collect forcing term in local array gsForc:
321     DO j=1-OLy,sNy+OLy
322     DO i=1-OLx,sNx+OLx
323     gsForc(i,j) = 0. _d 0
324     ENDDO
325     ENDDO
326     IF ( saltForcing ) THEN
327     CALL APPLY_FORCING_S(
328     U gsForc,
329     I iMin,iMax,jMin,jMax, k, bi,bj,
330     I myTime, myIter, myThid )
331     #ifdef ALLOW_DIAGNOSTICS
332     IF ( diagForcing ) THEN
333     CALL DIAGNOSTICS_FILL(gsForc,'gS_Forc ',k,1,2,bi,bj,myThid)
334     ENDIF
335     #endif /* ALLOW_DIAGNOSTICS */
336     ENDIF
337    
338 jmc 1.1 #ifdef ALLOW_ADAMSBASHFORTH_3
339 jmc 1.11 c m1 = 1 + MOD(iterNb+1,2)
340 jmc 1.1 m2 = 1 + MOD( iterNb ,2)
341     CALL GAD_CALC_RHS(
342     I bi, bj, iMin,iMax,jMin,jMax, k, kM1, kUp, kDown,
343     I xA, yA, maskUp, uFld(1-OLx,1-OLy,k),
344     I vFld(1-OLx,1-OLy,k), wFld(1-OLx,1-OLy,k),
345     I uTrans, vTrans, rTrans, rTransKp,
346     I diffKhS, diffK4S, KappaRk(1-OLx,1-OLy,k), diffKr4S,
347 jmc 1.15 I salt(1-OLx,1-OLy,1,bi,bj),
348     I gsNm(1-OLx,1-OLy,1,bi,bj,m2), dTtracerLev,
349 jmc 1.1 I GAD_SALINITY, saltAdvScheme, saltVertAdvScheme,
350     I calcAdvection, saltImplVertAdv, AdamsBashforth_S,
351     I saltVertDiff4, useGMRedi, useKPP,
352 jmc 1.5 O fZon, fMer,
353 jmc 1.14 U fVer, gS_loc,
354 jmc 1.1 I myTime, myIter, myThid )
355     #else /* ALLOW_ADAMSBASHFORTH_3 */
356     CALL GAD_CALC_RHS(
357     I bi, bj, iMin,iMax,jMin,jMax, k, kM1, kUp, kDown,
358     I xA, yA, maskUp, uFld(1-OLx,1-OLy,k),
359     I vFld(1-OLx,1-OLy,k), wFld(1-OLx,1-OLy,k),
360     I uTrans, vTrans, rTrans, rTransKp,
361     I diffKhS, diffK4S, KappaRk(1-OLx,1-OLy,k), diffKr4S,
362 jmc 1.15 I salt(1-OLx,1-OLy,1,bi,bj),
363     I gsNm1(1-OLx,1-OLy,1,bi,bj), dTtracerLev,
364 jmc 1.1 I GAD_SALINITY, saltAdvScheme, saltVertAdvScheme,
365     I calcAdvection, saltImplVertAdv, AdamsBashforth_S,
366     I saltVertDiff4, useGMRedi, useKPP,
367 jmc 1.5 O fZon, fMer,
368 jmc 1.14 U fVer, gS_loc,
369 jmc 1.1 I myTime, myIter, myThid )
370     #endif /* ALLOW_ADAMSBASHFORTH_3 */
371    
372     C-- External salinity forcing term(s) inside Adams-Bashforth:
373 jmc 1.10 IF ( saltForcing .AND. tracForcingOutAB.NE.1 ) THEN
374     DO j=1-OLy,sNy+OLy
375     DO i=1-OLx,sNx+OLx
376 jmc 1.14 gS_loc(i,j,k) = gS_loc(i,j,k) + gsForc(i,j)
377 jmc 1.10 ENDDO
378     ENDDO
379     ENDIF
380 jmc 1.1
381     IF ( AdamsBashforthGs ) THEN
382     #ifdef ALLOW_ADAMSBASHFORTH_3
383     CALL ADAMS_BASHFORTH3(
384     I bi, bj, k, Nr,
385 jmc 1.14 U gS_loc, gsNm,
386     O gs_AB,
387 jmc 1.1 I saltStartAB, iterNb, myThid )
388     #else
389     CALL ADAMS_BASHFORTH2(
390     I bi, bj, k, Nr,
391 jmc 1.14 U gS_loc, gsNm1(1-OLx,1-OLy,1,bi,bj),
392     O gs_AB,
393 jmc 1.1 I saltStartAB, iterNb, myThid )
394     #endif
395     #ifdef ALLOW_DIAGNOSTICS
396 jmc 1.10 IF ( diagAB_tend ) THEN
397 jmc 1.1 CALL DIAGNOSTICS_FILL(gs_AB,'AB_gS ',k,1,2,bi,bj,myThid)
398     ENDIF
399     #endif /* ALLOW_DIAGNOSTICS */
400     ENDIF
401    
402     C-- External salinity forcing term(s) outside Adams-Bashforth:
403 jmc 1.10 IF ( saltForcing .AND. tracForcingOutAB.EQ.1 ) THEN
404     DO j=1-OLy,sNy+OLy
405     DO i=1-OLx,sNx+OLx
406 jmc 1.14 gS_loc(i,j,k) = gS_loc(i,j,k) + gsForc(i,j)
407 jmc 1.10 ENDDO
408     ENDDO
409     ENDIF
410 jmc 1.1
411     #ifdef NONLIN_FRSURF
412     IF (nonlinFreeSurf.GT.0) THEN
413     CALL FREESURF_RESCALE_G(
414     I bi, bj, k,
415 jmc 1.14 U gS_loc,
416 jmc 1.1 I myThid )
417     IF ( AdamsBashforthGs ) THEN
418     #ifdef ALLOW_ADAMSBASHFORTH_3
419     # ifdef ALLOW_AUTODIFF_TAMC
420     CADJ STORE gsNm(:,:,k,bi,bj,1) = comlev1_bibj_k, key=kkey,
421     CADJ & byte=isbyte, kind = isbyte
422     CADJ STORE gsNm(:,:,k,bi,bj,2) = comlev1_bibj_k, key=kkey,
423     CADJ & kind = isbyte
424     # endif
425     CALL FREESURF_RESCALE_G(
426     I bi, bj, k,
427 jmc 1.14 U gsNm(1-OLx,1-OLy,1,bi,bj,1),
428 jmc 1.1 I myThid )
429     CALL FREESURF_RESCALE_G(
430     I bi, bj, k,
431 jmc 1.14 U gsNm(1-OLx,1-OLy,1,bi,bj,2),
432 jmc 1.1 I myThid )
433     #else
434     CALL FREESURF_RESCALE_G(
435     I bi, bj, k,
436 jmc 1.14 U gsNm1(1-OLx,1-OLy,1,bi,bj),
437 jmc 1.1 I myThid )
438     #endif
439     ENDIF
440     ENDIF
441     #endif /* NONLIN_FRSURF */
442    
443     C- end of vertical index (k) loop (Nr:1)
444     ENDDO
445    
446 jmc 1.5 #ifdef ALLOW_DOWN_SLOPE
447     IF ( useDOWN_SLOPE ) THEN
448     IF ( usingPCoords ) THEN
449     CALL DWNSLP_APPLY(
450     I GAD_SALINITY, bi, bj, kSurfC,
451 jmc 1.13 I salt(1-OLx,1-OLy,1,bi,bj),
452 jmc 1.14 U gS_loc,
453 jmc 1.13 I recip_hFac, recip_rA, recip_drF,
454     I dTtracerLev, myTime, myIter, myThid )
455 jmc 1.5 ELSE
456     CALL DWNSLP_APPLY(
457     I GAD_SALINITY, bi, bj, kLowC,
458 jmc 1.13 I salt(1-OLx,1-OLy,1,bi,bj),
459 jmc 1.14 U gS_loc,
460 jmc 1.13 I recip_hFac, recip_rA, recip_drF,
461     I dTtracerLev, myTime, myIter, myThid )
462 jmc 1.5 ENDIF
463     ENDIF
464     #endif /* ALLOW_DOWN_SLOPE */
465    
466 jmc 1.14 C- Integrate forward in time, storing in gS_loc: gS <= S + dt*gS
467 jmc 1.13 CALL TIMESTEP_TRACER(
468     I bi, bj, dTtracerLev,
469     I salt(1-OLx,1-OLy,1,bi,bj),
470 jmc 1.14 U gS_loc,
471 jmc 1.13 I myTime, myIter, myThid )
472    
473 jmc 1.5 C-- Implicit vertical advection & diffusion
474    
475     #ifdef INCLUDE_IMPLVERTADV_CODE
476 jmc 1.18 IF ( saltImplVertAdv .OR. implicitDiffusion ) THEN
477     C to recover older (prior to 2016-10-05) results:
478     c IF ( saltImplVertAdv ) THEN
479 jmc 1.5 #ifdef ALLOW_AUTODIFF_TAMC
480 jmc 1.7 CADJ STORE kappaRk(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
481 jmc 1.14 CADJ STORE gS_loc(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
482 jmc 1.7 CADJ STORE wFld(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
483     CADJ STORE salt(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
484     CADJ STORE recip_hFac(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
485 jmc 1.5 #endif /* ALLOW_AUTODIFF_TAMC */
486     CALL GAD_IMPLICIT_R(
487     I saltImplVertAdv, saltVertAdvScheme, GAD_SALINITY,
488 jmc 1.14 I dTtracerLev, kappaRk, recip_hFac, wFld,
489     I salt(1-OLx,1-OLy,1,bi,bj),
490     U gS_loc,
491 jmc 1.5 I bi, bj, myTime, myIter, myThid )
492     ELSEIF ( implicitDiffusion ) THEN
493     #else /* INCLUDE_IMPLVERTADV_CODE */
494     IF ( implicitDiffusion ) THEN
495     #endif /* INCLUDE_IMPLVERTADV_CODE */
496     #ifdef ALLOW_AUTODIFF_TAMC
497 jmc 1.7 CADJ STORE kappaRk(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
498 jmc 1.14 CADJ STORE gS_loc(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
499 jmc 1.5 #endif /* ALLOW_AUTODIFF_TAMC */
500     CALL IMPLDIFF(
501     I bi, bj, iMin, iMax, jMin, jMax,
502     I GAD_SALINITY, kappaRk, recip_hFac,
503 jmc 1.14 U gS_loc,
504 jmc 1.5 I myThid )
505     ENDIF
506    
507 jmc 1.12 IF ( AdamsBashforth_S ) THEN
508     C- Save current tracer field (for AB on tracer) and then update tracer
509 jmc 1.16 #ifdef ALLOW_ADAMSBASHFORTH_3
510 jmc 1.12 CALL CYCLE_AB_TRACER(
511 jmc 1.15 I bi, bj, gS_loc,
512     U salt(1-OLx,1-OLy,1,bi,bj),
513     O gsNm(1-OLx,1-OLy,1,bi,bj,m2),
514 jmc 1.12 I myTime, myIter, myThid )
515     #else /* ALLOW_ADAMSBASHFORTH_3 */
516 jmc 1.16 CALL CYCLE_AB_TRACER(
517     I bi, bj, gS_loc,
518     U salt(1-OLx,1-OLy,1,bi,bj),
519     O gsNm1(1-OLx,1-OLy,1,bi,bj),
520     I myTime, myIter, myThid )
521 jmc 1.12 #endif /* ALLOW_ADAMSBASHFORTH_3 */
522 jmc 1.16 ELSE
523 jmc 1.12 C- Update tracer fields: S(n) = S**
524     CALL CYCLE_TRACER(
525     I bi, bj,
526 jmc 1.15 O salt(1-OLx,1-OLy,1,bi,bj),
527     I gS_loc, myTime, myIter, myThid )
528 jmc 1.12 ENDIF
529    
530 jmc 1.1 #endif /* ALLOW_GENERIC_ADVDIFF */
531    
532     RETURN
533     END

  ViewVC Help
Powered by ViewVC 1.1.22