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

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

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


Revision 1.23.2.5 - (show annotations) (download)
Thu Nov 7 16:48:46 2002 UTC (21 years, 7 months ago) by heimbach
Branch: release1
CVS Tags: release1_p7
Changes since 1.23.2.4: +15 -10 lines
Added external_fields_load routine to TAF list.
Update corresponding checkpointing lists
(analog to exf handling of swapping).

1 C $Header: /u/u0/gcmpack/MITgcm/model/src/forward_step.F,v 1.34 2002/07/13 04:59:42 heimbach Exp $
2 C $Name: checkpoint46 $
3
4 #include "CPP_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: FORWARD_STEP
8 C !INTERFACE:
9 SUBROUTINE FORWARD_STEP( iloop, myTime, myIter, myThid )
10
11 C !DESCRIPTION: \bv
12 C *==================================================================
13 C | SUBROUTINE forward_step
14 C | o Run the ocean model and, optionally, evaluate a cost function.
15 C *==================================================================
16 C |
17 C | THE_MAIN_LOOP is the toplevel routine for the Tangent Linear and
18 C | Adjoint Model Compiler (TAMC). For this purpose the initialization
19 C | of the model was split into two parts. Those parameters that do
20 C | not depend on a specific model run are set in INITIALISE_FIXED,
21 C | whereas those that do depend on the specific realization are
22 C | initialized in INITIALISE_VARIA.
23 C |
24 C *==================================================================
25 C \ev
26
27 C !USES:
28 IMPLICIT NONE
29 C == Global variables ==
30 #include "SIZE.h"
31 #include "EEPARAMS.h"
32 #include "PARAMS.h"
33 #include "DYNVARS.h"
34 #include "FFIELDS.h"
35
36 #ifdef ALLOW_NONHYDROSTATIC
37 #include "CG3D.h"
38 #endif
39
40 #ifdef ALLOW_SHAP_FILT
41 #include "SHAP_FILT.h"
42 #endif
43 #ifdef ALLOW_ZONAL_FILT
44 #include "ZONAL_FILT.h"
45 #endif
46
47 #ifdef ALLOW_AUTODIFF_TAMC
48 # include "tamc.h"
49 # include "ctrl.h"
50 # include "ctrl_dummy.h"
51 # include "cost.h"
52 # ifdef INCLUDE_EXTERNAL_FORCING_PACKAGE
53 # include "exf_fields.h"
54 # ifdef ALLOW_BULKFORMULAE
55 # include "exf_constants.h"
56 # endif
57 # endif
58 # ifdef ALLOW_OBCS
59 # include "OBCS.h"
60 # endif
61 #endif
62
63 C !LOCAL VARIABLES:
64 C == Routine arguments ==
65 C note: under the multi-threaded model myiter and
66 C mytime are local variables passed around as routine
67 C arguments. Although this is fiddly it saves the need to
68 C impose additional synchronisation points when they are
69 C updated.
70 C myiter - iteration counter for this thread
71 C mytime - time counter for this thread
72 C mythid - thread number for this instance of the routine.
73 integer iloop
74 integer mythid
75 integer myiter
76 _RL mytime
77 INTEGER bi,bj
78
79 CEOP
80
81 #ifdef ALLOW_AUTODIFF_TAMC
82 C-- Reset the model iteration counter and the model time.
83 myiter = nIter0 + (iloop-1)
84 mytime = startTime + float(iloop-1)*deltaTclock
85 #endif
86
87 #if (defined (ALLOW_AUTODIFF_TAMC) && defined (ALLOW_AUTODIFF_MONITOR))
88 C Include call to a dummy routine. Its adjoint will be
89 C called at the proper place in the adjoint code.
90 C The adjoint routine will print out adjoint values
91 C if requested. The location of the call is important,
92 C it has to be after the adjoint of the exchanges
93 C (DO_GTERM_BLOCKING_EXCHANGES).
94 CALL DUMMY_IN_STEPPING( myTime, myIter, myThid )
95 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
96 #endif
97
98 #ifdef EXACT_CONSERV
99 IF (exactConserv) THEN
100 C-- Update etaH(n+1) :
101 DO bj=myByLo(myThid),myByHi(myThid)
102 DO bi=myBxLo(myThid),myBxHi(myThid)
103 CALL CALC_EXACT_ETA( .FALSE., bi,bj, uVel,vVel,
104 I myTime, myIter, myThid )
105 ENDDO
106 ENDDO
107 IF (implicDiv2Dflow .NE. 1. _d 0 .OR. useOBCS )
108 & _EXCH_XY_R8(etaH, myThid )
109 ENDIF
110 #endif /* EXACT_CONSERV */
111
112 #ifdef NONLIN_FRSURF
113 C-- compute the future surface level thickness
114 C according to etaH(n+1)
115 IF ( nonlinFreeSurf.GT.0) THEN
116 CALL CALC_SURF_DR(etaH, myTime, myIter, myThid )
117 ENDIF
118 #endif /* NONLIN_FRSURF */
119
120 C-- Load forcing/external data fields.
121 #ifdef ALLOW_AUTODIFF_TAMC
122 c**************************************
123 #include "checkpoint_lev1_directives.h"
124 c**************************************
125 #endif
126 #ifdef INCLUDE_EXTERNAL_FORCING_PACKAGE
127 C NOTE, that although the exf package is part of the
128 C distribution, it is not currently maintained, i.e.
129 C exf is disabled by default in genmake.
130 CALL TIMER_START('EXF_GETFORCING [THE_MAIN_LOOP]',mythid)
131 CALL EXF_GETFORCING( mytime, myiter, mythid )
132 CALL TIMER_STOP ('EXF_GETFORCING [THE_MAIN_LOOP]',mythid)
133 #else
134 CALL TIMER_START('EXTERNAL_FIELDS_LOAD[THE_MAIN_LOOP]',mythid)
135 CALL EXTERNAL_FIELDS_LOAD( mytime, myiter, mythid )
136 CALL TIMER_STOP ('EXTERNAL_FIELDS_LOAD[THE_MAIN_LOOP]',mythid)
137 #endif
138
139 C-- Step forward fields and calculate time tendency terms.
140 CALL TIMER_START('THERMODYNAMICS [THE_MAIN_LOOP]',mythid)
141 CALL THERMODYNAMICS( myTime, myIter, myThid )
142 CALL TIMER_STOP ('THERMODYNAMICS [THE_MAIN_LOOP]',mythid)
143
144 #ifdef ALLOW_SHAP_FILT
145 IF (useSHAP_FILT .AND.
146 & staggerTimeStep .AND. shap_filt_TrStagg ) THEN
147 CALL TIMER_START('SHAP_FILT [THE_MAIN_LOOP]',myThid)
148 CALL SHAP_FILT_APPLY_TS( gT, gS, myTime, myIter, myThid )
149 CALL TIMER_STOP ('SHAP_FILT [THE_MAIN_LOOP]',myThid)
150 ENDIF
151 #endif
152 #ifdef ALLOW_ZONAL_FILT
153 IF (useZONAL_FILT .AND.
154 & staggerTimeStep .AND. zonal_filt_TrStagg ) THEN
155 CALL TIMER_START('ZONAL_FILT_APPLY [THE_MAIN_LOOP]',myThid)
156 CALL ZONAL_FILT_APPLY_TS( gT, gS, myThid )
157 CALL TIMER_STOP ('ZONAL_FILT_APPLY [THE_MAIN_LOOP]',myThid)
158 ENDIF
159 #endif
160
161 C-- Step forward fields and calculate time tendency terms.
162 IF ( momStepping ) THEN
163 CALL TIMER_START('DYNAMICS [THE_MAIN_LOOP]',mythid)
164 CALL DYNAMICS( myTime, myIter, myThid )
165 CALL TIMER_STOP ('DYNAMICS [THE_MAIN_LOOP]',mythid)
166 ENDIF
167
168 #ifdef ALLOW_NONHYDROSTATIC
169 C-- Step forward W field in N-H algorithm
170 IF ( momStepping .AND. nonHydrostatic ) THEN
171 CALL TIMER_START('CALC_GW [THE_MAIN_LOOP]',myThid)
172 CALL CALC_GW(myThid)
173 CALL TIMER_STOP ('CALC_GW [THE_MAIN_LOOP]',myThid)
174 ENDIF
175 #endif
176
177 #ifdef NONLIN_FRSURF
178 C-- update hfacC,W,S and recip_hFac according to etaH(n+1) :
179 IF ( nonlinFreeSurf.GT.0) THEN
180 CALL TIMER_START('UPDATE_SURF_DR [THE_MAIN_LOOP]',myThid)
181 CALL UPDATE_SURF_DR( myTime, myIter, myThid )
182 CALL TIMER_STOP ('UPDATE_SURF_DR [THE_MAIN_LOOP]',myThid)
183 ENDIF
184 C- update also CG2D matrix (and preconditioner)
185 IF ( momStepping .AND. nonlinFreeSurf.GT.2 ) THEN
186 CALL TIMER_START('UPDATE_CG2D [THE_MAIN_LOOP]',myThid)
187 CALL UPDATE_CG2D( myTime, myIter, myThid )
188 CALL TIMER_STOP ('UPDATE_CG2D [THE_MAIN_LOOP]',myThid)
189 ENDIF
190 #endif
191
192 C-- Apply Filters to u*,v* before SOLVE_FOR_PRESSURE
193 #ifdef ALLOW_SHAP_FILT
194 IF (useSHAP_FILT .AND. shap_filt_uvStar) THEN
195 CALL TIMER_START('SHAP_FILT [THE_MAIN_LOOP]',myThid)
196 CALL SHAP_FILT_APPLY_UV( gUnm1,gVnm1, myTime,myIter,myThid )
197 IF (implicDiv2Dflow.LT.1.) THEN
198 C-- Explicit+Implicit part of the Barotropic Flow Divergence
199 C => Filtering of uVel,vVel is necessary
200 CALL SHAP_FILT_APPLY_UV( uVel,vVel, myTime,myIter,myThid )
201 ENDIF
202 CALL TIMER_STOP ('SHAP_FILT [THE_MAIN_LOOP]',myThid)
203 ENDIF
204 #endif
205 #ifdef ALLOW_ZONAL_FILT
206 IF (useZONAL_FILT .AND. zonal_filt_uvStar) THEN
207 CALL TIMER_START('ZONAL_FILT_APPLY [THE_MAIN_LOOP]',myThid)
208 CALL ZONAL_FILT_APPLY_UV( gUnm1, gVnm1, myThid )
209 IF (implicDiv2Dflow.LT.1.) THEN
210 C-- Explicit+Implicit part of the Barotropic Flow Divergence
211 C => Filtering of uVel,vVel is necessary
212 CALL ZONAL_FILT_APPLY_UV( uVel, vVel, myThid )
213 ENDIF
214 CALL TIMER_STOP ('ZONAL_FILT_APPLY [THE_MAIN_LOOP]',myThid)
215 ENDIF
216 #endif
217
218 C-- Solve elliptic equation(s).
219 C Two-dimensional only for conventional hydrostatic or
220 C three-dimensional for non-hydrostatic and/or IGW scheme.
221 IF ( momStepping ) THEN
222 CALL TIMER_START('SOLVE_FOR_PRESSURE [THE_MAIN_LOOP]',myThid)
223 CALL SOLVE_FOR_PRESSURE(myTime, myIter, myThid)
224 CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [THE_MAIN_LOOP]',myThid)
225 ENDIF
226
227 C-- Correct divergence in flow field and cycle time-stepping
228 C arrays (for all fields) ; update time-counter
229 myIter = nIter0 + iLoop
230 myTime = startTime + deltaTClock * float(iLoop)
231 CALL TIMER_START('THE_CORRECTION_STEP [THE_MAIN_LOOP]',myThid)
232 CALL THE_CORRECTION_STEP(myTime, myIter, myThid)
233 CALL TIMER_STOP ('THE_CORRECTION_STEP [THE_MAIN_LOOP]',myThid)
234
235 C-- Do "blocking" sends and receives for tendency "overlap" terms
236 c CALL TIMER_START('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
237 c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )
238 c CALL TIMER_STOP ('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
239
240 C-- Do "blocking" sends and receives for field "overlap" terms
241 CALL TIMER_START('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
242 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
243 CALL TIMER_STOP ('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
244
245 #ifdef ALLOW_FLT
246 C-- Calculate float trajectories
247 IF (useFLT) THEN
248 CALL TIMER_START('FLOATS [THE_MAIN_LOOP]',myThid)
249 CALL FLT_MAIN(myIter,myTime, myThid)
250 CALL TIMER_STOP ('FLOATS [THE_MAIN_LOOP]',myThid)
251 ENDIF
252 #endif
253
254 #ifndef EXCLUDE_MONITOR
255 C-- Check status of solution (statistics, cfl, etc...)
256 CALL TIMER_START('MONITOR [THE_MAIN_LOOP]',myThid)
257 CALL MONITOR( myIter, myTime, myThid )
258 CALL TIMER_STOP ('MONITOR [THE_MAIN_LOOP]',myThid)
259 #endif /* EXCLUDE_MONITOR */
260
261 C-- Do IO if needed.
262 CALL TIMER_START('DO_THE_MODEL_IO [THE_MAIN_LOOP]',myThid)
263 CALL DO_THE_MODEL_IO( myTime, myIter, myThid )
264 CALL TIMER_STOP ('DO_THE_MODEL_IO [THE_MAIN_LOOP]',myThid)
265
266 C-- Save state for restarts
267 C Note: (jmc: is it still the case after ckp35 ?)
268 C =====
269 C Because of the ordering of the timestepping code and
270 C tendency term code at end of loop model arrays hold
271 C U,V,T,S at "time-level" N but gu, gv, gs, gt, guNM1,...
272 C at "time-level" N+1/2 (guNM1 at "time-level" N+1/2 is
273 C gu at "time-level" N-1/2) and etaN at "time-level" N+1/2.
274 C where N = I+timeLevBase-1
275 C Thus a checkpoint contains U.0000000000, GU.0000000001 and
276 C etaN.0000000001 in the indexing scheme used for the model
277 C "state" files. This example is referred to as a checkpoint
278 C at time level 1
279 CALL TIMER_START('WRITE_CHECKPOINT [THE_MAIN_LOOP]',myThid)
280 CALL WRITE_CHECKPOINT(
281 & .FALSE., myTime, myIter, myThid )
282 CALL TIMER_STOP ('WRITE_CHECKPOINT [THE_MAIN_LOOP]',myThid)
283
284 END

  ViewVC Help
Powered by ViewVC 1.1.22