/[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.8 by jmc, Tue Feb 20 15:10:15 2001 UTC revision 1.127 by heimbach, Thu Dec 8 15:44:34 2005 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_OFFLINE
8       I                    iLoop,  # include "OFFLINE_OPTIONS.h"
9       U                    myCurrentTime, myCurrentIter,  #endif
10       &                    myThid)  #ifdef ALLOW_GMREDI
11  C     /==========================================================\  # include "GMREDI_OPTIONS.h"
12  C     | SUBROUTINE FORWARD_STEP                                  |  #endif
13  C     | o Does one instance of the model time stepping           |  
14  C     |   The time stepping loop in THE_MAIN_LOOP() calls        |  CBOP
15  C     |   this routine                                           |  C     !ROUTINE: FORWARD_STEP
16  C     |==========================================================|  C     !INTERFACE:
17  C     \==========================================================/        SUBROUTINE FORWARD_STEP( iloop, myTime, myIter, myThid )
18        IMPLICIT NONE  
19  C  C     !DESCRIPTION: \bv
20  C     Call Tree  C     *==================================================================
21  C     =========  C     | SUBROUTINE forward_step
22  C      C     | o Run the ocean model and, optionally, evaluate a cost function.
23  C      THE_MAIN_LOOP()  C     *==================================================================
24  C       |  C     |
25  C  ==>  | ** Time stepping loop starts here **  C     | THE_MAIN_LOOP is the toplevel routine for the Tangent Linear and
26  C  |    |  C     | Adjoint Model Compiler (TAMC). For this purpose the initialization
27  C /|\   |-FORWARD_STEP  C     | of the model was split into two parts. Those parameters that do
28  C  |    |  |  C     | not depend on a specific model run are set in INITIALISE_FIXED,  
29  C /|\   |  |--LOAD_EXTERNAL_DATA  C     | whereas those that do depend on the specific realization are
30  C  |    |  |   o Load and/or set time dependent forcing fields  C     | initialized in INITIALISE_VARIA.  
31  C /|\   |  |  C     |
32  C  |    |  |--DYNAMICS  C     *==================================================================
33  C /|\   |  |   o Evaluate "forward" terms  C     \ev
 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  |    |  |--DO_THE_MODEL_IO  
 C /|\   |  |   o Write model state  
 C  |    |  |  
 C /|\   |  |  
 C  |    |  |--WRITE_CHECKPOINT  
 C /|\   |  |   o Write restart file(s)  
 C  |    |  
 C /|\   |  
 C  |<== | ** Time stepping loop finishes here **  
 C  
34    
35  C     == Global variables ===  C     !USES:
36          IMPLICIT NONE
37    C     == Global variables ==
38  #include "SIZE.h"  #include "SIZE.h"
39  #include "EEPARAMS.h"  #include "EEPARAMS.h"
40  #include "PARAMS.h"  #include "PARAMS.h"
41  #include "DYNVARS.h"  #include "DYNVARS.h"
42  #include "CG2D.h"  
43  #ifdef ALLOW_NONHYDROSTATIC  #ifdef ALLOW_MNC
44  #include "CG3D.h"  #include "MNC_PARAMS.h"
45          EXTERNAL DIFFERENT_MULTIPLE
46          LOGICAL  DIFFERENT_MULTIPLE
47    #endif
48    
49    #ifdef HAVE_SIGREG
50    #include "SIGREG.h"
51    #endif
52    
53    #ifdef ALLOW_SHAP_FILT
54    # include "SHAP_FILT.h"
55    #endif
56    #ifdef ALLOW_ZONAL_FILT
57    # include "ZONAL_FILT.h"
58    #endif
59    #ifdef COMPONENT_MODULE
60    # include "CPL_PARAMS.h"
61  #endif  #endif
62    
63    #ifdef ALLOW_AUTODIFF_TAMC
64    # include "FFIELDS.h"
65    
66    # include "tamc.h"
67    # include "ctrl.h"
68    # include "ctrl_dummy.h"
69    # include "cost.h"
70    # include "EOS.h"
71    # ifdef NONLIN_FRSURF
72    #  include "GRID.h"
73    # endif
74    # ifdef ALLOW_EXF
75    #  include "exf_fields.h"
76    #  include "exf_clim_fields.h"
77    #  ifdef ALLOW_BULKFORMULAE
78    #   include "exf_constants.h"
79    #  endif
80    # endif
81    # ifdef ALLOW_OBCS
82    #  include "OBCS.h"
83    # endif
84    # ifdef ALLOW_PTRACERS
85    #  include "PTRACERS_SIZE.h"
86    #  include "PTRACERS.h"
87    # endif
88    # ifdef ALLOW_CD_CODE
89    #  include "CD_CODE_VARS.h"
90    # endif
91    # ifdef ALLOW_EBM
92    #  include "EBM.h"
93    # endif
94    # ifdef EXACT_CONSERV
95    #  include "SURFACE.h"
96    # endif
97    # ifdef ALLOW_KPP
98    #  include "KPP.h"
99    # endif
100    # ifdef ALLOW_GMREDI
101    #  include "GMREDI.h"
102    # endif
103    #endif /* ALLOW_AUTODIFF_TAMC */
104    
105    C     !LOCAL VARIABLES:
106  C     == Routine arguments ==  C     == Routine arguments ==
107  C     iLoop         - Invocation count (counter in THE_MAIN_LOOP)  C     note: under the multi-threaded model myiter and
108  C     myCurrentIter - Iteration counter for this thread  C           mytime are local variables passed around as routine
109  C     myCurrentTime - Time counter for this thread  C           arguments. Although this is fiddly it saves the need to
110  C     myThid - Thread number for this instance of the routine.  C           impose additional synchronisation points when they are
111        INTEGER iLoop  C           updated.
112        INTEGER myCurrentIter  C     myIter - iteration counter for this thread
113        _RL     myCurrentTime  C     myTime - time counter for this thread
114        INTEGER myThid    C     myThid - thread number for this instance of the routine.
115          INTEGER iloop
116          INTEGER myThid
117          INTEGER myIter
118          _RL     myTime
119    
120  C     == Local variables ==  C     == Local variables ==
121    #ifdef COMPONENT_MODULE
122          INTEGER myItP1
123    #endif
124    CEOP
125    
126    #ifdef ALLOW_DEBUG
127          IF ( debugLevel .GE. debLevB )
128         &    CALL DEBUG_ENTER('FORWARD_STEP',myThid)
129    #endif
130    
131    #ifdef ALLOW_AUTODIFF_TAMC
132    C--   Reset the model iteration counter and the model time.
133          myIter = nIter0 + (iloop-1)
134          myTime = startTime + float(iloop-1)*deltaTclock
135    #endif
136    
137  C--   Load forcing/external data fields  #ifdef ALLOW_AUTODIFF_TAMC
138        CALL TIMER_START('I/O (READ)         [FORWARD_STEP]',myThid)  c**************************************
139        CALL EXTERNAL_FIELDS_LOAD( myCurrentTime, myCurrentIter, myThid )  #include "checkpoint_lev1_directives.h"
140        CALL TIMER_STOP ('I/O (READ)         [FORWARD_STEP]',myThid)  c**************************************
141    #endif
142  C--   Step forward fields and calculate time tendency terms  
143        CALL TIMER_START('DYNAMICS           [FORWARD_STEP]',myThid)  C--   Switch on/off diagnostics for snap-shot output:
144        CALL DYNAMICS( myCurrentTime, myCurrentIter, myThid )  #ifdef ALLOW_DIAGNOSTICS
145        CALL TIMER_STOP ('DYNAMICS           [FORWARD_STEP]',myThid)        IF ( useDiagnostics ) THEN
146            CALL DIAGNOSTICS_SWITCH_ONOFF( myTime, myIter, myThid )
 #ifdef ALLOW_NONHYDROSTATIC  
 C--   Step forward W field in N-H algorithm  
       IF ( nonHydrostatic ) THEN  
        CALL TIMER_START('CALC_GW            [FORWARD_STEP]',myThid)  
        CALL CALC_GW( myThid)  
        CALL TIMER_STOP ('CALC_GW            [FORWARD_STEP]',myThid)  
