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

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

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


Revision 1.110 - (show annotations) (download)
Sun Nov 28 23:49:03 2004 UTC (19 years, 7 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint56b_post, checkpoint57, checkpoint57a_post, checkpoint56a_post, checkpoint56c_post, checkpoint57a_pre
Changes since 1.109: +30 -24 lines
o GCHEM: finish reorganizating the package
  - forward_step calls GCHEM_CALC_TENDENDY, which computes gchemTendency
    (introduces another 3D-array for each passive tracer, but only if
    GCHEM_SEPARATE_FORCING is undefined. For GCHEM_SEPARATE_FORCING
    gchemTendency is not needed because the timestep is done separately)

1 C $Header: /u/gcmpack/MITgcm/model/src/forward_step.F,v 1.109 2004/11/23 20:20:19 mlosch Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 #ifdef ALLOW_OFFLINE
8 # include "OFFLINE_OPTIONS.h"
9 #endif
10 #ifdef ALLOW_GMREDI
11 # include "GMREDI_OPTIONS.h"
12 #endif
13
14 CBOP
15 C !ROUTINE: FORWARD_STEP
16 C !INTERFACE:
17 SUBROUTINE FORWARD_STEP( iloop, myTime, myIter, myThid )
18
19 C !DESCRIPTION: \bv
20 C *==================================================================
21 C | SUBROUTINE forward_step
22 C | o Run the ocean model and, optionally, evaluate a cost function.
23 C *==================================================================
24 C |
25 C | THE_MAIN_LOOP is the toplevel routine for the Tangent Linear and
26 C | Adjoint Model Compiler (TAMC). For this purpose the initialization
27 C | of the model was split into two parts. Those parameters that do
28 C | not depend on a specific model run are set in INITIALISE_FIXED,
29 C | whereas those that do depend on the specific realization are
30 C | initialized in INITIALISE_VARIA.
31 C |
32 C *==================================================================
33 C \ev
34
35 C !USES:
36 IMPLICIT NONE
37 C == Global variables ==
38 #include "SIZE.h"
39 #include "EEPARAMS.h"
40 #include "PARAMS.h"
41 #include "DYNVARS.h"
42
43 #ifdef ALLOW_SHAP_FILT
44 # include "SHAP_FILT.h"
45 #endif
46 #ifdef ALLOW_ZONAL_FILT
47 # include "ZONAL_FILT.h"
48 #endif
49 #ifdef COMPONENT_MODULE
50 # include "CPL_PARAMS.h"
51 #endif
52
53 #ifdef ALLOW_AUTODIFF_TAMC
54 # include "FFIELDS.h"
55
56 # ifdef ALLOW_NONHYDROSTATIC
57 # include "CG3D.h"
58 # endif
59
60 # include "tamc.h"
61 # include "ctrl.h"
62 # include "ctrl_dummy.h"
63 # include "cost.h"
64 # include "EOS.h"
65 # ifdef ALLOW_EXF
66 # include "exf_fields.h"
67 # include "exf_clim_fields.h"
68 # ifdef ALLOW_BULKFORMULAE
69 # include "exf_constants.h"
70 # endif
71 # endif
72 # ifdef ALLOW_OBCS
73 # include "OBCS.h"
74 # endif
75 # ifdef ALLOW_PTRACERS
76 # include "PTRACERS_SIZE.h"
77 # include "PTRACERS.h"
78 # endif
79 # ifdef ALLOW_CD_CODE
80 # include "CD_CODE_VARS.h"
81 # endif
82 # ifdef ALLOW_EBM
83 # include "EBM.h"
84 # endif
85 # ifdef EXACT_CONSERV
86 # include "SURFACE.h"
87 # endif
88 # ifdef ALLOW_KPP
89 # include "KPP.h"
90 # endif
91 # ifdef ALLOW_GMREDI
92 # include "GMREDI.h"
93 # endif
94 #endif /* ALLOW_AUTODIFF_TAMC */
95
96 C !LOCAL VARIABLES:
97 C == Routine arguments ==
98 C note: under the multi-threaded model myiter and
99 C mytime are local variables passed around as routine
100 C arguments. Although this is fiddly it saves the need to
101 C impose additional synchronisation points when they are
102 C updated.
103 C myIter - iteration counter for this thread
104 C myTime - time counter for this thread
105 C myThid - thread number for this instance of the routine.
106 INTEGER iloop
107 INTEGER myThid
108 INTEGER myIter
109 _RL myTime
110
111 C == Local variables ==
112 INTEGER myItP1
113 CEOP
114
115 #ifdef ALLOW_DEBUG
116 IF ( debugLevel .GE. debLevB )
117 & CALL DEBUG_ENTER('FORWARD_STEP',myThid)
118 #endif
119
120 #ifdef ALLOW_AUTODIFF_TAMC
121 C-- Reset the model iteration counter and the model time.
122 myiter = nIter0 + (iloop-1)
123 mytime = startTime + float(iloop-1)*deltaTclock
124 #endif
125
126 #ifdef ALLOW_AUTODIFF_TAMC
127 c**************************************
128 #include "checkpoint_lev1_directives.h"
129 c**************************************
130 #endif
131
132 C-- Call external forcing package
133 #ifdef ALLOW_BULK_FORCE
134 IF ( useBulkForce ) THEN
135 #ifdef ALLOW_DEBUG
136 IF ( debugLevel .GE. debLevB )
137 & CALL DEBUG_CALL('BULKF_FIELDS_LOAD',myThid)
138 #endif
139 CALL TIMER_START('BULKF_FORCING [FORWARD_STEP]',mythid)
140 C- load all forcing fields at current time
141 CALL BULKF_FIELDS_LOAD( myTime, myIter, myThid )
142 C- calculate qnet and empmr (and wind stress)
143 CALL BULKF_FORCING( myTime, myIter, myThid )
144 CALL TIMER_STOP ('BULKF_FORCING [FORWARD_STEP]',mythid)
145 ELSE
146 #endif /* ALLOW_BULK_FORCE */
147
148 # ifdef ALLOW_EXF
149 # ifdef ALLOW_DEBUG
150 IF ( debugLevel .GE. debLevB )
151 & CALL DEBUG_CALL('EXF_GETFORCING',myThid)
152 # endif
153 CALL TIMER_START('EXF_GETFORCING [FORWARD_STEP]',mythid)
154 CALL EXF_GETFORCING( mytime, myiter, mythid )
155 CALL TIMER_STOP ('EXF_GETFORCING [FORWARD_STEP]',mythid)
156 # else /* ALLOW_EXF undef */
157 cph The following IF-statement creates an additional dependency
158 cph for the forcing fields requiring additional storing.
159 cph Therefore, the IF-statement will be put between CPP-OPTIONS,
160 cph assuming that ALLOW_SEAICE has not yet been differentiated.
161 # if (defined (ALLOW_SEAICE) || defined (ALLOW_EBM))
162 IF ( .NOT. useSEAICE .AND. .NOT. useEBM ) THEN
163 # endif
164 #ifdef ALLOW_DEBUG
165 IF ( debugLevel .GE. debLevB )
166 & CALL DEBUG_CALL('EXTERNAL_FIELDS_LOAD',myThid)
167 #endif
168 CALL TIMER_START('EXTERNAL_FIELDS_LOAD[FORWARD_STEP]',mythid)
169 CALL EXTERNAL_FIELDS_LOAD( mytime, myiter, mythid )
170 CALL TIMER_STOP ('EXTERNAL_FIELDS_LOAD[FORWARD_STEP]',mythid)
171 # if (defined (ALLOW_SEAICE) || defined (ALLOW_EBM))
172 ENDIF
173 # endif
174 # endif /* ALLOW_EXF */
175 #ifdef ALLOW_BULK_FORCE
176 C-- end of if/else block useBulfforce --
177 ENDIF
178 #endif /* ALLOW_BULK_FORCE */
179
180 #ifdef ALLOW_AUTODIFF
181 c-- Add control vector for forcing and parameter fields
182 if ( myiter .EQ. nIter0 )
183 & CALL CTRL_MAP_FORCING (mythid)
184 #endif
185
186 #if (defined (ALLOW_AUTODIFF_TAMC) && defined (ALLOW_AUTODIFF_MONITOR))
187 C Include call to a dummy routine. Its adjoint will be
188 C called at the proper place in the adjoint code.
189 C The adjoint routine will print out adjoint values
190 C if requested. The location of the call is important,
191 C it has to be after the adjoint of the exchanges
192 C (DO_GTERM_BLOCKING_EXCHANGES).
193 CALL DUMMY_IN_STEPPING( myTime, myIter, myThid )
194 cph I've commented this line since it may conflict with MITgcm's adjoint
195 cph However, need to check whether that's still consistent
196 cph with the ecco-branch (it should).
197 cph CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
198 #endif
199
200 # ifdef ALLOW_SEAICE
201 C-- Call sea ice model to compute forcing/external data fields. In
202 C addition to computing prognostic sea-ice variables and diagnosing the
203 C forcing/external data fields that drive the ocean model, SEAICE_MODEL
204 C also sets theta to the freezing point under sea-ice. The implied
205 C surface heat flux is then stored in variable surfaceTendencyTice,
206 C which is needed by KPP package (kpp_calc.F and kpp_transport_t.F)
207 C to diagnose surface buoyancy fluxes and for the non-local transport
208 C term. Because this call precedes model thermodynamics, temperature
209 C under sea-ice may not be "exactly" at the freezing point by the time
210 C theta is dumped or time-averaged.
211 IF ( useSEAICE ) THEN
212 #ifdef ALLOW_DEBUG
213 IF ( debugLevel .GE. debLevB )
214 & CALL DEBUG_CALL('SEAICE_MODEL',myThid)
215 #endif
216 CALL TIMER_START('SEAICE_MODEL [FORWARD_STEP]',myThid)
217 CALL SEAICE_MODEL( myTime, myIter, myThid )
218 CALL TIMER_STOP ('SEAICE_MODEL [FORWARD_STEP]',myThid)
219 ENDIF
220 # endif /* ALLOW_SEAICE */
221
222 #ifdef ALLOW_AUTODIFF_TAMC
223 # ifdef ALLOW_PTRACERS
224 cph this replaces _bibj storing of ptracer within thermodynamics
225 CADJ STORE ptracer = comlev1, key = ikey_dynamics
226 # endif
227 #endif
228
229 #ifdef ALLOW_OFFLINE
230 call OFFLINE_FIELDS_LOAD( myTime, myIter, myThid )
231 #endif
232
233 #ifdef ALLOW_PTRACERS
234 # ifdef ALLOW_GCHEM
235 IF ( useGCHEM ) THEN
236 CALL GCHEM_FIELDS_LOAD( mytime, myiter, mythid )
237 ENDIF
238 # endif
239 #endif
240
241 #ifdef COMPONENT_MODULE
242 IF ( useCoupler .AND. cpl_earlyExpImpCall ) THEN
243 C Post coupling data that I export.
244 C Read in coupling data that I import.
245 CALL TIMER_START('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid)
246 CALL CPL_EXPORT_MY_DATA( myIter, myTime, myThid )
247 CALL CPL_IMPORT_EXTERNAL_DATA( myIter, myTime, myThid )
248 CALL TIMER_STOP ('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid)
249 ENDIF
250 #endif /* COMPONENT_MODULE */
251
252 #ifdef ALLOW_EBM
253 IF ( useEBM ) THEN
254 # ifdef ALLOW_DEBUG
255 IF ( debugLevel .GE. debLevB )
256 & CALL DEBUG_CALL('EBM',myThid)
257 # endif
258 CALL TIMER_START('EBM [FORWARD_STEP]',mythid)
259 CALL EBM_DRIVER ( myTime, myIter, myThid )
260 CALL TIMER_STOP ('EBM [FORWARD_STEP]',mythid)
261 ENDIF
262 #endif
263
264 C-- Step forward fields and calculate time tendency terms.
265
266 #ifdef ALLOW_DEBUG
267 IF ( debugLevel .GE. debLevB )
268 & CALL DEBUG_CALL('DO_ATMOSPHERIC_PHYS',myThid)
269 #endif
270 CALL TIMER_START('DO_ATMOSPHERIC_PHYS [FORWARD_STEP]',mythid)
271 CALL DO_ATMOSPHERIC_PHYS( myTime, myIter, myThid )
272 CALL TIMER_STOP ('DO_ATMOSPHERIC_PHYS [FORWARD_STEP]',mythid)
273
274 #ifdef ALLOW_AUTODIFF_TAMC
275 CADJ STORE theta = comlev1, key = ikey_dynamics
276 CADJ STORE salt = comlev1, key = ikey_dynamics
277 CADJ STORE totphihyd = comlev1, key = ikey_dynamics
278 CADJ STORE surfaceforcingtice = comlev1, key = ikey_dynamics
279 # ifdef ALLOW_KPP
280 CADJ STORE uvel = comlev1, key = ikey_dynamics
281 CADJ STORE vvel = comlev1, key = ikey_dynamics
282 # endif
283 # ifdef EXACT_CONSERV
284 CADJ STORE empmr = comlev1, key = ikey_dynamics
285 CADJ STORE pmepr = comlev1, key = ikey_dynamics
286 # endif
287 #endif /* ALLOW_AUTODIFF_TAMC */
288
289 #ifndef ALLOW_OFFLINE
290 #ifdef ALLOW_DEBUG
291 IF ( debugLevel .GE. debLevB )
292 & CALL DEBUG_CALL('DO_OCEANIC_PHYS',myThid)
293 #endif
294 CALL TIMER_START('DO_OCEANIC_PHYS [FORWARD_STEP]',mythid)
295 CALL DO_OCEANIC_PHYS( myTime, myIter, myThid )
296 CALL TIMER_STOP ('DO_OCEANIC_PHYS [FORWARD_STEP]',mythid)
297 #endif
298
299 #ifdef ALLOW_GCHEM
300 C GCHEM package is an interface for any bio-geochemical or
301 C ecosystem model you would like to include.
302 C If GCHEM_SEPARATE_FORCING is not defined, you are
303 C responsible for computing tendency terms for passive
304 C tracers and storing them on a 3DxNumPtracers-array called
305 C gchemTendency in GCHEM_CALC_TENDENCY. This tendency is then added
306 C to gPtr in ptracers_forcing later-on.
307 C If GCHEM_SEPARATE_FORCING is defined, you are reponsible for
308 C UPDATING ptracers directly in GCHEM_FORCING_SEP. This amounts
309 C to a completely separate time step that you have to implement
310 C yourself (Eulerian seems to be fine in most cases).
311 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
312 C CAVEAT: Up to now, when GCHEM is turned on the field ptracerForcingSurf,
313 C which is needed for KPP is not set properly. ptracerForcingSurf must
314 C be treated differently depending on whether GCHEM_SEPARATE_FORCING
315 C is define or not. TBD.
316 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
317 #ifdef ALLOW_DEBUG
318 IF ( debugLevel .GE. debLevB )
319 & CALL DEBUG_CALL('GCHEM_CALC_TENDENCY',myThid)
320 #endif
321 IF ( useGCHEM ) THEN
322 CALL TIMER_START('GCHEM_CALC_TENDENCY [FORWARD_STEP]',myThid)
323 CALL GCHEM_CALC_TENDENCY( myTime, myIter, myThid )
324 CALL TIMER_STOP ('GCHEM_CALC_TENDENCY [FORWARD_STEP]',myThid)
325 ENDIF
326 #endif /* ALLOW_GCHEM */
327
328 #ifdef ALLOW_AUTODIFF_TAMC
329 cph needed to be moved here from do_oceanic_physics
330 cph to be visible down the road
331 c
332 CADJ STORE surfaceForcingS = comlev1, key = ikey_dynamics
333 CADJ STORE surfaceForcingT = comlev1, key = ikey_dynamics
334 CADJ STORE surfaceForcingTice = comlev1, key = ikey_dynamics
335 ctest(
336 CADJ STORE IVDConvCount = comlev1, key = ikey_dynamics
337 ctest)
338 # ifdef ALLOW_PTRACERS
339 CADJ STORE surfaceForcingPtr = comlev1, key = ikey_dynamics
340 # endif
341 c
342 # ifdef ALLOW_GMREDI
343 CADJ STORE Kwx = comlev1, key = ikey_dynamics
344 CADJ STORE Kwy = comlev1, key = ikey_dynamics
345 CADJ STORE Kwz = comlev1, key = ikey_dynamics
346 # ifdef GM_BOLUS_ADVEC
347 CADJ STORE GM_PsiX = comlev1, key = ikey_dynamics
348 CADJ STORE GM_PsiY = comlev1, key = ikey_dynamics
349 # endif
350 # endif
351 c
352 # ifdef ALLOW_KPP
353 CADJ STORE KPPghat = comlev1, key = ikey_dynamics
354 CADJ STORE KPPfrac = comlev1, key = ikey_dynamics
355 CADJ STORE KPPdiffKzS = comlev1, key = ikey_dynamics
356 CADJ STORE KPPdiffKzT = comlev1, key = ikey_dynamics
357 # endif
358 #endif /* ALLOW_AUTODIFF_TAMC */
359
360 IF ( .NOT.staggerTimeStep ) THEN
361 #ifdef ALLOW_DEBUG
362 IF ( debugLevel .GE. debLevB )
363 & CALL DEBUG_CALL('THERMODYNAMICS',myThid)
364 #endif
365 CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',mythid)
366 CALL THERMODYNAMICS( myTime, myIter, myThid )
367 CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',mythid)
368 C-- if not staggerTimeStep: end
369 ENDIF
370
371 #ifdef COMPONENT_MODULE
372 IF ( useCoupler .AND. .NOT.cpl_earlyExpImpCall ) THEN
373 C Post coupling data that I export.
374 C Read in coupling data that I import.
375 myItP1 = myIter + 1
376 CALL TIMER_START('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid)
377 CALL CPL_EXPORT_MY_DATA( myItP1, myTime, myThid )
378 CALL CPL_IMPORT_EXTERNAL_DATA( myItP1, myTime, myThid )
379 CALL TIMER_STOP ('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid)
380 # ifndef ALLOW_AIM
381 IF ( useRealFreshWaterFlux ) THEN
382 CALL OCN_APPLY_IMPORT( .FALSE., myTime, myIter, myThid )
383 ENDIF
384 # endif
385 ENDIF
386 #endif /* COMPONENT_MODULE */
387
388 C-- Step forward fields and calculate time tendency terms.
389 #ifndef ALLOW_OFFLINE
390 #ifndef ALLOW_AUTODIFF_TAMC
391 IF ( momStepping ) THEN
392 #endif
393 #ifdef ALLOW_DEBUG
394 IF ( debugLevel .GE. debLevB )
395 & CALL DEBUG_CALL('DYNAMICS',myThid)
396 #endif
397 CALL TIMER_START('DYNAMICS [FORWARD_STEP]',mythid)
398 CALL DYNAMICS( myTime, myIter, myThid )
399 CALL TIMER_STOP ('DYNAMICS [FORWARD_STEP]',mythid)
400 #ifndef ALLOW_AUTODIFF_TAMC
401 ENDIF
402 #endif
403 #endif
404
405 #ifdef ALLOW_NONHYDROSTATIC
406 C-- Step forward W field in N-H algorithm
407 IF ( momStepping .AND. nonHydrostatic ) THEN
408 #ifdef ALLOW_DEBUG
409 IF ( debugLevel .GE. debLevB )
410 & CALL DEBUG_CALL('CALC_GW',myThid)
411 #endif
412 CALL TIMER_START('CALC_GW [FORWARD_STEP]',myThid)
413 CALL CALC_GW(myThid)
414 CALL TIMER_STOP ('CALC_GW [FORWARD_STEP]',myThid)
415 ENDIF
416 #endif
417
418 C-- Update time-counter
419 myIter = nIter0 + iLoop
420 myTime = startTime + deltaTClock * float(iLoop)
421
422 C-- Update geometric factors:
423 #ifdef NONLIN_FRSURF
424 C- update hfacC,W,S and recip_hFac according to etaH(n+1) :
425 IF ( nonlinFreeSurf.GT.0) THEN
426 IF ( select_rStar.GT.0 ) THEN
427 CALL TIMER_START('UPDATE_R_STAR [FORWARD_STEP]',myThid)
428 CALL UPDATE_R_STAR( myTime, myIter, myThid )
429 CALL TIMER_STOP ('UPDATE_R_STAR [FORWARD_STEP]',myThid)
430 ELSE
431 CALL TIMER_START('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
432 CALL UPDATE_SURF_DR( myTime, myIter, myThid )
433 CALL TIMER_STOP ('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
434 ENDIF
435 ENDIF
436 C- update also CG2D matrix (and preconditioner)
437 IF ( momStepping .AND. nonlinFreeSurf.GT.2 ) THEN
438 CALL TIMER_START('UPDATE_CG2D [FORWARD_STEP]',myThid)
439 CALL UPDATE_CG2D( myTime, myIter, myThid )
440 CALL TIMER_STOP ('UPDATE_CG2D [FORWARD_STEP]',myThid)
441 ENDIF
442 #endif
443
444 C-- Apply Filters to u*,v* before SOLVE_FOR_PRESSURE
445 #ifdef ALLOW_SHAP_FILT
446 IF (useSHAP_FILT .AND. shap_filt_uvStar) THEN
447 CALL TIMER_START('SHAP_FILT [FORWARD_STEP]',myThid)
448 IF (implicDiv2Dflow.LT.1.) THEN
449 C-- Explicit+Implicit part of the Barotropic Flow Divergence
450 C => Filtering of uVel,vVel is necessary
451 CALL SHAP_FILT_APPLY_UV( uVel,vVel,
452 & myTime, myIter, myThid )
453 ENDIF
454 CALL SHAP_FILT_APPLY_UV( gU,gV,myTime,myIter,myThid)
455 CALL TIMER_STOP ('SHAP_FILT [FORWARD_STEP]',myThid)
456 ENDIF
457 #endif
458 #ifdef ALLOW_ZONAL_FILT
459 IF (useZONAL_FILT .AND. zonal_filt_uvStar) THEN
460 CALL TIMER_START('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid)
461 IF (implicDiv2Dflow.LT.1.) THEN
462 C-- Explicit+Implicit part of the Barotropic Flow Divergence
463 C => Filtering of uVel,vVel is necessary
464 CALL ZONAL_FILT_APPLY_UV( uVel, vVel, myThid )
465 ENDIF
466 CALL ZONAL_FILT_APPLY_UV( gU, gV, myThid )
467 CALL TIMER_STOP ('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid)
468 ENDIF
469 #endif
470
471 C-- Solve elliptic equation(s).
472 C Two-dimensional only for conventional hydrostatic or
473 C three-dimensional for non-hydrostatic and/or IGW scheme.
474 #ifndef ALLOW_OFFLINE
475 IF ( momStepping ) THEN
476 CALL TIMER_START('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
477 CALL SOLVE_FOR_PRESSURE(myTime, myIter, myThid)
478 CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
479 ENDIF
480 #endif
481
482 C-- Correct divergence in flow field and cycle time-stepping momentum
483 c IF ( momStepping ) THEN
484 #ifndef ALLOW_OFFLINE
485 CALL TIMER_START('UV_CORRECTION_STEP [FORWARD_STEP]',myThid)
486 CALL MOMENTUM_CORRECTION_STEP(myTime, myIter, myThid)
487 CALL TIMER_STOP ('UV_CORRECTION_STEP [FORWARD_STEP]',myThid)
488 #endif
489 c ENDIF
490
491 #ifdef EXACT_CONSERV
492 IF (exactConserv) THEN
493 C-- Update etaH(n+1) :
494 CALL TIMER_START('UPDATE_ETAH [FORWARD_STEP]',mythid)
495 CALL UPDATE_ETAH( myTime, myIter, myThid )
496 CALL TIMER_STOP ('UPDATE_ETAH [FORWARD_STEP]',mythid)
497 ENDIF
498 #endif /* EXACT_CONSERV */
499
500 #ifdef NONLIN_FRSURF
501 IF ( select_rStar.NE.0 ) THEN
502 C-- r* : compute the future level thickness according to etaH(n+1)
503 CALL TIMER_START('CALC_R_STAR [FORWARD_STEP]',mythid)
504 CALL CALC_R_STAR(etaH, myTime, myIter, myThid )
505 CALL TIMER_STOP ('CALC_R_STAR [FORWARD_STEP]',mythid)
506 ELSEIF ( nonlinFreeSurf.GT.0) THEN
507 C-- compute the future surface level thickness according to etaH(n+1)
508 CALL TIMER_START('CALC_SURF_DR [FORWARD_STEP]',mythid)
509 CALL CALC_SURF_DR(etaH, myTime, myIter, myThid )
510 CALL TIMER_STOP ('CALC_SURF_DR [FORWARD_STEP]',mythid)
511 ENDIF
512 #endif /* NONLIN_FRSURF */
513
514 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
515 IF ( staggerTimeStep ) THEN
516 C-- do exchanges of U,V (needed for multiDim) when using stagger time-step :
517 #ifdef ALLOW_DEBUG
518 IF ( debugLevel .GE. debLevB )
519 & CALL DEBUG_CALL('DO_STAGGER_FIELDS_EXCH.',myThid)
520 #endif
521 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
522 CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid )
523 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
524
525 #ifdef ALLOW_DEBUG
526 IF ( debugLevel .GE. debLevB )
527 & CALL DEBUG_CALL('THERMODYNAMICS',myThid)
528 #endif
529 CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',mythid)
530 CALL THERMODYNAMICS( myTime, myIter, myThid )
531 CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',mythid)
532
533 C-- if staggerTimeStep: end
534 ENDIF
535 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
536
537 #ifdef ALLOW_AUTODIFF_TAMC
538 cph This is needed because convective_adjustment calls
539 cph find_rho which may use pressure()
540 CADJ STORE totphihyd = comlev1, key = ikey_dynamics
541 #endif
542 C-- Cycle time-stepping Tracers arrays (T,S,+pTracers)
543 CALL TIMER_START('TS_CORRECTION_STEP [FORWARD_STEP]',myThid)
544 CALL TRACERS_CORRECTION_STEP(myTime, myIter, myThid)
545 CALL TIMER_STOP ('TS_CORRECTION_STEP [FORWARD_STEP]',myThid)
546
547 #ifdef ALLOW_GCHEM
548 C Add separate timestepping of chemical/biological/forcing
549 C of ptracers here in GCHEM_FORCING_SEP
550 #ifdef ALLOW_DEBUG
551 IF ( debugLevel .GE. debLevB )
552 & CALL DEBUG_CALL('GCHEM_FORCING_SEP',myThid)
553 #endif /* ALLOW_DEBUG */
554 IF ( useGCHEM ) THEN
555 CALL TIMER_START('GCHEM_FORCING_SEP [FORWARD_STEP]',myThid)
556 CALL GCHEM_FORCING_SEP( myTime,myIter,myThid )
557 CALL TIMER_STOP ('GCHEM_FORCING_SEP [FORWARD_STEP]',myThid)
558 ENDIF
559 #endif /* ALLOW_GCHEM */
560
561 C-- Do "blocking" sends and receives for tendency "overlap" terms
562 c CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
563 c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )
564 c CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
565
566 C-- Do "blocking" sends and receives for field "overlap" terms
567 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
568 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
569 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
570
571
572 C AMM
573 #ifdef ALLOW_GRIDALT
574 if (useGRIDALT) then
575 CALL GRIDALT_UPDATE(myThid)
576 endif
577 #endif
578 C AMM
579
580 C AMM
581 #ifdef ALLOW_FIZHI
582 if( useFIZHI) then
583 CALL TIMER_START('FIZHI [FORWARD_STEP]',mythid)
584 CALL STEP_FIZHI_CORR ( myTime, myIter, myThid )
585 CALL TIMER_STOP('FIZHI [FORWARD_STEP]',mythid)
586 endif
587 #endif
588 C AMM
589
590 #ifdef ALLOW_FLT
591 C-- Calculate float trajectories
592 IF (useFLT) THEN
593 CALL TIMER_START('FLOATS [FORWARD_STEP]',myThid)
594 CALL FLT_MAIN(myIter,myTime, myThid)
595 CALL TIMER_STOP ('FLOATS [FORWARD_STEP]',myThid)
596 ENDIF
597 #endif
598
599 C-- State-variables statistics (time-aver, diagnostics ...)
600 CALL TIMER_START('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
601 CALL DO_STATEVARS_DIAGS( myTime, myIter, myThid )
602 CALL TIMER_STOP ('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
603
604 #ifndef ALLOW_OFFLINE
605 #ifdef ALLOW_MONITOR
606 C-- Check status of solution (statistics, cfl, etc...)
607 CALL TIMER_START('MONITOR [FORWARD_STEP]',myThid)
608 CALL MONITOR( myIter, myTime, myThid )
609 CALL TIMER_STOP ('MONITOR [FORWARD_STEP]',myThid)
610 #endif /* ALLOW_MONITOR */
611 #endif
612
613 #ifdef ALLOW_COST
614 C-- compare model with data and compute cost function
615 C-- this is done after exchanges to allow interpolation
616 CALL TIMER_START('COST_TILE [FORWARD_STEP]',myThid)
617 CALL COST_TILE ( mytime, myiter, myThid )
618 CALL TIMER_STOP ('COST_TILE [FORWARD_STEP]',myThid)
619 #endif
620
621 C-- Do IO if needed.
622 #ifdef ALLOW_OFFLINE
623 CALL TIMER_START('OFFLINE_MODEL_IO [FORWARD_STEP]',myThid)
624 CALL OFFLINE_MODEL_IO( myTime, myIter, myThid )
625 CALL TIMER_STOP ('OFFLINE_MODEL_IO [FORWARD_STEP]',myThid)
626 #else
627 CALL TIMER_START('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
628 CALL DO_THE_MODEL_IO( myTime, myIter, myThid )
629 CALL TIMER_STOP ('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
630 #endif
631
632 C-- Save state for restarts
633 CALL TIMER_START('WRITE_CHECKPOINT [FORWARD_STEP]',myThid)
634 CALL PACKAGES_WRITE_PICKUP(
635 I .FALSE., myTime, myIter, myThid )
636 #ifndef ALLOW_OFFLINE
637 CALL WRITE_CHECKPOINT(
638 I .FALSE., myTime, myIter, myThid )
639 #endif
640 CALL TIMER_STOP ('WRITE_CHECKPOINT [FORWARD_STEP]',myThid)
641
642 #ifdef ALLOW_DEBUG
643 IF ( debugLevel .GE. debLevB )
644 & CALL DEBUG_LEAVE('FORWARD_STEP',myThid)
645 #endif
646
647 RETURN
648 END

  ViewVC Help
Powered by ViewVC 1.1.22