1 |
C $Header: /u/gcmpack/MITgcm/model/src/temp_integrate.F,v 1.14 2014/08/08 19:55:52 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "PACKAGES_CONFIG.h" |
5 |
#include "CPP_OPTIONS.h" |
6 |
#ifdef ALLOW_AUTODIFF |
7 |
# include "AUTODIFF_OPTIONS.h" |
8 |
#endif |
9 |
#ifdef ALLOW_GENERIC_ADVDIFF |
10 |
# include "GAD_OPTIONS.h" |
11 |
#endif |
12 |
|
13 |
CBOP |
14 |
C !ROUTINE: TEMP_INTEGRATE |
15 |
C !INTERFACE: |
16 |
SUBROUTINE TEMP_INTEGRATE( |
17 |
I bi, bj, recip_hFac, |
18 |
I uFld, vFld, wFld, |
19 |
U KappaRk, |
20 |
I myTime, myIter, myThid ) |
21 |
C !DESCRIPTION: \bv |
22 |
C *==========================================================* |
23 |
C | SUBROUTINE TEMP_INTEGRATE |
24 |
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 |
C *==========================================================* |
29 |
C | A procedure called APPLY_FORCING_T is called from |
30 |
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 |
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 |
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 |
#include "GRID.h" |
59 |
#include "DYNVARS.h" |
60 |
#include "RESTART.h" |
61 |
#ifdef ALLOW_GENERIC_ADVDIFF |
62 |
# include "GAD.h" |
63 |
# include "GAD_SOM_VARS.h" |
64 |
#endif |
65 |
#ifdef ALLOW_TIMEAVE |
66 |
# include "TIMEAVE_STATV.h" |
67 |
#endif |
68 |
#ifdef ALLOW_AUTODIFF |
69 |
# include "tamc.h" |
70 |
# include "tamc_keys.h" |
71 |
#endif |
72 |
|
73 |
C !INPUT/OUTPUT PARAMETERS: |
74 |
C == Routine arguments == |
75 |
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 |
INTEGER bi, bj |
84 |
_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 |
_RL myTime |
90 |
INTEGER myIter |
91 |
INTEGER myThid |
92 |
CEOP |
93 |
|
94 |
#ifdef ALLOW_GENERIC_ADVDIFF |
95 |
#ifdef ALLOW_DIAGNOSTICS |
96 |
C !FUNCTIONS: |
97 |
LOGICAL DIAGNOSTICS_IS_ON |
98 |
EXTERNAL DIAGNOSTICS_IS_ON |
99 |
#endif /* ALLOW_DIAGNOSTICS */ |
100 |
|
101 |
C !LOCAL VARIABLES: |
102 |
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 |
C gT_loc :: Temperature tendency (local to this S/R) |
120 |
C gtForc :: Temperature forcing tendency |
121 |
C gt_AB :: Adams-Bashforth temperature tendency increment |
122 |
C useVariableK :: T when vertical diffusion is not constant |
123 |
INTEGER iMin, iMax, jMin, jMax |
124 |
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 |
_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 |
_RL gT_loc (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
137 |
_RL gtForc (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
138 |
_RL gt_AB (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
139 |
#ifdef ALLOW_DIAGNOSTICS |
140 |
LOGICAL diagForcing, diagAB_tend |
141 |
#endif |
142 |
LOGICAL calcAdvection |
143 |
INTEGER iterNb |
144 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
145 |
INTEGER m2 |
146 |
#endif |
147 |
#ifdef ALLOW_TIMEAVE |
148 |
LOGICAL useVariableK |
149 |
#endif |
150 |
|
151 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
152 |
|
153 |
C- Loop ranges for daughter routines |
154 |
iMin = 1-OLx+2 |
155 |
iMax = sNx+OLx-1 |
156 |
jMin = 1-OLy+2 |
157 |
jMax = sNy+OLy-1 |
158 |
|
159 |
iterNb = myIter |
160 |
IF (staggerTimeStep) iterNb = myIter - 1 |
161 |
|
162 |
#ifdef ALLOW_DIAGNOSTICS |
163 |
diagForcing = .FALSE. |
164 |
diagAB_tend = .FALSE. |
165 |
IF ( useDiagnostics .AND. tempForcing ) |
166 |
& diagForcing = DIAGNOSTICS_IS_ON( 'gT_Forc ', myThid ) |
167 |
IF ( useDiagnostics .AND. AdamsBashforthGt ) |
168 |
& diagAB_tend = DIAGNOSTICS_IS_ON( 'AB_gT ', myThid ) |
169 |
#endif |
170 |
|
171 |
#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 |
C- Tracer tendency needs to be set to zero (moved here from gad_calc_rhs): |
185 |
DO k=1,Nr |
186 |
DO j=1-OLy,sNy+OLy |
187 |
DO i=1-OLx,sNx+OLx |
188 |
gT_loc(i,j,k) = 0. _d 0 |
189 |
ENDDO |
190 |
ENDDO |
191 |
ENDDO |
192 |
DO j=1-OLy,sNy+OLy |
193 |
DO i=1-OLx,sNx+OLx |
194 |
fVer(i,j,1) = 0. _d 0 |
195 |
fVer(i,j,2) = 0. _d 0 |
196 |
ENDDO |
197 |
ENDDO |
198 |
#ifdef ALLOW_AUTODIFF |
199 |
DO k=1,Nr |
200 |
DO j=1-OLy,sNy+OLy |
201 |
DO i=1-OLx,sNx+OLx |
202 |
kappaRk(i,j,k) = 0. _d 0 |
203 |
ENDDO |
204 |
ENDDO |
205 |
ENDDO |
206 |
CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte |
207 |
CADJ STORE wFld(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte |
208 |
#endif /* ALLOW_AUTODIFF */ |
209 |
|
210 |
#ifdef INCLUDE_CALC_DIFFUSIVITY_CALL |
211 |
CALL CALC_3D_DIFFUSIVITY( |
212 |
I bi, bj, iMin, iMax, jMin, jMax, |
213 |
I GAD_TEMPERATURE, useGMredi, useKPP, |
214 |
O kappaRk, |
215 |
I myThid ) |
216 |
#endif /* INCLUDE_CALC_DIFFUSIVITY_CALL */ |
217 |
|
218 |
#ifndef DISABLE_MULTIDIM_ADVECTION |
219 |
C-- Some advection schemes are better calculated using a multi-dimensional |
220 |
C method in the absence of any other terms and, if used, is done here. |
221 |
C |
222 |
C The CPP flag DISABLE_MULTIDIM_ADVECTION is currently unset in GAD_OPTIONS.h |
223 |
C The default is to use multi-dimensinal advection for non-linear advection |
224 |
C schemes. However, for the sake of efficiency of the adjoint it is necessary |
225 |
C to be able to exclude this scheme to avoid excessive storage and |
226 |
C recomputation. It *is* differentiable, if you need it. |
227 |
C Edit GAD_OPTIONS.h and #define DISABLE_MULTIDIM_ADVECTION to |
228 |
C disable this section of code. |
229 |
#ifdef GAD_ALLOW_TS_SOM_ADV |
230 |
# ifdef ALLOW_AUTODIFF_TAMC |
231 |
CADJ STORE som_T = comlev1_bibj, key=itdkey, byte=isbyte |
232 |
# endif |
233 |
IF ( tempSOM_Advection ) THEN |
234 |
# ifdef ALLOW_DEBUG |
235 |
IF (debugMode) CALL DEBUG_CALL('GAD_SOM_ADVECT',myThid) |
236 |
# endif |
237 |
CALL GAD_SOM_ADVECT( |
238 |
I tempImplVertAdv, |
239 |
I tempAdvScheme, tempVertAdvScheme, GAD_TEMPERATURE, |
240 |
I dTtracerLev, uFld, vFld, wFld, theta, |
241 |
U som_T, |
242 |
O gT_loc, |
243 |
I bi, bj, myTime, myIter, myThid ) |
244 |
ELSEIF (tempMultiDimAdvec) THEN |
245 |
#else /* GAD_ALLOW_TS_SOM_ADV */ |
246 |
IF (tempMultiDimAdvec) THEN |
247 |
#endif /* GAD_ALLOW_TS_SOM_ADV */ |
248 |
# ifdef ALLOW_DEBUG |
249 |
IF (debugMode) CALL DEBUG_CALL('GAD_ADVECTION',myThid) |
250 |
# endif |
251 |
CALL GAD_ADVECTION( |
252 |
I tempImplVertAdv, |
253 |
I tempAdvScheme, tempVertAdvScheme, GAD_TEMPERATURE, |
254 |
I dTtracerLev, uFld, vFld, wFld, theta, |
255 |
O gT_loc, |
256 |
I bi, bj, myTime, myIter, myThid ) |
257 |
ENDIF |
258 |
#endif /* DISABLE_MULTIDIM_ADVECTION */ |
259 |
|
260 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
261 |
|
262 |
C- Start vertical index (k) loop (Nr:1) |
263 |
calcAdvection = tempAdvection .AND. .NOT.tempMultiDimAdvec |
264 |
DO k=Nr,1,-1 |
265 |
#ifdef ALLOW_AUTODIFF_TAMC |
266 |
kkey = (itdkey-1)*Nr + k |
267 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
268 |
kM1 = MAX(1,k-1) |
269 |
kUp = 1+MOD(k+1,2) |
270 |
kDown= 1+MOD(k,2) |
271 |
|
272 |
#ifdef ALLOW_AUTODIFF_TAMC |
273 |
CADJ STORE fVer(:,:,:) = comlev1_bibj_k, key=kkey, |
274 |
CADJ & byte=isbyte, kind = isbyte |
275 |
CADJ STORE gT_loc(:,:,k) = comlev1_bibj_k, key=kkey, |
276 |
CADJ & byte=isbyte, kind = isbyte |
277 |
# ifdef ALLOW_ADAMSBASHFORTH_3 |
278 |
CADJ STORE gtNm(:,:,k,bi,bj,1) = comlev1_bibj_k, key=kkey, |
279 |
CADJ & byte=isbyte, kind = isbyte |
280 |
CADJ STORE gtNm(:,:,k,bi,bj,2) = comlev1_bibj_k, key=kkey, |
281 |
CADJ & byte=isbyte, kind = isbyte |
282 |
# else |
283 |
CADJ STORE gtNm1(:,:,k,bi,bj) = comlev1_bibj_k, key=kkey, |
284 |
CADJ & byte=isbyte, kind = isbyte |
285 |
# endif |
286 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
287 |
CALL CALC_ADV_FLOW( |
288 |
I uFld, vFld, wFld, |
289 |
U rTrans, |
290 |
O uTrans, vTrans, rTransKp, |
291 |
O maskUp, xA, yA, |
292 |
I k, bi, bj, myThid ) |
293 |
|
294 |
C-- Collect forcing term in local array gtForc: |
295 |
DO j=1-OLy,sNy+OLy |
296 |
DO i=1-OLx,sNx+OLx |
297 |
gtForc(i,j) = 0. _d 0 |
298 |
ENDDO |
299 |
ENDDO |
300 |
IF ( tempForcing ) THEN |
301 |
CALL APPLY_FORCING_T( |
302 |
U gtForc, |
303 |
I iMin,iMax,jMin,jMax, k, bi,bj, |
304 |
I myTime, myIter, myThid ) |
305 |
#ifdef ALLOW_DIAGNOSTICS |
306 |
IF ( diagForcing ) THEN |
307 |
CALL DIAGNOSTICS_FILL(gtForc,'gT_Forc ',k,1,2,bi,bj,myThid) |
308 |
ENDIF |
309 |
#endif /* ALLOW_DIAGNOSTICS */ |
310 |
ENDIF |
311 |
|
312 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
313 |
c m1 = 1 + MOD(iterNb+1,2) |
314 |
m2 = 1 + MOD( iterNb ,2) |
315 |
CALL GAD_CALC_RHS( |
316 |
I bi, bj, iMin,iMax,jMin,jMax, k, kM1, kUp, kDown, |
317 |
I xA, yA, maskUp, uFld(1-OLx,1-OLy,k), |
318 |
I vFld(1-OLx,1-OLy,k), wFld(1-OLx,1-OLy,k), |
319 |
I uTrans, vTrans, rTrans, rTransKp, |
320 |
I diffKhT, diffK4T, KappaRk(1-OLx,1-OLy,k), diffKr4T, |
321 |
I theta, gtNm(1-OLx,1-OLy,1,1,1,m2), dTtracerLev, |
322 |
I GAD_TEMPERATURE, tempAdvScheme, tempVertAdvScheme, |
323 |
I calcAdvection, tempImplVertAdv, AdamsBashforth_T, |
324 |
I tempVertDiff4, useGMRedi, useKPP, |
325 |
O fZon, fMer, |
326 |
U fVer, gT_loc, |
327 |
I myTime, myIter, myThid ) |
328 |
#else /* ALLOW_ADAMSBASHFORTH_3 */ |
329 |
CALL GAD_CALC_RHS( |
330 |
I bi, bj, iMin,iMax,jMin,jMax, k, kM1, kUp, kDown, |
331 |
I xA, yA, maskUp, uFld(1-OLx,1-OLy,k), |
332 |
I vFld(1-OLx,1-OLy,k), wFld(1-OLx,1-OLy,k), |
333 |
I uTrans, vTrans, rTrans, rTransKp, |
334 |
I diffKhT, diffK4T, KappaRk(1-OLx,1-OLy,k), diffKr4T, |
335 |
I theta, gtNm1, dTtracerLev, |
336 |
I GAD_TEMPERATURE, tempAdvScheme, tempVertAdvScheme, |
337 |
I calcAdvection, tempImplVertAdv, AdamsBashforth_T, |
338 |
I tempVertDiff4, useGMRedi, useKPP, |
339 |
O fZon, fMer, |
340 |
U fVer, gT_loc, |
341 |
I myTime, myIter, myThid ) |
342 |
#endif |
343 |
|
344 |
C-- External thermal forcing term(s) inside Adams-Bashforth: |
345 |
IF ( tempForcing .AND. tracForcingOutAB.NE.1 ) THEN |
346 |
DO j=1-OLy,sNy+OLy |
347 |
DO i=1-OLx,sNx+OLx |
348 |
gT_loc(i,j,k) = gT_loc(i,j,k) + gtForc(i,j) |
349 |
ENDDO |
350 |
ENDDO |
351 |
ENDIF |
352 |
|
353 |
IF ( AdamsBashforthGt ) THEN |
354 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
355 |
CALL ADAMS_BASHFORTH3( |
356 |
I bi, bj, k, Nr, |
357 |
U gT_loc, gtNm, |
358 |
O gt_AB, |
359 |
I tempStartAB, iterNb, myThid ) |
360 |
#else |
361 |
CALL ADAMS_BASHFORTH2( |
362 |
I bi, bj, k, Nr, |
363 |
U gT_loc, gtNm1(1-OLx,1-OLy,1,bi,bj), |
364 |
O gt_AB, |
365 |
I tempStartAB, iterNb, myThid ) |
366 |
#endif |
367 |
#ifdef ALLOW_DIAGNOSTICS |
368 |
IF ( diagAB_tend ) THEN |
369 |
CALL DIAGNOSTICS_FILL(gt_AB,'AB_gT ',k,1,2,bi,bj,myThid) |
370 |
ENDIF |
371 |
#endif /* ALLOW_DIAGNOSTICS */ |
372 |
ENDIF |
373 |
|
374 |
C-- External thermal forcing term(s) outside Adams-Bashforth: |
375 |
IF ( tempForcing .AND. tracForcingOutAB.EQ.1 ) THEN |
376 |
DO j=1-OLy,sNy+OLy |
377 |
DO i=1-OLx,sNx+OLx |
378 |
gT_loc(i,j,k) = gT_loc(i,j,k) + gtForc(i,j) |
379 |
ENDDO |
380 |
ENDDO |
381 |
ENDIF |
382 |
|
383 |
#ifdef NONLIN_FRSURF |
384 |
IF (nonlinFreeSurf.GT.0) THEN |
385 |
CALL FREESURF_RESCALE_G( |
386 |
I bi, bj, k, |
387 |
U gT_loc, |
388 |
I myThid ) |
389 |
IF ( AdamsBashforthGt ) THEN |
390 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
391 |
# ifdef ALLOW_AUTODIFF_TAMC |
392 |
CADJ STORE gtNm(:,:,k,bi,bj,1) = comlev1_bibj_k, key=kkey, |
393 |
CADJ & byte=isbyte, kind = isbyte |
394 |
CADJ STORE gtNm(:,:,k,bi,bj,2) = comlev1_bibj_k, key=kkey, |
395 |
CADJ & byte=isbyte, kind = isbyte |
396 |
# endif |
397 |
CALL FREESURF_RESCALE_G( |
398 |
I bi, bj, k, |
399 |
U gtNm(1-OLx,1-OLy,1,bi,bj,1), |
400 |
I myThid ) |
401 |
CALL FREESURF_RESCALE_G( |
402 |
I bi, bj, k, |
403 |
U gtNm(1-OLx,1-OLy,1,bi,bj,2), |
404 |
I myThid ) |
405 |
#else |
406 |
CALL FREESURF_RESCALE_G( |
407 |
I bi, bj, k, |
408 |
U gtNm1(1-OLx,1-OLy,1,bi,bj), |
409 |
I myThid ) |
410 |
#endif |
411 |
ENDIF |
412 |
ENDIF |
413 |
#endif /* NONLIN_FRSURF */ |
414 |
|
415 |
C- end of vertical index (k) loop (Nr:1) |
416 |
ENDDO |
417 |
|
418 |
#ifdef ALLOW_DOWN_SLOPE |
419 |
IF ( useDOWN_SLOPE ) THEN |
420 |
IF ( usingPCoords ) THEN |
421 |
CALL DWNSLP_APPLY( |
422 |
I GAD_TEMPERATURE, bi, bj, kSurfC, |
423 |
I theta(1-OLx,1-OLy,1,bi,bj), |
424 |
U gT_loc, |
425 |
I recip_hFac, recip_rA, recip_drF, |
426 |
I dTtracerLev, myTime, myIter, myThid ) |
427 |
ELSE |
428 |
CALL DWNSLP_APPLY( |
429 |
I GAD_TEMPERATURE, bi, bj, kLowC, |
430 |
I theta(1-OLx,1-OLy,1,bi,bj), |
431 |
U gT_loc, |
432 |
I recip_hFac, recip_rA, recip_drF, |
433 |
I dTtracerLev, myTime, myIter, myThid ) |
434 |
ENDIF |
435 |
ENDIF |
436 |
#endif /* ALLOW_DOWN_SLOPE */ |
437 |
|
438 |
C- Integrate forward in time, storing in gT_loc: gT <= T + dt*gT |
439 |
CALL TIMESTEP_TRACER( |
440 |
I bi, bj, dTtracerLev, |
441 |
I theta(1-OLx,1-OLy,1,bi,bj), |
442 |
U gT_loc, |
443 |
I myTime, myIter, myThid ) |
444 |
|
445 |
iMin = 0 |
446 |
iMax = sNx+1 |
447 |
jMin = 0 |
448 |
jMax = sNy+1 |
449 |
|
450 |
C-- Implicit vertical advection & diffusion |
451 |
|
452 |
#ifdef INCLUDE_IMPLVERTADV_CODE |
453 |
IF ( tempImplVertAdv ) THEN |
454 |
#ifdef ALLOW_AUTODIFF_TAMC |
455 |
CADJ STORE kappaRk(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte |
456 |
CADJ STORE gT_loc(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte |
457 |
CADJ STORE wFld(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte |
458 |
CADJ STORE theta(:,:,:,bi,bj) = comlev1_bibj , key=itdkey, byte=isbyte |
459 |
CADJ STORE recip_hFac(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte |
460 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
461 |
CALL GAD_IMPLICIT_R( |
462 |
I tempImplVertAdv, tempVertAdvScheme, GAD_TEMPERATURE, |
463 |
I dTtracerLev, kappaRk, recip_hFac, wFld, |
464 |
I theta(1-OLx,1-OLy,1,bi,bj), |
465 |
U gT_loc, |
466 |
I bi, bj, myTime, myIter, myThid ) |
467 |
ELSEIF ( implicitDiffusion ) THEN |
468 |
#else /* INCLUDE_IMPLVERTADV_CODE */ |
469 |
IF ( implicitDiffusion ) THEN |
470 |
#endif /* INCLUDE_IMPLVERTADV_CODE */ |
471 |
#ifdef ALLOW_AUTODIFF_TAMC |
472 |
CADJ STORE kappaRk(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte |
473 |
CADJ STORE gT_loc(:,:,:) = comlev1_bibj , key=itdkey, byte=isbyte |
474 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
475 |
CALL IMPLDIFF( |
476 |
I bi, bj, iMin, iMax, jMin, jMax, |
477 |
I GAD_TEMPERATURE, kappaRk, recip_hFac, |
478 |
U gT_loc, |
479 |
I myThid ) |
480 |
ENDIF |
481 |
|
482 |
#ifdef ALLOW_TIMEAVE |
483 |
useVariableK = useKPP .OR. usePP81 .OR. useKL10 .OR. useMY82 |
484 |
& .OR. useGGL90 .OR. useGMredi .OR. ivdc_kappa.NE.0. |
485 |
IF ( taveFreq.GT.0. .AND. useVariableK |
486 |
& .AND.implicitDiffusion ) THEN |
487 |
CALL TIMEAVE_CUMUL_DIF_1T( TdiffRtave, |
488 |
I gT_loc, kappaRk, |
489 |
I Nr, 3, deltaTClock, bi, bj, myThid ) |
490 |
ENDIF |
491 |
#endif /* ALLOW_TIMEAVE */ |
492 |
|
493 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
494 |
IF ( AdamsBashforth_T ) THEN |
495 |
C- Save current tracer field (for AB on tracer) and then update tracer |
496 |
CALL CYCLE_AB_TRACER( |
497 |
I bi, bj, |
498 |
U theta(1-OLx,1-OLy,1,bi,bj), gT_loc, |
499 |
U gtNm(1-OLx,1-OLy,1,bi,bj,m2), |
500 |
I myTime, myIter, myThid ) |
501 |
ELSEIF ( tempStepping ) THEN |
502 |
#else /* ALLOW_ADAMSBASHFORTH_3 */ |
503 |
IF ( tempStepping ) THEN |
504 |
#endif /* ALLOW_ADAMSBASHFORTH_3 */ |
505 |
C- Update tracer fields: T(n) = T** |
506 |
CALL CYCLE_TRACER( |
507 |
I bi, bj, |
508 |
U theta(1-OLx,1-OLy,1,bi,bj), gT_loc, |
509 |
I myTime, myIter, myThid ) |
510 |
ENDIF |
511 |
|
512 |
#endif /* ALLOW_GENERIC_ADVDIFF */ |
513 |
|
514 |
RETURN |
515 |
END |