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

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

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


Revision 1.12 - (hide annotations) (download)
Fri Jul 13 20:14:08 2001 UTC (22 years, 11 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint40pre3
Changes since 1.11: +126 -102 lines
o Extracted core part of the_main_loop and re-created forward_step
  N.B.: Time-dependent part of cost function should remain in
        the_main_loop (or contributions must be stored)
o Added some parameter recomputations of nIter0

1 adcroft 1.1
2     #include "CPP_OPTIONS.h"
3    
4 heimbach 1.12 subroutine forward_step( iloop, mytime, myiter, mythid )
5    
6     c ==================================================================
7     c SUBROUTINE forward_step
8     c ==================================================================
9     c
10     c o Run the ocean model and evaluate the specified cost function.
11     c
12     c *the_main_loop* is the toplevel routine for the Tangent Linear and
13     c Adjoint Model Compiler (TAMC). For this purpose the initialization
14     c of the model was split into two parts. Those parameters that do
15     c not depend on a specific model run are set in *initialise_fixed*,
16     c whereas those that do depend on the specific realization are
17     c initialized in *initialise_varia*.
18     c This routine is to be used in conjuction with the MITgcmuv
19     c checkpoint 37.
20     c
21     c ==================================================================
22     c SUBROUTINE forward_step
23     c ==================================================================
24    
25     implicit none
26    
27     c == global variables ==
28 adcroft 1.1
29     #include "SIZE.h"
30     #include "EEPARAMS.h"
31     #include "PARAMS.h"
32     #include "DYNVARS.h"
33 heimbach 1.12 #include "FFIELDS.h"
34     #include "TR1.h"
35    
36 adcroft 1.1 #ifdef ALLOW_NONHYDROSTATIC
37     #include "CG3D.h"
38     #endif
39    
40 heimbach 1.12 #ifdef ALLOW_AUTODIFF_TAMC
41     #include "tamc.h"
42     #include "ctrl.h"
43     #include "ctrl_dummy.h"
44     #include "cost.h"
45     #endif
46    
47     c == routine arguments ==
48     c note: under the multi-threaded model myiter and
49     c mytime are local variables passed around as routine
50     c arguments. Although this is fiddly it saves the need to
51     c impose additional synchronisation points when they are
52     c updated.
53     c myiter - iteration counter for this thread
54     c mytime - time counter for this thread
55     c mythid - thread number for this instance of the routine.
56     integer iloop
57     integer mythid
58     integer myiter
59     _RL mytime
60    
61     c == local variables ==
62    
63     c-- == end of interface ==
64    
65     c-- >>> Loop body start <<<
66    
67     #ifdef ALLOW_AUTODIFF_TAMC
68     c-- Set the model iteration counter and the model time.
69     myiter = nIter0 + (iloop-1)
70     mytime = startTime + float(iloop-1)*deltaTclock
71    
72     c Include call to a dummy routine. Its adjoint will be
73     c called at the proper place in the adjoint code.
74     c The adjoint routine will print out adjoint values
75     c if requested. The location of the call is important,
76     c it has to be after the adjoint of the exchanges
77     c (DO_GTERM_BLOCKING_EXCHANGES).
78     call dummy_in_stepping( myTime, myIter, myThid )
79     #endif
80    
81     c-- Load forcing/external data fields.
82     #ifdef INCLUDE_EXTERNAL_FORCING_PACKAGE
83     c NOTE, that although the exf package is part of the
84     c distribution, it is not currently maintained, i.e.
85     c exf is disabled by default in genmake.
86     CALL EXF_GETFORCING( mytime, myiter, mythid )
87     #else
88     CALL TIMER_START('EXTERNAL_FIELDS_LOAD[THE_MAIN_LOOP]',mythid)
89     CALL EXTERNAL_FIELDS_LOAD( mytime, myiter, mythid )
90     CALL TIMER_STOP ('EXTERNAL_FIELDS_LOAD[THE_MAIN_LOOP]',mythid)
91     #endif
92    
93     c-- Step forward fields and calculate time tendency terms.
94     CALL TIMER_START('DYNAMICS [THE_MAIN_LOOP]',mythid)
95     CALL DYNAMICS( myTime, myIter, myThid )
96     CALL TIMER_STOP ('DYNAMICS [THE_MAIN_LOOP]',mythid)
97    
98     #ifndef ALLOW_AUTODIFF_TAMC
99 adcroft 1.1
100     #ifdef ALLOW_NONHYDROSTATIC
101     C-- Step forward W field in N-H algorithm
102 heimbach 1.12 IF ( nonHydrostatic ) THEN
103     CALL TIMER_START('CALC_GW [THE_MAIN_LOOP]',myThid)
104     CALL CALC_GW(myThid)
105     CALL TIMER_STOP ('CALC_GW [THE_MAIN_LOOP]',myThid)
106     ENDIF
107 adcroft 1.1 #endif
108    
109 adcroft 1.5 #ifdef ALLOW_SHAP_FILT
110     C-- Step forward all tiles, filter and exchange.
111 heimbach 1.12 CALL TIMER_START('SHAP_FILT [THE_MAIN_LOOP]',myThid)
112 adcroft 1.5 CALL SHAP_FILT_APPLY(
113     I gUnm1, gVnm1, gTnm1, gSnm1,
114 heimbach 1.12 I myTime, myIter, myThid )
115 jmc 1.9 IF (implicDiv2Dflow.LT.1.) THEN
116     C-- Explicit+Implicit part of the Barotropic Flow Divergence
117     C => Filtering of uVel,vVel is necessary
118 heimbach 1.12 CALL SHAP_FILT_UV( uVel, vVel, myTime, myThid )
119 jmc 1.9 ENDIF
120 heimbach 1.12 CALL TIMER_STOP ('SHAP_FILT [THE_MAIN_LOOP]',myThid)
121 adcroft 1.5 #endif
122    
123     #ifdef ALLOW_ZONAL_FILT
124 jmc 1.8 IF (zonal_filt_lat.LT.90.) THEN
125 heimbach 1.12 CALL TIMER_START('ZONAL_FILT_APPLY [THE_MAIN_LOOP]',myThid)
126     CALL ZONAL_FILT_APPLY(
127 adcroft 1.5 U gUnm1, gVnm1, gTnm1, gSnm1,
128     I myThid )
129 heimbach 1.12 CALL TIMER_STOP ('ZONAL_FILT_APPLY [THE_MAIN_LOOP]',myThid)
130 jmc 1.8 ENDIF
131 adcroft 1.5 #endif
132    
133 heimbach 1.12 #endif /* ALLOW_AUTODIFF_TAMC */
134    
135 adcroft 1.1 C-- Solve elliptic equation(s).
136     C Two-dimensional only for conventional hydrostatic or
137     C three-dimensional for non-hydrostatic and/or IGW scheme.
138 heimbach 1.12 CALL TIMER_START('SOLVE_FOR_PRESSURE [THE_MAIN_LOOP]',myThid)
139 adcroft 1.1 CALL SOLVE_FOR_PRESSURE( myThid )
140 heimbach 1.12 CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [THE_MAIN_LOOP]',myThid)
141 adcroft 1.1
142 adcroft 1.5 C-- Correct divergence in flow field and cycle time-stepping
143 jmc 1.7 C arrays (for all fields) ; update time-counter
144 heimbach 1.12 myIter = nIter0 + iLoop
145     myTime = startTime + deltaTClock * float(iLoop)
146     CALL TIMER_START('THE_CORRECTION_STEP [THE_MAIN_LOOP]',myThid)
147     CALL THE_CORRECTION_STEP(myTime, myIter, myThid)
148     CALL TIMER_STOP ('THE_CORRECTION_STEP [THE_MAIN_LOOP]',myThid)
149 adcroft 1.5
150 adcroft 1.1 C-- Do "blocking" sends and receives for tendency "overlap" terms
151 heimbach 1.12 c CALL TIMER_START('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
152 jmc 1.7 c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )
153 heimbach 1.12 c CALL TIMER_STOP ('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
154 adcroft 1.5
155     C-- Do "blocking" sends and receives for field "overlap" terms
156 heimbach 1.12 CALL TIMER_START('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
157 adcroft 1.5 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
158 heimbach 1.12 CALL TIMER_STOP ('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
159    
160     #ifndef EXCLUDE_MONITOR
161     C-- Check status of solution (statistics, cfl, etc...)
162     CALL MONITOR( myIter, myTime, myThid )
163     #endif /* EXCLUDE_MONITOR */
164 adcroft 1.1
165 heimbach 1.12 #ifndef ALLOW_AUTODIFF_TAMC
166 jmc 1.7 C-- Do IO if needed.
167 heimbach 1.12 CALL TIMER_START('DO_THE_MODEL_IO [THE_MAIN_LOOP]',myThid)
168     CALL DO_THE_MODEL_IO( myTime, myIter, myThid )
169     CALL TIMER_STOP ('DO_THE_MODEL_IO [THE_MAIN_LOOP]',myThid)
170 adcroft 1.1
171     C-- Save state for restarts
172 jmc 1.7 C Note: (jmc: is it still the case after ckp35 ?)
173 adcroft 1.1 C =====
174     C Because of the ordering of the timestepping code and
175     C tendency term code at end of loop model arrays hold
176     C U,V,T,S at "time-level" N but gu, gv, gs, gt, guNM1,...
177     C at "time-level" N+1/2 (guNM1 at "time-level" N+1/2 is
178 jmc 1.10 C gu at "time-level" N-1/2) and etaN at "time-level" N+1/2.
179 adcroft 1.1 C where N = I+timeLevBase-1
180     C Thus a checkpoint contains U.0000000000, GU.0000000001 and
181 jmc 1.10 C etaN.0000000001 in the indexing scheme used for the model
182 adcroft 1.1 C "state" files. This example is referred to as a checkpoint
183     C at time level 1
184 heimbach 1.12 CALL TIMER_START('WRITE_CHECKPOINT [THE_MAIN_LOOP]',myThid)
185 adcroft 1.1 CALL WRITE_CHECKPOINT(
186 heimbach 1.12 & .FALSE., myTime, myIter, myThid )
187     CALL TIMER_STOP ('WRITE_CHECKPOINT [THE_MAIN_LOOP]',myThid)
188    
189     #endif /* ALLOW_AUTODIFF_TAMC */
190 adcroft 1.1
191 heimbach 1.12 c-- >>> Loop body end <<<
192 adcroft 1.1
193     END

  ViewVC Help
Powered by ViewVC 1.1.22