147        ENDIF        ENDIF
148  #endif  #endif
149    
150  #ifdef ALLOW_SHAP_FILT  C--   State-variables diagnostics
151  C--   Step forward all tiles, filter and exchange.        IF ( usediagnostics ) THEN
152        CALL TIMER_START('SHAP_FILT          [FORWARD_STEP]',myThid)          CALL TIMER_START('DO_STATEVARS_DIAGS  [FORWARD_STEP]',myThid)
153        CALL SHAP_FILT_APPLY(          CALL DO_STATEVARS_DIAGS( myTime, 0, myIter, myThid )
154       I                     gUnm1, gVnm1, gTnm1, gSnm1,          CALL TIMER_STOP ('DO_STATEVARS_DIAGS  [FORWARD_STEP]',myThid)
155       I                     myCurrentTime, myCurrentIter, myThid )        ENDIF
156        CALL TIMER_STOP ('SHAP_FILT          [FORWARD_STEP]',myThid)  
157    C--   Call external forcing package
158    #ifdef ALLOW_BULK_FORCE
159          IF ( useBulkForce ) THEN
160    #ifdef ALLOW_DEBUG
161           IF ( debugLevel .GE. debLevB )
162         &    CALL DEBUG_CALL('BULKF_FIELDS_LOAD',myThid)
163    #endif
164           CALL TIMER_START('BULKF_FORCING      [FORWARD_STEP]',mythid)
165    C-    load all forcing fields at current time
166           CALL BULKF_FIELDS_LOAD( myTime, myIter, myThid )
167    C-    calculate qnet and empmr (and wind stress)
168           CALL BULKF_FORCING( myTime, myIter, myThid )
169           CALL TIMER_STOP ('BULKF_FORCING      [FORWARD_STEP]',mythid)
170          ELSE
171    #endif /* ALLOW_BULK_FORCE */
172    
173    # ifdef ALLOW_EXF
174    #  ifdef ALLOW_DEBUG
175          IF ( debugLevel .GE. debLevB )
176         &    CALL DEBUG_CALL('EXF_GETFORCING',myThid)
177    #  endif
178          CALL TIMER_START('EXF_GETFORCING     [FORWARD_STEP]',mythid)
179          CALL EXF_GETFORCING( mytime, myiter, mythid )
180          CALL TIMER_STOP ('EXF_GETFORCING     [FORWARD_STEP]',mythid)
181    # else /* ALLOW_EXF undef */
182    cph The following IF-statement creates an additional dependency
183    cph for the forcing fields requiring additional storing.
184    cph Therefore, the IF-statement will be put between CPP-OPTIONS,
185    cph assuming that ALLOW_SEAICE has not yet been differentiated.
186    #  if (defined (ALLOW_SEAICE) || defined (ALLOW_EBM))
187          IF ( .NOT. useSEAICE .AND. .NOT. useEBM ) THEN
188    #  endif
189    #ifdef ALLOW_DEBUG
190           IF ( debugLevel .GE. debLevB )
191         &    CALL DEBUG_CALL('EXTERNAL_FIELDS_LOAD',myThid)
192  #endif  #endif
193           CALL TIMER_START('EXTERNAL_FIELDS_LOAD[FORWARD_STEP]',mythid)
194           CALL EXTERNAL_FIELDS_LOAD( mytime, myiter, mythid )
195           CALL TIMER_STOP ('EXTERNAL_FIELDS_LOAD[FORWARD_STEP]',mythid)
196    #  if (defined (ALLOW_SEAICE) || defined (ALLOW_EBM))
197          ENDIF
198    #  endif
199    # endif /* ALLOW_EXF */
200    #ifdef ALLOW_BULK_FORCE
201    C--   end of if/else block useBulfforce --
202          ENDIF
203    #endif /* ALLOW_BULK_FORCE */
204    
205  #ifdef ALLOW_ZONAL_FILT  #ifdef ALLOW_AUTODIFF
206        IF (zonal_filt_lat.LT.90.) THEN  c--   Add control vector for forcing and parameter fields
207        CALL ZONAL_FILT_APPLY(        if ( myiter .EQ. nIter0 )
208       U           gUnm1, gVnm1, gTnm1, gSnm1,       &     CALL CTRL_MAP_FORCING (mythid)
209       I           myThid )  #endif
210    
211    #if (defined (ALLOW_AUTODIFF_TAMC) && defined (ALLOW_AUTODIFF_MONITOR))
212    C     Include call to a dummy routine. Its adjoint will be
213    C     called at the proper place in the adjoint code.
214    C     The adjoint routine will print out adjoint values
215    C     if requested. The location of the call is important,
216    C     it has to be after the adjoint of the exchanges
217    C     (DO_GTERM_BLOCKING_EXCHANGES).
218          CALL DUMMY_IN_STEPPING( myTime, myIter, myThid )
219    cph   I've commented this line since it may conflict with MITgcm's adjoint
220    cph   However, need to check whether that's still consistent
221    cph   with the ecco-branch (it should).
222    cph      CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
223    #endif
224    
225    # ifdef ALLOW_SEAICE
226    C--   Call sea ice model to compute forcing/external data fields.  In
227    C     addition to computing prognostic sea-ice variables and diagnosing the
228    C     forcing/external data fields that drive the ocean model, SEAICE_MODEL
229    C     also sets theta to the freezing point under sea-ice.  The implied
230    C     surface heat flux is then stored in variable surfaceTendencyTice,
231    C     which is needed by KPP package (kpp_calc.F and kpp_transport_t.F)
232    C     to diagnose surface buoyancy fluxes and for the non-local transport
233    C     term.  Because this call precedes model thermodynamics, temperature
234    C     under sea-ice may not be "exactly" at the freezing point by the time
235    C     theta is dumped or time-averaged.
236          IF ( useSEAICE ) THEN
237    #ifdef ALLOW_DEBUG
238             IF ( debugLevel .GE. debLevB )
239         &    CALL DEBUG_CALL('SEAICE_MODEL',myThid)
240    #endif
241             CALL TIMER_START('SEAICE_MODEL       [FORWARD_STEP]',myThid)
242             CALL SEAICE_MODEL( myTime, myIter, myThid )
243             CALL TIMER_STOP ('SEAICE_MODEL       [FORWARD_STEP]',myThid)
244          ENDIF
245    # endif /* ALLOW_SEAICE */
246    
247    #ifdef ALLOW_AUTODIFF_TAMC
248    # ifdef ALLOW_PTRACERS
249    cph this replaces _bibj storing of ptracer within thermodynamics
250    CADJ STORE ptracer  = comlev1, key = ikey_dynamics
251    # endif
252    #endif
253    
254    #ifdef ALLOW_OFFLINE
255            call OFFLINE_FIELDS_LOAD( myTime, myIter, myThid )
256    #endif
257    
258    #ifdef ALLOW_PTRACERS
259    # ifdef ALLOW_GCHEM
260            IF ( useGCHEM ) THEN
261    #ifdef ALLOW_DEBUG
262             IF ( debugLevel .GE. debLevB )
263         &        CALL DEBUG_CALL('GCHEM_FIELDS_LOAD',myThid)
264    #endif /* ALLOW_DEBUG */
265             CALL GCHEM_FIELDS_LOAD( mytime, myiter, mythid )
266            ENDIF
267    # endif
268    #endif
269    
270    #ifdef ALLOW_RBCS
271            IF ( useRBCS ) THEN
272             CALL RBCS_FIELDS_LOAD( mytime, myiter, mythid )
273            ENDIF
274    #endif
275    
276    #ifdef COMPONENT_MODULE
277           IF ( useCoupler .AND. cpl_earlyExpImpCall ) THEN
278    C      Post coupling data that I export.
279    C      Read in coupling data that I import.
280             CALL TIMER_START('CPL_EXPORT-IMPORT  [FORWARD_STEP]',myThid)
281             CALL CPL_EXPORT_MY_DATA(       myIter, myTime, myThid )
282             CALL CPL_IMPORT_EXTERNAL_DATA( myIter, myTime, myThid )
283             CALL TIMER_STOP ('CPL_EXPORT-IMPORT  [FORWARD_STEP]',myThid)
284           ENDIF
285    #endif /* COMPONENT_MODULE */
286    
287    #ifdef ALLOW_EBM
288            IF ( useEBM ) THEN
289    # ifdef ALLOW_DEBUG
290             IF ( debugLevel .GE. debLevB )
291         &    CALL DEBUG_CALL('EBM',myThid)
292    # endif
293             CALL TIMER_START('EBM                [FORWARD_STEP]',mythid)
294             CALL EBM_DRIVER ( myTime, myIter, myThid )
295             CALL TIMER_STOP ('EBM                [FORWARD_STEP]',mythid)
296            ENDIF
297    #endif
298    
299    C--     Step forward fields and calculate time tendency terms.
300    
301    #ifdef ALLOW_DEBUG
302           IF ( debugLevel .GE. debLevB )
303         &    CALL DEBUG_CALL('DO_ATMOSPHERIC_PHYS',myThid)
304    #endif
305           CALL TIMER_START('DO_ATMOSPHERIC_PHYS [FORWARD_STEP]',mythid)
306           CALL DO_ATMOSPHERIC_PHYS( myTime, myIter, myThid )
307           CALL TIMER_STOP ('DO_ATMOSPHERIC_PHYS [FORWARD_STEP]',mythid)
308    
309    #ifdef ALLOW_AUTODIFF_TAMC
310    CADJ STORE theta              = comlev1, key = ikey_dynamics
311    CADJ STORE salt               = comlev1, key = ikey_dynamics
312    CADJ STORE totphihyd          = comlev1, key = ikey_dynamics
313    CADJ STORE surfaceforcingtice = comlev1, key = ikey_dynamics
314    # ifdef ALLOW_KPP
315    CADJ STORE uvel               = comlev1, key = ikey_dynamics
316    CADJ STORE vvel               = comlev1, key = ikey_dynamics
317    # endif
318    # ifdef EXACT_CONSERV
319    CADJ STORE empmr              = comlev1, key = ikey_dynamics
320    CADJ STORE pmepr              = comlev1, key = ikey_dynamics
321    # endif
322    # ifdef NONLIN_FRSURF
323    cph-test
324    CADJ STORE hFac_surfC         = comlev1, key = ikey_dynamics
325    CADJ STORE hfac_surfs         = comlev1, key = ikey_dynamics
326    CADJ STORE hfac_surfw         = comlev1, key = ikey_dynamics
327    CADJ STORE hFacC, hFacS, hFacW
328    CADJ &     = comlev1, key = ikey_dynamics
329    CADJ STORE recip_hFacC, recip_hFacS, recip_hFacW
330    CADJ &     = comlev1, key = ikey_dynamics
331    # endif
332    #endif /* ALLOW_AUTODIFF_TAMC */
333    
334    #ifndef ALLOW_OFFLINE
335    #ifdef ALLOW_DEBUG
336           IF ( debugLevel .GE. debLevB )
337         &    CALL DEBUG_CALL('DO_OCEANIC_PHYS',myThid)
338    #endif
339           CALL TIMER_START('DO_OCEANIC_PHYS     [FORWARD_STEP]',mythid)
340           CALL DO_OCEANIC_PHYS( myTime, myIter, myThid )
341           CALL TIMER_STOP ('DO_OCEANIC_PHYS     [FORWARD_STEP]',mythid)
342    #endif
343    
344    #ifdef ALLOW_GCHEM
345    C     GCHEM package is an interface for any bio-geochemical or
346    C     ecosystem model you would like to include.
347    C     If GCHEM_SEPARATE_FORCING is not defined, you are
348    C     responsible for computing tendency terms for passive
349    C     tracers and storing them on a 3DxNumPtracers-array called
350    C     gchemTendency in GCHEM_CALC_TENDENCY. This tendency is then added
351    C     to gPtr in ptracers_forcing later-on.
352    C     If GCHEM_SEPARATE_FORCING is defined, you are reponsible for
353    C     UPDATING ptracers directly in GCHEM_FORCING_SEP. This amounts
354    C     to a completely separate time step that you have to implement
355    C     yourself (Eulerian seems to be fine in most cases).
356    CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
357    C     CAVEAT: Up to now, when GCHEM is turned on the field ptracerForcingSurf,
358    C     which is needed for KPP is not set properly. ptracerForcingSurf must
359    C     be treated differently depending on whether GCHEM_SEPARATE_FORCING
360    C     is define or not. TBD.
361    CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
362           IF ( useGCHEM ) THEN
363    #ifdef ALLOW_DEBUG
364            IF ( debugLevel .GE. debLevB )
365         &       CALL DEBUG_CALL('GCHEM_CALC_TENDENCY',myThid)
366    #endif
367            CALL TIMER_START('GCHEM_CALC_TENDENCY [FORWARD_STEP]',myThid)
368            CALL GCHEM_CALC_TENDENCY( myTime, myIter, myThid )
369            CALL TIMER_STOP ('GCHEM_CALC_TENDENCY [FORWARD_STEP]',myThid)
370           ENDIF
371    #endif /* ALLOW_GCHEM */
372    
373    #ifdef ALLOW_AUTODIFF_TAMC
374    cph needed to be moved here from do_oceanic_physics
375    cph to be visible down the road
376    c
377    CADJ STORE surfaceForcingS    = comlev1, key = ikey_dynamics
378    CADJ STORE surfaceForcingT    = comlev1, key = ikey_dynamics
379    CADJ STORE surfaceForcingTice = comlev1, key = ikey_dynamics
380    ctest(
381    CADJ STORE IVDConvCount       = comlev1, key = ikey_dynamics
382    ctest)
383    # ifdef ALLOW_PTRACERS
384    CADJ STORE surfaceForcingPtr  = comlev1, key = ikey_dynamics
385    # endif
386    c
387    # ifdef ALLOW_GMREDI
388    CADJ STORE Kwx                = comlev1, key = ikey_dynamics
389    CADJ STORE Kwy                = comlev1, key = ikey_dynamics
390    CADJ STORE Kwz                = comlev1, key = ikey_dynamics
391    #  ifdef GM_BOLUS_ADVEC
392    CADJ STORE GM_PsiX            = comlev1, key = ikey_dynamics
393    CADJ STORE GM_PsiY            = comlev1, key = ikey_dynamics
394    #  endif
395    # endif
396    c
397    # ifdef ALLOW_KPP
398    CADJ STORE KPPghat            = comlev1, key = ikey_dynamics
399    CADJ STORE KPPfrac            = comlev1, key = ikey_dynamics
400    CADJ STORE KPPdiffKzS         = comlev1, key = ikey_dynamics
401    CADJ STORE KPPdiffKzT         = comlev1, key = ikey_dynamics
402    # endif
403    c
404    # ifdef NONLIN_FRSURF
405    cph-test
406    CADJ STORE etaN               = comlev1, key = ikey_dynamics
407    CADJ STORE etaH               = comlev1, key = ikey_dynamics
408    CADJ STORE uVel               = comlev1, key = ikey_dynamics
409    CADJ STORE vVel               = comlev1, key = ikey_dynamics
410    #  ifdef ALLOW_CD_CODE
411    CADJ STORE etanm1             = comlev1, key = ikey_dynamics
412    #  endif
413    # endif
414    #endif /* ALLOW_AUTODIFF_TAMC */
415    
416          IF ( .NOT.staggerTimeStep ) THEN
417    #ifdef ALLOW_DEBUG
418            IF ( debugLevel .GE. debLevB )
419         &    CALL DEBUG_CALL('THERMODYNAMICS',myThid)
420    #endif
421            CALL TIMER_START('THERMODYNAMICS      [FORWARD_STEP]',mythid)
422            CALL THERMODYNAMICS( myTime, myIter, myThid )
423            CALL TIMER_STOP ('THERMODYNAMICS      [FORWARD_STEP]',mythid)
424    C--    if not staggerTimeStep: end
425        ENDIF        ENDIF
426    
427    #ifdef COMPONENT_MODULE
428           IF ( useCoupler .AND. .NOT.cpl_earlyExpImpCall ) THEN
429    C      Post coupling data that I export.
430    C      Read in coupling data that I import.
431             myItP1 = myIter + 1
432             CALL TIMER_START('CPL_EXPORT-IMPORT  [FORWARD_STEP]',myThid)
433             CALL CPL_EXPORT_MY_DATA(       myItP1, myTime, myThid )
434             CALL CPL_IMPORT_EXTERNAL_DATA( myItP1, myTime, myThid )
435             CALL TIMER_STOP ('CPL_EXPORT-IMPORT  [FORWARD_STEP]',myThid)
436    # ifndef ALLOW_AIM
437            IF ( useRealFreshWaterFlux ) THEN
438             CALL OCN_APPLY_IMPORT( .FALSE., myTime, myIter, myThid )
439            ENDIF
440    # endif
441           ENDIF
442    #endif /* COMPONENT_MODULE */
443    
444    C--   Step forward fields and calculate time tendency terms.
445    #ifndef ALLOW_OFFLINE
446    #ifndef ALLOW_AUTODIFF_TAMC
447          IF ( momStepping ) THEN
448    #endif
449    #ifdef ALLOW_DEBUG
450            IF ( debugLevel .GE. debLevB )
451         &    CALL DEBUG_CALL('DYNAMICS',myThid)
452    #endif
453            CALL TIMER_START('DYNAMICS            [FORWARD_STEP]',mythid)
454            CALL DYNAMICS( myTime, myIter, myThid )
455            CALL TIMER_STOP ('DYNAMICS            [FORWARD_STEP]',mythid)
456    #ifndef ALLOW_AUTODIFF_TAMC
457          ENDIF
458    #endif
459    #endif
460    
461    C--   Update time-counter
462          myIter = nIter0 + iLoop
463          myTime = startTime + deltaTClock * float(iLoop)
464    
465    #ifdef ALLOW_MNC
466    C     Update the default next iter for MNC
467          IF ( useMNC ) THEN
468             CALL MNC_CW_CITER_SETG( 1, 1, -1, myIter , myThid )
469    
470    C        TODO: Logic should be added here so that users can specify, on
471    C        a per-citer-group basis, when it is time to update the
472    C        "current" (and not just the "next") iteration
473    
474    C        TODO: the following is just a temporary band-aid (mostly, for
475    C        Baylor) until someone writes a routine that better handles time
476    C        boundaries such as weeks, months, years, etc.
477             IF ( mnc_filefreq .GT. 0 ) THEN
478               IF (DIFFERENT_MULTIPLE(mnc_filefreq,myTime,deltaTClock))
479         &          THEN
480                 CALL MNC_CW_CITER_SETG( 1, 1, myIter, -1 , myThid )
481               ENDIF
482             ENDIF
483           ENDIF
484  #endif  #endif
485    
486    C--   Update geometric factors:
487    #ifdef NONLIN_FRSURF
488    C-    update hfacC,W,S and recip_hFac according to etaH(n+1) :
489          IF ( nonlinFreeSurf.GT.0) THEN
490           IF ( select_rStar.GT.0 ) THEN
491    # ifndef DISABLE_RSTAR_CODE
492            CALL TIMER_START('UPDATE_R_STAR      [FORWARD_STEP]',myThid)
493            CALL UPDATE_R_STAR( myTime, myIter, myThid )
494            CALL TIMER_STOP ('UPDATE_R_STAR      [FORWARD_STEP]',myThid)
495    # endif /* DISABLE_RSTAR_CODE */
496           ELSE
497    #ifdef ALLOW_AUTODIFF_TAMC
498    cph-test
499    CADJ STORE hFac_surfC, hFac_surfS, hFac_surfW
500    CADJ &     = comlev1, key = ikey_dynamics
501    #endif
502            CALL TIMER_START('UPDATE_SURF_DR     [FORWARD_STEP]',myThid)
503            CALL UPDATE_SURF_DR( myTime, myIter, myThid )
504            CALL TIMER_STOP ('UPDATE_SURF_DR     [FORWARD_STEP]',myThid)
505           ENDIF
506          ENDIF
507    C-    update also CG2D matrix (and preconditioner)
508          IF ( momStepping .AND. nonlinFreeSurf.GT.2 ) THEN
509            CALL TIMER_START('UPDATE_CG2D        [FORWARD_STEP]',myThid)
510            CALL UPDATE_CG2D( myTime, myIter, myThid )
511            CALL TIMER_STOP ('UPDATE_CG2D        [FORWARD_STEP]',myThid)
512          ENDIF
513    #endif /* NONLIN_FRSURF */
514    
515    C--   Apply Filters to u*,v* before SOLVE_FOR_PRESSURE
516    #ifdef ALLOW_SHAP_FILT
517          IF (useSHAP_FILT .AND. shap_filt_uvStar) THEN
518            CALL TIMER_START('SHAP_FILT_UV        [FORWARD_STEP]',myThid)
519            IF (implicDiv2Dflow.LT.1.) THEN
520    C--   Explicit+Implicit part of the Barotropic Flow Divergence
521    C      => Filtering of uVel,vVel is necessary
522              CALL SHAP_FILT_APPLY_UV( uVel,vVel,
523         &                             myTime, myIter, myThid )
524            ENDIF
525            CALL SHAP_FILT_APPLY_UV( gU,gV,myTime,myIter,myThid)
526            CALL TIMER_STOP ('SHAP_FILT_UV        [FORWARD_STEP]',myThid)
527          ENDIF
528    #endif
529    #ifdef ALLOW_ZONAL_FILT
530          IF (useZONAL_FILT .AND. zonal_filt_uvStar) THEN
531            CALL TIMER_START('ZONAL_FILT_UV       [FORWARD_STEP]',myThid)
532            IF (implicDiv2Dflow.LT.1.) THEN
533    C--   Explicit+Implicit part of the Barotropic Flow Divergence
534    C      => Filtering of uVel,vVel is necessary
535              CALL ZONAL_FILT_APPLY_UV( uVel, vVel, myThid )
536            ENDIF
537            CALL ZONAL_FILT_APPLY_UV( gU, gV, myThid )
538            CALL TIMER_STOP ('ZONAL_FILT_UV       [FORWARD_STEP]',myThid)
539          ENDIF
540    #endif  
541    
542  C--   Solve elliptic equation(s).  C--   Solve elliptic equation(s).
543  C     Two-dimensional only for conventional hydrostatic or  C     Two-dimensional only for conventional hydrostatic or
544  C     three-dimensional for non-hydrostatic and/or IGW scheme.  C     three-dimensional for non-hydrostatic and/or IGW scheme.
545        CALL TIMER_START('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)  #ifndef ALLOW_OFFLINE
546        CALL SOLVE_FOR_PRESSURE( myThid )        IF ( momStepping ) THEN
547        CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)          CALL TIMER_START('SOLVE_FOR_PRESSURE  [FORWARD_STEP]',myThid)
548            CALL SOLVE_FOR_PRESSURE(myTime, myIter, myThid)
549  C--   Correct divergence in flow field and cycle time-stepping          CALL TIMER_STOP ('SOLVE_FOR_PRESSURE  [FORWARD_STEP]',myThid)
550  C     arrays (for all fields) ; update time-counter        ENDIF
551        myCurrentIter = nIter0 + iLoop  #endif
552        myCurrentTime = startTime + deltaTClock * float(iLoop)  
553        CALL THE_CORRECTION_STEP(myCurrentTime, myCurrentIter, myThid)  C--   Correct divergence in flow field and cycle time-stepping momentum
554    c     IF ( momStepping ) THEN
555    #ifndef ALLOW_OFFLINE
556            CALL TIMER_START('MOM_CORRECTION_STEP [FORWARD_STEP]',myThid)
557            CALL MOMENTUM_CORRECTION_STEP(myTime, myIter, myThid)
558            CALL TIMER_STOP ('MOM_CORRECTION_STEP [FORWARD_STEP]',myThid)
559    #endif
560    c     ENDIF
561    
562    #ifdef EXACT_CONSERV
563          IF (exactConserv) THEN
564    C--   Update etaH(n+1) :
565            CALL TIMER_START('UPDATE_ETAH         [FORWARD_STEP]',mythid)
566            CALL UPDATE_ETAH( myTime, myIter, myThid )
567            CALL TIMER_STOP ('UPDATE_ETAH         [FORWARD_STEP]',mythid)
568          ENDIF
569    #endif /* EXACT_CONSERV */
570    
571    #ifdef NONLIN_FRSURF
572          IF ( select_rStar.NE.0 ) THEN
573    # ifndef DISABLE_RSTAR_CODE
574    C--   r* : compute the future level thickness according to etaH(n+1)
575            CALL TIMER_START('CALC_R_STAR       [FORWARD_STEP]',mythid)
576            CALL CALC_R_STAR(etaH, myTime, myIter, myThid )
577            CALL TIMER_STOP ('CALC_R_STAR       [FORWARD_STEP]',mythid)
578    # endif /* DISABLE_RSTAR_CODE */
579          ELSEIF ( nonlinFreeSurf.GT.0) THEN
580    C--   compute the future surface level thickness according to etaH(n+1)
581    #ifdef ALLOW_AUTODIFF_TAMC
582    CADJ STORE etaH          = comlev1, key = ikey_dynamics
583    #endif
584            CALL TIMER_START('CALC_SURF_DR      [FORWARD_STEP]',mythid)
585            CALL CALC_SURF_DR(etaH, myTime, myIter, myThid )
586            CALL TIMER_STOP ('CALC_SURF_DR      [FORWARD_STEP]',mythid)
587          ENDIF
588    #endif /* NONLIN_FRSURF */
589    
590    C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
591          IF ( staggerTimeStep ) THEN
592    C--   do exchanges of U,V (needed for multiDim) when using stagger time-step :
593    #ifdef ALLOW_DEBUG
594            IF ( debugLevel .GE. debLevB )
595         &    CALL DEBUG_CALL('DO_STAGGER_FIELDS_EXCH.',myThid)
596    #endif
597            CALL TIMER_START('BLOCKING_EXCHANGES  [FORWARD_STEP]',myThid)
598            CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid )
599            CALL TIMER_STOP ('BLOCKING_EXCHANGES  [FORWARD_STEP]',myThid)
600    
601    C--   State-variables diagnostics
602           IF ( usediagnostics ) THEN
603            CALL TIMER_START('DO_STATEVARS_DIAGS  [FORWARD_STEP]',myThid)
604            CALL DO_STATEVARS_DIAGS( myTime, 1, myIter, myThid )
605            CALL TIMER_STOP ('DO_STATEVARS_DIAGS  [FORWARD_STEP]',myThid)
606           ENDIF
607    
608    #ifdef ALLOW_DEBUG
609            IF ( debugLevel .GE. debLevB )
610         &    CALL DEBUG_CALL('THERMODYNAMICS',myThid)
611    #endif
612            CALL TIMER_START('THERMODYNAMICS      [FORWARD_STEP]',mythid)
613            CALL THERMODYNAMICS( myTime, myIter, myThid )
614            CALL TIMER_STOP ('THERMODYNAMICS      [FORWARD_STEP]',mythid)
615    
616    C--    if staggerTimeStep: end
617          ENDIF
618    C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
619    
620    #ifdef ALLOW_AUTODIFF_TAMC
621    cph This is needed because convective_adjustment calls
622    cph find_rho which may use pressure()
623    CADJ STORE totphihyd  = comlev1, key = ikey_dynamics
624    #endif
625    C--   Cycle time-stepping Tracers arrays (T,S,+pTracers)
626            CALL TIMER_START('TRC_CORRECTION_STEP [FORWARD_STEP]',myThid)
627            CALL TRACERS_CORRECTION_STEP(myTime, myIter, myThid)
628            CALL TIMER_STOP ('TRC_CORRECTION_STEP [FORWARD_STEP]',myThid)
629    
630    #ifdef ALLOW_GCHEM
631    C     Add separate timestepping of chemical/biological/forcing
632    C     of ptracers here in GCHEM_FORCING_SEP
633            IF ( useGCHEM ) THEN
634    #ifdef ALLOW_DEBUG
635             IF ( debugLevel .GE. debLevB )
636         &        CALL DEBUG_CALL('GCHEM_FORCING_SEP',myThid)
637    #endif /* ALLOW_DEBUG */
638             CALL TIMER_START('GCHEM_FORCING_SEP  [FORWARD_STEP]',myThid)
639             CALL GCHEM_FORCING_SEP( myTime,myIter,myThid )
640             CALL TIMER_STOP ('GCHEM_FORCING_SEP  [FORWARD_STEP]',myThid)
641            ENDIF  
642    #endif /* ALLOW_GCHEM */
643    
644  C--   Do "blocking" sends and receives for tendency "overlap" terms  C--   Do "blocking" sends and receives for tendency "overlap" terms
645  c     CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)  c     CALL TIMER_START('BLOCKING_EXCHANGES  [FORWARD_STEP]',myThid)
646  c     CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )  c     CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )
647  c     CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)  c     CALL TIMER_STOP ('BLOCKING_EXCHANGES  [FORWARD_STEP]',myThid)
648    
649  C--   Do "blocking" sends and receives for field "overlap" terms  C--   Do "blocking" sends and receives for field "overlap" terms
650        CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)        CALL TIMER_START('BLOCKING_EXCHANGES  [FORWARD_STEP]',myThid)
651        CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )        CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
652        CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)        CALL TIMER_STOP ('BLOCKING_EXCHANGES  [FORWARD_STEP]',myThid)
653    
654    
655    C AMM
656    #ifdef ALLOW_GRIDALT
657            if (useGRIDALT) then
658             CALL GRIDALT_UPDATE(myThid)
659            endif
660    #endif
661    C AMM
662    
663    C AMM
664    #ifdef ALLOW_FIZHI
665            if( useFIZHI) then
666             CALL TIMER_START('FIZHI               [FORWARD_STEP]',mythid)
667             CALL STEP_FIZHI_CORR ( myTime, myIter, myThid )
668             CALL TIMER_STOP('FIZHI               [FORWARD_STEP]',mythid)
669            endif
670    #endif
671    C AMM
672    
673    #ifdef ALLOW_FLT
674    C--   Calculate float trajectories
675          IF (useFLT) THEN
676            CALL TIMER_START('FLOATS            [FORWARD_STEP]',myThid)
677            CALL FLT_MAIN(myIter,myTime, myThid)
678            CALL TIMER_STOP ('FLOATS            [FORWARD_STEP]',myThid)
679          ENDIF
680    #endif
681    
682    C--   State-variables time-averaging
683          CALL TIMER_START('DO_STATEVARS_TAVE   [FORWARD_STEP]',myThid)
684          CALL DO_STATEVARS_TAVE( myTime, myIter, myThid )
685          CALL TIMER_STOP ('DO_STATEVARS_TAVE   [FORWARD_STEP]',myThid)
686    
687    #ifndef ALLOW_OFFLINE
688    #ifdef ALLOW_MONITOR
689    C--   Check status of solution (statistics, cfl, etc...)
690          CALL TIMER_START('MONITOR             [FORWARD_STEP]',myThid)
691          CALL MONITOR( myIter, myTime, myThid )
692          CALL TIMER_STOP ('MONITOR             [FORWARD_STEP]',myThid)
693    #endif /* ALLOW_MONITOR */
694    #endif
695    
696    #ifdef ALLOW_COST
697    C--     compare model with data and compute cost function
698    C--     this is done after exchanges to allow interpolation
699          CALL TIMER_START('COST_TILE           [FORWARD_STEP]',myThid)
700          CALL COST_TILE  ( mytime, myiter, myThid )
701          CALL TIMER_STOP ('COST_TILE           [FORWARD_STEP]',myThid)
702    #endif
703    
704  C--   Do IO if needed.  C--   Do IO if needed.
705        CALL TIMER_START('I/O (WRITE)        [FORWARD_STEP]',myThid)  #ifdef ALLOW_OFFLINE
706        CALL DO_THE_MODEL_IO(.FALSE.,        CALL TIMER_START('OFFLINE_MODEL_IO    [FORWARD_STEP]',myThid)
707       &                     myCurrentTime, myCurrentIter, myThid )        CALL OFFLINE_MODEL_IO( myTime, myIter, myThid )
708        CALL TIMER_STOP ('I/O (WRITE)        [FORWARD_STEP]',myThid)        CALL TIMER_STOP ('OFFLINE_MODEL_IO    [FORWARD_STEP]',myThid)
709    #else
710          CALL TIMER_START('DO_THE_MODEL_IO     [FORWARD_STEP]',myThid)
711          CALL DO_THE_MODEL_IO( myTime, myIter, myThid )
712          CALL TIMER_STOP ('DO_THE_MODEL_IO     [FORWARD_STEP]',myThid)
713    #endif
714    
715    #ifdef HAVE_SIGREG
716          IF ( useSIGREG ) THEN
717            IF ( i_got_signal .GT. 0 ) THEN
718              CALL PACKAGES_WRITE_PICKUP(
719         I         .TRUE., myTime, myIter, myThid )
720    #ifndef ALLOW_OFFLINE
721              CALL WRITE_CHECKPOINT(
722         I         .TRUE., myTime, myIter, myThid )  
723    #endif
724              STOP 'Checkpoint completed -- killed by signal handler'
725            ENDIF
726          ENDIF
727    #endif
728    
729  C--   Save state for restarts  C--   Save state for restarts
730  C     Note:    (jmc: is it still the case after ckp35 ?)        CALL TIMER_START('WRITE_CHECKPOINT    [FORWARD_STEP]',myThid)
731  C     =====        CALL PACKAGES_WRITE_PICKUP(
732  C     Because of the ordering of the timestepping code and       I               .FALSE., myTime, myIter, myThid )
733  C     tendency term code at end of loop model arrays hold  #ifndef ALLOW_OFFLINE
 C     U,V,T,S  at "time-level" N but gu, gv, gs, gt, guNM1,...  
 C     at "time-level" N+1/2 (guNM1 at "time-level" N+1/2 is  
 C     gu at "time-level" N-1/2) and cg2d_x at "time-level" N+1/2.  
 C      where N = I+timeLevBase-1  
 C     Thus a checkpoint contains U.0000000000, GU.0000000001 and  
 C     cg2d_x.0000000001 in the indexing scheme used for the model  
 C     "state" files. This example is referred to as a checkpoint  
 C     at time level 1  
       CALL TIMER_START('I/O (WRITE)        [FORWARD_STEP]',myThid)  
734        CALL WRITE_CHECKPOINT(        CALL WRITE_CHECKPOINT(
735       &        .FALSE., myCurrentTime, myCurrentIter, myThid )       I               .FALSE., myTime, myIter, myThid )  
736        CALL TIMER_STOP ('I/O (WRITE)        [FORWARD_STEP]',myThid)  #endif
737          CALL TIMER_STOP ('WRITE_CHECKPOINT    [FORWARD_STEP]',myThid)
738    
739    #ifdef ALLOW_DEBUG
740          IF ( debugLevel .GE. debLevB )
741         &    CALL DEBUG_LEAVE('FORWARD_STEP',myThid)
742    #endif
743    
744        RETURN        RETURN
745        END        END

Legend:
Removed from v.1.8  
changed lines
  Added in v.1.127

  ViewVC Help
Powered by ViewVC 1.1.22