/[MITgcm]/MITgcm_contrib/dgoldberg/code_ad_flowline/forward_step.F
ViewVC logotype

Contents of /MITgcm_contrib/dgoldberg/code_ad_flowline/forward_step.F

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


Revision 1.1 - (show annotations) (download)
Fri Jul 20 13:42:08 2012 UTC (13 years ago) by heimbach
Branch: MAIN
CVS Tags: HEAD
Add a complete config. directory

1 C $Header: /u/gcmpack/MITgcm/model/src/forward_step.F,v 1.196 2011/10/04 00:33:51 heimbach Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 #ifdef ALLOW_GMREDI
8 # include "GMREDI_OPTIONS.h"
9 #endif
10 #ifdef ALLOW_OBCS
11 # include "OBCS_OPTIONS.h"
12 #endif
13 #ifdef ALLOW_SEAICE
14 # include "SEAICE_OPTIONS.h"
15 #endif
16 #ifdef ALLOW_PTRACERS
17 # include "PTRACERS_OPTIONS.h"
18 #endif
19 #ifdef ALLOW_STREAMICE
20 # include "STREAMICE_OPTIONS.h"
21 #endif
22
23 CBOP
24 C !ROUTINE: FORWARD_STEP
25 C !INTERFACE:
26 SUBROUTINE FORWARD_STEP( iloop, myTime, myIter, myThid )
27
28 C !DESCRIPTION: \bv
29 C *==================================================================
30 C | SUBROUTINE forward_step
31 C | o Run the ocean model and, optionally, evaluate a cost function.
32 C *==================================================================
33 C |
34 C | THE_MAIN_LOOP is the toplevel routine for the Tangent Linear and
35 C | Adjoint Model Compiler (TAMC). For this purpose the initialization
36 C | of the model was split into two parts. Those parameters that do
37 C | not depend on a specific model run are set in INITIALISE_FIXED,
38 C | whereas those that do depend on the specific realization are
39 ! C | initialized in INITIALISE_VARIA.
40 C |
41 C *==================================================================
42 C \ev
43
44 C !USES:
45 IMPLICIT NONE
46 C == Global variables ==
47 #include "SIZE.h"
48 #include "EEPARAMS.h"
49 #include "PARAMS.h"
50 #include "DYNVARS.h"
51
52 #ifdef HAVE_SIGREG
53 #include "SIGREG.h"
54 #endif
55
56 #ifdef ALLOW_SHAP_FILT
57 # include "SHAP_FILT.h"
58 #endif
59 #ifdef ALLOW_ZONAL_FILT
60 # include "ZONAL_FILT.h"
61 #endif
62 #ifdef COMPONENT_MODULE
63 # include "CPL_PARAMS.h"
64 #endif
65
66 #ifdef ALLOW_LONGSTEP
67 # include "LONGSTEP_PARAMS.h"
68 # include "LONGSTEP.h"
69 #endif
70
71 #ifdef ALLOW_AUTODIFF_TAMC
72 # include "AUTODIFF_MYFIELDS.h"
73 # include "FFIELDS.h"
74 # include "SURFACE.h"
75
76 # include "tamc.h"
77 # include "ctrl.h"
78 # include "ctrl_dummy.h"
79 # include "cost.h"
80 # ifdef ALLOW_ECCO
81 # include "ecco_cost.h"
82 # endif
83 # include "EOS.h"
84 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
85 # include "GRID.h"
86 # endif
87 # ifdef ALLOW_EXF
88 # include "EXF_FIELDS.h"
89 # ifdef ALLOW_BULKFORMULAE
90 # include "EXF_CONSTANTS.h"
91 # endif
92 # endif
93 # ifdef ALLOW_PTRACERS
94 # include "PTRACERS_SIZE.h"
95 # include "PTRACERS_FIELDS.h"
96 # endif
97 # ifdef ALLOW_GCHEM
98 # include "GCHEM_FIELDS.h"
99 # endif
100 # ifdef ALLOW_CFC
101 # include "CFC.h"
102 # endif
103 # ifdef ALLOW_DIC
104 # include "DIC_VARS.h"
105 # include "DIC_LOAD.h"
106 # include "DIC_ATMOS.h"
107 # include "DIC_COST.h"
108 # endif
109 # ifdef ALLOW_OBCS
110 # include "OBCS_FIELDS.h"
111 # include "OBCS_SEAICE.h"
112 # ifdef ALLOW_PTRACERS
113 # include "OBCS_PTRACERS.h"
114 # endif
115 # endif
116 # ifdef ALLOW_CD_CODE
117 # include "CD_CODE_VARS.h"
118 # endif
119 # ifdef ALLOW_THSICE
120 # include "THSICE_VARS.h"
121 # endif
122 # ifdef ALLOW_SEAICE
123 # include "SEAICE.h"
124 # include "SEAICE_COST.h"
125 # endif
126 # ifdef ALLOW_SALT_PLUME
127 # include "SALT_PLUME.h"
128 # endif
129 # ifdef ALLOW_SHELFICE
130 # include "SHELFICE.h"
131 # include "SHELFICE_COST.h"
132 # endif
133 # ifdef ALLOW_STREAMICE
134 # include "STREAMICE.h"
135 # include "STREAMICE_ADV.h"
136 # include "STREAMICE_BDRY.h"
137 # include "STREAMICE_CG.h"
138 # endif
139 # ifdef ALLOW_EBM
140 # include "EBM.h"
141 # endif
142 # ifdef ALLOW_KPP
143 # include "KPP.h"
144 # endif
145 # ifdef ALLOW_GGL90
146 # include "GGL90.h"
147 # endif
148 # ifdef ALLOW_GMREDI
149 # include "GMREDI.h"
150 # endif
151 # ifdef ALLOW_RBCS
152 # include "RBCS_SIZE.h"
153 # include "RBCS_FIELDS.h"
154 # endif
155 # ifdef ALLOW_OFFLINE
156 # include "OFFLINE.h"
157 # endif
158 #endif /* ALLOW_AUTODIFF_TAMC */
159
160 C !INPUT/OUTPUT PARAMETERS:
161 C == Routine arguments ==
162 C note: under the multi-threaded model myIter and
163 C myTime are local variables passed around as routine
164 C arguments. Although this is fiddly it saves the need to
165 C impose additional synchronisation points when they are
166 C updated.
167 C myTime :: time counter for this thread
168 C myIter :: iteration counter for this thread
169 C myThid :: thread number for this instance of the routine.
170 INTEGER iloop
171 _RL myTime
172 INTEGER myIter
173 INTEGER myThid
174
175 C !LOCAL VARIABLES:
176 C == Local variables ==
177 C modelEnd :: true if reaching the end of the run
178 C myTimeBeg :: time at beginning of time step (needed by longstep)
179 C myIterBeg :: iteration number at beginning of time step
180 LOGICAL modelEnd
181 #ifdef COMPONENT_MODULE
182 INTEGER myItP1
183 #endif
184 #ifdef ALLOW_LONGSTEP
185 INTEGER myIterBeg
186 _RL myTimeBeg
187 #endif /* ALLOW_LONGSTEP */
188 CEOP
189
190 #ifdef ALLOW_DEBUG
191 IF (debugMode) CALL DEBUG_ENTER('FORWARD_STEP',myThid)
192 #endif
193
194 #ifdef ALLOW_AUTODIFF_TAMC
195 CALL AUTODIFF_INADMODE_UNSET( myThid )
196 #endif
197
198 #ifdef ALLOW_AUTODIFF_TAMC
199 C-- Reset the model iteration counter and the model time.
200 myIter = nIter0 + (iloop-1)
201 myTime = startTime + float(iloop-1)*deltaTclock
202 #endif
203
204 #ifdef ALLOW_LONGSTEP
205 C store this for longstep_average with staggerTimeStep
206 C which is called after myIter and myTime are incremented
207 C but needs iter/time at beginning of time step
208 myIterBeg = myIter
209 myTimeBeg = myTime
210 #endif /* ALLOW_LONGSTEP */
211
212 #ifdef ALLOW_AUTODIFF_TAMC
213 c**************************************
214 #include "checkpoint_lev1_directives.h"
215 #include "checkpoint_lev1_template.h"
216 c**************************************
217 #endif
218
219 C-- Reset geometric factors to their current values
220 Cgf (only has an impact for the adjoint)
221 #ifdef NONLIN_FRSURF
222 C- update hfacC,W,S and recip_hFac according to etaH(n+1) :
223 IF ( select_rStar.GT.0 ) THEN
224 # ifndef DISABLE_RSTAR_CODE
225 # ifdef ALLOW_AUTODIFF_TAMC
226 CADJ STORE rstarFacC, rstarFacS, rstarFacW =
227 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
228 # endif
229 CALL TIMER_START('UPDATE_R_STAR [FORWARD_STEP]',myThid)
230 CALL UPDATE_R_STAR( .FALSE., myTime, myIter, myThid )
231 CALL TIMER_STOP ('UPDATE_R_STAR [FORWARD_STEP]',myThid)
232 # endif /* DISABLE_RSTAR_CODE */
233 ELSE
234 #ifdef ALLOW_AUTODIFF_TAMC
235 CADJ STORE hFac_surfC, hFac_surfS, hFac_surfW
236 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
237 #endif
238 CALL TIMER_START('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
239 CALL UPDATE_SURF_DR( .FALSE., myTime, myIter, myThid )
240 CALL TIMER_STOP ('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
241 ENDIF
242 # ifdef ALLOW_AUTODIFF_TAMC
243 CADJ STORE hFacC, hFacS, hFacW =
244 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
245 CADJ STORE recip_hFacC, recip_hFacS, recip_hFacW =
246 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
247 # endif
248 #endif /* NONLIN_FRSURF */
249
250 C-- Switch on/off diagnostics for snap-shot output:
251 #ifdef ALLOW_DIAGNOSTICS
252 IF ( useDiagnostics ) THEN
253 CALL DIAGNOSTICS_SWITCH_ONOFF( myTime, myIter, myThid )
254 C-- State-variables diagnostics
255 CALL TIMER_START('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
256 CALL DO_STATEVARS_DIAGS( myTime, 0, myIter, myThid )
257 CALL TIMER_STOP ('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
258 ENDIF
259 #endif
260
261 #ifdef ALLOW_NEST_CHILD
262 IF ( useNEST_CHILD) THEN
263 CALL NEST_CHILD_SETMEMO( myTime, myIter, myThid )
264 ENDIF
265 #endif /* ALLOW_NEST_CHILD */
266
267 #ifdef ALLOW_NEST_PARENT
268 IF ( useNEST_PARENT) THEN
269 CALL NEST_PARENT_IO_1( myTime, myIter, myThid )
270 ENDIF
271 #endif /* ALLOW_NEST_PARENT */
272
273 #ifdef ALLOW_PROFILES
274 #ifdef ALLOW_DEBUG
275 IF (debugMode) CALL DEBUG_CALL('',myThid)
276 #endif
277 c-- Accumulate in-situ time averages of theta, salt, and SSH.
278 CALL TIMER_START('PROFILES_INLOOP [FORWARD_STEP]', mythid)
279 CALL PROFILES_INLOOP( mytime, mythid )
280 CALL TIMER_STOP ('PROFILES_INLOOP [FORWARD_STEP]', mythid)
281 #endif
282
283 C-- Call driver to load external forcing fields from file
284 #ifdef ALLOW_DEBUG
285 IF (debugMode) CALL DEBUG_CALL('LOAD_FIELDS_DRIVER',myThid)
286 #endif
287 #ifdef ALLOW_AUTODIFF_TAMC
288 cph Important STORE that avoids hidden recomp. of load_fields_driver
289 CADJ STORE theta = comlev1, key = ikey_dynamics,
290 CADJ & kind = isbyte
291 CADJ STORE uvel, vvel = comlev1, key = ikey_dynamics,
292 CADJ & kind = isbyte
293 #endif
294 CALL TIMER_START('LOAD_FIELDS_DRIVER [FORWARD_STEP]',myThid)
295 CALL LOAD_FIELDS_DRIVER( myTime, myIter, myThid )
296 CALL TIMER_STOP ('LOAD_FIELDS_DRIVER [FORWARD_STEP]',myThid)
297
298 C-- Call Bulk-Formulae forcing package
299 #ifdef ALLOW_BULK_FORCE
300 IF ( useBulkForce ) THEN
301 #ifdef ALLOW_DEBUG
302 IF (debugMode) CALL DEBUG_CALL('BULKF_FORCING',myThid)
303 #endif
304 CALL TIMER_START('BULKF_FORCING [FORWARD_STEP]',myThid)
305 C- calculate qnet and empmr (and wind stress)
306 CALL BULKF_FORCING( myTime, myIter, myThid )
307 CALL TIMER_STOP ('BULKF_FORCING [FORWARD_STEP]',myThid)
308 ENDIF
309 #endif /* ALLOW_BULK_FORCE */
310
311 C-- Call external chepaml forcing package
312 #ifdef ALLOW_CHEAPAML
313 IF ( useCheapAML ) THEN
314 #ifdef ALLOW_DEBUG
315 IF (debugMode) CALL DEBUG_CALL('CHEAPAML',myThid)
316 #endif
317 CALL TIMER_START('CHEAPAML [FORWARD_STEP]',mythid)
318 C- calculate qnet (and wind stress)
319 CALL CHEAPAML( myTime, myIter,myThid )
320 CALL TIMER_STOP ('CHEAPAML [FORWARD_STEP]',mythid)
321 ENDIF
322 #endif /*ALLOW_CHEAPAML */
323
324
325 #ifdef ALLOW_AUTODIFF
326 c-- Add control vector for forcing and parameter fields
327 IF ( myIter .EQ. nIter0 )
328 & CALL CTRL_MAP_FORCING (myThid)
329 #endif
330
331 #if (defined (ALLOW_AUTODIFF_TAMC) && defined (ALLOW_AUTODIFF_MONITOR))
332 CALL DUMMY_IN_STEPPING( myTime, myIter, myThid )
333 #endif
334
335 #ifdef COMPONENT_MODULE
336 IF ( useCoupler .AND. cpl_earlyExpImpCall ) THEN
337 C Post coupling data that I export.
338 C Read in coupling data that I import.
339 CALL TIMER_START('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid)
340 CALL CPL_EXPORT_MY_DATA( myTime, myIter, myThid )
341 CALL CPL_IMPORT_EXTERNAL_DATA( myTime, myIter, myThid )
342 CALL TIMER_STOP ('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid)
343 ENDIF
344 #endif /* COMPONENT_MODULE */
345 #ifdef ALLOW_OASIS
346 IF ( useOASIS ) THEN
347 CALL TIMER_START('OASIS_PUT-GET [FORWARD_STEP]',myThid)
348 C Post coupling data that I export.
349 CALL OASIS_PUT( myTime, myIter, myThid )
350 C Read in coupling data that I import.
351 CALL OASIS_GET( myTime, myIter, myThid )
352 CALL TIMER_STOP ('OASIS_PUT-GET [FORWARD_STEP]',myThid)
353 ENDIF
354 #endif /* ALLOW_OASIS */
355
356
357 #ifdef ALLOW_EBM
358 IF ( useEBM ) THEN
359 # ifdef ALLOW_DEBUG
360 IF (debugMode) CALL DEBUG_CALL('EBM',myThid)
361 # endif
362 CALL TIMER_START('EBM [FORWARD_STEP]',myThid)
363 CALL EBM_DRIVER ( myTime, myIter, myThid )
364 CALL TIMER_STOP ('EBM [FORWARD_STEP]',myThid)
365 ENDIF
366 #endif /* ALLOW_EBM */
367
368 C-- Step forward fields and calculate time tendency terms.
369
370 #ifdef ALLOW_DEBUG
371 IF (debugMode) CALL DEBUG_CALL('DO_ATMOSPHERIC_PHYS',myThid)
372 #endif
373 CALL TIMER_START('DO_ATMOSPHERIC_PHYS [FORWARD_STEP]',myThid)
374 CALL DO_ATMOSPHERIC_PHYS( myTime, myIter, myThid )
375 CALL TIMER_STOP ('DO_ATMOSPHERIC_PHYS [FORWARD_STEP]',myThid)
376
377 #ifdef ALLOW_AUTODIFF_TAMC
378 CADJ STORE surfaceforcingtice = comlev1, key = ikey_dynamics,
379 CADJ & kind = isbyte
380 # ifdef ALLOW_KPP
381 CADJ STORE uvel = comlev1, key = ikey_dynamics,
382 CADJ & kind = isbyte
383 CADJ STORE vvel = comlev1, key = ikey_dynamics,
384 CADJ & kind = isbyte
385 # endif /* ALLOW_KPP */
386 # ifdef EXACT_CONSERV
387 CADJ STORE empmr = comlev1, key=ikey_dynamics, kind=isbyte
388 CADJ STORE pmepr = comlev1, key=ikey_dynamics, kind=isbyte
389 # endif
390 # ifdef ALLOW_OBCS
391 CADJ STORE salt = comlev1, key=ikey_dynamics, kind=isbyte
392 CADJ STORE totphihyd = comlev1, key=ikey_dynamics, kind=isbyte
393 # ifdef ALLOW_OBCS_STEVENS
394 CADJ STORE gsnm1 = comlev1, key=ikey_dynamics, kind=isbyte
395 CADJ STORE gtnm1 = comlev1, key=ikey_dynamics, kind=isbyte
396 # endif
397 # endif /* ALLOW_OBCS */
398 # ifdef ALLOW_PTRACERS
399 CADJ STORE ptracer = comlev1, key = ikey_dynamics,
400 CADJ & kind = isbyte
401 # endif /* ALLOW_PTRACERS */
402 # if (defined ALLOW_DEPTH_CONTROL)
403 CADJ STORE hFacC = comlev1, key = ikey_dynamics, kind = isbyte
404 # ifndef DISABLE_RSTAR_CODE
405 CADJ STORE rstarexpc = comlev1, key = ikey_dynamics, kind = isbyte
406 # endif
407 # endif
408 #endif /* ALLOW_AUTODIFF_TAMC */
409
410 #ifdef ALLOW_OFFLINE
411 IF ( .NOT. useOffLine ) THEN
412 #endif
413 #ifdef ALLOW_DEBUG
414 IF (debugMode) CALL DEBUG_CALL('DO_OCEANIC_PHYS',myThid)
415 #endif
416 CALL TIMER_START('DO_OCEANIC_PHYS [FORWARD_STEP]',myThid)
417 CALL DO_OCEANIC_PHYS( myTime, myIter, myThid )
418 CALL TIMER_STOP ('DO_OCEANIC_PHYS [FORWARD_STEP]',myThid)
419
420
421 # ifdef ALLOW_STREAMICE
422 CALL STREAMICE_SOLO_TIMESTEP ( myThid, myiter,
423 & iLoop, mytime )
424 #endif
425
426 #ifdef ALLOW_AUTODIFF_TAMC
427 CADJ STORE EmPmR = comlev1, key = ikey_dynamics,
428 CADJ & kind = isbyte
429 # ifdef EXACT_CONSERV
430 CADJ STORE pmepr = comlev1, key = ikey_dynamics,
431 CADJ & kind = isbyte
432 # endif
433 #endif
434 #ifdef ALLOW_OFFLINE
435 ENDIF
436 #endif
437
438 #ifdef ALLOW_AUTODIFF_TAMC
439 # if (defined ALLOW_DEPTH_CONTROL)
440 CADJ STORE hFacC, hFacS, hFacW
441 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
442 CADJ STORE recip_hFacC, recip_hFacS, recip_hFacW
443 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
444 c
445 CADJ STORE surfaceforcingu, surfaceforcingv =
446 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
447 # endif
448 #endif /* ALLOW_AUTODIFF_TAMC */
449
450 #ifdef ALLOW_GCHEM
451 #ifdef ALLOW_AUTODIFF_TAMC
452 CADJ STORE ptracer = comlev1, key = ikey_dynamics,
453 CADJ & kind = isbyte
454 CADJ STORE theta = comlev1, key = ikey_dynamics,
455 CADJ & kind = isbyte
456 CADJ STORE salt = comlev1, key = ikey_dynamics,
457 CADJ & kind = isbyte
458 #endif
459 IF ( useGCHEM ) THEN
460 #ifdef ALLOW_DEBUG
461 IF (debugMode) CALL DEBUG_CALL('GCHEM_CALC_TENDENCY',myThid)
462 #endif
463 CALL TIMER_START('GCHEM_CALC_TENDENCY [FORWARD_STEP]',myThid)
464 CALL GCHEM_CALC_TENDENCY( myTime, myIter, myThid )
465 CALL TIMER_STOP ('GCHEM_CALC_TENDENCY [FORWARD_STEP]',myThid)
466 ENDIF
467 #endif /* ALLOW_GCHEM */
468
469 #ifdef ALLOW_AUTODIFF_TAMC
470 cph needed to be moved here from do_oceanic_physics
471 cph to be visible down the road
472 c
473 CADJ STORE rhoInSitu = comlev1, key = ikey_dynamics,
474 CADJ & kind = isbyte
475 CADJ STORE surfaceForcingS = comlev1, key = ikey_dynamics,
476 CADJ & kind = isbyte
477 CADJ STORE surfaceForcingT = comlev1, key = ikey_dynamics,
478 CADJ & kind = isbyte
479 CADJ STORE surfaceForcingTice = comlev1, key = ikey_dynamics,
480 CADJ & kind = isbyte
481 CADJ STORE IVDConvCount = comlev1, key = ikey_dynamics,
482 CADJ & kind = isbyte
483 # ifdef ALLOW_PTRACERS
484 CADJ STORE surfaceForcingPTr = comlev1, key = ikey_dynamics,
485 CADJ & kind = isbyte
486 # endif
487 c
488 # ifdef ALLOW_GMREDI
489 CADJ STORE Kwx = comlev1, key = ikey_dynamics,
490 CADJ & kind = isbyte
491 CADJ STORE Kwy = comlev1, key = ikey_dynamics,
492 CADJ & kind = isbyte
493 CADJ STORE Kwz = comlev1, key = ikey_dynamics,
494 CADJ & kind = isbyte
495 # ifdef GM_BOLUS_ADVEC
496 CADJ STORE GM_PsiX = comlev1, key = ikey_dynamics,
497 CADJ & kind = isbyte
498 CADJ STORE GM_PsiY = comlev1, key = ikey_dynamics,
499 CADJ & kind = isbyte
500 # endif
501 # endif
502 c
503 # ifdef ALLOW_KPP
504 CADJ STORE KPPghat = comlev1, key = ikey_dynamics,
505 CADJ & kind = isbyte
506 CADJ STORE KPPfrac = comlev1, key = ikey_dynamics,
507 CADJ & kind = isbyte
508 CADJ STORE KPPdiffKzS = comlev1, key = ikey_dynamics,
509 CADJ & kind = isbyte
510 CADJ STORE KPPdiffKzT = comlev1, key = ikey_dynamics,
511 CADJ & kind = isbyte
512 # endif
513 c
514 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
515 CADJ STORE theta,salt,wvel = comlev1, key = ikey_dynamics, kind = isbyte
516 CADJ STORE etaH = comlev1, key = ikey_dynamics, kind = isbyte
517 # ifdef ALLOW_CD_CODE
518 CADJ STORE etanm1 = comlev1, key = ikey_dynamics, kind = isbyte
519 # endif
520 # ifndef DISABLE_RSTAR_CODE
521 CADJ STORE rstarexpc = comlev1, key = ikey_dynamics, kind = isbyte
522 # endif
523 # endif
524 #endif /* ALLOW_AUTODIFF_TAMC */
525
526 #ifdef ALLOW_LONGSTEP
527 IF ( usePTRACERS ) THEN
528 IF ( LS_whenToSample .EQ. 0 ) THEN
529 C Average all variables before advection (but after do_oceanic_phys
530 C where Qsw, KPP and GMRedi stuff is computed).
531 C This is like diagnostics package and will reproduce offline
532 C results.
533 #ifdef ALLOW_DEBUG
534 IF (debugMode) CALL DEBUG_CALL('LONGSTEP_AVERAGE',myThid)
535 #endif
536 CALL TIMER_START('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid)
537 CALL LONGSTEP_AVERAGE( myTime, myIter, myThid )
538 CALL TIMER_STOP ('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid)
539
540 #ifdef ALLOW_DEBUG
541 IF (debugMode)
542 & CALL DEBUG_CALL('LONGSTEP_THERMODYNAMICS',myThid)
543 #endif
544 CALL TIMER_START('LONGSTEP_THERMODYNAMICS [FORWARD_STEP]',
545 & myThid)
546 CALL LONGSTEP_THERMODYNAMICS( myTime, myIter, myThid )
547 CALL TIMER_STOP ('LONGSTEP_THERMODYNAMICS [FORWARD_STEP]',
548 & myThid)
549 ENDIF
550 ENDIF
551 #endif /* ALLOW_LONGSTEP */
552
553 IF ( .NOT.staggerTimeStep ) THEN
554 #ifdef ALLOW_DEBUG
555 IF (debugMode) CALL DEBUG_CALL('THERMODYNAMICS',myThid)
556 #endif
557 #ifdef ALLOW_AUTODIFF_TAMC
558 CADJ STORE salt = comlev1, key = ikey_dynamics,
559 CADJ & kind = isbyte
560 #endif
561 CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',myThid)
562 CALL THERMODYNAMICS( myTime, myIter, myThid )
563 CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',myThid)
564 C-- if not staggerTimeStep: end
565 ENDIF
566
567 #ifdef ALLOW_LONGSTEP
568 IF ( usePTRACERS ) THEN
569 IF ( LS_whenToSample .EQ. 1 ) THEN
570 C Average T and S after thermodynamics, but U,V,W before dynamics.
571 C This will reproduce online results with staggerTimeStep=.FALSE.
572 C for LS_nIter=1
573 #ifdef ALLOW_DEBUG
574 IF (debugMode) CALL DEBUG_CALL('LONGSTEP_AVERAGE',myThid)
575 #endif
576 CALL TIMER_START('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid)
577 CALL LONGSTEP_AVERAGE( myTime, myIter, myThid )
578 CALL TIMER_STOP ('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid)
579
580 #ifdef ALLOW_DEBUG
581 IF (debugMode)
582 & CALL DEBUG_CALL('LONGSTEP_THERMODYNAMICS',myThid)
583 #endif
584 CALL TIMER_START('LONGSTEP_THERMODYNAMICS [FORWARD_STEP]',
585 & myThid)
586 CALL LONGSTEP_THERMODYNAMICS( myTime, myIter, myThid )
587 CALL TIMER_STOP ('LONGSTEP_THERMODYNAMICS [FORWARD_STEP]',
588 & myThid)
589 ENDIF
590 ENDIF
591 #endif /* ALLOW_LONGSTEP */
592
593 c #ifdef ALLOW_NONHYDROSTATIC
594 IF ( implicitIntGravWave ) THEN
595 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
596 CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid )
597 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
598 ENDIF
599 c #endif
600
601 #ifdef COMPONENT_MODULE
602 IF ( useCoupler .AND. .NOT.cpl_earlyExpImpCall ) THEN
603 C Post coupling data that I export.
604 C Read in coupling data that I import.
605 myItP1 = myIter + 1
606 CALL TIMER_START('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid)
607 CALL CPL_EXPORT_MY_DATA( myTime, myItP1, myThid )
608 CALL CPL_IMPORT_EXTERNAL_DATA( myTime, myItP1, myThid )
609 CALL TIMER_STOP ('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid)
610 # ifdef ALLOW_OCN_COMPON_INTERF
611 IF ( useRealFreshWaterFlux ) THEN
612 CALL OCN_APPLY_IMPORT( .FALSE., myTime, myIter, myThid )
613 ENDIF
614 # endif /* ALLOW_OCN_COMPON_INTERF */
615 ENDIF
616 #endif /* COMPONENT_MODULE */
617
618 #ifdef ALLOW_AUTODIFF_TAMC
619 CADJ STORE etaN = comlev1, key = ikey_dynamics, kind = isbyte
620 # if (defined ALLOW_DEPTH_CONTROL)
621 CADJ STORE hFacC, hFacS, hFacW
622 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
623 CADJ STORE recip_hFacC, recip_hFacS, recip_hFacW
624 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
625 c
626 # ifndef DISABLE_RSTAR_CODE
627 CADJ STORE rstarFacC, rstarFacS, rstarFacW =
628 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
629 c
630 CADJ STORE h0facc,h0facs,h0facw
631 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
632 CADJ STORE rstardhcdt,rstardhsdt,rstardhwdt
633 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
634 CADJ STORE rstarexpc,rstarexps,rstarexpw
635 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
636 # endif
637 # endif
638 #endif
639
640 #ifndef ALLOW_OFFLINE
641 C-- Step forward fields and calculate time tendency terms.
642 #ifndef ALLOW_AUTODIFF_TAMC
643 IF ( momStepping ) THEN
644 #endif
645 #ifdef ALLOW_DEBUG
646 IF (debugMode) CALL DEBUG_CALL('DYNAMICS',myThid)
647 #endif
648 CALL TIMER_START('DYNAMICS [FORWARD_STEP]',myThid)
649 CALL DYNAMICS( myTime, myIter, myThid )
650 CALL TIMER_STOP ('DYNAMICS [FORWARD_STEP]',myThid)
651 #ifndef ALLOW_AUTODIFF_TAMC
652 ENDIF
653 #endif
654 #endif /* ndfef ALLOW_OFFLINE */
655
656 #ifdef ALLOW_AUTODIFF_TAMC
657 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
658 CADJ STORE gU, gV = comlev1, key = ikey_dynamics,
659 CADJ & kind = isbyte
660 # endif
661 #endif
662
663 C-- Update time-counter
664
665 myIter = nIter0 + iLoop
666 myTime = startTime + deltaTClock * float(iLoop)
667
668 #ifdef ALLOW_MNC
669 C Update MNC time information
670 IF ( useMNC ) THEN
671 CALL MNC_UPDATE_TIME( myTime, myIter, myThid )
672 ENDIF
673 #endif /* ALLOW_MNC */
674
675 C-- Update geometric factors:
676 #ifdef NONLIN_FRSURF
677 C- update hfacC,W,S and recip_hFac according to etaH(n+1) :
678 IF ( select_rStar.GT.0 ) THEN
679 # ifndef DISABLE_RSTAR_CODE
680 # ifdef ALLOW_AUTODIFF_TAMC
681 CADJ STORE rstarFacC, rstarFacS, rstarFacW =
682 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
683 # endif
684 CALL TIMER_START('UPDATE_R_STAR [FORWARD_STEP]',myThid)
685 CALL UPDATE_R_STAR( .TRUE., myTime, myIter, myThid )
686 CALL TIMER_STOP ('UPDATE_R_STAR [FORWARD_STEP]',myThid)
687 # endif /* DISABLE_RSTAR_CODE */
688 ELSEIF ( selectSigmaCoord.NE.0 ) THEN
689 # ifndef DISABLE_SIGMA_CODE
690 CALL UPDATE_SIGMA( etaH, myTime, myIter, myThid )
691 # endif /* DISABLE_RSTAR_CODE */
692 ELSE
693 #ifdef ALLOW_AUTODIFF_TAMC
694 CADJ STORE hFac_surfC, hFac_surfS, hFac_surfW
695 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
696 #endif
697 CALL TIMER_START('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
698 CALL UPDATE_SURF_DR( .TRUE., myTime, myIter, myThid )
699 CALL TIMER_STOP ('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
700 ENDIF
701 # ifdef ALLOW_AUTODIFF_TAMC
702 CADJ STORE hFacC, hFacS, hFacW =
703 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
704 CADJ STORE recip_hFacC, recip_hFacS, recip_hFacW =
705 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
706 # endif
707 C- update also CG2D matrix (and preconditioner)
708 IF ( momStepping .AND. nonlinFreeSurf.GT.2 ) THEN
709 CALL TIMER_START('UPDATE_CG2D [FORWARD_STEP]',myThid)
710 CALL UPDATE_CG2D( myTime, myIter, myThid )
711 CALL TIMER_STOP ('UPDATE_CG2D [FORWARD_STEP]',myThid)
712 ENDIF
713 #endif /* NONLIN_FRSURF */
714
715 C-- Apply Filters to u*,v* before SOLVE_FOR_PRESSURE
716 #ifdef ALLOW_SHAP_FILT
717 IF (useSHAP_FILT .AND. shap_filt_uvStar) THEN
718 CALL TIMER_START('SHAP_FILT_UV [FORWARD_STEP]',myThid)
719 IF (implicDiv2Dflow.LT.1.) THEN
720 C-- Explicit+Implicit part of the Barotropic Flow Divergence
721 C => Filtering of uVel,vVel is necessary
722 CALL SHAP_FILT_APPLY_UV( uVel,vVel,
723 & myTime, myIter, myThid )
724 ENDIF
725 CALL SHAP_FILT_APPLY_UV( gU,gV,myTime,myIter,myThid)
726 CALL TIMER_STOP ('SHAP_FILT_UV [FORWARD_STEP]',myThid)
727 ENDIF
728 #endif
729 #ifdef ALLOW_ZONAL_FILT
730 IF (useZONAL_FILT .AND. zonal_filt_uvStar) THEN
731 CALL TIMER_START('ZONAL_FILT_UV [FORWARD_STEP]',myThid)
732 IF (implicDiv2Dflow.LT.1.) THEN
733 C-- Explicit+Implicit part of the Barotropic Flow Divergence
734 C => Filtering of uVel,vVel is necessary
735 CALL ZONAL_FILT_APPLY_UV( uVel, vVel, myThid )
736 ENDIF
737 CALL ZONAL_FILT_APPLY_UV( gU, gV, myThid )
738 CALL TIMER_STOP ('ZONAL_FILT_UV [FORWARD_STEP]',myThid)
739 ENDIF
740 #endif
741
742 #ifndef ALLOW_OFFLINE
743
744 C-- Solve elliptic equation(s).
745 C Two-dimensional only for conventional hydrostatic or
746 C three-dimensional for non-hydrostatic and/or IGW scheme.
747 IF ( momStepping ) THEN
748 #ifdef ALLOW_AUTODIFF_TAMC
749 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
750 CADJ STORE uvel, vvel
751 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
752 CADJ STORE empmr,hfacs,hfacw
753 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
754 # endif
755 #endif
756 CALL TIMER_START('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
757 CALL SOLVE_FOR_PRESSURE(myTime, myIter, myThid)
758 CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
759 ENDIF
760
761 C-- Correct divergence in flow field and cycle time-stepping momentum
762 #ifndef ALLOW_AUTODIFF_TAMC
763 IF ( momStepping ) THEN
764 #endif
765 #ifdef ALLOW_AUTODIFF_TAMC
766 # if (defined ALLOW_DEPTH_CONTROL)
767 # ifndef DISABLE_RSTAR_CODE
768 cph-test
769 cph not clear, why this one
770 CADJ STORE h0facc = comlev1, key = ikey_dynamics,
771 CADJ & kind = isbyte
772 # endif
773 CADJ STORE etan, uvel,vvel
774 CADJ & = comlev1, key = ikey_dynamics
775 # endif
776 #endif
777 CALL TIMER_START('MOM_CORRECTION_STEP [FORWARD_STEP]',myThid)
778 CALL MOMENTUM_CORRECTION_STEP(myTime, myIter, myThid)
779 CALL TIMER_STOP ('MOM_CORRECTION_STEP [FORWARD_STEP]',myThid)
780 #ifndef ALLOW_AUTODIFF_TAMC
781 ENDIF
782 #endif
783 #endif /* ndfef ALLOW_OFFLINE */
784
785 #ifdef EXACT_CONSERV
786 IF (exactConserv) THEN
787 C-- Update etaH(n+1) :
788 CALL TIMER_START('UPDATE_ETAH [FORWARD_STEP]',myThid)
789 CALL UPDATE_ETAH( myTime, myIter, myThid )
790 CALL TIMER_STOP ('UPDATE_ETAH [FORWARD_STEP]',myThid)
791 ENDIF
792 #endif /* EXACT_CONSERV */
793
794 #ifdef NONLIN_FRSURF
795 IF ( select_rStar.NE.0 ) THEN
796 # ifndef DISABLE_RSTAR_CODE
797 # ifdef ALLOW_AUTODIFF_TAMC
798 CADJ STORE rstarfacc,rstarfacs,rstarfacw =
799 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
800 # endif
801 C-- r* : compute the future level thickness according to etaH(n+1)
802 CALL TIMER_START('CALC_R_STAR [FORWARD_STEP]',myThid)
803 CALL CALC_R_STAR(etaH, myTime, myIter, myThid )
804 CALL TIMER_STOP ('CALC_R_STAR [FORWARD_STEP]',myThid)
805 # endif /* DISABLE_RSTAR_CODE */
806 ELSEIF ( nonlinFreeSurf.GT.0 .AND. selectSigmaCoord.EQ.0 ) THEN
807 C-- compute the future surface level thickness according to etaH(n+1)
808 # ifdef ALLOW_AUTODIFF_TAMC
809 CADJ STORE etaH = comlev1, key = ikey_dynamics,
810 CADJ & kind = isbyte
811 # endif
812 CALL TIMER_START('CALC_SURF_DR [FORWARD_STEP]',myThid)
813 CALL CALC_SURF_DR(etaH, myTime, myIter, myThid )
814 CALL TIMER_STOP ('CALC_SURF_DR [FORWARD_STEP]',myThid)
815 ENDIF
816 # ifdef ALLOW_AUTODIFF_TAMC
817 CADJ STORE hFac_surfC = comlev1, key = ikey_dynamics,
818 CADJ & kind = isbyte
819 CADJ STORE salt,theta,vvel = comlev1, key = ikey_dynamics,
820 CADJ & kind = isbyte
821 # endif
822 #endif /* NONLIN_FRSURF */
823
824 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
825 IF ( staggerTimeStep ) THEN
826 C-- do exchanges of U,V (needed for multiDim) when using stagger time-step :
827 #ifdef ALLOW_DEBUG
828 IF (debugMode)
829 & CALL DEBUG_CALL('DO_STAGGER_FIELDS_EXCH.',myThid)
830 #endif
831 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
832 CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid )
833 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
834
835 #ifdef ALLOW_DIAGNOSTICS
836 C-- State-variables diagnostics
837 IF ( useDiagnostics ) THEN
838 CALL TIMER_START('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
839 CALL DO_STATEVARS_DIAGS( myTime, 1, myIter, myThid )
840 CALL TIMER_STOP ('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
841 ENDIF
842 #endif
843
844 #ifdef ALLOW_DEBUG
845 IF (debugMode) CALL DEBUG_CALL('THERMODYNAMICS',myThid)
846 #endif
847 CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',myThid)
848 CALL THERMODYNAMICS( myTime, myIter, myThid )
849 CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',myThid)
850
851 C-- if staggerTimeStep: end
852 ENDIF
853 C---+--------+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
854
855 #ifdef ALLOW_AUTODIFF_TAMC
856 cph This is needed because convective_adjustment calls
857 cph find_rho which may use pressure()
858 CADJ STORE totphihyd = comlev1, key = ikey_dynamics,
859 CADJ & kind = isbyte
860 #endif
861 C-- Cycle time-stepping Tracers arrays (T,S,+pTracers)
862 CALL TIMER_START('TRC_CORRECTION_STEP [FORWARD_STEP]',myThid)
863 CALL TRACERS_CORRECTION_STEP(myTime, myIter, myThid)
864 CALL TIMER_STOP ('TRC_CORRECTION_STEP [FORWARD_STEP]',myThid)
865
866 #ifdef ALLOW_LONGSTEP
867 IF ( usePTRACERS ) THEN
868 IF ( LS_whenToSample .EQ. 2 ) THEN
869 C Average everything at the end of the timestep. This will
870 C reproduce online results with staggerTimeStep=.TRUE.
871 C when LS_nIter=1
872 #ifdef ALLOW_DEBUG
873 IF (debugMode) CALL DEBUG_CALL('LONGSTEP_AVERAGE',myThid)
874 #endif
875 CALL TIMER_START('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid)
876 C myIter has been update after dynamics, but the averaging window
877 C should be determined by myIter at beginning of timestep
878 CALL LONGSTEP_AVERAGE( myTimeBeg, myIterBeg, myThid )
879 CALL TIMER_STOP ('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid)
880
881 #ifdef ALLOW_DEBUG
882 IF (debugMode)
883 & CALL DEBUG_CALL('LONGSTEP_THERMODYNAMICS',myThid)
884 #endif
885 CALL TIMER_START('LONGSTEP_THERMODYNAMICS [FORWARD_STEP]',
886 & myThid)
887 CALL LONGSTEP_THERMODYNAMICS( myTime, myIter, myThid )
888 CALL TIMER_STOP ('LONGSTEP_THERMODYNAMICS [FORWARD_STEP]',
889 & myThid)
890 C-- if LS_whenToSample.EQ.2: end
891 ENDIF
892
893 C-- Cycle time-stepping Tracers arrays (pTracers)
894 CALL TIMER_START('LS_CORRECTION_STEP [FORWARD_STEP]',myThid)
895 CALL LONGSTEP_CORRECTION_STEP(myTime, myIter, myThid)
896 CALL TIMER_STOP ('LS_CORRECTION_STEP [FORWARD_STEP]',myThid)
897 C-- if usePTRACERS: end
898 ENDIF
899 #endif /* ALLOW_LONGSTEP */
900
901 #ifdef ALLOW_GCHEM
902 C Add separate timestepping of chemical/biological/forcing
903 C of ptracers here in GCHEM_FORCING_SEP
904 #ifdef ALLOW_AUTODIFF_TAMC
905 CADJ STORE ptracer = comlev1, key = ikey_dynamics,
906 CADJ & kind = isbyte
907 CADJ STORE theta = comlev1, key = ikey_dynamics,
908 CADJ & kind = isbyte
909 CADJ STORE salt = comlev1, key = ikey_dynamics,
910 CADJ & kind = isbyte
911 #endif
912
913 #ifdef ALLOW_LONGSTEP
914 IF ( LS_doTimeStep ) THEN
915 #else
916 IF ( .TRUE. ) THEN
917 #endif
918 IF ( useGCHEM ) THEN
919 #ifdef ALLOW_DEBUG
920 IF (debugMode) CALL DEBUG_CALL('GCHEM_FORCING_SEP',myThid)
921 #endif /* ALLOW_DEBUG */
922 CALL TIMER_START('GCHEM_FORCING_SEP [FORWARD_STEP]',myThid)
923 CALL GCHEM_FORCING_SEP( myTime,myIter,myThid )
924 CALL TIMER_STOP ('GCHEM_FORCING_SEP [FORWARD_STEP]',myThid)
925 ENDIF
926 C endif LS_doTimeStep
927 ENDIF
928 #endif /* ALLOW_GCHEM */
929
930 C-- Do "blocking" sends and receives for tendency "overlap" terms
931 c CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
932 c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )
933 c CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
934
935 C-- Do "blocking" sends and receives for field "overlap" terms
936 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
937 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
938 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
939
940 #ifdef ALLOW_DIAGNOSTICS
941 IF ( useDiagnostics ) THEN
942 CALL TIMER_START('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
943 CALL DO_STATEVARS_DIAGS( myTime, 2, myIter, myThid )
944 CALL TIMER_STOP ('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
945 ENDIF
946 #endif
947
948 #ifdef ALLOW_GRIDALT
949 IF (useGRIDALT) THEN
950 CALL GRIDALT_UPDATE(myThid)
951 ENDIF
952 #endif
953
954 #ifdef ALLOW_FIZHI
955 IF (useFIZHI) THEN
956 CALL TIMER_START('FIZHI [FORWARD_STEP]',myThid)
957 CALL STEP_FIZHI_CORR ( myTime, myIter, myThid, dTtracerLev(1) )
958 CALL TIMER_STOP ('FIZHI [FORWARD_STEP]',myThid)
959 ENDIF
960 #endif
961
962 #ifdef ALLOW_FLT
963 C-- Calculate float trajectories
964 IF (useFLT) THEN
965 CALL TIMER_START('FLOATS [FORWARD_STEP]',myThid)
966 CALL FLT_MAIN( myTime, myIter, myThid )
967 CALL TIMER_STOP ('FLOATS [FORWARD_STEP]',myThid)
968 ENDIF
969 #endif
970
971 #ifdef ALLOW_TIMEAVE
972 C-- State-variables time-averaging
973 CALL TIMER_START('DO_STATEVARS_TAVE [FORWARD_STEP]',myThid)
974 CALL DO_STATEVARS_TAVE( myTime, myIter, myThid )
975 CALL TIMER_STOP ('DO_STATEVARS_TAVE [FORWARD_STEP]',myThid)
976 #endif
977
978 #ifdef ALLOW_NEST_PARENT
979 IF ( useNEST_PARENT) THEN
980 CALL NEST_PARENT_IO_2( myTime, myIter, myThid )
981 ENDIF
982 #endif /* ALLOW_NEST_PARENT */
983
984 #ifdef ALLOW_NEST_CHILD
985 IF ( useNEST_CHILD) THEN
986 CALL NEST_CHILD_TRANSP( myTime, myIter, myThid )
987 ENDIF
988 #endif /* ALLOW_NEST_CHILD */
989
990 #ifdef ALLOW_MONITOR
991 IF ( .NOT.useOffLine ) THEN
992 C-- Check status of solution (statistics, cfl, etc...)
993 CALL TIMER_START('MONITOR [FORWARD_STEP]',myThid)
994 CALL MONITOR( myTime, myIter, myThid )
995 CALL TIMER_STOP ('MONITOR [FORWARD_STEP]',myThid)
996 ENDIF
997 #endif /* ALLOW_MONITOR */
998
999 #ifdef ALLOW_COST
1000 C-- compare model with data and compute cost function
1001 C-- this is done after exchanges to allow interpolation
1002 CALL TIMER_START('COST_TILE [FORWARD_STEP]',myThid)
1003 CALL COST_TILE ( myTime, myIter, myThid )
1004 CALL TIMER_STOP ('COST_TILE [FORWARD_STEP]',myThid)
1005 #endif
1006
1007 C-- Check if it has reached the end of simulation
1008 modelEnd = myTime.EQ.endTime .OR. myIter.EQ.nEndIter
1009 #ifdef HAVE_SIGREG
1010 IF ( useSIGREG ) THEN
1011 modelEnd = modelEnd .OR. ( i_got_signal.GT.0 )
1012 ENDIF
1013 #endif /* HAVE_SIGREG */
1014
1015 C-- Do IO if needed.
1016 CALL TIMER_START('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
1017 CALL DO_THE_MODEL_IO( modelEnd, myTime, myIter, myThid )
1018 CALL TIMER_STOP ('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
1019
1020 C-- Save state for restarts
1021 CALL TIMER_START('DO_WRITE_PICKUP [FORWARD_STEP]',myThid)
1022 CALL DO_WRITE_PICKUP( modelEnd, myTime, myIter, myThid )
1023 CALL TIMER_STOP ('DO_WRITE_PICKUP [FORWARD_STEP]',myThid)
1024
1025 #ifdef HAVE_SIGREG
1026 IF ( useSIGREG ) THEN
1027 IF ( modelEnd .AND. i_got_signal.GT.0 ) THEN
1028 STOP 'Checkpoint completed -- killed by signal handler'
1029 ENDIF
1030 ENDIF
1031 #endif /* HAVE_SIGREG */
1032
1033 #ifdef ALLOW_AUTODIFF_TAMC
1034 CALL AUTODIFF_INADMODE_SET( myThid )
1035 #endif
1036
1037 #ifdef ALLOW_SHOWFLOPS
1038 CALL TIMER_START('SHOWFLOPS_INLOOP [THE_MAIN_LOOP]', mythid)
1039 CALL SHOWFLOPS_INLOOP( iloop, mythid )
1040 CALL TIMER_STOP ('SHOWFLOPS_INLOOP [THE_MAIN_LOOP]', mythid)
1041 #endif
1042
1043 #ifdef ALLOW_DEBUG
1044 IF (debugMode) CALL DEBUG_LEAVE('FORWARD_STEP',myThid)
1045 #endif
1046
1047 RETURN
1048 END

  ViewVC Help
Powered by ViewVC 1.1.22