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

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

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


Revision 1.19 - (hide annotations) (download)
Sun Oct 9 18:13:09 2016 UTC (7 years, 7 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.18: +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.19 C $Header: /u/gcmpack/MITgcm/model/src/temp_integrate.F,v 1.18 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: TEMP_INTEGRATE
15     C !INTERFACE:
16     SUBROUTINE TEMP_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 TEMP_INTEGRATE
24 jmc 1.13 C | o Calculate tendency for temperature and integrates
25     C | forward in time. The temperature 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_T is called from
30 jmc 1.1 C | here. These procedures can be used to add per problem
31     C | heat 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     C | simply setting the diffusion tensor terms appropriately.
49 jmc 1.1 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.5 #ifdef ALLOW_TIMEAVE
66     # include "TIMEAVE_STATV.h"
67     #endif
68 jmc 1.7 #ifdef ALLOW_AUTODIFF
69 jmc 1.1 # include "tamc.h"
70     # include "tamc_keys.h"
71     #endif
72    
73     C !INPUT/OUTPUT PARAMETERS:
74     C == Routine arguments ==
75 jmc 1.5 C bi, bj, :: tile indices
76     C recip_hFac :: reciprocal of cell open-depth factor (@ next iter)
77     C uFld,vFld :: Local copy of horizontal velocity field
78     C wFld :: Local copy of vertical velocity field
79     C KappaRk :: Vertical diffusion for Tempertature
80     C myTime :: current time
81     C myIter :: current iteration number
82     C myThid :: my Thread Id. number
83 jmc 1.4 INTEGER bi, bj
84 jmc 1.5 _RS recip_hFac(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
85     _RL uFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
86     _RL vFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
87     _RL wFld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
88     _RL KappaRk (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
89 jmc 1.1 _RL myTime
90     INTEGER myIter
91     INTEGER myThid
92     CEOP
93    
94     #ifdef ALLOW_GENERIC_ADVDIFF
95 jmc 1.10 #ifdef ALLOW_DIAGNOSTICS
96     C !FUNCTIONS:
97     LOGICAL DIAGNOSTICS_IS_ON
98     EXTERNAL DIAGNOSTICS_IS_ON
99     #endif /* ALLOW_DIAGNOSTICS */
100    
101 jmc 1.1 C !LOCAL VARIABLES:
102 jmc 1.5 C iMin, iMax :: 1rst index loop range
103     C jMin, jMax :: 2nd index loop range
104     C k :: vertical index
105     C kM1 :: =k-1 for k>1, =1 for k=1
106     C kUp :: index into 2 1/2D array, toggles between 1|2
107     C kDown :: index into 2 1/2D array, toggles between 2|1
108     C xA :: Tracer cell face area normal to X
109     C yA :: Tracer cell face area normal to X
110     C maskUp :: Land/water mask for Wvel points (interface k)
111     C uTrans :: Zonal volume transport through cell face
112     C vTrans :: Meridional volume transport through cell face
113     C rTrans :: Vertical volume transport at interface k
114     C rTransKp :: Vertical volume transport at inteface k+1
115     C fZon :: Flux of temperature (T) in the zonal direction
116     C fMer :: Flux of temperature (T) in the meridional direction
117     C fVer :: Flux of temperature (T) in the vertical direction
118     C at the upper(U) and lower(D) faces of a cell.
119 jmc 1.15 C gT_loc :: Temperature tendency (local to this S/R)
120 jmc 1.10 C gtForc :: Temperature forcing tendency
121 jmc 1.8 C gt_AB :: Adams-Bashforth temperature tendency increment
122 jmc 1.5 C useVariableK :: T when vertical diffusion is not constant
123 jmc 1.4 INTEGER iMin, iMax, jMin, jMax
124 jmc 1.1 INTEGER i, j, k
125     INTEGER kUp, kDown, kM1
126     _RS xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
127     _RS yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
128     _RS maskUp (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
129     _RL uTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
130     _RL vTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
131     _RL rTrans (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
132     _RL rTransKp(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
133 jmc 1.5 _RL fZon (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
134     _RL fMer (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
135     _RL fVer (1-OLx:sNx+OLx,1-OLy:sNy+OLy,2)
136 jmc 1.15 _RL gT_loc (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
137 jmc 1.10 _RL gtForc (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
138 jmc 1.1 _RL gt_AB (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
139 jmc 1.10 #ifdef ALLOW_DIAGNOSTICS
140     LOGICAL diagForcing, diagAB_tend
141     #endif
142 jmc 1.1 LOGICAL calcAdvection
143     INTEGER iterNb
144     #ifdef ALLOW_ADAMSBASHFORTH_3
145 jmc 1.11 INTEGER m2
146 jmc 1.1 #endif
147 jmc 1.5 #ifdef ALLOW_TIMEAVE
148     LOGICAL useVariableK
149     #endif
150    
151 jmc 1.4 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
152    
153     iterNb = myIter
154 jmc 1.5 IF (staggerTimeStep) iterNb = myIter - 1
155 jmc 1.1
156 jmc 1.18 C- Loop ranges for daughter routines:
157     c iMin = 1
158     c iMax = sNx
159     c jMin = 1
160     c jMax = sNy
161     C Regarding model dynamics, only needs to get correct tracer tendency
162     C (gT_loc) in tile interior (1:sNx,1:sNy);
163     C However, for some diagnostics, we may want to get valid tendency
164     C extended over 1 point in tile halo region (0:sNx+1,0:sNy=1).
165     iMin = 0
166     iMax = sNx+1
167     jMin = 0
168     jMax = sNy+1
169    
170 jmc 1.10 #ifdef ALLOW_DIAGNOSTICS
171     diagForcing = .FALSE.
172     diagAB_tend = .FALSE.
173     IF ( useDiagnostics .AND. tempForcing )
174     & diagForcing = DIAGNOSTICS_IS_ON( 'gT_Forc ', myThid )
175     IF ( useDiagnostics .AND. AdamsBashforthGt )
176     & diagAB_tend = DIAGNOSTICS_IS_ON( 'AB_gT ', myThid )
177     #endif
178    
179 jmc 1.1 #ifdef ALLOW_AUTODIFF_TAMC
180     act1 = bi - myBxLo(myThid)
181     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
182     act2 = bj - myByLo(myThid)
183     max2 = myByHi(myThid) - myByLo(myThid) + 1
184     act3 = myThid - 1
185     max3 = nTx*nTy
186     act4 = ikey_dynamics - 1
187     itdkey = (act1 + 1) + act2*max1
188     & + act3*max1*max2
189     & + act4*max1*max2*max3
190     #endif /* ALLOW_AUTODIFF_TAMC */
191    
192 jmc 1.17 C- Apply AB on T :
193     IF ( AdamsBashforth_T ) THEN
194     C compute T^n+1/2 (stored in gtNm) extrapolating T forward in time
195     #ifdef ALLOW_ADAMSBASHFORTH_3
196     c m1 = 1 + MOD(iterNb+1,2)
197     c m2 = 1 + MOD( iterNb ,2)
198     CALL ADAMS_BASHFORTH3(
199     I bi, bj, 0, Nr,
200     I theta(1-OLx,1-OLy,1,bi,bj),
201     U gtNm, gt_AB,
202     I tempStartAB, iterNb, myThid )
203     #else /* ALLOW_ADAMSBASHFORTH_3 */
204     CALL ADAMS_BASHFORTH2(
205     I bi, bj, 0, Nr,
206     I theta(1-OLx,1-OLy,1,bi,bj),
207     U gtNm1(1-OLx,1-OLy,1,bi,bj), gt_AB,
208     I tempStartAB, iterNb, myThid )
209     #endif /* ALLOW_ADAMSBASHFORTH_3 */
210     ENDIF
211    
212 jmc 1.4 C- Tracer tendency needs to be set to zero (moved here from gad_calc_rhs):
213     DO k=1,Nr
214     DO j=1-OLy,sNy+OLy
215     DO i=1-OLx,sNx+OLx
216 jmc 1.15 gT_loc(i,j,k) = 0. _d 0
217 jmc 1.4 ENDDO
218     ENDDO
219     ENDDO
220 jmc 1.1 DO j=1-OLy,sNy+OLy
221     DO i=1-OLx,sNx+OLx
222 jmc 1.5 fVer(i,j,1) = 0. _d 0
223     fVer(i,j,2) = 0. _d 0
224 jmc 1.1 ENDDO
225     ENDDO
226 jmc 1.4 #ifdef ALLOW_AUTODIFF
227     DO k=1,Nr
228     DO j=1-OLy,sNy+OLy
229     DO i=1-OLx,sNx+OLx
230     kappaRk(i,j,k) = 0. _d 0
231     ENDDO
232     ENDDO
233     ENDDO
234 jmc 1.17 CADJ STORE wFld(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
235     CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
236     # ifdef ALLOW_ADAMSBASHFORTH_3
237     CADJ STORE gtNm(:,:,:,bi,bj,1) = comlev1_bibj, key=itdkey, byte=isbyte
238     CADJ STORE gtNm(:,:,:,bi,bj,2) = comlev1_bibj, key=itdkey, byte=isbyte
239     # else
240     CADJ STORE gtNm1(:,:,:,bi,bj) = comlev1_bibj, key=itdkey, byte=isbyte
241     # endif
242 jmc 1.4 #endif /* ALLOW_AUTODIFF */
243    
244     #ifdef INCLUDE_CALC_DIFFUSIVITY_CALL
245 jmc 1.5 CALL CALC_3D_DIFFUSIVITY(
246 jmc 1.4 I bi, bj, iMin, iMax, jMin, jMax,
247     I GAD_TEMPERATURE, useGMredi, useKPP,
248     O kappaRk,
249     I myThid )
250     #endif /* INCLUDE_CALC_DIFFUSIVITY_CALL */
251    
252     #ifndef DISABLE_MULTIDIM_ADVECTION
253     C-- Some advection schemes are better calculated using a multi-dimensional
254     C method in the absence of any other terms and, if used, is done here.
255     C
256     C The CPP flag DISABLE_MULTIDIM_ADVECTION is currently unset in GAD_OPTIONS.h
257     C The default is to use multi-dimensinal advection for non-linear advection
258     C schemes. However, for the sake of efficiency of the adjoint it is necessary
259     C to be able to exclude this scheme to avoid excessive storage and
260     C recomputation. It *is* differentiable, if you need it.
261     C Edit GAD_OPTIONS.h and #define DISABLE_MULTIDIM_ADVECTION to
262     C disable this section of code.
263     #ifdef GAD_ALLOW_TS_SOM_ADV
264     # ifdef ALLOW_AUTODIFF_TAMC
265     CADJ STORE som_T = comlev1_bibj, key=itdkey, byte=isbyte
266     # endif
267     IF ( tempSOM_Advection ) THEN
268     # ifdef ALLOW_DEBUG
269     IF (debugMode) CALL DEBUG_CALL('GAD_SOM_ADVECT',myThid)
270     # endif
271     CALL GAD_SOM_ADVECT(
272 jmc 1.8 I tempImplVertAdv,
273     I tempAdvScheme, tempVertAdvScheme, GAD_TEMPERATURE,
274     I dTtracerLev, uFld, vFld, wFld, theta,
275 jmc 1.4 U som_T,
276 jmc 1.15 O gT_loc,
277 jmc 1.4 I bi, bj, myTime, myIter, myThid )
278     ELSEIF (tempMultiDimAdvec) THEN
279     #else /* GAD_ALLOW_TS_SOM_ADV */
280     IF (tempMultiDimAdvec) THEN
281     #endif /* GAD_ALLOW_TS_SOM_ADV */
282     # ifdef ALLOW_DEBUG
283     IF (debugMode) CALL DEBUG_CALL('GAD_ADVECTION',myThid)
284     # endif
285     CALL GAD_ADVECTION(
286 jmc 1.8 I tempImplVertAdv,
287     I tempAdvScheme, tempVertAdvScheme, GAD_TEMPERATURE,
288     I dTtracerLev, uFld, vFld, wFld, theta,
289 jmc 1.15 O gT_loc,
290 jmc 1.4 I bi, bj, myTime, myIter, myThid )
291     ENDIF
292     #endif /* DISABLE_MULTIDIM_ADVECTION */
293 jmc 1.1
294 jmc 1.4 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
295    
296     C- Start vertical index (k) loop (Nr:1)
297     calcAdvection = tempAdvection .AND. .NOT.tempMultiDimAdvec
298 jmc 1.1 DO k=Nr,1,-1
299     #ifdef ALLOW_AUTODIFF_TAMC
300     kkey = (itdkey-1)*Nr + k
301     #endif /* ALLOW_AUTODIFF_TAMC */
302     kM1 = MAX(1,k-1)
303     kUp = 1+MOD(k+1,2)
304     kDown= 1+MOD(k,2)
305    
306     #ifdef ALLOW_AUTODIFF_TAMC
307 jmc 1.5 CADJ STORE fVer(:,:,:) = comlev1_bibj_k, key=kkey,
308 jmc 1.1 CADJ & byte=isbyte, kind = isbyte
309 jmc 1.15 CADJ STORE gT_loc(:,:,k) = comlev1_bibj_k, key=kkey,
310 jmc 1.1 CADJ & byte=isbyte, kind = isbyte
311     # ifdef ALLOW_ADAMSBASHFORTH_3
312     CADJ STORE gtNm(:,:,k,bi,bj,1) = comlev1_bibj_k, key=kkey,
313     CADJ & byte=isbyte, kind = isbyte
314     CADJ STORE gtNm(:,:,k,bi,bj,2) = comlev1_bibj_k, key=kkey,
315     CADJ & byte=isbyte, kind = isbyte
316     # else
317     CADJ STORE gtNm1(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey,
318     CADJ & byte=isbyte, kind = isbyte
319     # endif
320     #endif /* ALLOW_AUTODIFF_TAMC */
321     CALL CALC_ADV_FLOW(
322     I uFld, vFld, wFld,
323     U rTrans,
324     O uTrans, vTrans, rTransKp,
325     O maskUp, xA, yA,
326     I k, bi, bj, myThid )
327    
328 jmc 1.10 C-- Collect forcing term in local array gtForc:
329     DO j=1-OLy,sNy+OLy
330     DO i=1-OLx,sNx+OLx
331     gtForc(i,j) = 0. _d 0
332     ENDDO
333     ENDDO
334     IF ( tempForcing ) THEN
335     CALL APPLY_FORCING_T(
336     U gtForc,
337     I iMin,iMax,jMin,jMax, k, bi,bj,
338     I myTime, myIter, myThid )
339     #ifdef ALLOW_DIAGNOSTICS
340     IF ( diagForcing ) THEN
341     CALL DIAGNOSTICS_FILL(gtForc,'gT_Forc ',k,1,2,bi,bj,myThid)
342     ENDIF
343     #endif /* ALLOW_DIAGNOSTICS */
344     ENDIF
345    
346 jmc 1.1 #ifdef ALLOW_ADAMSBASHFORTH_3
347 jmc 1.11 c m1 = 1 + MOD(iterNb+1,2)
348 jmc 1.1 m2 = 1 + MOD( iterNb ,2)
349     CALL GAD_CALC_RHS(
350     I bi, bj, iMin,iMax,jMin,jMax, k, kM1, kUp, kDown,
351     I xA, yA, maskUp, uFld(1-OLx,1-OLy,k),
352     I vFld(1-OLx,1-OLy,k), wFld(1-OLx,1-OLy,k),
353     I uTrans, vTrans, rTrans, rTransKp,
354     I diffKhT, diffK4T, KappaRk(1-OLx,1-OLy,k), diffKr4T,
355 jmc 1.16 I theta(1-OLx,1-OLy,1,bi,bj),
356     I gtNm(1-OLx,1-OLy,1,bi,bj,m2), dTtracerLev,
357 jmc 1.1 I GAD_TEMPERATURE, tempAdvScheme, tempVertAdvScheme,
358     I calcAdvection, tempImplVertAdv, AdamsBashforth_T,
359     I tempVertDiff4, useGMRedi, useKPP,
360 jmc 1.5 O fZon, fMer,
361 jmc 1.15 U fVer, gT_loc,
362 jmc 1.1 I myTime, myIter, myThid )
363     #else /* ALLOW_ADAMSBASHFORTH_3 */
364     CALL GAD_CALC_RHS(
365     I bi, bj, iMin,iMax,jMin,jMax, k, kM1, kUp, kDown,
366     I xA, yA, maskUp, uFld(1-OLx,1-OLy,k),
367     I vFld(1-OLx,1-OLy,k), wFld(1-OLx,1-OLy,k),
368     I uTrans, vTrans, rTrans, rTransKp,
369     I diffKhT, diffK4T, KappaRk(1-OLx,1-OLy,k), diffKr4T,
370 jmc 1.16 I theta(1-OLx,1-OLy,1,bi,bj),
371     I gtNm1(1-OLx,1-OLy,1,bi,bj), dTtracerLev,
372 jmc 1.1 I GAD_TEMPERATURE, tempAdvScheme, tempVertAdvScheme,
373     I calcAdvection, tempImplVertAdv, AdamsBashforth_T,
374     I tempVertDiff4, useGMRedi, useKPP,
375 jmc 1.5 O fZon, fMer,
376 jmc 1.15 U fVer, gT_loc,
377 jmc 1.1 I myTime, myIter, myThid )
378     #endif
379    
380     C-- External thermal forcing term(s) inside Adams-Bashforth:
381 jmc 1.10 IF ( tempForcing .AND. tracForcingOutAB.NE.1 ) THEN
382     DO j=1-OLy,sNy+OLy
383     DO i=1-OLx,sNx+OLx
384 jmc 1.15 gT_loc(i,j,k) = gT_loc(i,j,k) + gtForc(i,j)
385 jmc 1.10 ENDDO
386     ENDDO
387     ENDIF
388 jmc 1.1
389     IF ( AdamsBashforthGt ) THEN
390     #ifdef ALLOW_ADAMSBASHFORTH_3
391     CALL ADAMS_BASHFORTH3(
392     I bi, bj, k, Nr,
393 jmc 1.15 U gT_loc, gtNm,
394     O gt_AB,
395 jmc 1.1 I tempStartAB, iterNb, myThid )
396     #else
397     CALL ADAMS_BASHFORTH2(
398     I bi, bj, k, Nr,
399 jmc 1.15 U gT_loc, gtNm1(1-OLx,1-OLy,1,bi,bj),
400     O gt_AB,
401 jmc 1.1 I tempStartAB, iterNb, myThid )
402     #endif
403     #ifdef ALLOW_DIAGNOSTICS
404 jmc 1.10 IF ( diagAB_tend ) THEN
405 jmc 1.1 CALL DIAGNOSTICS_FILL(gt_AB,'AB_gT ',k,1,2,bi,bj,myThid)
406     ENDIF
407     #endif /* ALLOW_DIAGNOSTICS */
408     ENDIF
409    
410     C-- External thermal forcing term(s) outside Adams-Bashforth:
411 jmc 1.10 IF ( tempForcing .AND. tracForcingOutAB.EQ.1 ) THEN
412     DO j=1-OLy,sNy+OLy
413     DO i=1-OLx,sNx+OLx
414 jmc 1.15 gT_loc(i,j,k) = gT_loc(i,j,k) + gtForc(i,j)
415 jmc 1.10 ENDDO
416     ENDDO
417     ENDIF
418 jmc 1.1
419     #ifdef NONLIN_FRSURF
420     IF (nonlinFreeSurf.GT.0) THEN
421     CALL FREESURF_RESCALE_G(
422     I bi, bj, k,
423 jmc 1.15 U gT_loc,
424 jmc 1.1 I myThid )
425     IF ( AdamsBashforthGt ) THEN
426     #ifdef ALLOW_ADAMSBASHFORTH_3
427     # ifdef ALLOW_AUTODIFF_TAMC
428     CADJ STORE gtNm(:,:,k,bi,bj,1) = comlev1_bibj_k, key=kkey,
429     CADJ & byte=isbyte, kind = isbyte
430     CADJ STORE gtNm(:,:,k,bi,bj,2) = comlev1_bibj_k, key=kkey,
431     CADJ & byte=isbyte, kind = isbyte
432     # endif
433     CALL FREESURF_RESCALE_G(
434     I bi, bj, k,
435 jmc 1.15 U gtNm(1-OLx,1-OLy,1,bi,bj,1),
436 jmc 1.1 I myThid )
437     CALL FREESURF_RESCALE_G(
438     I bi, bj, k,
439 jmc 1.15 U gtNm(1-OLx,1-OLy,1,bi,bj,2),
440 jmc 1.1 I myThid )
441     #else
442     CALL FREESURF_RESCALE_G(
443     I bi, bj, k,
444 jmc 1.15 U gtNm1(1-OLx,1-OLy,1,bi,bj),
445 jmc 1.1 I myThid )
446     #endif
447     ENDIF
448     ENDIF
449     #endif /* NONLIN_FRSURF */
450    
451     C- end of vertical index (k) loop (Nr:1)
452     ENDDO
453    
454 jmc 1.5 #ifdef ALLOW_DOWN_SLOPE
455     IF ( useDOWN_SLOPE ) THEN
456     IF ( usingPCoords ) THEN
457     CALL DWNSLP_APPLY(
458     I GAD_TEMPERATURE, bi, bj, kSurfC,
459 jmc 1.14 I theta(1-OLx,1-OLy,1,bi,bj),
460 jmc 1.15 U gT_loc,
461 jmc 1.14 I recip_hFac, recip_rA, recip_drF,
462     I dTtracerLev, myTime, myIter, myThid )
463 jmc 1.5 ELSE
464     CALL DWNSLP_APPLY(
465     I GAD_TEMPERATURE, bi, bj, kLowC,
466 jmc 1.14 I theta(1-OLx,1-OLy,1,bi,bj),
467 jmc 1.15 U gT_loc,
468 jmc 1.14 I recip_hFac, recip_rA, recip_drF,
469     I dTtracerLev, myTime, myIter, myThid )
470 jmc 1.5 ENDIF
471     ENDIF
472     #endif /* ALLOW_DOWN_SLOPE */
473    
474 jmc 1.15 C- Integrate forward in time, storing in gT_loc: gT <= T + dt*gT
475 jmc 1.14 CALL TIMESTEP_TRACER(
476     I bi, bj, dTtracerLev,
477     I theta(1-OLx,1-OLy,1,bi,bj),
478 jmc 1.15 U gT_loc,
479 jmc 1.14 I myTime, myIter, myThid )
480    
481 jmc 1.5 C-- Implicit vertical advection & diffusion
482    
483     #ifdef INCLUDE_IMPLVERTADV_CODE
484 jmc 1.19 IF ( tempImplVertAdv .OR. implicitDiffusion ) THEN
485     C to recover older (prior to 2016-10-05) results:
486     c IF ( tempImplVertAdv ) THEN
487 jmc 1.5 #ifdef ALLOW_AUTODIFF_TAMC
488 jmc 1.7 CADJ STORE kappaRk(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
489 jmc 1.15 CADJ STORE gT_loc(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
490 jmc 1.7 CADJ STORE wFld(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
491     CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte
492     CADJ STORE recip_hFac(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
493 jmc 1.5 #endif /* ALLOW_AUTODIFF_TAMC */
494     CALL GAD_IMPLICIT_R(
495     I tempImplVertAdv, tempVertAdvScheme, GAD_TEMPERATURE,
496 jmc 1.15 I dTtracerLev, kappaRk, recip_hFac, wFld,
497     I theta(1-OLx,1-OLy,1,bi,bj),
498     U gT_loc,
499 jmc 1.5 I bi, bj, myTime, myIter, myThid )
500     ELSEIF ( implicitDiffusion ) THEN
501     #else /* INCLUDE_IMPLVERTADV_CODE */
502     IF ( implicitDiffusion ) THEN
503     #endif /* INCLUDE_IMPLVERTADV_CODE */
504     #ifdef ALLOW_AUTODIFF_TAMC
505 jmc 1.7 CADJ STORE kappaRk(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
506 jmc 1.15 CADJ STORE gT_loc(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte
507 jmc 1.5 #endif /* ALLOW_AUTODIFF_TAMC */
508     CALL IMPLDIFF(
509     I bi, bj, iMin, iMax, jMin, jMax,
510     I GAD_TEMPERATURE, kappaRk, recip_hFac,
511 jmc 1.15 U gT_loc,
512 jmc 1.5 I myThid )
513     ENDIF
514    
515     #ifdef ALLOW_TIMEAVE
516 jmc 1.12 useVariableK = useKPP .OR. usePP81 .OR. useKL10 .OR. useMY82
517     & .OR. useGGL90 .OR. useGMredi .OR. ivdc_kappa.NE.0.
518 jmc 1.5 IF ( taveFreq.GT.0. .AND. useVariableK
519     & .AND.implicitDiffusion ) THEN
520 jmc 1.15 CALL TIMEAVE_CUMUL_DIF_1T( TdiffRtave,
521     I gT_loc, kappaRk,
522     I Nr, 3, deltaTClock, bi, bj, myThid )
523 jmc 1.5 ENDIF
524     #endif /* ALLOW_TIMEAVE */
525    
526 jmc 1.13 IF ( AdamsBashforth_T ) THEN
527     C- Save current tracer field (for AB on tracer) and then update tracer
528 jmc 1.17 #ifdef ALLOW_ADAMSBASHFORTH_3
529 jmc 1.13 CALL CYCLE_AB_TRACER(
530 jmc 1.16 I bi, bj, gT_loc,
531     U theta(1-OLx,1-OLy,1,bi,bj),
532     O gtNm(1-OLx,1-OLy,1,bi,bj,m2),
533 jmc 1.13 I myTime, myIter, myThid )
534     #else /* ALLOW_ADAMSBASHFORTH_3 */
535 jmc 1.17 CALL CYCLE_AB_TRACER(
536     I bi, bj, gT_loc,
537     U theta(1-OLx,1-OLy,1,bi,bj),
538     O gtNm1(1-OLx,1-OLy,1,bi,bj),
539     I myTime, myIter, myThid )
540 jmc 1.13 #endif /* ALLOW_ADAMSBASHFORTH_3 */
541 jmc 1.17 ELSE
542 jmc 1.13 C- Update tracer fields: T(n) = T**
543     CALL CYCLE_TRACER(
544     I bi, bj,
545 jmc 1.16 O theta(1-OLx,1-OLy,1,bi,bj),
546     I gT_loc, myTime, myIter, myThid )
547 jmc 1.13 ENDIF
548    
549 jmc 1.1 #endif /* ALLOW_GENERIC_ADVDIFF */
550    
551     RETURN
552     END

  ViewVC Help
Powered by ViewVC 1.1.22