/[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.37 - (show annotations) (download)
Fri Nov 15 03:01:21 2002 UTC (21 years, 6 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint47
Changes since 1.36: +9 -14 lines
differentiable version of checkpoint46n_post
o external_fields_load now part of differentiation list
o pressure needs multiple storing;
  would be nice to have store_pressure at beginning or
  end of forward_step, e.g. by having phiHyd global (5-dim.)
  (NB: pressure is needed for certain cases in find_rho,
  which is also invoked through convective_adjustment).
o recomputations in find_rho for cases
 'JMD95'/'UNESCO' or 'MDJWF' are OK.
o #define ATMOSPHERIC_LOADING should be differentiable
o ini_forcing shifted to begining of initialise_varia

1 C $Header: /u/gcmpack/MITgcm/model/src/forward_step.F,v 1.36 2002/11/12 20:39:46 heimbach Exp $
2 C $Name: $
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 # include "EOS.h"
53 # ifdef INCLUDE_EXTERNAL_FORCING_PACKAGE
54 # include "exf_fields.h"
55 # ifdef ALLOW_BULKFORMULAE
56 # include "exf_constants.h"
57 # endif
58 # endif
59 # ifdef ALLOW_OBCS
60 # include "OBCS.h"
61 # endif
62 #endif
63
64 C !LOCAL VARIABLES:
65 C == Routine arguments ==
66 C note: under the multi-threaded model myiter and
67 C mytime are local variables passed around as routine
68 C arguments. Although this is fiddly it saves the need to
69 C impose additional synchronisation points when they are
70 C updated.
71 C myiter - iteration counter for this thread
72 C mytime - time counter for this thread
73 C mythid - thread number for this instance of the routine.
74 integer iloop
75 integer mythid
76 integer myiter
77 _RL mytime
78 INTEGER bi,bj
79
80 CEOP
81
82
83 #ifdef ALLOW_AUTODIFF_TAMC
84 C-- Reset the model iteration counter and the model time.
85 myiter = nIter0 + (iloop-1)
86 mytime = startTime + float(iloop-1)*deltaTclock
87 #endif
88
89 #if (defined (ALLOW_AUTODIFF_TAMC) && defined (ALLOW_AUTODIFF_MONITOR))
90 C Include call to a dummy routine. Its adjoint will be
91 C called at the proper place in the adjoint code.
92 C The adjoint routine will print out adjoint values
93 C if requested. The location of the call is important,
94 C it has to be after the adjoint of the exchanges
95 C (DO_GTERM_BLOCKING_EXCHANGES).
96 CALL DUMMY_IN_STEPPING( myTime, myIter, myThid )
97 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
98 #endif
99
100 #ifdef EXACT_CONSERV
101 IF (exactConserv) THEN
102 C-- Update etaH(n+1) :
103 CALL TIMER_START('UPDATE_ETAH [FORWARD_STEP]',mythid)
104 CALL UPDATE_ETAH( myTime, myIter, myThid )
105 CALL TIMER_STOP ('UPDATE_ETAH [FORWARD_STEP]',mythid)
106 ENDIF
107 #endif /* EXACT_CONSERV */
108
109 #ifdef NONLIN_FRSURF
110 C-- compute the future surface level thickness
111 C according to etaH(n+1)
112 IF ( nonlinFreeSurf.GT.0) THEN
113 CALL TIMER_START('CALC_SURF_DR [FORWARD_STEP]',mythid)
114 CALL CALC_SURF_DR(etaH, myTime, myIter, myThid )
115 CALL TIMER_STOP ('CALC_SURF_DR [FORWARD_STEP]',mythid)
116 ENDIF
117 #endif /* NONLIN_FRSURF */
118
119 C-- Load forcing/external data fields.
120 #ifdef ALLOW_AUTODIFF_TAMC
121 c**************************************
122 #include "checkpoint_lev1_directives.h"
123 c**************************************
124 #endif
125 #ifdef INCLUDE_EXTERNAL_FORCING_PACKAGE
126 C-- Call external forcing package
127 IF ( .not. useSEAICE ) THEN
128 CALL TIMER_START('EXF_GETFORCING [FORWARD_STEP]',mythid)
129 CALL EXF_GETFORCING( mytime, myiter, mythid )
130 CALL TIMER_STOP ('EXF_GETFORCING [FORWARD_STEP]',mythid)
131 ENDIF
132 #else
133 IF ( .not. useSEAICE ) THEN
134 CALL TIMER_START('EXTERNAL_FIELDS_LOAD[FORWARD_STEP]',mythid)
135 CALL EXTERNAL_FIELDS_LOAD( mytime, myiter, mythid )
136 CALL TIMER_STOP ('EXTERNAL_FIELDS_LOAD[FORWARD_STEP]',mythid)
137 ENDIF
138 #endif /* INCLUDE_EXTERNAL_FORCING_PACKAGE */
139
140 #ifdef ALLOW_SEAICE
141 C-- Call sea ice model to compute forcing/external data fields.
142 IF ( useSEAICE ) THEN
143 CALL TIMER_START('SEAICE_MODEL [FORWARD_STEP]',myThid)
144 CALL SEAICE_MODEL( myTime, myIter, myThid )
145 CALL TIMER_STOP ('SEAICE_MODEL [FORWARD_STEP]',myThid)
146 ENDIF
147 #endif ALLOW_SEAICE
148
149 C-- Step forward fields and calculate time tendency terms.
150 CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',mythid)
151 CALL THERMODYNAMICS( myTime, myIter, myThid )
152 CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',mythid)
153
154 C-- do exchanges (needed for DYNAMICS) when using stagger time-step :
155 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
156 CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid )
157 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
158
159 #ifdef ALLOW_SHAP_FILT
160 IF (useSHAP_FILT .AND.
161 & staggerTimeStep .AND. shap_filt_TrStagg ) THEN
162 CALL TIMER_START('SHAP_FILT [FORWARD_STEP]',myThid)
163 CALL SHAP_FILT_APPLY_TS( gT, gS, myTime, myIter, myThid )
164 CALL TIMER_STOP ('SHAP_FILT [FORWARD_STEP]',myThid)
165 ENDIF
166 #endif
167 #ifdef ALLOW_ZONAL_FILT
168 IF (useZONAL_FILT .AND.
169 & staggerTimeStep .AND. zonal_filt_TrStagg ) THEN
170 CALL TIMER_START('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid)
171 CALL ZONAL_FILT_APPLY_TS( gT, gS, myThid )
172 CALL TIMER_STOP ('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid)
173 ENDIF
174 #endif
175
176 C-- Step forward fields and calculate time tendency terms.
177 IF ( momStepping ) THEN
178 CALL TIMER_START('DYNAMICS [FORWARD_STEP]',mythid)
179 CALL DYNAMICS( myTime, myIter, myThid )
180 CALL TIMER_STOP ('DYNAMICS [FORWARD_STEP]',mythid)
181 ENDIF
182
183 #ifdef ALLOW_NONHYDROSTATIC
184 C-- Step forward W field in N-H algorithm
185 IF ( momStepping .AND. nonHydrostatic ) THEN
186 CALL TIMER_START('CALC_GW [FORWARD_STEP]',myThid)
187 CALL CALC_GW(myThid)
188 CALL TIMER_STOP ('CALC_GW [FORWARD_STEP]',myThid)
189 ENDIF
190 #endif
191
192 #ifdef NONLIN_FRSURF
193 C-- update hfacC,W,S and recip_hFac according to etaH(n+1) :
194 IF ( nonlinFreeSurf.GT.0) THEN
195 CALL TIMER_START('UPDATE_SURF_DR [FORWARD_STEP]',mythid)
196 CALL UPDATE_SURF_DR( myTime, myIter, myThid )
197 CALL TIMER_STOP ('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
198 ENDIF
199 C- update also CG2D matrix (and preconditioner)
200 IF ( momStepping .AND. nonlinFreeSurf.GT.2 ) THEN
201 CALL TIMER_START('UPDATE_CG2D [FORWARD_STEP]',mythid)
202 CALL UPDATE_CG2D( myTime, myIter, myThid )
203 CALL TIMER_START('UPDATE_CG2D [FORWARD_STEP]',mythid)
204 ENDIF
205 #endif
206
207 C-- Apply Filters to u*,v* before SOLVE_FOR_PRESSURE
208 #ifdef ALLOW_SHAP_FILT
209 IF (useSHAP_FILT .AND. shap_filt_uvStar) THEN
210 CALL TIMER_START('SHAP_FILT [FORWARD_STEP]',myThid)
211 CALL SHAP_FILT_APPLY_UV( gUnm1,gVnm1, myTime,myIter,myThid )
212 IF (implicDiv2Dflow.LT.1.) THEN
213 C-- Explicit+Implicit part of the Barotropic Flow Divergence
214 C => Filtering of uVel,vVel is necessary
215 CALL SHAP_FILT_APPLY_UV( uVel,vVel, myTime,myIter,myThid )
216 ENDIF
217 CALL TIMER_STOP ('SHAP_FILT [FORWARD_STEP]',myThid)
218 ENDIF
219 #endif
220 #ifdef ALLOW_ZONAL_FILT
221 IF (useZONAL_FILT .AND. zonal_filt_uvStar) THEN
222 CALL TIMER_START('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid)
223 CALL ZONAL_FILT_APPLY_UV( gUnm1, gVnm1, myThid )
224 IF (implicDiv2Dflow.LT.1.) THEN
225 C-- Explicit+Implicit part of the Barotropic Flow Divergence
226 C => Filtering of uVel,vVel is necessary
227 CALL ZONAL_FILT_APPLY_UV( uVel, vVel, myThid )
228 ENDIF
229 CALL TIMER_STOP ('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid)
230 ENDIF
231 #endif
232
233 C-- Solve elliptic equation(s).
234 C Two-dimensional only for conventional hydrostatic or
235 C three-dimensional for non-hydrostatic and/or IGW scheme.
236 IF ( momStepping ) THEN
237 CALL TIMER_START('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
238 CALL SOLVE_FOR_PRESSURE(myTime, myIter, myThid)
239 CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
240 ENDIF
241
242 #ifdef ALLOW_AUTODIFF_TAMC
243 cph This is needed because convective_adjustment calls
244 cph find_rho which may use pressure()
245 CADJ STORE pressure = comlev1, key = ikey_dynamics
246 #endif
247 C-- Correct divergence in flow field and cycle time-stepping
248 C arrays (for all fields) ; update time-counter
249 myIter = nIter0 + iLoop
250 myTime = startTime + deltaTClock * float(iLoop)
251 CALL TIMER_START('THE_CORRECTION_STEP [FORWARD_STEP]',myThid)
252 CALL THE_CORRECTION_STEP(myTime, myIter, myThid)
253 CALL TIMER_STOP ('THE_CORRECTION_STEP [FORWARD_STEP]',myThid)
254
255 C-- Do "blocking" sends and receives for tendency "overlap" terms
256 c CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
257 c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )
258 c CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
259
260 C-- Do "blocking" sends and receives for field "overlap" terms
261 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
262 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
263 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
264
265 #ifdef ALLOW_FLT
266 C-- Calculate float trajectories
267 IF (useFLT) THEN
268 CALL TIMER_START('FLOATS [FORWARD_STEP]',myThid)
269 CALL FLT_MAIN(myIter,myTime, myThid)
270 CALL TIMER_STOP ('FLOATS [FORWARD_STEP]',myThid)
271 ENDIF
272 #endif
273
274 #ifndef EXCLUDE_MONITOR
275 C-- Check status of solution (statistics, cfl, etc...)
276 CALL TIMER_START('MONITOR [FORWARD_STEP]',mythid)
277 CALL MONITOR( myIter, myTime, myThid )
278 CALL TIMER_STOP ('MONITOR [FORWARD_STEP]',myThid)
279 #endif /* EXCLUDE_MONITOR */
280
281 C-- Do IO if needed.
282 CALL TIMER_START('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
283 CALL DO_THE_MODEL_IO( myTime, myIter, myThid )
284 CALL TIMER_STOP ('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
285
286 C-- Save state for restarts
287 C Note: (jmc: is it still the case after ckp35 ?)
288 C =====
289 C Because of the ordering of the timestepping code and
290 C tendency term code at end of loop model arrays hold
291 C U,V,T,S at "time-level" N but gu, gv, gs, gt, guNM1,...
292 C at "time-level" N+1/2 (guNM1 at "time-level" N+1/2 is
293 C gu at "time-level" N-1/2) and etaN at "time-level" N+1/2.
294 C where N = I+timeLevBase-1
295 C Thus a checkpoint contains U.0000000000, GU.0000000001 and
296 C etaN.0000000001 in the indexing scheme used for the model
297 C "state" files. This example is referred to as a checkpoint
298 C at time level 1
299 CALL TIMER_START('WRITE_CHECKPOINT [FORWARD_STEP]',myThid)
300 CALL WRITE_CHECKPOINT(
301 & .FALSE., myTime, myIter, myThid )
302 CALL TIMER_STOP ('WRITE_CHECKPOINT [FORWARD_STEP]',myThid)
303
304 END

  ViewVC Help
Powered by ViewVC 1.1.22