/[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.3 - (hide annotations) (download)
Thu May 5 18:16:04 2016 UTC (9 years, 2 months ago) by dgoldberg
Branch: MAIN
CVS Tags: checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65w
Changes since 1.2: +6 -1 lines
renaming files, merged  shelfice_thermodynamics, added CPP directives

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

  ViewVC Help
Powered by ViewVC 1.1.22