/[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.1 - (hide annotations) (download)
Tue Jul 21 14:52:40 2015 UTC (10 years ago) by dgoldberg
Branch: MAIN
*** empty log message ***

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

  ViewVC Help
Powered by ViewVC 1.1.22