/[MITgcm]/MITgcm/model/src/forward_step.F
ViewVC logotype

Diff of /MITgcm/model/src/forward_step.F

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

revision 1.6 by cnh, Sun Feb 4 14:38:47 2001 UTC revision 1.157 by heimbach, Thu May 1 23:52:24 2008 UTC
# Line 1  Line 1 
1  C $Header$  C $Header$
2  C $Name$  C $Name$
3    
4    #include "PACKAGES_CONFIG.h"
5  #include "CPP_OPTIONS.h"  #include "CPP_OPTIONS.h"
6    
7        SUBROUTINE FORWARD_STEP(  #ifdef ALLOW_GMREDI
8       I                    doHalfStep, iLoop,  # include "GMREDI_OPTIONS.h"
9       U                    myCurrentTime, myCurrentIter,  #endif
10       &                    myThid)  #ifdef ALLOW_OBCS
11  C     /==========================================================\  # include "OBCS_OPTIONS.h"
12  C     | SUBROUTINE FORWARD_STEP                                  |  #endif
13  C     | o Does one instance of the model time stepping           |  #ifdef ALLOW_SEAICE
14  C     |   The time stepping loop in THE_MAIN_LOOP() calls        |  # include "SEAICE_OPTIONS.h"
15  C     |   this routine                                           |  #endif
 C     |==========================================================|  
 C     \==========================================================/  
       IMPLICIT NONE  
 C  
 C     Call Tree  
 C     =========  
 C      
 C      THE_MAIN_LOOP()  
 C       |  
 C  ==>  | ** Time stepping loop starts here **  
 C  |    |  
 C /|\   |-FORWARD_STEP  
 C  |    |  |  
 C /|\   |  |--LOAD_EXTERNAL_DATA  
 C  |    |  |   o Load and/or set time dependent forcing fields  
 C /|\   |  |  
 C  |    |  |--DYNAMICS  
 C /|\   |  |   o Evaluate "forward" terms  
 C  |    |  |  
 C /|\   |  |  
 C  |    |  |--DO_THE_MODEL_IO  
 C /|\   |  |   o Write model state  
 C  |    |  |  
 C /|\   |  |  
 C  |    |  |--SOLVE_FOR_PRESSURE  
 C /|\   |  |   o Find pressure field to keep flow non-divergent  
 C  |    |  |  
 C /|\   |  |  
 C  |    |  |--THE_CORRECTION_STEP  
 C /|\   |  |   o Correct flow field with pressure gradient  
 C  |    |  |     and cycle time-stepping arrays (for all fields)  
 C /|\   |  |  
 C  |    |  |--DO_GTERM_BLOCKING_EXCHANGES  
 C /|\   |  |   o Update overlap regions  
 C  |    |  |  
 C /|\   |  |  
 C  |    |  |--WRITE_CHECKPOINT  
 C /|\   |  |   o Write restart file(s)  
 C  |    |  
 C /|\   |  
 C  |<== | ** Time stepping loop finishes here **  
 C  
