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

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

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


Revision 1.13 - (show annotations) (download)
Thu Aug 7 17:32:03 2014 UTC (9 years, 10 months ago) by jmc
Branch: MAIN
Changes since 1.12: +25 -3 lines
move CYCLE_TRACER calls from tracers_correction_step.F to temp/salt/ptracer_integrate.F
 so that theta,salt and pTracers arrays are already updated when leaving
 S/R THERMODYANMICS while adjustments (filters, conv.adjustment) are still
 applied later, in S/R TRACERS_CORRECTION_STEP.

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

  ViewVC Help
Powered by ViewVC 1.1.22