/[MITgcm]/MITgcm_contrib/shelfice_remeshing/TEST_TAF/code/forward_step.F
ViewVC logotype

Annotation of /MITgcm_contrib/shelfice_remeshing/TEST_TAF/code/forward_step.F

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


Revision 1.2 - (hide annotations) (download)
Thu Dec 10 21:31:28 2015 UTC (9 years, 7 months ago) by dgoldberg
Branch: MAIN
CVS Tags: HEAD
Changes since 1.1: +12 -1 lines
test taf

1 dgoldberg 1.2 C $Header: /u/gcmpack/MITgcm_contrib/shelfice_remeshing/TEST_TAF/code/forward_step.F,v 1.1 2015/12/10 21:22:54 dgoldberg Exp $
2 dgoldberg 1.1 C $Name: $
3    
4     #include "PACKAGES_CONFIG.h"
5     #include "CPP_OPTIONS.h"
6    
7     #ifdef ALLOW_AUTODIFF
8     # include "AUTODIFF_OPTIONS.h"
9     #endif
10     #ifdef ALLOW_GENERIC_ADVDIFF
11     # include "GAD_OPTIONS.h"
12     #endif
13     #ifdef ALLOW_GGL90
14     # include "GGL90_OPTIONS.h"
15     #endif
16     #ifdef ALLOW_GMREDI
17     # include "GMREDI_OPTIONS.h"
18     #endif
19     #ifdef ALLOW_OBCS
20     # include "OBCS_OPTIONS.h"
21     #endif
22     #ifdef ALLOW_THSICE
23     # include "THSICE_OPTIONS.h"
24     #endif
25     #ifdef ALLOW_SEAICE
26     # include "SEAICE_OPTIONS.h"
27     #endif
28     #ifdef ALLOW_PTRACERS
29     # include "PTRACERS_OPTIONS.h"
30     #endif
31     #ifdef ALLOW_EXF
32     # include "EXF_OPTIONS.h"
33     #endif
34     #ifdef ALLOW_STREAMICE
35     # include "STREAMICE_OPTIONS.h"
36     #endif
37     #ifdef ALLOW_COST
38     # include "COST_OPTIONS.h"
39     #endif
40     #ifdef ALLOW_CTRL
41     # include "CTRL_OPTIONS.h"
42     #endif
43     #ifdef ALLOW_ECCO
44     # include "ECCO_OPTIONS.h"
45     #endif
46    
47     #define ALLOW_MOM_STEPPING
48     #if ( defined (ALLOW_AUTODIFF) && defined (ALLOW_OFFLINE) )
49     # undef ALLOW_MOM_STEPPING
50     #endif
51    
52     CBOP
53     C !ROUTINE: FORWARD_STEP
54     C !INTERFACE:
55     SUBROUTINE FORWARD_STEP( iloop, myTime, myIter, myThid )
56    
57     C !DESCRIPTION: \bv
58     C *=================================================================
59     C | SUBROUTINE forward_step
60     C | o Step forward in time the model variables for one time-step
61     C *=================================================================
62     C | The algorithm...
63     C |
64     C | "Calculation of Gs"
65     C | ===================
66     C | This is where all the accelerations and tendencies (ie.
67     C | physics, parameterizations etc...) are calculated
68     C | rho = rho ( theta[n], salt[n] )
69     C | b = b(rho, theta)
70     C | K31 = K31 ( rho )
71     C | Gu[n] = Gu( u[n], v[n], wVel, b, ... )
72     C | Gv[n] = Gv( u[n], v[n], wVel, b, ... )
73     C | Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
74     C | Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
75     C |
76     C | "Time-stepping" or "Prediction"
77     C | ================================
78     C | The models variables are stepped forward with the appropriate
79     C | time-stepping scheme (currently we use Adams-Bashforth II)
80     C | - For momentum, the result is always *only* a "prediction"
81     C | in that the flow may be divergent and will be "corrected"
82     C | later with a surface pressure gradient.
83     C | - Normally for tracers the result is the new field at time
84     C | level [n+1} *BUT* in the case of implicit diffusion the result
85     C | is also *only* a prediction.
86     C | - We denote "predictors" with an asterisk (*).
87     C | U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
88     C | V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
89     C | theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
90     C | salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
91     C | With implicit diffusion:
92     C | theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
93     C | salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
94     C | (1 + dt * K * d_zz) theta[n+1] = theta*
95     C | (1 + dt * K * d_zz) salt[n+1] = salt*
96     C |
97     C | "Correction Step"
98     C | =================
99     C | Here we update the horizontal velocities with the surface
100     C | pressure such that the resulting flow is either consistent
101     C | with the free-surface evolution or the rigid-lid:
102     C | U[n] = U* + dt x d/dx P
103     C | V[n] = V* + dt x d/dy P
104     C | W[n] = W* + dt x d/dz P (NH mode)
105     C *=================================================================
106     C \ev
107    
108     C !CALLING SEQUENCE:
109     C FORWARD_STEP
110     C |
111     C |-- SHELFICE_REMESHING
112     C |
113     C |-- AUTODIFF_INADMODE_UNSET
114     C |
115     C |-- RESET_NLFS_VARS
116     C |-- UPDATE_R_STAR
117     C |-- UPDATE_SURF_DR
118     C |
119     C |-- PTRACERS_SWITCH_ONOFF
120     C |
121     C |-- DIAGNOSTICS_SWITCH_ONOFF
122     C |-- DO_STATEVARS_DIAGS
123     C |
124     C |-- NEST_CHILD_SETMEMO
125     C |-- NEST_PARENT_IO_1
126     C |
127     C |-- LOAD_FIELDS_DRIVER
128     C |
129     C |-- BULKF_FORCING
130     C |
131     C |-- CHEAPAML
132     C |
133     C |-- CTRL_MAP_FORCING
134     C |-- DUMMY_IN_STEPPING
135     C |
136     C |-- CPL_EXPORT_IMPORT_DATA
137     C |
138     C |-- OASIS_PUT
139     C |-- OASIS_GET
140     C |
141     C |-- EBM_DRIVER
142     C |
143     C |-- DO_ATMOSPHERIC_PHYS
144     C |
145     C |-- DO_OCEANIC_PHYS
146     C |
147     C |-- STREAMICE_TIMESTEP
148     C |
149     C |-- GCHEM_CALC_TENDENCY
150     C |
151     C |-- LONGSTEP_AVERAGE
152     C |-- LONGSTEP_THERMODYNAMICS
153     C |
154     C |-- THERMODYNAMICS
155     C |
156     C |-- LONGSTEP_AVERAGE
157     C |-- LONGSTEP_THERMODYNAMICS
158     C |
159     C |-- DO_STAGGER_FIELDS_EXCHANGES
160     C |
161     C |-- DYNAMICS
162     C |
163     C |-- MNC_UPDATE_TIME
164     C |
165     C |-- OFFLINE_FIELDS_LOAD
166     C |
167     C |-- UPDATE_R_STAR
168     C |-- UPDATE_SIGMA
169     C |-- UPDATE_SURF_DR
170     C |-- UPDATE_CG2D
171     C |
172     C |-- SHAP_FILT_APPLY_UV
173     C |-- ZONAL_FILT_APPLY_UV
174     C |
175     C |-- SOLVE_FOR_PRESSURE
176     C |
177     C |-- MOMENTUM_CORRECTION_STEP
178     C |
179     C |-- INTEGR_CONTINUITY
180     C |
181     C |-- CALC_R_STAR
182     C |-- CALC_SURF_DR
183     C |
184     C |-- DO_STAGGER_FIELDS_EXCHANGES
185     C |
186     C |-- DO_STATEVARS_DIAGS
187     C |
188     C |-- THERMODYNAMICS
189     C |
190     C |-- TRACERS_CORRECTION_STEP
191     C |
192     C |-- LONGSTEP_AVERAGE
193     C |-- LONGSTEP_THERMODYNAMICS
194     C |
195     C |-- GCHEM_FORCING_SEP
196     C |
197     C |-- DO_FIELDS_BLOCKING_EXCHANGES
198     C |
199     C |-- DO_STATEVARS_DIAGS
200     C |
201     C |-- GRIDALT_UPDATE
202     C |-- STEP_FIZHI_CORR
203     C |
204     C |-- FLT_MAIN
205     C |
206     C |-- DO_STATEVARS_TAVE
207     C |
208     C |-- NEST_PARENT_IO_2
209     C |-- NEST_CHILD_TRANSP
210     C |
211     C |-- MONITOR
212     C |
213     C |-- COST_TILE
214     C |
215     C |-- DO_THE_MODEL_IO
216     C |
217     C |-- PTRACERS_RESET
218     C |
219     C |-- DO_WRITE_PICKUP
220     C |
221     C |-- AUTODIFF_INADMODE_SET
222     C |
223     C |-- SHOWFLOPS_INLOOP
224    
225     C !USES:
226     IMPLICIT NONE
227     C == Global variables ==
228     #include "SIZE.h"
229     #include "EEPARAMS.h"
230     #include "PARAMS.h"
231     #include "DYNVARS.h"
232    
233     #ifdef HAVE_SIGREG
234     #include "SIGREG.h"
235     #endif
236    
237     #ifdef ALLOW_SHAP_FILT
238     # include "SHAP_FILT.h"
239     #endif
240     #ifdef ALLOW_ZONAL_FILT
241     # include "ZONAL_FILT.h"
242     #endif
243    
244     #ifdef ALLOW_LONGSTEP
245     # include "LONGSTEP_PARAMS.h"
246     # include "LONGSTEP.h"
247     #endif
248    
249     #ifdef ALLOW_AUTODIFF
250     # include "AUTODIFF_MYFIELDS.h"
251     # include "FFIELDS.h"
252     # include "SURFACE.h"
253    
254     # include "tamc.h"
255     # ifdef ALLOW_CTRL
256     # include "CTRL_SIZE.h"
257     # include "ctrl.h"
258     # include "ctrl_dummy.h"
259     # include "CTRL_GENARR.h"
260     # include "CTRL_OBCS.h"
261     # endif
262     # ifdef ALLOW_COST
263     # include "cost.h"
264     # endif
265     # ifdef ALLOW_ECCO
266     # include "ecco_cost.h"
267     # endif
268     # include "EOS.h"
269     # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
270     # include "GRID.h"
271     # endif
272     # ifdef ALLOW_GMREDI
273     # include "GMREDI.h"
274     # endif
275     # ifdef ALLOW_EXF
276     # include "EXF_FIELDS.h"
277     # include "EXF_PARAM.h"
278     # ifdef ALLOW_BULKFORMULAE
279     # include "EXF_CONSTANTS.h"
280     # endif
281     # endif
282     # ifdef ALLOW_CD_CODE
283     # include "CD_CODE_VARS.h"
284     # endif
285     # ifdef ALLOW_GENERIC_ADVDIFF
286     # include "GAD.h"
287     # include "GAD_SOM_VARS.h"
288     # endif
289     # ifdef ALLOW_GGL90
290     # include "GGL90.h"
291     # endif
292     # ifdef ALLOW_PTRACERS
293     # include "PTRACERS_SIZE.h"
294     # include "PTRACERS_FIELDS.h"
295     # endif
296     # ifdef ALLOW_GCHEM
297     # include "GCHEM_FIELDS.h"
298     # endif
299     # ifdef ALLOW_CFC
300     # include "CFC.h"
301     # endif
302     # ifdef ALLOW_DIC
303     # include "DIC_VARS.h"
304     # include "DIC_LOAD.h"
305     # include "DIC_ATMOS.h"
306     # include "DIC_COST.h"
307     # endif
308     # ifdef ALLOW_OBCS
309     # include "OBCS_PARAMS.h"
310     # include "OBCS_FIELDS.h"
311     # include "OBCS_SEAICE.h"
312     # ifdef ALLOW_PTRACERS
313     # include "OBCS_PTRACERS.h"
314     # endif
315     # endif
316     # ifdef ALLOW_THSICE
317     # include "THSICE_PARAMS.h"
318     # include "THSICE_SIZE.h"
319     # include "THSICE_VARS.h"
320     # include "THSICE_COST.h"
321     # endif
322     # ifdef ALLOW_SEAICE
323     # include "SEAICE_SIZE.h"
324     # include "SEAICE.h"
325     # include "SEAICE_COST.h"
326     # endif
327     # ifdef ALLOW_SALT_PLUME
328     # include "SALT_PLUME.h"
329     # endif
330     # ifdef ALLOW_SHELFICE
331     # include "SHELFICE.h"
332     # include "SHELFICE_COST.h"
333     # endif
334     # ifdef ALLOW_STREAMICE
335     # include "STREAMICE.h"
336     # include "STREAMICE_ADV.h"
337     # include "STREAMICE_BDRY.h"
338     # include "STREAMICE_CG.h"
339     # endif
340     # ifdef ALLOW_EBM
341     # include "EBM.h"
342     # endif
343     # ifdef ALLOW_KPP
344     # include "KPP.h"
345     # endif
346     # ifdef ALLOW_RBCS
347     # include "RBCS_SIZE.h"
348     # include "RBCS_FIELDS.h"
349     # endif
350     # ifdef ALLOW_OFFLINE
351     # include "OFFLINE.h"
352     # endif
353     # ifdef ALLOW_CG2D_NSA
354     # include "CG2D.h"
355     # endif
356     #endif /* ALLOW_AUTODIFF */
357    
358     C !INPUT/OUTPUT PARAMETERS:
359     C == Routine arguments ==
360     C note: under the multi-threaded model myIter and
361     C myTime are local variables passed around as routine
362     C arguments. Although this is fiddly it saves the need to
363     C impose additional synchronisation points when they are
364     C updated.
365     C myTime :: time counter for this thread
366     C myIter :: iteration counter for this thread
367     C myThid :: thread number for this instance of the routine.
368     INTEGER iloop
369     _RL myTime
370     INTEGER myIter
371     INTEGER myThid
372    
373     C !LOCAL VARIABLES:
374     C == Local variables ==
375     C modelEnd :: true if reaching the end of the run
376     C myTimeBeg :: time at beginning of time step (needed by longstep)
377     C myIterBeg :: iteration number at beginning of time step
378     LOGICAL modelEnd
379     #ifdef ALLOW_LONGSTEP
380     INTEGER myIterBeg
381     _RL myTimeBeg
382     #endif /* ALLOW_LONGSTEP */
383     CEOP
384    
385     #ifdef ALLOW_SHELFICE
386     C-- Remesh shelfice
387     CALL TIMER_START('SHELFICE_REMESHING [FORWARD_STEP]',myThid)
388     CALL SHELFICE_REMESHING (myTime, myIter, myThid )
389     CALL TIMER_STOP('SHELFICE_REMESHING [FORWARD_STEP]',myThid)
390     #endif /* ALLOW_SHELFICE */
391    
392     #ifdef ALLOW_DEBUG
393     IF (debugMode) CALL DEBUG_ENTER('FORWARD_STEP',myThid)
394     #endif
395    
396     #ifdef ALLOW_AUTODIFF
397     CALL AUTODIFF_INADMODE_UNSET( myThid )
398     #endif
399    
400     #ifdef ALLOW_AUTODIFF
401     C-- Reset the model iteration counter and the model time.
402     myIter = nIter0 + (iloop-1)
403     myTime = startTime + deltaTClock*(iLoop-1)
404     #endif
405    
406     #ifdef ALLOW_LONGSTEP
407     C store this for longstep_average with staggerTimeStep
408     C which is called after myIter and myTime are incremented
409     C but needs iter/time at beginning of time step
410     myIterBeg = myIter
411     myTimeBeg = myTime
412     #endif /* ALLOW_LONGSTEP */
413    
414     #ifdef ALLOW_AUTODIFF_TAMC
415     c**************************************
416     #include "checkpoint_lev1_directives.h"
417     #include "checkpoint_lev1_template.h"
418     c**************************************
419     #endif
420    
421     C-- Reset geometric factors (hFacC,W,S & recip_hFac) to their current values:
422     C added to simplify adjoint derivation - no effect in forward run
423     #ifdef NONLIN_FRSURF
424     #ifndef ALLOW_AUTODIFF
425     IF ( doResetHFactors ) THEN
426     #endif
427     CALL RESET_NLFS_VARS( myTime, myIter, myThid )
428     IF ( select_rStar.GT.0 ) THEN
429     # ifndef DISABLE_RSTAR_CODE
430     # ifdef ALLOW_AUTODIFF_TAMC
431     CADJ STORE rStarFacC, rStarFacS, rStarFacW =
432     CADJ & comlev1, key = ikey_dynamics, kind = isbyte
433     # endif
434     CALL TIMER_START('UPDATE_R_STAR [FORWARD_STEP]',myThid)
435     CALL UPDATE_R_STAR( .FALSE., myTime, myIter, myThid )
436     CALL TIMER_STOP ('UPDATE_R_STAR [FORWARD_STEP]',myThid)
437     # endif /* DISABLE_RSTAR_CODE */
438     ELSE
439     #ifdef ALLOW_AUTODIFF_TAMC
440     CADJ STORE hFac_surfC, hFac_surfS, hFac_surfW
441     CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
442     #endif
443     CALL TIMER_START('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
444     CALL UPDATE_SURF_DR( .FALSE., myTime, myIter, myThid )
445     CALL TIMER_STOP ('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
446     ENDIF
447     #ifdef ALLOW_AUTODIFF_TAMC
448     CADJ STORE hFacC, hFacS, hFacW =
449     CADJ & comlev1, key = ikey_dynamics, kind = isbyte
450     CADJ STORE recip_hFacC, recip_hFacS, recip_hFacW =
451     CADJ & comlev1, key = ikey_dynamics, kind = isbyte
452     #endif
453     #ifndef ALLOW_AUTODIFF
454     ENDIF
455     #endif
456     #endif /* NONLIN_FRSURF */
457    
458     #ifdef ALLOW_PTRACERS
459     C-- Switch on/off individual tracer time-stepping
460     IF ( usePTRACERS ) THEN
461     CALL PTRACERS_SWITCH_ONOFF( myTime, myIter, myThid )
462     ENDIF
463     #endif /* ALLOW_PTRACERS */
464    
465     C-- Switch on/off diagnostics for snap-shot output:
466     #ifdef ALLOW_DIAGNOSTICS
467     IF ( useDiagnostics ) THEN
468     CALL DIAGNOSTICS_SWITCH_ONOFF( myTime, myIter, myThid )
469     C-- State-variables diagnostics
470     CALL TIMER_START('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
471     CALL DO_STATEVARS_DIAGS( myTime, 0, myIter, myThid )
472     CALL TIMER_STOP ('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
473     ENDIF
474     #endif /* ALLOW_DIAGNOSTICS */
475    
476     #ifdef ALLOW_NEST_CHILD
477     IF ( useNEST_CHILD) THEN
478     CALL NEST_CHILD_SETMEMO( myTime, myIter, myThid )
479     ENDIF
480     #endif /* ALLOW_NEST_CHILD */
481    
482     #ifdef ALLOW_NEST_PARENT
483     IF ( useNEST_PARENT) THEN
484     CALL NEST_PARENT_IO_1( myTime, myIter, myThid )
485     ENDIF
486     #endif /* ALLOW_NEST_PARENT */
487    
488     C-- Call driver to load external forcing fields from file
489     #ifdef ALLOW_DEBUG
490     IF (debugMode) CALL DEBUG_CALL('LOAD_FIELDS_DRIVER',myThid)
491     #endif
492     #ifdef ALLOW_AUTODIFF_TAMC
493     cph Important STORE that avoids hidden recomp. of load_fields_driver
494     CADJ STORE theta = comlev1, key = ikey_dynamics,
495     CADJ & kind = isbyte
496     CADJ STORE uVel, vVel = comlev1, key = ikey_dynamics,
497     CADJ & kind = isbyte
498     #endif
499     CALL TIMER_START('LOAD_FIELDS_DRIVER [FORWARD_STEP]',myThid)
500     CALL LOAD_FIELDS_DRIVER( myTime, myIter, myThid )
501     CALL TIMER_STOP ('LOAD_FIELDS_DRIVER [FORWARD_STEP]',myThid)
502    
503     C-- Call Bulk-Formulae forcing package
504     #ifdef ALLOW_BULK_FORCE
505     IF ( useBulkForce ) THEN
506     #ifdef ALLOW_DEBUG
507     IF (debugMode) CALL DEBUG_CALL('BULKF_FORCING',myThid)
508     #endif
509     CALL TIMER_START('BULKF_FORCING [FORWARD_STEP]',myThid)
510     C- calculate qnet and empmr (and wind stress)
511     CALL BULKF_FORCING( myTime, myIter, myThid )
512     CALL TIMER_STOP ('BULKF_FORCING [FORWARD_STEP]',myThid)
513     ENDIF
514     #endif /* ALLOW_BULK_FORCE */
515    
516     C-- Call external chepaml forcing package
517     #ifdef ALLOW_CHEAPAML
518     IF ( useCheapAML ) THEN
519     #ifdef ALLOW_DEBUG
520     IF (debugMode) CALL DEBUG_CALL('CHEAPAML',myThid)
521     #endif
522     CALL TIMER_START('CHEAPAML [FORWARD_STEP]',myThid)
523     C- calculate qnet (and wind stress)
524     CALL CHEAPAML( myTime, myIter,myThid )
525     CALL TIMER_STOP ('CHEAPAML [FORWARD_STEP]',myThid)
526     ENDIF
527     #endif /*ALLOW_CHEAPAML */
528    
529     #ifdef ALLOW_CTRL
530     C-- Add control vector for forcing and parameter fields
531     IF ( ( useCTRL ).AND.( myIter .EQ. nIter0 ) )
532     & CALL CTRL_MAP_FORCING (myThid)
533     #endif
534    
535     #ifdef ALLOW_AUTODIFF_MONITOR
536     CALL DUMMY_IN_STEPPING( myTime, myIter, myThid )
537     #endif
538    
539     #ifdef COMPONENT_MODULE
540     IF ( useCoupler ) THEN
541     C Post coupling data that I export.
542     C Read in coupling data that I import.
543     CALL TIMER_START('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid)
544     CALL CPL_EXPORT_IMPORT_DATA( myTime, myIter, myThid )
545     CALL TIMER_STOP ('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid)
546     ENDIF
547     #endif /* COMPONENT_MODULE */
548     #ifdef ALLOW_OASIS
549     IF ( useOASIS ) THEN
550     CALL TIMER_START('OASIS_PUT-GET [FORWARD_STEP]',myThid)
551     C Post coupling data that I export.
552     CALL OASIS_PUT( myTime, myIter, myThid )
553     C Read in coupling data that I import.
554     CALL OASIS_GET( myTime, myIter, myThid )
555     CALL TIMER_STOP ('OASIS_PUT-GET [FORWARD_STEP]',myThid)
556     ENDIF
557     #endif /* ALLOW_OASIS */
558    
559     #ifdef ALLOW_EBM
560     IF ( useEBM ) THEN
561     # ifdef ALLOW_DEBUG
562     IF (debugMode) CALL DEBUG_CALL('EBM',myThid)
563     # endif
564     CALL TIMER_START('EBM [FORWARD_STEP]',myThid)
565     CALL EBM_DRIVER ( myTime, myIter, myThid )
566     CALL TIMER_STOP ('EBM [FORWARD_STEP]',myThid)
567     ENDIF
568     #endif /* ALLOW_EBM */
569    
570     C-- Step forward fields and calculate time tendency terms.
571    
572     #ifdef ALLOW_DEBUG
573     IF (debugMode) CALL DEBUG_CALL('DO_ATMOSPHERIC_PHYS',myThid)
574     #endif
575     CALL TIMER_START('DO_ATMOSPHERIC_PHYS [FORWARD_STEP]',myThid)
576     CALL DO_ATMOSPHERIC_PHYS( myTime, myIter, myThid )
577     CALL TIMER_STOP ('DO_ATMOSPHERIC_PHYS [FORWARD_STEP]',myThid)
578    
579     #ifdef ALLOW_AUTODIFF_TAMC
580     CADJ STORE surfaceForcingTice = comlev1, key = ikey_dynamics,
581     CADJ & kind = isbyte
582     # ifdef ALLOW_KPP
583     CADJ STORE uVel, vVel = comlev1, key = ikey_dynamics,
584     CADJ & kind = isbyte
585     # endif /* ALLOW_KPP */
586     # ifdef EXACT_CONSERV
587     CADJ STORE EmPmR = comlev1, key=ikey_dynamics, kind=isbyte
588     CADJ STORE PmEpR = comlev1, key=ikey_dynamics, kind=isbyte
589     # endif
590     # ifdef ALLOW_OBCS
591     CADJ STORE salt = comlev1, key=ikey_dynamics, kind=isbyte
592     CADJ STORE totphihyd = comlev1, key=ikey_dynamics, kind=isbyte
593     # ifdef ALLOW_OBCS_STEVENS
594     CADJ STORE gsNm1 = comlev1, key=ikey_dynamics, kind=isbyte
595     CADJ STORE gtNm1 = comlev1, key=ikey_dynamics, kind=isbyte
596     # endif
597     # endif /* ALLOW_OBCS */
598     # ifdef ALLOW_PTRACERS
599     CADJ STORE pTracer = comlev1, key = ikey_dynamics,
600     CADJ & kind = isbyte
601     # endif /* ALLOW_PTRACERS */
602     # ifdef ALLOW_DEPTH_CONTROL
603     CADJ STORE hFacC = comlev1, key = ikey_dynamics, kind = isbyte
604     # endif
605     #endif /* ALLOW_AUTODIFF_TAMC */
606    
607     #ifdef ALLOW_DEBUG
608     IF (debugMode) CALL DEBUG_CALL('DO_OCEANIC_PHYS',myThid)
609     #endif
610     CALL TIMER_START('DO_OCEANIC_PHYS [FORWARD_STEP]',myThid)
611     CALL DO_OCEANIC_PHYS( myTime, myIter, myThid )
612     CALL TIMER_STOP ('DO_OCEANIC_PHYS [FORWARD_STEP]',myThid)
613    
614     #ifdef ALLOW_STREAMICE
615     IF (useStreamIce) THEN
616     CALL STREAMICE_TIMESTEP ( myThid, myIter,
617     & iLoop, myTime )
618     ENDIF
619     #endif
620    
621     #ifdef ALLOW_AUTODIFF_TAMC
622     CADJ STORE EmPmR = comlev1, key = ikey_dynamics,
623     CADJ & kind = isbyte
624     # ifdef EXACT_CONSERV
625     CADJ STORE PmEpR = comlev1, key = ikey_dynamics,
626     CADJ & kind = isbyte
627     # endif
628     cph-test(
629     CADJ STORE surfaceForcingU = comlev1, key=ikey_dynamics, kind=isbyte
630     CADJ STORE surfaceForcingV = comlev1, key=ikey_dynamics, kind=isbyte
631     CADJ STORE qsw = comlev1, key = ikey_dynamics, kind = isbyte
632     # ifdef ATMOSPHERIC_LOADING
633     CADJ STORE phi0surf = comlev1, key = ikey_dynamics, kind = isbyte
634     # endif
635     cph-test)
636    
637     # ifdef ALLOW_DEPTH_CONTROL
638     CADJ STORE hFacC, hFacS, hFacW
639     CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
640     CADJ STORE recip_hFacC, recip_hFacS, recip_hFacW
641     CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
642     CADJ STORE surfaceForcingU, surfaceForcingV =
643     CADJ & comlev1, key = ikey_dynamics, kind = isbyte
644     # endif
645     #endif /* ALLOW_AUTODIFF_TAMC */
646    
647     #ifdef ALLOW_GCHEM
648     IF ( useGCHEM ) THEN
649     #ifdef ALLOW_AUTODIFF_TAMC
650     CADJ STORE pTracer = comlev1, key = ikey_dynamics,
651     CADJ & kind = isbyte
652     CADJ STORE theta, salt = comlev1, key = ikey_dynamics,
653     CADJ & kind = isbyte
654     #endif
655     #ifdef ALLOW_DEBUG
656     IF (debugMode) CALL DEBUG_CALL('GCHEM_CALC_TENDENCY',myThid)
657     #endif
658     CALL TIMER_START('GCHEM_CALC_TENDENCY [FORWARD_STEP]',myThid)
659     CALL GCHEM_CALC_TENDENCY( myTime, myIter, myThid )
660     CALL TIMER_STOP ('GCHEM_CALC_TENDENCY [FORWARD_STEP]',myThid)
661     ENDIF
662     #endif /* ALLOW_GCHEM */
663    
664     #ifdef ALLOW_AUTODIFF_TAMC
665     cph needed to be moved here from do_oceanic_physics
666     cph to be visible down the road
667    
668     CADJ STORE rhoInSitu = comlev1, key = ikey_dynamics,
669     CADJ & kind = isbyte
670     CADJ STORE surfaceForcingS = comlev1, key = ikey_dynamics,
671     CADJ & kind = isbyte
672     CADJ STORE surfaceForcingT = comlev1, key = ikey_dynamics,
673     CADJ & kind = isbyte
674     CADJ STORE surfaceForcingTice = comlev1, key = ikey_dynamics,
675     CADJ & kind = isbyte
676     CADJ STORE IVDConvCount = comlev1, key = ikey_dynamics,
677     CADJ & kind = isbyte
678     # ifdef ALLOW_PTRACERS
679     CADJ STORE surfaceForcingPTr = comlev1, key = ikey_dynamics,
680     CADJ & kind = isbyte
681     # endif
682    
683     # ifdef ALLOW_GMREDI
684     CADJ STORE Kwx = comlev1, key = ikey_dynamics,
685     CADJ & kind = isbyte
686     CADJ STORE Kwy = comlev1, key = ikey_dynamics,
687     CADJ & kind = isbyte
688     CADJ STORE Kwz = comlev1, key = ikey_dynamics,
689     CADJ & kind = isbyte
690     # ifdef GM_BOLUS_ADVEC
691     CADJ STORE GM_PsiX = comlev1, key = ikey_dynamics,
692     CADJ & kind = isbyte
693     CADJ STORE GM_PsiY = comlev1, key = ikey_dynamics,
694     CADJ & kind = isbyte
695     # endif
696     # endif
697    
698     # ifdef ALLOW_KPP
699     CADJ STORE KPPghat = comlev1, key = ikey_dynamics,
700     CADJ & kind = isbyte
701     CADJ STORE KPPfrac = comlev1, key = ikey_dynamics,
702     CADJ & kind = isbyte
703     CADJ STORE KPPdiffKzS = comlev1, key = ikey_dynamics,
704     CADJ & kind = isbyte
705     CADJ STORE KPPdiffKzT = comlev1, key = ikey_dynamics,
706     CADJ & kind = isbyte
707     # endif
708    
709     # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
710     CADJ STORE theta,salt = comlev1, key = ikey_dynamics, kind = isbyte
711     CADJ STORE etaH = comlev1, key = ikey_dynamics, kind = isbyte
712     # ifdef ALLOW_CD_CODE
713     CADJ STORE etaNm1 = comlev1, key = ikey_dynamics, kind = isbyte
714     # endif
715     # ifndef DISABLE_RSTAR_CODE
716     CADJ STORE rStarExpC = comlev1, key = ikey_dynamics, kind = isbyte
717     # endif
718     # endif
719     #endif /* ALLOW_AUTODIFF_TAMC */
720    
721     #ifdef ALLOW_LONGSTEP
722     IF ( usePTRACERS .AND. LS_whenToSample .EQ. 0 ) THEN
723     C Average all variables before advection (but after do_oceanic_phys
724     C where Qsw, KPP and GMRedi stuff is computed).
725     C This is like diagnostics package and will reproduce offline
726     C results.
727     #ifdef ALLOW_DEBUG
728     IF (debugMode) CALL DEBUG_CALL('LONGSTEP_AVERAGE',myThid)
729     #endif
730     CALL TIMER_START('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid)
731     CALL LONGSTEP_AVERAGE( myTime, myIter, myThid )
732     CALL TIMER_STOP ('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid)
733    
734     #ifdef ALLOW_DEBUG
735     IF (debugMode)
736     & CALL DEBUG_CALL('LONGSTEP_THERMODYNAMICS',myThid)
737     #endif
738     CALL TIMER_START('LONGSTEP_THERMODYNAMICS [FORWARD_STEP]',
739     & myThid)
740     CALL LONGSTEP_THERMODYNAMICS( myTime, myIter, myThid )
741     CALL TIMER_STOP ('LONGSTEP_THERMODYNAMICS [FORWARD_STEP]',
742     & myThid)
743     ENDIF
744     #endif /* ALLOW_LONGSTEP */
745    
746     IF ( .NOT.staggerTimeStep ) THEN
747     #ifdef ALLOW_AUTODIFF_TAMC
748     CADJ STORE wVel = comlev1, key = ikey_dynamics, kind = isbyte
749     #endif
750     #ifdef ALLOW_DEBUG
751     IF (debugMode) CALL DEBUG_CALL('THERMODYNAMICS',myThid)
752     #endif
753     CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',myThid)
754     CALL THERMODYNAMICS( myTime, myIter, myThid )
755     CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',myThid)
756     C-- if not staggerTimeStep: end
757     ENDIF
758    
759     #ifdef ALLOW_LONGSTEP
760     IF ( usePTRACERS .AND. LS_whenToSample .EQ. 1 ) THEN
761     C Average T and S after thermodynamics, but U,V,W before dynamics.
762     C This will reproduce online results with staggerTimeStep=.FALSE.
763     C for LS_nIter=1
764     #ifdef ALLOW_DEBUG
765     IF (debugMode) CALL DEBUG_CALL('LONGSTEP_AVERAGE',myThid)
766     #endif
767     CALL TIMER_START('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid)
768     CALL LONGSTEP_AVERAGE( myTime, myIter, myThid )
769     CALL TIMER_STOP ('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid)
770    
771     #ifdef ALLOW_DEBUG
772     IF (debugMode)
773     & CALL DEBUG_CALL('LONGSTEP_THERMODYNAMICS',myThid)
774     #endif
775     CALL TIMER_START('LONGSTEP_THERMODYNAMICS [FORWARD_STEP]',
776     & myThid)
777     CALL LONGSTEP_THERMODYNAMICS( myTime, myIter, myThid )
778     CALL TIMER_STOP ('LONGSTEP_THERMODYNAMICS [FORWARD_STEP]',
779     & myThid)
780     ENDIF
781     #endif /* ALLOW_LONGSTEP */
782    
783     c #ifdef ALLOW_NONHYDROSTATIC
784     IF ( implicitIntGravWave ) THEN
785     CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
786     CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid )
787     CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
788     ENDIF
789     c #endif
790    
791     #ifdef ALLOW_AUTODIFF_TAMC
792     CADJ STORE etaN = comlev1, key = ikey_dynamics, kind = isbyte
793     # ifdef ALLOW_DEPTH_CONTROL
794     CADJ STORE hFacC, hFacS, hFacW
795     CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
796     CADJ STORE recip_hFacC, recip_hFacS, recip_hFacW
797     CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
798     # endif /* ALLOW_DEPTH_CONTROL */
799     #endif /* ALLOW_AUTODIFF_TAMC */
800    
801     C-- Step forward fields and calculate time tendency terms.
802     #ifdef ALLOW_MOM_STEPPING
803     #ifndef ALLOW_AUTODIFF
804     IF ( momStepping ) THEN
805     #endif
806     #ifdef ALLOW_DEBUG
807     IF (debugMode) CALL DEBUG_CALL('DYNAMICS',myThid)
808     #endif
809     CALL TIMER_START('DYNAMICS [FORWARD_STEP]',myThid)
810     CALL DYNAMICS( myTime, myIter, myThid )
811     CALL TIMER_STOP ('DYNAMICS [FORWARD_STEP]',myThid)
812     #ifndef ALLOW_AUTODIFF
813     ENDIF
814     #endif
815     #endif /* ALLOW_MOM_STEPPING */
816    
817     #ifdef ALLOW_AUTODIFF_TAMC
818     # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
819     CADJ STORE gU, gV = comlev1, key = ikey_dynamics,
820     CADJ & kind = isbyte
821     # endif
822     #endif
823    
824     C-- Update time-counter
825     myIter = nIter0 + iLoop
826     myTime = startTime + deltaTClock*iLoop
827    
828     #ifdef ALLOW_MNC
829     C Update MNC time information
830     IF ( useMNC ) THEN
831     CALL MNC_UPDATE_TIME( myTime, myIter, myThid )
832     ENDIF
833     #endif /* ALLOW_MNC */
834    
835     #ifdef ALLOW_OFFLINE
836     C Load new Offline fields and update state-variable
837     IF ( useOffLine ) THEN
838     #ifdef ALLOW_DEBUG
839     IF (debugMode) CALL DEBUG_CALL('OFFLINE_FIELDS_LOAD',myThid)
840     #endif /* ALLOW_DEBUG */
841     CALL TIMER_START('OFFLINE_FLDS_LOAD [FORWARD_STEP]',myThid)
842     CALL OFFLINE_FIELDS_LOAD( myTime, myIter, myThid )
843     CALL TIMER_STOP ('OFFLINE_FLDS_LOAD [FORWARD_STEP]',myThid)
844     ENDIF
845     #endif /* ALLOW_OFFLINE */
846    
847     C-- Update geometric factors:
848     #ifdef NONLIN_FRSURF
849     C- update hfacC,W,S and recip_hFac according to etaH(n+1) :
850     IF ( select_rStar.GT.0 ) THEN
851     # ifndef DISABLE_RSTAR_CODE
852     # ifdef ALLOW_AUTODIFF_TAMC
853     CADJ STORE rStarFacC, rStarFacS, rStarFacW =
854     CADJ & comlev1, key = ikey_dynamics, kind = isbyte
855     # endif
856     CALL TIMER_START('UPDATE_R_STAR [FORWARD_STEP]',myThid)
857     CALL UPDATE_R_STAR( .TRUE., myTime, myIter, myThid )
858     CALL TIMER_STOP ('UPDATE_R_STAR [FORWARD_STEP]',myThid)
859     # endif /* DISABLE_RSTAR_CODE */
860     ELSEIF ( selectSigmaCoord.NE.0 ) THEN
861     # ifndef DISABLE_SIGMA_CODE
862     CALL UPDATE_SIGMA( etaH, myTime, myIter, myThid )
863     # endif /* DISABLE_RSTAR_CODE */
864     ELSE
865     # ifdef ALLOW_AUTODIFF_TAMC
866     CADJ STORE hFac_surfC, hFac_surfS, hFac_surfW
867     CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
868     # endif
869     CALL TIMER_START('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
870     CALL UPDATE_SURF_DR( .TRUE., myTime, myIter, myThid )
871     CALL TIMER_STOP ('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
872     ENDIF
873     # ifdef ALLOW_AUTODIFF_TAMC
874     CADJ STORE hFacC, hFacS, hFacW =
875     CADJ & comlev1, key = ikey_dynamics, kind = isbyte
876     CADJ STORE recip_hFacC, recip_hFacS, recip_hFacW =
877     CADJ & comlev1, key = ikey_dynamics, kind = isbyte
878     # ifdef ALLOW_CG2D_NSA
879     CADJ STORE aW2d, aS2d, aC2d =
880     CADJ & comlev1, key = ikey_dynamics, kind = isbyte
881     CADJ STORE pC, pS, pW =
882     CADJ & comlev1, key = ikey_dynamics, kind = isbyte
883     # endif
884     # endif
885     C- update also CG2D matrix (and preconditioner)
886     IF ( momStepping .AND. nonlinFreeSurf.GT.2 ) THEN
887     CALL TIMER_START('UPDATE_CG2D [FORWARD_STEP]',myThid)
888     CALL UPDATE_CG2D( myTime, myIter, myThid )
889     CALL TIMER_STOP ('UPDATE_CG2D [FORWARD_STEP]',myThid)
890     ENDIF
891     #endif /* NONLIN_FRSURF */
892    
893     C-- Apply Filters to u*,v* before SOLVE_FOR_PRESSURE
894     #ifdef ALLOW_SHAP_FILT
895     IF (useSHAP_FILT .AND. shap_filt_uvStar) THEN
896     CALL TIMER_START('SHAP_FILT_UV [FORWARD_STEP]',myThid)
897     IF (implicDiv2Dflow.LT.1.) THEN
898     C-- Explicit+Implicit part of the Barotropic Flow Divergence
899     C => Filtering of uVel,vVel is necessary
900     CALL SHAP_FILT_APPLY_UV( uVel,vVel,
901     & myTime, myIter, myThid )
902     ENDIF
903     CALL SHAP_FILT_APPLY_UV( gU,gV,myTime,myIter,myThid)
904     CALL TIMER_STOP ('SHAP_FILT_UV [FORWARD_STEP]',myThid)
905     ENDIF
906     #endif
907     #ifdef ALLOW_ZONAL_FILT
908     IF (useZONAL_FILT .AND. zonal_filt_uvStar) THEN
909     CALL TIMER_START('ZONAL_FILT_UV [FORWARD_STEP]',myThid)
910     IF (implicDiv2Dflow.LT.1.) THEN
911     C-- Explicit+Implicit part of the Barotropic Flow Divergence
912     C => Filtering of uVel,vVel is necessary
913     CALL ZONAL_FILT_APPLY_UV( uVel, vVel, myThid )
914     ENDIF
915     CALL ZONAL_FILT_APPLY_UV( gU, gV, myThid )
916     CALL TIMER_STOP ('ZONAL_FILT_UV [FORWARD_STEP]',myThid)
917     ENDIF
918     #endif
919    
920     C-- Solve elliptic equation(s).
921     C Two-dimensional only for conventional hydrostatic or
922     C three-dimensional for non-hydrostatic and/or IGW scheme.
923     IF ( momStepping ) THEN
924     #ifdef ALLOW_AUTODIFF_TAMC
925     # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
926     CADJ STORE uVel, vVel
927     CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
928     CADJ STORE EmPmR,hFacS,hFacW
929     CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
930     # endif
931     #endif
932     CALL TIMER_START('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
933     CALL SOLVE_FOR_PRESSURE(myTime, myIter, myThid)
934     CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
935     ENDIF
936    
937     C-- Correct divergence in flow field and cycle time-stepping momentum
938     #ifdef ALLOW_MOM_STEPPING
939     #ifndef ALLOW_AUTODIFF
940     IF ( momStepping ) THEN
941     #endif
942     #ifdef ALLOW_AUTODIFF_TAMC
943     # ifdef ALLOW_DEPTH_CONTROL
944     CADJ STORE etaN, uVel,vVel
945     CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
946     # endif /* ALLOW_DEPTH_CONTROL */
947     #endif /* ALLOW_AUTODIFF_TAMC */
948     CALL TIMER_START('MOM_CORRECTION_STEP [FORWARD_STEP]',myThid)
949     CALL MOMENTUM_CORRECTION_STEP(myTime, myIter, myThid)
950     CALL TIMER_STOP ('MOM_CORRECTION_STEP [FORWARD_STEP]',myThid)
951     #ifndef ALLOW_AUTODIFF
952     ENDIF
953     #endif
954     #endif /* ALLOW_MOM_STEPPING */
955     #ifdef ALLOW_AUTODIFF_TAMC
956     CADJ STORE uVel, vVel = comlev1, key = ikey_dynamics, kind = isbyte
957     #endif
958    
959     IF ( calc_wVelocity ) THEN
960     C-- Integrate continuity vertically for vertical velocity
961     C (+ update "etaN" & "etaH", exact volume conservation):
962     CALL TIMER_START('INTEGR_CONTINUITY [FORWARD_STEP]',myThid)
963     CALL INTEGR_CONTINUITY( uVel, vVel, myTime, myIter, myThid)
964     CALL TIMER_STOP ('INTEGR_CONTINUITY [FORWARD_STEP]',myThid)
965     ENDIF
966    
967     #ifdef NONLIN_FRSURF
968     IF ( select_rStar.NE.0 ) THEN
969     # ifndef DISABLE_RSTAR_CODE
970     # ifdef ALLOW_AUTODIFF_TAMC
971     CADJ STORE etaH
972     CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
973     CADJ STORE rStarFacC,rStarFacS,rStarFacW
974     CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
975     # endif
976     C-- r* : compute the future level thickness according to etaH(n+1)
977     CALL TIMER_START('CALC_R_STAR [FORWARD_STEP]',myThid)
978     CALL CALC_R_STAR(etaH, myTime, myIter, myThid )
979     CALL TIMER_STOP ('CALC_R_STAR [FORWARD_STEP]',myThid)
980     # endif /* DISABLE_RSTAR_CODE */
981     ELSEIF ( nonlinFreeSurf.GT.0 .AND. selectSigmaCoord.EQ.0 ) THEN
982     C-- compute the future surface level thickness according to etaH(n+1)
983     # ifdef ALLOW_AUTODIFF_TAMC
984     CADJ STORE etaH = comlev1, key = ikey_dynamics,
985     CADJ & kind = isbyte
986     # endif
987     CALL TIMER_START('CALC_SURF_DR [FORWARD_STEP]',myThid)
988     CALL CALC_SURF_DR(etaH, myTime, myIter, myThid )
989     CALL TIMER_STOP ('CALC_SURF_DR [FORWARD_STEP]',myThid)
990     ENDIF
991     # ifdef ALLOW_AUTODIFF_TAMC
992     CADJ STORE rStarExpC
993     CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
994     CADJ STORE hFac_surfC = comlev1, key = ikey_dynamics,
995     CADJ & kind = isbyte
996     CADJ STORE salt,theta = comlev1, key = ikey_dynamics,
997     CADJ & kind = isbyte
998     # endif
999     #endif /* NONLIN_FRSURF */
1000    
1001     C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
1002     IF ( staggerTimeStep ) THEN
1003     C-- do exchanges of U,V (needed for multiDim) when using stagger time-step :
1004     #ifdef ALLOW_DEBUG
1005     IF (debugMode)
1006     & CALL DEBUG_CALL('DO_STAGGER_FIELDS_EXCH.',myThid)
1007     #endif
1008     CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
1009     CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid )
1010     CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
1011    
1012     #ifdef ALLOW_DIAGNOSTICS
1013     C-- State-variables diagnostics
1014     IF ( useDiagnostics ) THEN
1015     CALL TIMER_START('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
1016     CALL DO_STATEVARS_DIAGS( myTime, 1, myIter, myThid )
1017     CALL TIMER_STOP ('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
1018     ENDIF
1019     #endif
1020    
1021    
1022    
1023    
1024    
1025    
1026     #ifdef ALLOW_AUTODIFF_TAMC
1027     CADJ STORE wVel = comlev1, key = ikey_dynamics, kind = isbyte
1028     #endif
1029     #ifdef ALLOW_DEBUG
1030     IF (debugMode) CALL DEBUG_CALL('THERMODYNAMICS',myThid)
1031     #endif
1032     CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',myThid)
1033     CALL THERMODYNAMICS( myTime, myIter, myThid )
1034     CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',myThid)
1035    
1036     C-- if staggerTimeStep: end
1037     ENDIF
1038     C---+--------+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
1039    
1040 dgoldberg 1.2 #ifdef ALLOW_SHELFICE
1041     IF (SHELFICEMassStepping) THEN
1042     CALL SHELFICE_STEP_ICEMASS( myTime, myIter, myThid )
1043     ENDIF
1044     C-- Calculate new loading anomaly (in case the ice-shelf mass was updated)
1045     !#ifndef ALLOW_AUTODIFF
1046     c IF ( SHELFICEloadAnomalyFile .EQ. ' ' ) THEN
1047     c ENDIF
1048     !#endif /* ndef ALLOW_AUTODIFF */
1049     #endif
1050    
1051 dgoldberg 1.1 #ifdef ALLOW_AUTODIFF_TAMC
1052     cph This is needed because convective_adjustment calls
1053     cph find_rho which may use pressure()
1054     CADJ STORE totPhiHyd = comlev1, key = ikey_dynamics,
1055     CADJ & kind = isbyte
1056     #endif
1057     C-- Apply adjustments to Tracers arrays (T,S,+pTracers)
1058     CALL TIMER_START('TRC_CORRECTION_STEP [FORWARD_STEP]',myThid)
1059     CALL TRACERS_CORRECTION_STEP(myTime, myIter, myThid)
1060     CALL TIMER_STOP ('TRC_CORRECTION_STEP [FORWARD_STEP]',myThid)
1061    
1062     #ifdef ALLOW_LONGSTEP
1063     IF ( usePTRACERS ) THEN
1064     IF ( LS_whenToSample .EQ. 2 ) THEN
1065     C Average everything at the end of the timestep. This will
1066     C reproduce online results with staggerTimeStep=.TRUE.
1067     C when LS_nIter=1
1068     #ifdef ALLOW_DEBUG
1069     IF (debugMode) CALL DEBUG_CALL('LONGSTEP_AVERAGE',myThid)
1070     #endif
1071     CALL TIMER_START('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid)
1072     C myIter has been update after dynamics, but the averaging window
1073     C should be determined by myIter at beginning of timestep
1074     CALL LONGSTEP_AVERAGE( myTimeBeg, myIterBeg, myThid )
1075     CALL TIMER_STOP ('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid)
1076    
1077     #ifdef ALLOW_DEBUG
1078     IF (debugMode)
1079     & CALL DEBUG_CALL('LONGSTEP_THERMODYNAMICS',myThid)
1080     #endif
1081     CALL TIMER_START('LONGSTEP_THERMODYNAMICS [FORWARD_STEP]',
1082     & myThid)
1083     CALL LONGSTEP_THERMODYNAMICS( myTime, myIter, myThid )
1084     CALL TIMER_STOP ('LONGSTEP_THERMODYNAMICS [FORWARD_STEP]',
1085     & myThid)
1086     C-- if LS_whenToSample.EQ.2: end
1087     ENDIF
1088    
1089     C-- Apply adjustments to passive Tracers arrays (pTracers)
1090     c CALL TIMER_START('LS_CORRECTION_STEP [FORWARD_STEP]',myThid)
1091     c CALL LONGSTEP_CORRECTION_STEP(myTime, myIter, myThid)
1092     c CALL TIMER_STOP ('LS_CORRECTION_STEP [FORWARD_STEP]',myThid)
1093     C-- if usePTRACERS: end
1094     ENDIF
1095     #endif /* ALLOW_LONGSTEP */
1096    
1097     #ifdef ALLOW_GCHEM
1098     C Add separate timestepping of chemical/biological/forcing
1099     C of ptracers here in GCHEM_FORCING_SEP
1100     #ifdef ALLOW_LONGSTEP
1101     IF ( useGCHEM .AND. LS_doTimeStep ) THEN
1102     #else
1103     IF ( useGCHEM ) THEN
1104     #endif
1105     #ifdef ALLOW_AUTODIFF_TAMC
1106     CADJ STORE pTracer = comlev1, key = ikey_dynamics,
1107     CADJ & kind = isbyte
1108     CADJ STORE theta, salt = comlev1, key = ikey_dynamics,
1109     CADJ & kind = isbyte
1110     #endif
1111     #ifdef ALLOW_DEBUG
1112     IF (debugMode) CALL DEBUG_CALL('GCHEM_FORCING_SEP',myThid)
1113     #endif /* ALLOW_DEBUG */
1114     CALL TIMER_START('GCHEM_FORCING_SEP [FORWARD_STEP]',myThid)
1115     CALL GCHEM_FORCING_SEP( myTime,myIter,myThid )
1116     CALL TIMER_STOP ('GCHEM_FORCING_SEP [FORWARD_STEP]',myThid)
1117     ENDIF
1118     #endif /* ALLOW_GCHEM */
1119    
1120     C-- Do "blocking" sends and receives for tendency "overlap" terms
1121     c CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
1122     c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )
1123     c CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
1124    
1125     C-- Do "blocking" sends and receives for field "overlap" terms
1126     CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
1127     CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
1128     CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
1129    
1130     #ifdef ALLOW_DIAGNOSTICS
1131     IF ( useDiagnostics ) THEN
1132     CALL TIMER_START('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
1133     CALL DO_STATEVARS_DIAGS( myTime, 2, myIter, myThid )
1134     CALL TIMER_STOP ('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
1135     ENDIF
1136     #endif
1137    
1138     #ifdef ALLOW_GRIDALT
1139     IF (useGRIDALT) THEN
1140     CALL GRIDALT_UPDATE(myThid)
1141     ENDIF
1142     #endif
1143    
1144     #ifdef ALLOW_FIZHI
1145     IF (useFIZHI) THEN
1146     CALL TIMER_START('FIZHI [FORWARD_STEP]',myThid)
1147     CALL STEP_FIZHI_CORR ( myTime, myIter, myThid, dTtracerLev(1) )
1148     CALL TIMER_STOP ('FIZHI [FORWARD_STEP]',myThid)
1149     ENDIF
1150     #endif
1151    
1152     #ifdef ALLOW_FLT
1153     C-- Calculate float trajectories
1154     IF (useFLT) THEN
1155     CALL TIMER_START('FLOATS [FORWARD_STEP]',myThid)
1156     CALL FLT_MAIN( myTime, myIter, myThid )
1157     CALL TIMER_STOP ('FLOATS [FORWARD_STEP]',myThid)
1158     ENDIF
1159     #endif
1160    
1161     #ifdef ALLOW_TIMEAVE
1162     C-- State-variables time-averaging
1163     CALL TIMER_START('DO_STATEVARS_TAVE [FORWARD_STEP]',myThid)
1164     CALL DO_STATEVARS_TAVE( myTime, myIter, myThid )
1165     CALL TIMER_STOP ('DO_STATEVARS_TAVE [FORWARD_STEP]',myThid)
1166     #endif
1167    
1168     #ifdef ALLOW_NEST_PARENT
1169     IF ( useNEST_PARENT) THEN
1170     CALL NEST_PARENT_IO_2( myTime, myIter, myThid )
1171     ENDIF
1172     #endif /* ALLOW_NEST_PARENT */
1173    
1174     #ifdef ALLOW_NEST_CHILD
1175     IF ( useNEST_CHILD) THEN
1176     CALL NEST_CHILD_TRANSP( myTime, myIter, myThid )
1177     ENDIF
1178     #endif /* ALLOW_NEST_CHILD */
1179    
1180     #ifdef ALLOW_MONITOR
1181     IF ( monitorFreq.GT.0. .OR. adjMonitorFreq.GT.0. ) THEN
1182     C-- Check status of solution (statistics, cfl, etc...)
1183     CALL TIMER_START('MONITOR [FORWARD_STEP]',myThid)
1184     CALL MONITOR( myTime, myIter, myThid )
1185     CALL TIMER_STOP ('MONITOR [FORWARD_STEP]',myThid)
1186     ENDIF
1187     #endif /* ALLOW_MONITOR */
1188    
1189     #ifdef ALLOW_COST
1190     C-- compare model with data and compute cost function
1191     C-- this is done after exchanges to allow interpolation
1192     CALL TIMER_START('COST_TILE [FORWARD_STEP]',myThid)
1193     CALL COST_TILE ( myTime, myIter, myThid )
1194     CALL TIMER_STOP ('COST_TILE [FORWARD_STEP]',myThid)
1195     #endif
1196    
1197     C-- Check if it has reached the end of simulation
1198     modelEnd = myTime.EQ.endTime .OR. myIter.EQ.nEndIter
1199     #ifdef HAVE_SIGREG
1200     IF ( useSIGREG ) THEN
1201     modelEnd = modelEnd .OR. ( i_got_signal.GT.0 )
1202     ENDIF
1203     #endif /* HAVE_SIGREG */
1204    
1205     C-- Do IO if needed.
1206     CALL TIMER_START('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
1207     CALL DO_THE_MODEL_IO( modelEnd, myTime, myIter, myThid )
1208     CALL TIMER_STOP ('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
1209    
1210     #ifdef ALLOW_PTRACERS
1211     C Reset the ptracers (but after the io is done)
1212     IF ( usePTRACERS ) THEN
1213     CALL TIMER_START('PTRACERS_RESET [FORWARD_STEP]',myThid)
1214     CALL PTRACERS_RESET( myTime, myIter, myThid )
1215     CALL TIMER_STOP ('PTRACERS_RESET [FORWARD_STEP]',myThid)
1216     ENDIF
1217     #endif /* ALLOW_PTRACERS */
1218    
1219     C-- Save state for restarts
1220     CALL TIMER_START('DO_WRITE_PICKUP [FORWARD_STEP]',myThid)
1221     CALL DO_WRITE_PICKUP( modelEnd, myTime, myIter, myThid )
1222     CALL TIMER_STOP ('DO_WRITE_PICKUP [FORWARD_STEP]',myThid)
1223    
1224     #ifdef HAVE_SIGREG
1225     IF ( useSIGREG ) THEN
1226     IF ( modelEnd .AND. i_got_signal.GT.0 ) THEN
1227     STOP 'Checkpoint completed -- killed by signal handler'
1228     ENDIF
1229     ENDIF
1230     #endif /* HAVE_SIGREG */
1231    
1232     #ifdef ALLOW_AUTODIFF
1233     CALL AUTODIFF_INADMODE_SET( myThid )
1234     #endif
1235    
1236     #ifdef ALLOW_SHOWFLOPS
1237     CALL TIMER_START('SHOWFLOPS_INLOOP [THE_MAIN_LOOP]', myThid)
1238     CALL SHOWFLOPS_INLOOP( iloop, myThid )
1239     CALL TIMER_STOP ('SHOWFLOPS_INLOOP [THE_MAIN_LOOP]', myThid)
1240     #endif
1241    
1242     #ifdef ALLOW_DEBUG
1243     IF (debugMode) CALL DEBUG_LEAVE('FORWARD_STEP',myThid)
1244     #endif
1245    
1246     RETURN
1247     END

  ViewVC Help
Powered by ViewVC 1.1.22