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

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

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


Revision 1.1 - (hide annotations) (download)
Fri Dec 11 19:48:31 2015 UTC (9 years, 7 months ago) by dgoldberg
Branch: MAIN
CVS Tags: checkpoint65r, checkpoint65s
shelfice_remeshing verification

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

  ViewVC Help
Powered by ViewVC 1.1.22