/[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.181 - (show annotations) (download)
Sat Sep 25 23:09:54 2010 UTC (13 years, 9 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint62l
Changes since 1.180: +12 -1 lines
add hooks for new packages OASIS, the package will follow
(hooks exclude the seaice pkg for now)

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

  ViewVC Help
Powered by ViewVC 1.1.22