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

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

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


Revision 1.3 - (hide annotations) (download)
Wed Aug 19 10:28:51 2015 UTC (9 years, 11 months ago) by dgoldberg
Branch: MAIN
Changes since 1.2: +2 -2 lines
*** empty log message ***

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

  ViewVC Help
Powered by ViewVC 1.1.22