/[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.45 - (show annotations) (download)
Sat Dec 28 10:11:10 2002 UTC (21 years, 5 months ago) by dimitri
Branch: MAIN
CVS Tags: checkpoint47i_post, checkpoint47g_post, checkpoint47f_post, checkpoint47h_post
Changes since 1.44: +55 -49 lines
checkpoint47f_post
Merging from release1_p10:
o modifications for using pkg/exf with pkg/seaice
  - pkg/seaice CPP options SEAICE_EXTERNAL_FORCING
    and SEAICE_EXTERNAL_FLUXES
  - pkg/exf CPP options EXF_READ_EVAP and
    EXF_NO_BULK_COMPUTATIONS
  - usage examples are Experiments 8 and 9 in
    verification/lab_sea/README
  - verification/lab_sea default experiment now uses
    pkg/gmredi, pkg/kpp, pkg/seaice, and pkg/exf

1 C $Header: /u/gcmpack/MITgcm/model/src/forward_step.F,v 1.44 2002/12/21 00:33:54 jmc 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 # if (defined (ALLOW_BULKFORMULAE) || defined (ALLOW_BULK_FORCE))
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 #ifdef ALLOW_BULK_FORCE
79 INTEGER bi,bj
80 #endif ALLOW_BULK_FORCE
81
82 CEOP
83
84 #ifdef ALLOW_AUTODIFF_TAMC
85 C-- Reset the model iteration counter and the model time.
86 myiter = nIter0 + (iloop-1)
87 mytime = startTime + float(iloop-1)*deltaTclock
88 #endif
89
90 #if (defined (ALLOW_AUTODIFF_TAMC) && defined (ALLOW_AUTODIFF_MONITOR))
91 C Include call to a dummy routine. Its adjoint will be
92 C called at the proper place in the adjoint code.
93 C The adjoint routine will print out adjoint values
94 C if requested. The location of the call is important,
95 C it has to be after the adjoint of the exchanges
96 C (DO_GTERM_BLOCKING_EXCHANGES).
97 CALL DUMMY_IN_STEPPING( myTime, myIter, myThid )
98 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
99 #endif
100
101 #ifdef EXACT_CONSERV
102 IF (exactConserv) THEN
103 C-- Update etaH(n+1) :
104 CALL TIMER_START('UPDATE_ETAH [FORWARD_STEP]',mythid)
105 CALL UPDATE_ETAH( myTime, myIter, myThid )
106 CALL TIMER_STOP ('UPDATE_ETAH [FORWARD_STEP]',mythid)
107 ENDIF
108 #endif /* EXACT_CONSERV */
109
110 #ifdef NONLIN_FRSURF
111 C-- compute the future surface level thickness
112 C according to etaH(n+1)
113 IF ( nonlinFreeSurf.GT.0) THEN
114 CALL TIMER_START('CALC_SURF_DR [FORWARD_STEP]',mythid)
115 CALL CALC_SURF_DR(etaH, myTime, myIter, myThid )
116 CALL TIMER_STOP ('CALC_SURF_DR [FORWARD_STEP]',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
127
128 C-- Call external forcing package
129 cswdblk -- add ---
130 #ifdef ALLOW_BULK_FORCE
131 CALL TIMER_START('BULKF_FIELDS_LOAD[THE_MAIN_LOOP]',mythid)
132 CALL BULKF_FIELDS_LOAD( mytime, myiter, mythid )
133 CALL TIMER_STOP ('BULKF_FIELDS_LOAD[THE_MAIN_LOOP]',mythid)
134 c calculate qnet and empmr (and wind stress)
135 DO bj=myByLo(myThid),myByHi(myThid)
136 DO bi=myBxLo(myThid),myBxHi(myThid)
137 CALL BULKF_FORCING( bi,bj, mytime, myiter, mythid )
138 ENDDO
139 ENDDO
140 c Update the tile edges.
141 _EXCH_XY_R8(Qnet, mythid)
142 _EXCH_XY_R8(EmPmR, mythid)
143 CALL EXCH_UV_XY_RS(fu, fv, .TRUE., myThid)
144 C _EXCH_XY_R8(fu , mythid)
145 C _EXCH_XY_R8(fv , mythid)
146 cswdblk -- end add ---
147 #else
148 #ifdef INCLUDE_EXTERNAL_FORCING_PACKAGE
149 C NOTE, that although the exf package is part of the
150 C distribution, it is not currently maintained, i.e.
151 C exf is disabled by default in genmake.
152 CALL TIMER_START('EXF_GETFORCING [FORWARD_STEP]',mythid)
153 CALL EXF_GETFORCING( mytime, myiter, mythid )
154 CALL TIMER_STOP ('EXF_GETFORCING [FORWARD_STEP]',mythid)
155 #else
156 IF ( .NOT. useSEAICE ) THEN
157 CALL TIMER_START('EXTERNAL_FIELDS_LOAD[FORWARD_STEP]',mythid)
158 CALL EXTERNAL_FIELDS_LOAD( mytime, myiter, mythid )
159 CALL TIMER_STOP ('EXTERNAL_FIELDS_LOAD[FORWARD_STEP]',mythid)
160 ENDIF
161 #endif /* INCLUDE_EXTERNAL_FORCING_PACKAGE */
162
163 #ifdef ALLOW_SEAICE
164 C-- Call sea ice model to compute forcing/external data fields. In
165 C addition to computing prognostic sea-ice variables and diagnosing the
166 C forcing/external data fields that drive the ocean model, SEAICE_MODEL
167 C also sets theta to the freezing point under sea-ice. The implied
168 C surface heat flux is then stored in variable surfaceTendencyTice,
169 C which is needed by KPP package (kpp_calc.F and kpp_transport_t.F)
170 C to diagnose surface buoyancy fluxes and for the non-local transport
171 C term. Because this call precedes model thermodynamics, temperature
172 C under sea-ice may not be "exactly" at the freezing point by the time
173 C theta is dumped or time-averaged.
174 IF ( useSEAICE ) THEN
175 CALL TIMER_START('SEAICE_MODEL [FORWARD_STEP]',myThid)
176 CALL SEAICE_MODEL( myTime, myIter, myThid )
177 CALL TIMER_STOP ('SEAICE_MODEL [FORWARD_STEP]',myThid)
178 ENDIF
179 #endif ALLOW_SEAICE
180 #endif ALLOW_BULK_FORCE
181
182 C-- Step forward fields and calculate time tendency terms.
183 CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',mythid)
184 CALL THERMODYNAMICS( myTime, myIter, myThid )
185 CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',mythid)
186
187 C-- do exchanges (needed for DYNAMICS) when using stagger time-step :
188 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
189 CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid )
190 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
191
192 #ifdef ALLOW_SHAP_FILT
193 IF (useSHAP_FILT .AND.
194 & staggerTimeStep .AND. shap_filt_TrStagg ) THEN
195 CALL TIMER_START('SHAP_FILT [FORWARD_STEP]',myThid)
196 CALL SHAP_FILT_APPLY_TS( gT, gS, myTime, myIter, myThid )
197 CALL TIMER_STOP ('SHAP_FILT [FORWARD_STEP]',myThid)
198 ENDIF
199 #endif
200 #ifdef ALLOW_ZONAL_FILT
201 IF (useZONAL_FILT .AND.
202 & staggerTimeStep .AND. zonal_filt_TrStagg ) THEN
203 CALL TIMER_START('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid)
204 CALL ZONAL_FILT_APPLY_TS( gT, gS, myThid )
205 CALL TIMER_STOP ('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid)
206 ENDIF
207 #endif
208
209 C-- Step forward fields and calculate time tendency terms.
210 IF ( momStepping ) THEN
211 CALL TIMER_START('DYNAMICS [FORWARD_STEP]',mythid)
212 CALL DYNAMICS( myTime, myIter, myThid )
213 CALL TIMER_STOP ('DYNAMICS [FORWARD_STEP]',mythid)
214 ENDIF
215
216 #ifdef ALLOW_NONHYDROSTATIC
217 C-- Step forward W field in N-H algorithm
218 IF ( momStepping .AND. nonHydrostatic ) THEN
219 CALL TIMER_START('CALC_GW [FORWARD_STEP]',myThid)
220 CALL CALC_GW(myThid)
221 CALL TIMER_STOP ('CALC_GW [FORWARD_STEP]',myThid)
222 ENDIF
223 #endif
224
225 #ifdef NONLIN_FRSURF
226 C-- update hfacC,W,S and recip_hFac according to etaH(n+1) :
227 IF ( nonlinFreeSurf.GT.0) THEN
228 CALL TIMER_START('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
229 CALL UPDATE_SURF_DR( myTime, myIter, myThid )
230 CALL TIMER_STOP ('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
231 ENDIF
232 C- update also CG2D matrix (and preconditioner)
233 IF ( momStepping .AND. nonlinFreeSurf.GT.2 ) THEN
234 CALL TIMER_START('UPDATE_CG2D [FORWARD_STEP]',myThid)
235 CALL UPDATE_CG2D( myTime, myIter, myThid )
236 CALL TIMER_STOP ('UPDATE_CG2D [FORWARD_STEP]',myThid)
237 ENDIF
238 #endif
239
240 C-- Apply Filters to u*,v* before SOLVE_FOR_PRESSURE
241 #ifdef ALLOW_SHAP_FILT
242 IF (useSHAP_FILT .AND. shap_filt_uvStar) THEN
243 CALL TIMER_START('SHAP_FILT [FORWARD_STEP]',myThid)
244 CALL SHAP_FILT_APPLY_UV( gUnm1,gVnm1, myTime,myIter,myThid )
245 IF (implicDiv2Dflow.LT.1.) THEN
246 C-- Explicit+Implicit part of the Barotropic Flow Divergence
247 C => Filtering of uVel,vVel is necessary
248 CALL SHAP_FILT_APPLY_UV( uVel,vVel, myTime,myIter,myThid )
249 ENDIF
250 CALL TIMER_STOP ('SHAP_FILT [FORWARD_STEP]',myThid)
251 ENDIF
252 #endif
253 #ifdef ALLOW_ZONAL_FILT
254 IF (useZONAL_FILT .AND. zonal_filt_uvStar) THEN
255 CALL TIMER_START('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid)
256 CALL ZONAL_FILT_APPLY_UV( gUnm1, gVnm1, myThid )
257 IF (implicDiv2Dflow.LT.1.) THEN
258 C-- Explicit+Implicit part of the Barotropic Flow Divergence
259 C => Filtering of uVel,vVel is necessary
260 CALL ZONAL_FILT_APPLY_UV( uVel, vVel, myThid )
261 ENDIF
262 CALL TIMER_STOP ('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid)
263 ENDIF
264 #endif
265
266 C-- Solve elliptic equation(s).
267 C Two-dimensional only for conventional hydrostatic or
268 C three-dimensional for non-hydrostatic and/or IGW scheme.
269 IF ( momStepping ) THEN
270 CALL TIMER_START('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
271 CALL SOLVE_FOR_PRESSURE(myTime, myIter, myThid)
272 CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
273 ENDIF
274
275 #ifdef ALLOW_AUTODIFF_TAMC
276 cph This is needed because convective_adjustment calls
277 cph find_rho which may use pressure()
278 CADJ STORE pressure = comlev1, key = ikey_dynamics
279 #endif
280 C-- Correct divergence in flow field and cycle time-stepping
281 C arrays (for all fields) ; update time-counter
282 myIter = nIter0 + iLoop
283 myTime = startTime + deltaTClock * float(iLoop)
284 CALL TIMER_START('THE_CORRECTION_STEP [FORWARD_STEP]',myThid)
285 CALL THE_CORRECTION_STEP(myTime, myIter, myThid)
286 CALL TIMER_STOP ('THE_CORRECTION_STEP [FORWARD_STEP]',myThid)
287
288 C-- Do "blocking" sends and receives for tendency "overlap" terms
289 c CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
290 c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )
291 c CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
292
293 C-- Do "blocking" sends and receives for field "overlap" terms
294 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
295 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
296 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
297
298 #ifdef ALLOW_FLT
299 C-- Calculate float trajectories
300 IF (useFLT) THEN
301 CALL TIMER_START('FLOATS [FORWARD_STEP]',myThid)
302 CALL FLT_MAIN(myIter,myTime, myThid)
303 CALL TIMER_STOP ('FLOATS [FORWARD_STEP]',myThid)
304 ENDIF
305 #endif
306
307 #ifndef EXCLUDE_MONITOR
308 C-- Check status of solution (statistics, cfl, etc...)
309 CALL TIMER_START('MONITOR [FORWARD_STEP]',myThid)
310 CALL MONITOR( myIter, myTime, myThid )
311 CALL TIMER_STOP ('MONITOR [FORWARD_STEP]',myThid)
312 #endif /* EXCLUDE_MONITOR */
313
314 C-- Do IO if needed.
315 CALL TIMER_START('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
316 CALL DO_THE_MODEL_IO( myTime, myIter, myThid )
317 CALL TIMER_STOP ('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
318
319 C-- Save state for restarts
320 C Note: (jmc: is it still the case after ckp35 ?)
321 C =====
322 C Because of the ordering of the timestepping code and
323 C tendency term code at end of loop model arrays hold
324 C U,V,T,S at "time-level" N but gu, gv, gs, gt, guNM1,...
325 C at "time-level" N+1/2 (guNM1 at "time-level" N+1/2 is
326 C gu at "time-level" N-1/2) and etaN at "time-level" N+1/2.
327 C where N = I+timeLevBase-1
328 C Thus a checkpoint contains U.0000000000, GU.0000000001 and
329 C etaN.0000000001 in the indexing scheme used for the model
330 C "state" files. This example is referred to as a checkpoint
331 C at time level 1
332 CALL TIMER_START('WRITE_CHECKPOINT [FORWARD_STEP]',myThid)
333 CALL WRITE_CHECKPOINT(
334 & .FALSE., myTime, myIter, myThid )
335 CALL TIMER_STOP ('WRITE_CHECKPOINT [FORWARD_STEP]',myThid)
336
337 END

  ViewVC Help
Powered by ViewVC 1.1.22