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

Annotation of /MITgcm_contrib/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)
Thu Aug 27 12:40:19 2015 UTC (9 years, 10 months ago) by dgoldberg
Branch: MAIN
Changes since 1.3: +5 -4 lines
*** empty log message ***

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

  ViewVC Help
Powered by ViewVC 1.1.22