16    
17  C     == Global variables ===  CBOP
18    C     !ROUTINE: FORWARD_STEP
19    C     !INTERFACE:
20          SUBROUTINE FORWARD_STEP( iloop, myTime, myIter, myThid )
21    
22    C     !DESCRIPTION: \bv
23    C     *==================================================================
24    C     | SUBROUTINE forward_step
25    C     | o Run the ocean model and, optionally, evaluate a cost function.
26    C     *==================================================================
27    C     |
28    C     | THE_MAIN_LOOP is the toplevel routine for the Tangent Linear and
29    C     | Adjoint Model Compiler (TAMC). For this purpose the initialization
30    C     | of the model was split into two parts. Those parameters that do
31    C     | not depend on a specific model run are set in INITIALISE_FIXED,
32    C     | whereas those that do depend on the specific realization are
33    C     | initialized in INITIALISE_VARIA.
34    C     |
35    C     *==================================================================
36    C     \ev
37    
38    C     !USES:
39          IMPLICIT NONE
40    C     == Global variables ==
41  #include "SIZE.h"  #include "SIZE.h"
42  #include "EEPARAMS.h"  #include "EEPARAMS.h"
43  #include "PARAMS.h"  #include "PARAMS.h"
44  #include "DYNVARS.h"  #include "DYNVARS.h"
45  #include "CG2D.h"  
46  #ifdef ALLOW_NONHYDROSTATIC  #ifdef ALLOW_MNC
47  #include "CG3D.h"  #include "MNC_PARAMS.h"
48          EXTERNAL DIFFERENT_MULTIPLE
49          LOGICAL  DIFFERENT_MULTIPLE
50    #endif
51    
52    #ifdef HAVE_SIGREG
53    #include "SIGREG.h"
54    #endif
55    
56    #ifdef ALLOW_SHAP_FILT
57    # include "SHAP_FILT.h"
58  #endif  #endif
59    #ifdef ALLOW_ZONAL_FILT
60    # include "ZONAL_FILT.h"
61    #endif
62    #ifdef COMPONENT_MODULE
63    # include "CPL_PARAMS.h"
64    #endif
65    
66    #ifdef ALLOW_AUTODIFF_TAMC
67    # include "FFIELDS.h"
68    # include "SURFACE.h"
69    
70    # include "tamc.h"
71    # include "ctrl.h"
72    # include "ctrl_dummy.h"
73    # include "cost.h"
74    # include "EOS.h"
75    # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
76    #  include "GRID.h"
77    # endif
78    # ifdef ALLOW_EXF
79    #  include "EXF_FIELDS.h"
80    #  ifdef ALLOW_BULKFORMULAE
81    #   include "EXF_CONSTANTS.h"
82    #  endif
83    # endif
84    # ifdef ALLOW_PTRACERS
85    #  include "PTRACERS_SIZE.h"
86    #  include "PTRACERS_FIELDS.h"
87    # endif
88    # ifdef ALLOW_GCHEM
89    #  include "GCHEM_FIELDS.h"
90    # endif
91    # ifdef ALLOW_CFC
92    #  include "CFC.h"
93    # endif
94    # ifdef ALLOW_DIC
95    #  include "DIC_VARS.h"
96    #  include "DIC_LOAD.h"
97    #  include "DIC_ATMOS.h"
98    # endif
99    # ifdef ALLOW_OBCS
100    #  include "OBCS.h"
101    #  ifdef ALLOW_PTRACERS
102    #   include "OBCS_PTRACERS.h"
103    #  endif
104    # endif
105    # ifdef ALLOW_CD_CODE
106    #  include "CD_CODE_VARS.h"
107    # endif
108    # ifdef ALLOW_THSICE
109    #  include "THSICE_VARS.h"
110    # endif
111    # ifdef ALLOW_SEAICE
112    #  include "SEAICE.h"
113    # endif
114    # ifdef ALLOW_EBM
115    #  include "EBM.h"
116    # endif
117    # ifdef ALLOW_KPP
118    #  include "KPP.h"
119    # endif
120    # ifdef ALLOW_GMREDI
121    #  include "GMREDI.h"
122    # endif
123    # ifdef ALLOW_RBCS
124    #  include "RBCS.h"
125    # endif
126    #endif /* ALLOW_AUTODIFF_TAMC */
127    
128    C     !INPUT/OUTPUT PARAMETERS:
129  C     == Routine arguments ==  C     == Routine arguments ==
130  C     doHalfStep - If .TRUE. then this is the last half step  C     note: under the multi-threaded model myIter and
131  C     iLoop         - Invocation count (counter in THE_MAIN_LOOP)  C           myTime are local variables passed around as routine
132  C     myCurrentIter - Iteration counter for this thread  C           arguments. Although this is fiddly it saves the need to
133  C     myCurrentTime - Time counter for this thread  C           impose additional synchronisation points when they are
134  C     myThid - Thread number for this instance of the routine.  C           updated.
135        LOGICAL doHalfStep  C     myTime :: time counter for this thread
136        INTEGER iLoop  C     myIter :: iteration counter for this thread
137        INTEGER myCurrentIter  C     myThid :: thread number for this instance of the routine.
138        _RL     myCurrentTime        INTEGER iloop
139        INTEGER myThid          _RL     myTime
140          INTEGER myIter
141          INTEGER myThid
142    
143    C     !LOCAL VARIABLES:
144  C     == Local variables ==  C     == Local variables ==
145    C     modelEnd  :: true if reaching the end of the run
146          LOGICAL modelEnd
147    #ifdef COMPONENT_MODULE
148          INTEGER myItP1
149    #endif
150    CEOP
151    
152    #ifdef ALLOW_DEBUG
153          IF ( debugLevel .GE. debLevB )
154         &    CALL DEBUG_ENTER('FORWARD_STEP',myThid)
155    #endif
156    
157    #ifdef ALLOW_AUTODIFF_TAMC
158          CALL AUTODIFF_INADMODE_UNSET( myThid )
159    #endif
160    
161    #ifdef ALLOW_AUTODIFF_TAMC
162    C--   Reset the model iteration counter and the model time.
163          myIter = nIter0 + (iloop-1)
164          myTime = startTime + float(iloop-1)*deltaTclock
165    #endif
166    
167  C--   Load forcing/external data fields  #ifdef ALLOW_AUTODIFF_TAMC
168        CALL TIMER_START('I/O (READ)         [FORWARD_STEP]',myThid)  c**************************************
169        CALL EXTERNAL_FIELDS_LOAD( myCurrentTime, myCurrentIter, myThid )  #include "checkpoint_lev1_directives.h"
170        CALL TIMER_STOP ('I/O (READ)         [FORWARD_STEP]',myThid)  #include "checkpoint_lev1_template.h"
171    c**************************************
172  C--   Step forward fields and calculate time tendency terms  #endif
173        CALL TIMER_START('DYNAMICS           [FORWARD_STEP]',myThid)  
174        CALL DYNAMICS( myCurrentTime, myCurrentIter, myThid )  C--   Switch on/off diagnostics for snap-shot output:
175        CALL TIMER_STOP ('DYNAMICS           [FORWARD_STEP]',myThid)  #ifdef ALLOW_DIAGNOSTICS
176          IF ( useDiagnostics ) THEN
177  #ifdef ALLOW_NONHYDROSTATIC          CALL DIAGNOSTICS_SWITCH_ONOFF( myTime, myIter, myThid )
178  C--   Step forward W field in N-H algorithm  C--   State-variables diagnostics
179        IF ( nonHydrostatic ) THEN          CALL TIMER_START('DO_STATEVARS_DIAGS  [FORWARD_STEP]',myThid)
180         CALL TIMER_START('CALC_GW            [FORWARD_STEP]',myThid)          CALL DO_STATEVARS_DIAGS( myTime, 0, myIter, myThid )
181         CALL CALC_GW( myThid)          CALL TIMER_STOP ('DO_STATEVARS_DIAGS  [FORWARD_STEP]',myThid)
        CALL TIMER_STOP ('CALC_GW            [FORWARD_STEP]',myThid)  
