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

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

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


Revision 1.2 - (hide annotations) (download)
Mon Sep 28 12:59:59 2015 UTC (9 years, 10 months ago) by dgoldberg
Branch: MAIN
Changes since 1.1: +1 -376 lines
*** empty log message ***

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

  ViewVC Help
Powered by ViewVC 1.1.22