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

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

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


Revision 1.4 - (hide annotations) (download)
Mon Aug 21 16:06:48 2017 UTC (7 years, 10 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint67a, checkpoint67b, checkpoint67d, HEAD
Changes since 1.3: +28 -31 lines
bring customized version up to date with standard version (in model/src)

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

  ViewVC Help
Powered by ViewVC 1.1.22