182        ENDIF        ENDIF
183  #endif  #endif
184    
185  #ifdef ALLOW_SHAP_FILT  #ifdef ALLOW_PROFILES
186  C--   Step forward all tiles, filter and exchange.  #ifdef ALLOW_DEBUG
187        CALL TIMER_START('SHAP_FILT          [FORWARD_STEP]',myThid)        IF (debugMode) CALL DEBUG_CALL('',myThid)
188        CALL SHAP_FILT_APPLY(  #endif
189       I                     gUnm1, gVnm1, gTnm1, gSnm1,  c--     Accumulate in-situ time averages of theta, salt, and SSH.
190       I                     myCurrentTime, myCurrentIter, myThid )          CALL TIMER_START('PROFILES_INLOOP    [THE_MAIN_LOOP]', mythid)
191        CALL TIMER_STOP ('SHAP_FILT          [FORWARD_STEP]',myThid)          CALL PROFILES_INLOOP( mytime, mythid )
192            CALL TIMER_STOP ('PROFILES_INLOOP    [THE_MAIN_LOOP]', mythid)
193  #endif  #endif
194    
195  #ifdef ALLOW_ZONAL_FILT  C--   Call driver to load external forcing fields from file
196        CALL ZONAL_FILT_APPLY(  #ifdef ALLOW_DEBUG
197       U           gUnm1, gVnm1, gTnm1, gSnm1,        IF ( debugLevel .GE. debLevB )
198       I           myThid )       & CALL DEBUG_CALL('LOAD_FIELDS_DRIVER',myThid)
199    #endif
200    #ifdef ALLOW_AUTODIFF_TAMC
201    cph Important STORE that avoids hidden recomp. of load_fields_driver
202    CADJ STORE theta      = comlev1, key = ikey_dynamics
203    #endif
204          CALL TIMER_START('LOAD_FIELDS_DRIVER  [FORWARD_STEP]',myThid)
205          CALL LOAD_FIELDS_DRIVER( myTime, myIter, myThid )
206          CALL TIMER_STOP ('LOAD_FIELDS_DRIVER  [FORWARD_STEP]',myThid)
207    
208    C--   Call Bulk-Formulae forcing package
209    #ifdef ALLOW_BULK_FORCE
210          IF ( useBulkForce ) THEN
211    #ifdef ALLOW_DEBUG
212            IF ( debugLevel .GE. debLevB )
213         &   CALL DEBUG_CALL('BULKF_FORCING',myThid)
214    #endif
215            CALL TIMER_START('BULKF_FORCING       [FORWARD_STEP]',myThid)
216    C-    calculate qnet and empmr (and wind stress)
217            CALL BULKF_FORCING( myTime, myIter, myThid )
218            CALL TIMER_STOP ('BULKF_FORCING       [FORWARD_STEP]',myThid)
219          ENDIF
220    #endif /* ALLOW_BULK_FORCE */
221    
222    #ifdef ALLOW_AUTODIFF
223    c--   Add control vector for forcing and parameter fields
224          IF ( myIter .EQ. nIter0 )
225         &     CALL CTRL_MAP_FORCING (myThid)
226  #endif  #endif
227    
228  C--   Do IO if needed.  #if (defined (ALLOW_AUTODIFF_TAMC) && defined (ALLOW_AUTODIFF_MONITOR))
229  C     Note:        CALL DUMMY_IN_STEPPING( myTime, myIter, myThid )
230  C     =====  #endif
231  C     At this point model arrays hold U,V,T,S  at "time-level" N  
232  C     and cg2d_x at "time-level" N-1/2 where N = I+timeLevBase-1.  #ifdef COMPONENT_MODULE
233  C     By convention this is taken to be the model "state".        IF ( useCoupler .AND. cpl_earlyExpImpCall ) THEN
234        CALL TIMER_START('I/O (WRITE)        [FORWARD_STEP]',myThid)  C      Post coupling data that I export.
235        CALL DO_THE_MODEL_IO(  C      Read in coupling data that I import.
236       &        doHalfStep, myCurrentTime, myCurrentIter, myThid )           CALL TIMER_START('CPL_EXPORT-IMPORT  [FORWARD_STEP]',myThid)
237        CALL TIMER_STOP ('I/O (WRITE)        [FORWARD_STEP]',myThid)           CALL CPL_EXPORT_MY_DATA(       myIter, myTime, myThid )
238             CALL CPL_IMPORT_EXTERNAL_DATA( myIter, myTime, myThid )
239             CALL TIMER_STOP ('CPL_EXPORT-IMPORT  [FORWARD_STEP]',myThid)
240          ENDIF
241    #endif /* COMPONENT_MODULE */
242    
243    #ifdef ALLOW_EBM
244          IF ( useEBM ) THEN
245    # ifdef ALLOW_DEBUG
246             IF ( debugLevel .GE. debLevB )
247         &    CALL DEBUG_CALL('EBM',myThid)
248    # endif
249             CALL TIMER_START('EBM                [FORWARD_STEP]',myThid)
250             CALL EBM_DRIVER ( myTime, myIter, myThid )
251             CALL TIMER_STOP ('EBM                [FORWARD_STEP]',myThid)
252          ENDIF
253    #endif /* ALLOW_EBM */
254    
255    C--     Step forward fields and calculate time tendency terms.
256    
257    #ifdef ALLOW_DEBUG
258          IF ( debugLevel .GE. debLevB )
259         & CALL DEBUG_CALL('DO_ATMOSPHERIC_PHYS',myThid)
260    #endif
261          CALL TIMER_START('DO_ATMOSPHERIC_PHYS [FORWARD_STEP]',myThid)
262          CALL DO_ATMOSPHERIC_PHYS( myTime, myIter, myThid )
263          CALL TIMER_STOP ('DO_ATMOSPHERIC_PHYS [FORWARD_STEP]',myThid)
264    
265    #ifdef ALLOW_AUTODIFF_TAMC
266    CADJ STORE surfaceforcingtice    = comlev1, key = ikey_dynamics
267    # ifdef ALLOW_KPP
268    CADJ STORE uvel               = comlev1, key = ikey_dynamics
269    CADJ STORE vvel               = comlev1, key = ikey_dynamics
270    # endif
271    # ifdef EXACT_CONSERV
272    cphCADJ STORE empmr              = comlev1, key = ikey_dynamics
273    cphCADJ STORE pmepr              = comlev1, key = ikey_dynamics
274    # endif
275    # ifdef ALLOW_PTRACERS
276    CADJ STORE ptracer               = comlev1, key = ikey_dynamics
277    # endif
278    # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
279    cph-test
280    CADJ STORE hFacC                 = comlev1, key = ikey_dynamics
281    #  ifndef DISABLE_RSTAR_CODE
282    CADJ STORE rstarexpc             = comlev1, key = ikey_dynamics
283    #  endif
284    # endif
285    #endif /* ALLOW_AUTODIFF_TAMC */
286    
287    #ifndef ALLOW_AUTODIFF_TAMC
288          IF ( .NOT. useOffLine ) THEN
289    #endif
290    #ifdef ALLOW_DEBUG
291           IF ( debugLevel .GE. debLevB )
292         &    CALL DEBUG_CALL('DO_OCEANIC_PHYS',myThid)
293    #endif
294           CALL TIMER_START('DO_OCEANIC_PHYS     [FORWARD_STEP]',myThid)
295           CALL DO_OCEANIC_PHYS( myTime, myIter, myThid )
296           CALL TIMER_STOP ('DO_OCEANIC_PHYS     [FORWARD_STEP]',myThid)
297    #ifdef ALLOW_AUTODIFF_TAMC
298    CADJ STORE EmPmR    = comlev1, key = ikey_dynamics
299    # ifdef EXACT_CONSERV
300    CADJ STORE pmepr    = comlev1, key = ikey_dynamics
301    # endif
302    #else
303          ENDIF
304    #endif
305    
306    #ifdef ALLOW_AUTODIFF_TAMC
307    # ifdef NONLIN_FRSURF
308    cph-test
309    CADJ STORE hFac_surfC         = comlev1, key = ikey_dynamics
310    CADJ STORE hfac_surfs         = comlev1, key = ikey_dynamics
311    CADJ STORE hfac_surfw         = comlev1, key = ikey_dynamics
312    # endif
313    # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
314    CADJ STORE hFacC, hFacS, hFacW
315    CADJ &     = comlev1, key = ikey_dynamics
316    CADJ STORE recip_hFacC, recip_hFacS, recip_hFacW
317    CADJ &     = comlev1, key = ikey_dynamics
318    c
319    CADJ STORE surfaceforcingu = comlev1, key = ikey_dynamics
320    CADJ STORE surfaceforcingv = comlev1, key = ikey_dynamics
321    # endif
322    #endif /* ALLOW_AUTODIFF_TAMC */
323    
324    #ifdef ALLOW_GCHEM
325    #ifdef ALLOW_AUTODIFF_TAMC
326    CADJ STORE ptracer  = comlev1, key = ikey_dynamics
327    CADJ STORE theta  = comlev1, key = ikey_dynamics
328    CADJ STORE salt  = comlev1, key = ikey_dynamics
329    #endif
330          IF ( useGCHEM ) THEN
331    #ifdef ALLOW_DEBUG
332            IF ( debugLevel .GE. debLevB )
333         &       CALL DEBUG_CALL('GCHEM_CALC_TENDENCY',myThid)
334    #endif
335            CALL TIMER_START('GCHEM_CALC_TENDENCY [FORWARD_STEP]',myThid)
336            CALL GCHEM_CALC_TENDENCY( myTime, myIter, myThid )
337            CALL TIMER_STOP ('GCHEM_CALC_TENDENCY [FORWARD_STEP]',myThid)
338          ENDIF
339    #endif /* ALLOW_GCHEM */
340    
341    #ifdef ALLOW_AUTODIFF_TAMC
342    cph needed to be moved here from do_oceanic_physics
343    cph to be visible down the road
344    c
345    CADJ STORE surfaceForcingS    = comlev1, key = ikey_dynamics
346    CADJ STORE surfaceForcingT    = comlev1, key = ikey_dynamics
347    CADJ STORE surfaceForcingTice = comlev1, key = ikey_dynamics
348    ctest(
349    CADJ STORE IVDConvCount       = comlev1, key = ikey_dynamics
350    ctest)
351    # ifdef ALLOW_PTRACERS
352    CADJ STORE surfaceForcingPTr  = comlev1, key = ikey_dynamics
353    # endif
354    c
355    # ifdef ALLOW_GMREDI
356    CADJ STORE Kwx                = comlev1, key = ikey_dynamics
357    CADJ STORE Kwy                = comlev1, key = ikey_dynamics
358    CADJ STORE Kwz                = comlev1, key = ikey_dynamics
359    #  ifdef GM_BOLUS_ADVEC
360    CADJ STORE GM_PsiX            = comlev1, key = ikey_dynamics
361    CADJ STORE GM_PsiY            = comlev1, key = ikey_dynamics
362    #  endif
363    # endif
364    c
365    # ifdef ALLOW_KPP
366    CADJ STORE KPPghat            = comlev1, key = ikey_dynamics
367    CADJ STORE KPPfrac            = comlev1, key = ikey_dynamics
368    CADJ STORE KPPdiffKzS         = comlev1, key = ikey_dynamics
369    CADJ STORE KPPdiffKzT         = comlev1, key = ikey_dynamics
370    # endif
371    c
372    # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
373    cph-adj-test(
374    CADJ STORE theta,salt,wvel          = comlev1, key = ikey_dynamics
375    cph-adj-test)
376    CADJ STORE etaH               = comlev1, key = ikey_dynamics
377    #  ifdef ALLOW_CD_CODE
378    CADJ STORE etanm1             = comlev1, key = ikey_dynamics
379    #  endif
380    #  ifndef DISABLE_RSTAR_CODE
381    cph-test
382    CADJ STORE rstarexpc = comlev1, key = ikey_dynamics
383    #  endif
384    # endif
385    #endif /* ALLOW_AUTODIFF_TAMC */
386    
387          IF ( .NOT.staggerTimeStep ) THEN
388    #ifdef ALLOW_DEBUG
389            IF ( debugLevel .GE. debLevB )
390         &    CALL DEBUG_CALL('THERMODYNAMICS',myThid)
391    #endif
392    cph-adj-test(
393    CADJ STORE salt,wvel          = comlev1, key = ikey_dynamics
394    cph-adj-test)
395            CALL TIMER_START('THERMODYNAMICS      [FORWARD_STEP]',myThid)
396            CALL THERMODYNAMICS( myTime, myIter, myThid )
397            CALL TIMER_STOP ('THERMODYNAMICS      [FORWARD_STEP]',myThid)
398    C--    if not staggerTimeStep: end
399          ENDIF
400    c #ifdef ALLOW_NONHYDROSTATIC
401          IF ( implicitIntGravWave ) THEN
402            CALL TIMER_START('BLOCKING_EXCHANGES  [FORWARD_STEP]',myThid)
403            CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid )
404            CALL TIMER_STOP ('BLOCKING_EXCHANGES  [FORWARD_STEP]',myThid)
405          ENDIF
406    c #endif
407    
408    #ifdef COMPONENT_MODULE
409          IF ( useCoupler .AND. .NOT.cpl_earlyExpImpCall ) THEN
410    C      Post coupling data that I export.
411    C      Read in coupling data that I import.
412             myItP1 = myIter + 1
413             CALL TIMER_START('CPL_EXPORT-IMPORT  [FORWARD_STEP]',myThid)
414             CALL CPL_EXPORT_MY_DATA(       myItP1, myTime, myThid )
415             CALL CPL_IMPORT_EXTERNAL_DATA( myItP1, myTime, myThid )
416             CALL TIMER_STOP ('CPL_EXPORT-IMPORT  [FORWARD_STEP]',myThid)
417    # ifdef ALLOW_OCN_COMPON_INTERF
418            IF ( useRealFreshWaterFlux ) THEN
419             CALL OCN_APPLY_IMPORT( .FALSE., myTime, myIter, myThid )
420            ENDIF
421    # endif /* ALLOW_OCN_COMPON_INTERF */
422          ENDIF
423    #endif /* COMPONENT_MODULE */
424    
425    #ifdef ALLOW_AUTODIFF_TAMC
426    # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
427    CADJ STORE hFacC          = comlev1, key = ikey_dynamics
428    CADJ STORE hFacS          = comlev1, key = ikey_dynamics
429    CADJ STORE hFacW          = comlev1, key = ikey_dynamics
430    CADJ STORE recip_hFacC    = comlev1, key = ikey_dynamics
431    CADJ STORE recip_hFacS    = comlev1, key = ikey_dynamics
432    CADJ STORE recip_hFacW    = comlev1, key = ikey_dynamics
433    CADJ STORE etaN           = comlev1, key = ikey_dynamics
434    c
435    #  ifndef DISABLE_RSTAR_CODE
436    CADJ STORE rstarFacC    = comlev1, key = ikey_dynamics
437    CADJ STORE rstarFacS    = comlev1, key = ikey_dynamics
438    CADJ STORE rstarFacW    = comlev1, key = ikey_dynamics
439    c
440    CADJ STORE h0facc,h0facs,h0facw  
441    CADJ &     = comlev1, key = ikey_dynamics
442    CADJ STORE rstardhcdt,rstardhsdt,rstardhwdt
443    CADJ &     = comlev1, key = ikey_dynamics
444    CADJ STORE rstarexpc,rstarexps,rstarexpw
445    CADJ &     = comlev1, key = ikey_dynamics
446    #  endif
447    # endif
448    #endif
449    C--   Step forward fields and calculate time tendency terms.
450    #ifndef ALLOW_AUTODIFF_TAMC
451          IF ( momStepping ) THEN
452    #endif
453    #ifdef ALLOW_DEBUG
454            IF ( debugLevel .GE. debLevB )
455         &    CALL DEBUG_CALL('DYNAMICS',myThid)
456    #endif
457            CALL TIMER_START('DYNAMICS            [FORWARD_STEP]',myThid)
458            CALL DYNAMICS( myTime, myIter, myThid )
459            CALL TIMER_STOP ('DYNAMICS            [FORWARD_STEP]',myThid)
460    #ifndef ALLOW_AUTODIFF_TAMC
461          ENDIF
462    #endif
463    
464    #ifdef ALLOW_AUTODIFF_TAMC
465    # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
466    cph-test
467    CADJ STORE gU, gV  = comlev1, key = ikey_dynamics
468    # endif
469    #endif
470    
471    C--   Update time-counter
472          myIter = nIter0 + iLoop
473          myTime = startTime + deltaTClock * float(iLoop)
474    
475    #ifdef ALLOW_MNC
476    C     Update the default next iter for MNC
477          IF ( useMNC ) THEN
478             CALL MNC_CW_CITER_SETG( 1, 1, -1, myIter , myThid )
479    
480    C        TODO: Logic should be added here so that users can specify, on
481    C        a per-citer-group basis, when it is time to update the
482    C        "current" (and not just the "next") iteration
483    
484    C        TODO: the following is just a temporary band-aid (mostly, for
485    C        Baylor) until someone writes a routine that better handles time
486    C        boundaries such as weeks, months, years, etc.
487             IF ( mnc_filefreq .GT. 0 ) THEN
488               IF (DIFFERENT_MULTIPLE(mnc_filefreq,myTime,deltaTClock))
489         &          THEN
490                 CALL MNC_CW_CITER_SETG( 1, 1, myIter, -1 , myThid )
491               ENDIF
492             ENDIF
493          ENDIF
494    #endif /* ALLOW_MNC */
495    
496    C--   Update geometric factors:
497    #ifdef NONLIN_FRSURF
498    C-    update hfacC,W,S and recip_hFac according to etaH(n+1) :
499          IF ( nonlinFreeSurf.GT.0) THEN
500           IF ( select_rStar.GT.0 ) THEN
501    # ifndef DISABLE_RSTAR_CODE
502    # ifdef ALLOW_AUTODIFF_TAMC
503    cph-test
504    CADJ STORE hFacC    = comlev1, key = ikey_dynamics
505    CADJ STORE hFacS    = comlev1, key = ikey_dynamics
506    CADJ STORE hFacW    = comlev1, key = ikey_dynamics
507    CADJ STORE recip_hFacC    = comlev1, key = ikey_dynamics
508    CADJ STORE recip_hFacS    = comlev1, key = ikey_dynamics
509    CADJ STORE recip_hFacW    = comlev1, key = ikey_dynamics
510    c
511    CADJ STORE rstarFacC    = comlev1, key = ikey_dynamics
512    CADJ STORE rstarFacS    = comlev1, key = ikey_dynamics
513    CADJ STORE rstarFacW    = comlev1, key = ikey_dynamics
514    c
515    CADJ STORE h0facc,h0facs,h0facw  = comlev1, key = ikey_dynamics
516    # endif
517            CALL TIMER_START('UPDATE_R_STAR       [FORWARD_STEP]',myThid)
518            CALL UPDATE_R_STAR( myTime, myIter, myThid )
519            CALL TIMER_STOP ('UPDATE_R_STAR       [FORWARD_STEP]',myThid)
520    # ifdef ALLOW_AUTODIFF_TAMC
521    cph-test
522    CADJ STORE hFacC    = comlev1, key = ikey_dynamics
523    CADJ STORE hFacS    = comlev1, key = ikey_dynamics
524    CADJ STORE hFacW    = comlev1, key = ikey_dynamics
525    CADJ STORE recip_hFacC    = comlev1, key = ikey_dynamics
526    CADJ STORE recip_hFacS    = comlev1, key = ikey_dynamics
527    CADJ STORE recip_hFacW    = comlev1, key = ikey_dynamics
528    # endif
529    # endif /* DISABLE_RSTAR_CODE */
530           ELSE
531    #ifdef ALLOW_AUTODIFF_TAMC
532    CADJ STORE hFac_surfC, hFac_surfS, hFac_surfW
533    CADJ &     = comlev1, key = ikey_dynamics
534    #endif
535            CALL TIMER_START('UPDATE_SURF_DR      [FORWARD_STEP]',myThid)
536            CALL UPDATE_SURF_DR( myTime, myIter, myThid )
537            CALL TIMER_STOP ('UPDATE_SURF_DR      [FORWARD_STEP]',myThid)
538           ENDIF
539          ENDIF
540    # ifdef ALLOW_AUTODIFF_TAMC
541    cph-test
542    CADJ STORE hFacC    = comlev1, key = ikey_dynamics
543    CADJ STORE hFacS    = comlev1, key = ikey_dynamics
544    CADJ STORE hFacW    = comlev1, key = ikey_dynamics
545    CADJ STORE recip_hFacC    = comlev1, key = ikey_dynamics
546    CADJ STORE recip_hFacS    = comlev1, key = ikey_dynamics
547    CADJ STORE recip_hFacW    = comlev1, key = ikey_dynamics
548    # endif
549    C-    update also CG2D matrix (and preconditioner)
550          IF ( momStepping .AND. nonlinFreeSurf.GT.2 ) THEN
551            CALL TIMER_START('UPDATE_CG2D         [FORWARD_STEP]',myThid)
552            CALL UPDATE_CG2D( myTime, myIter, myThid )
553            CALL TIMER_STOP ('UPDATE_CG2D         [FORWARD_STEP]',myThid)
554          ENDIF
555    #endif /* NONLIN_FRSURF */
556    
557        IF (.NOT. doHalfStep) THEN  C--   Apply Filters to u*,v* before SOLVE_FOR_PRESSURE
558    #ifdef ALLOW_SHAP_FILT
559          IF (useSHAP_FILT .AND. shap_filt_uvStar) THEN
560            CALL TIMER_START('SHAP_FILT_UV        [FORWARD_STEP]',myThid)
561            IF (implicDiv2Dflow.LT.1.) THEN
562    C--   Explicit+Implicit part of the Barotropic Flow Divergence
563    C      => Filtering of uVel,vVel is necessary
564              CALL SHAP_FILT_APPLY_UV( uVel,vVel,
565         &                             myTime, myIter, myThid )
566            ENDIF
567            CALL SHAP_FILT_APPLY_UV( gU,gV,myTime,myIter,myThid)
568            CALL TIMER_STOP ('SHAP_FILT_UV        [FORWARD_STEP]',myThid)
569          ENDIF
570    #endif
571    #ifdef ALLOW_ZONAL_FILT
572          IF (useZONAL_FILT .AND. zonal_filt_uvStar) THEN
573            CALL TIMER_START('ZONAL_FILT_UV       [FORWARD_STEP]',myThid)
574            IF (implicDiv2Dflow.LT.1.) THEN
575    C--   Explicit+Implicit part of the Barotropic Flow Divergence
576    C      => Filtering of uVel,vVel is necessary
577              CALL ZONAL_FILT_APPLY_UV( uVel, vVel, myThid )
578            ENDIF
579            CALL ZONAL_FILT_APPLY_UV( gU, gV, myThid )
580            CALL TIMER_STOP ('ZONAL_FILT_UV       [FORWARD_STEP]',myThid)
581          ENDIF
582    #endif
583    
584  C--   Solve elliptic equation(s).  C--   Solve elliptic equation(s).
585  C     Two-dimensional only for conventional hydrostatic or  C     Two-dimensional only for conventional hydrostatic or
586  C     three-dimensional for non-hydrostatic and/or IGW scheme.  C     three-dimensional for non-hydrostatic and/or IGW scheme.
587        CALL TIMER_START('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)        IF ( momStepping ) THEN
588        CALL SOLVE_FOR_PRESSURE( myThid )  #ifdef ALLOW_AUTODIFF_TAMC
589        CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)  # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
590    CADJ STORE uvel, vvel
591  C--   Correct divergence in flow field and cycle time-stepping  CADJ &     = comlev1, key = ikey_dynamics
592  C     arrays (for all fields)  CADJ STORE empmr,hfacs,hfacw
593        CALL THE_CORRECTION_STEP(myCurrentTime, myCurrentIter, myThid)  CADJ &     = comlev1, key = ikey_dynamics
594    # endif
595    #endif
596            CALL TIMER_START('SOLVE_FOR_PRESSURE  [FORWARD_STEP]',myThid)
597            CALL SOLVE_FOR_PRESSURE(myTime, myIter, myThid)
598            CALL TIMER_STOP ('SOLVE_FOR_PRESSURE  [FORWARD_STEP]',myThid)
599          ENDIF
600    
601    C--   Correct divergence in flow field and cycle time-stepping momentum
602    #ifndef ALLOW_AUTODIFF_TAMC
603          IF ( momStepping ) THEN
604    #endif
605    #ifdef ALLOW_AUTODIFF_TAMC
606    # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
607    #  ifndef DISABLE_RSTAR_CODE
608    cph-test
609    cph not clear, why this one
610    CADJ STORE h0facc = comlev1, key = ikey_dynamics
611    #  endif
612    # endif
613    # ifdef ALLOW_DEPTH_CONTROL
614    CADJ STORE etan, uvel,vvel
615    CADJ &     = comlev1, key = ikey_dynamics
616    # endif
617    #endif
618            CALL TIMER_START('MOM_CORRECTION_STEP [FORWARD_STEP]',myThid)
619            CALL MOMENTUM_CORRECTION_STEP(myTime, myIter, myThid)
620            CALL TIMER_STOP ('MOM_CORRECTION_STEP [FORWARD_STEP]',myThid)
621    #ifndef ALLOW_AUTODIFF_TAMC
622          ENDIF
623    #endif
624    
625    #ifdef EXACT_CONSERV
626          IF (exactConserv) THEN
627    C--   Update etaH(n+1) :
628            CALL TIMER_START('UPDATE_ETAH         [FORWARD_STEP]',myThid)
629            CALL UPDATE_ETAH( myTime, myIter, myThid )
630            CALL TIMER_STOP ('UPDATE_ETAH         [FORWARD_STEP]',myThid)
631          ENDIF
632    #endif /* EXACT_CONSERV */
633    
634    #ifdef NONLIN_FRSURF
635          IF ( select_rStar.NE.0 ) THEN
636    # ifndef DISABLE_RSTAR_CODE
637    #  ifdef ALLOW_AUTODIFF_TAMC
638    cph-test
639    CADJ STORE rstarfacc,rstarfacs,rstarfacw = comlev1, key = ikey_dynamics
640    #  endif
641    C--   r* : compute the future level thickness according to etaH(n+1)
642            CALL TIMER_START('CALC_R_STAR         [FORWARD_STEP]',myThid)
643            CALL CALC_R_STAR(etaH, myTime, myIter, myThid )
644            CALL TIMER_STOP ('CALC_R_STAR         [FORWARD_STEP]',myThid)
645    # endif /* DISABLE_RSTAR_CODE */
646          ELSEIF ( nonlinFreeSurf.GT.0) THEN
647    C--   compute the future surface level thickness according to etaH(n+1)
648    # ifdef ALLOW_AUTODIFF_TAMC
649    CADJ STORE etaH          = comlev1, key = ikey_dynamics
650    # endif
651            CALL TIMER_START('CALC_SURF_DR      [FORWARD_STEP]',myThid)
652            CALL CALC_SURF_DR(etaH, myTime, myIter, myThid )
653            CALL TIMER_STOP ('CALC_SURF_DR      [FORWARD_STEP]',myThid)
654          ENDIF
655    # ifdef ALLOW_AUTODIFF_TAMC
656    CADJ STORE hFac_surfC    = comlev1, key = ikey_dynamics
657    cph-adj-test(
658    CADJ STORE salt,theta,vvel  = comlev1, key = ikey_dynamics
659    cph-adj-test)
660    # endif
661    #endif /* NONLIN_FRSURF */
662    
663    C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
664          IF ( staggerTimeStep ) THEN
665    C--   do exchanges of U,V (needed for multiDim) when using stagger time-step :
666    #ifdef ALLOW_DEBUG
667            IF ( debugLevel .GE. debLevB )
668         &   CALL DEBUG_CALL('DO_STAGGER_FIELDS_EXCH.',myThid)
669    #endif
670            CALL TIMER_START('BLOCKING_EXCHANGES  [FORWARD_STEP]',myThid)
671            CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid )
672            CALL TIMER_STOP ('BLOCKING_EXCHANGES  [FORWARD_STEP]',myThid)
673    
674    #ifdef ALLOW_DIAGNOSTICS
675    C--   State-variables diagnostics
676            IF ( useDiagnostics ) THEN
677              CALL TIMER_START('DO_STATEVARS_DIAGS  [FORWARD_STEP]',myThid)
678              CALL DO_STATEVARS_DIAGS( myTime, 1, myIter, myThid )
679              CALL TIMER_STOP ('DO_STATEVARS_DIAGS  [FORWARD_STEP]',myThid)
680            ENDIF
681    #endif
682    
683    #ifdef ALLOW_DEBUG
684            IF ( debugLevel .GE. debLevB )
685         &   CALL DEBUG_CALL('THERMODYNAMICS',myThid)
686    #endif
687            CALL TIMER_START('THERMODYNAMICS      [FORWARD_STEP]',myThid)
688            CALL THERMODYNAMICS( myTime, myIter, myThid )
689            CALL TIMER_STOP ('THERMODYNAMICS      [FORWARD_STEP]',myThid)
690    
691    C--    if staggerTimeStep: end
692          ENDIF
693    C---+--------+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
694    
695    #ifdef ALLOW_AUTODIFF_TAMC
696    cph This is needed because convective_adjustment calls
697    cph find_rho which may use pressure()
698    CADJ STORE totphihyd  = comlev1, key = ikey_dynamics
699    #endif
700    C--   Cycle time-stepping Tracers arrays (T,S,+pTracers)
701          CALL TIMER_START('TRC_CORRECTION_STEP [FORWARD_STEP]',myThid)
702          CALL TRACERS_CORRECTION_STEP(myTime, myIter, myThid)
703          CALL TIMER_STOP ('TRC_CORRECTION_STEP [FORWARD_STEP]',myThid)
704    
705    #ifdef ALLOW_GCHEM
706    C     Add separate timestepping of chemical/biological/forcing
707    C     of ptracers here in GCHEM_FORCING_SEP
708    #ifdef ALLOW_AUTODIFF_TAMC
709    CADJ STORE ptracer  = comlev1, key = ikey_dynamics
710    CADJ STORE theta  = comlev1, key = ikey_dynamics
711    CADJ STORE salt  = comlev1, key = ikey_dynamics
712    #endif
713          IF ( useGCHEM ) THEN
714    #ifdef ALLOW_DEBUG
715             IF ( debugLevel .GE. debLevB )
716         &    CALL DEBUG_CALL('GCHEM_FORCING_SEP',myThid)
717    #endif /* ALLOW_DEBUG */
718             CALL TIMER_START('GCHEM_FORCING_SEP  [FORWARD_STEP]',myThid)
719             CALL GCHEM_FORCING_SEP( myTime,myIter,myThid )
720             CALL TIMER_STOP ('GCHEM_FORCING_SEP  [FORWARD_STEP]',myThid)
721          ENDIF
722    #endif /* ALLOW_GCHEM */
723    
724  C--   Do "blocking" sends and receives for tendency "overlap" terms  C--   Do "blocking" sends and receives for tendency "overlap" terms
725  C     CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)  c     CALL TIMER_START('BLOCKING_EXCHANGES  [FORWARD_STEP]',myThid)
726  C     CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )  c     CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )
727  C     CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)  c     CALL TIMER_STOP ('BLOCKING_EXCHANGES  [FORWARD_STEP]',myThid)
728    
729  C--   Do "blocking" sends and receives for field "overlap" terms  C--   Do "blocking" sends and receives for field "overlap" terms
730        CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)        CALL TIMER_START('BLOCKING_EXCHANGES  [FORWARD_STEP]',myThid)
731        CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )        CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
732        CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)        CALL TIMER_STOP ('BLOCKING_EXCHANGES  [FORWARD_STEP]',myThid)
733    
734  c     myCurrentIter = myCurrentIter + 1  #ifdef ALLOW_DIAGNOSTICS
735  c     myCurrentTime = myCurrentTime + deltaTClock        IF ( useDiagnostics ) THEN
736        myCurrentIter = nIter0 + iLoop         CALL TIMER_START('DO_STATEVARS_DIAGS  [FORWARD_STEP]',myThid)
737        myCurrentTime = startTime + deltaTClock * float(iLoop)         CALL DO_STATEVARS_DIAGS( myTime, 2, myIter, myThid )
738           CALL TIMER_STOP ('DO_STATEVARS_DIAGS  [FORWARD_STEP]',myThid)
739          ENDIF
740    #endif
741    
742  C--   Save state for restarts  #ifdef ALLOW_GRIDALT
743  C     Note:        IF (useGRIDALT) THEN
744  C     =====           CALL GRIDALT_UPDATE(myThid)
745  C     Because of the ordering of the timestepping code and        ENDIF
746  C     tendency term code at end of loop model arrays hold  #endif
747  C     U,V,T,S  at "time-level" N but gu, gv, gs, gt, guNM1,...  
748  C     at "time-level" N+1/2 (guNM1 at "time-level" N+1/2 is  #ifdef ALLOW_FIZHI
749  C     gu at "time-level" N-1/2) and cg2d_x at "time-level" N+1/2.        IF (useFIZHI) THEN
750  C      where N = I+timeLevBase-1           CALL TIMER_START('FIZHI               [FORWARD_STEP]',myThid)
751  C     Thus a checkpoint contains U.0000000000, GU.0000000001 and           CALL STEP_FIZHI_CORR ( myTime, myIter, myThid, dTtracerLev(1) )
752  C     cg2d_x.0000000001 in the indexing scheme used for the model           CALL TIMER_STOP ('FIZHI               [FORWARD_STEP]',myThid)
753  C     "state" files. This example is referred to as a checkpoint        ENDIF
754  C     at time level 1  #endif
       CALL TIMER_START('I/O (WRITE)        [FORWARD_STEP]',myThid)  
       CALL WRITE_CHECKPOINT(  
      &        .FALSE., myCurrentTime, myCurrentIter, myThid )  
       CALL TIMER_STOP ('I/O (WRITE)        [FORWARD_STEP]',myThid)  
755    
756    #ifdef ALLOW_FLT
757    C--   Calculate float trajectories
758          IF (useFLT) THEN
759            CALL TIMER_START('FLOATS            [FORWARD_STEP]',myThid)
760            CALL FLT_MAIN(myIter,myTime, myThid)
761            CALL TIMER_STOP ('FLOATS            [FORWARD_STEP]',myThid)
762        ENDIF        ENDIF
763    #endif
764    
765    #ifdef ALLOW_TIMEAVE
766    C--   State-variables time-averaging
767          CALL TIMER_START('DO_STATEVARS_TAVE   [FORWARD_STEP]',myThid)
768          CALL DO_STATEVARS_TAVE( myTime, myIter, myThid )
769          CALL TIMER_STOP ('DO_STATEVARS_TAVE   [FORWARD_STEP]',myThid)
770    #endif
771    
772    #ifdef ALLOW_MONITOR
773          IF ( .NOT.useOffLine ) THEN
774    C--   Check status of solution (statistics, cfl, etc...)
775            CALL TIMER_START('MONITOR             [FORWARD_STEP]',myThid)
776            CALL MONITOR( myIter, myTime, myThid )
777            CALL TIMER_STOP ('MONITOR             [FORWARD_STEP]',myThid)
778          ENDIF
779    #endif /* ALLOW_MONITOR */
780    
781    #ifdef ALLOW_COST
782    C--     compare model with data and compute cost function
783    C--     this is done after exchanges to allow interpolation
784          CALL TIMER_START('COST_TILE           [FORWARD_STEP]',myThid)
785          CALL COST_TILE  ( myTime, myIter, myThid )
786          CALL TIMER_STOP ('COST_TILE           [FORWARD_STEP]',myThid)
787    #endif
788    
789    C--   Do IO if needed.
790          CALL TIMER_START('DO_THE_MODEL_IO     [FORWARD_STEP]',myThid)
791          CALL DO_THE_MODEL_IO( myTime, myIter, myThid )
792          CALL TIMER_STOP ('DO_THE_MODEL_IO     [FORWARD_STEP]',myThid)
793    
794          modelEnd = myTime.EQ.endTime .OR. myIter.EQ.nEndIter
795    #ifdef HAVE_SIGREG
796          IF ( useSIGREG ) THEN
797            modelEnd = modelEnd .OR. ( i_got_signal.GT.0 )
798          ENDIF
799    #endif /* HAVE_SIGREG */
800    
801    C--   Save state for restarts
802          CALL TIMER_START('DO_WRITE_PICKUP     [FORWARD_STEP]',myThid)
803          CALL DO_WRITE_PICKUP(
804         I                      modelEnd, myTime, myIter, myThid )
805          CALL TIMER_STOP ('DO_WRITE_PICKUP     [FORWARD_STEP]',myThid)
806    
807    #ifdef HAVE_SIGREG
808          IF ( useSIGREG ) THEN
809            IF ( modelEnd .AND. i_got_signal.GT.0 ) THEN
810              STOP 'Checkpoint completed -- killed by signal handler'
811            ENDIF
812          ENDIF
813    #endif /* HAVE_SIGREG */
814    
815    #ifdef ALLOW_AUTODIFF_TAMC
816          CALL AUTODIFF_INADMODE_SET( myThid )
817    #endif
818    
819    #ifdef ALLOW_SHOWFLOPS
820          CALL TIMER_START('SHOWFLOPS_INLOOP   [THE_MAIN_LOOP]', mythid)
821          CALL SHOWFLOPS_INLOOP( iloop, mythid )
822          CALL TIMER_STOP ('SHOWFLOPS_INLOOP   [THE_MAIN_LOOP]', mythid)
823    #endif
824    
825    #ifdef ALLOW_DEBUG
826          IF ( debugLevel .GE. debLevB )
827         &    CALL DEBUG_LEAVE('FORWARD_STEP',myThid)
828    #endif
829    
830        RETURN        RETURN
831        END        END

Legend:
Removed from v.1.6  
changed lines
  Added in v.1.157

  ViewVC Help
Powered by ViewVC 1.1.22