/[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.97 - (show annotations) (download)
Tue Jul 6 21:12:51 2004 UTC (19 years, 10 months ago) by jmc
Branch: MAIN
Changes since 1.96: +16 -46 lines
put atmospheric physics & state-vars diagnostics calls in 2 dedicated S/R.

1 C $Header: /u/gcmpack/MITgcm/model/src/forward_step.F,v 1.96 2004/07/06 01:05:53 jmc Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 cswdptr -- add --
8 #ifdef ALLOW_GCHEM
9 # include "GCHEM_OPTIONS.h"
10 #endif
11 cswdptr -- end add ---
12
13 CBOP
14 C !ROUTINE: FORWARD_STEP
15 C !INTERFACE:
16 SUBROUTINE FORWARD_STEP( iloop, myTime, myIter, myThid )
17
18 C !DESCRIPTION: \bv
19 C *==================================================================
20 C | SUBROUTINE forward_step
21 C | o Run the ocean model and, optionally, evaluate a cost function.
22 C *==================================================================
23 C |
24 C | THE_MAIN_LOOP is the toplevel routine for the Tangent Linear and
25 C | Adjoint Model Compiler (TAMC). For this purpose the initialization
26 C | of the model was split into two parts. Those parameters that do
27 C | not depend on a specific model run are set in INITIALISE_FIXED,
28 C | whereas those that do depend on the specific realization are
29 C | initialized in INITIALISE_VARIA.
30 C |
31 C *==================================================================
32 C \ev
33
34 C !USES:
35 IMPLICIT NONE
36 C == Global variables ==
37 #include "SIZE.h"
38 #include "EEPARAMS.h"
39 #include "PARAMS.h"
40 #include "DYNVARS.h"
41
42 #ifdef ALLOW_SHAP_FILT
43 # include "SHAP_FILT.h"
44 #endif
45 #ifdef ALLOW_ZONAL_FILT
46 # include "ZONAL_FILT.h"
47 #endif
48 #ifdef COMPONENT_MODULE
49 # include "CPL_PARAMS.h"
50 #endif
51
52 #ifdef ALLOW_AUTODIFF_TAMC
53 # include "FFIELDS.h"
54
55 # ifdef ALLOW_NONHYDROSTATIC
56 # include "CG3D.h"
57 # endif
58
59 # include "tamc.h"
60 # include "ctrl.h"
61 # include "ctrl_dummy.h"
62 # include "cost.h"
63 # include "EOS.h"
64 # ifdef ALLOW_EXF
65 # include "exf_fields.h"
66 # include "exf_clim_fields.h"
67 # ifdef ALLOW_BULKFORMULAE
68 # include "exf_constants.h"
69 # endif
70 # endif
71 # ifdef ALLOW_OBCS
72 # include "OBCS.h"
73 # endif
74 # ifdef ALLOW_PTRACERS
75 # include "PTRACERS.h"
76 # endif
77 # ifdef ALLOW_CD_CODE
78 # include "CD_CODE_VARS.h"
79 # endif
80 # ifdef ALLOW_EBM
81 # include "EBM.h"
82 # endif
83 #endif /* ALLOW_AUTODIFF_TAMC */
84
85 C !LOCAL VARIABLES:
86 C == Routine arguments ==
87 C note: under the multi-threaded model myiter and
88 C mytime are local variables passed around as routine
89 C arguments. Although this is fiddly it saves the need to
90 C impose additional synchronisation points when they are
91 C updated.
92 C myIter - iteration counter for this thread
93 C myTime - time counter for this thread
94 C myThid - thread number for this instance of the routine.
95 INTEGER iloop
96 INTEGER myThid
97 INTEGER myIter
98 _RL myTime
99
100 C == Local variables ==
101 INTEGER myItP1
102 CEOP
103
104 #ifdef ALLOW_DEBUG
105 IF ( debugLevel .GE. debLevB )
106 & CALL DEBUG_ENTER('FORWARD_STEP',myThid)
107 #endif
108
109 #ifdef ALLOW_AUTODIFF_TAMC
110 C-- Reset the model iteration counter and the model time.
111 myiter = nIter0 + (iloop-1)
112 mytime = startTime + float(iloop-1)*deltaTclock
113 #endif
114
115 #if (defined (ALLOW_AUTODIFF_TAMC) && defined (ALLOW_AUTODIFF_MONITOR))
116 C Include call to a dummy routine. Its adjoint will be
117 C called at the proper place in the adjoint code.
118 C The adjoint routine will print out adjoint values
119 C if requested. The location of the call is important,
120 C it has to be after the adjoint of the exchanges
121 C (DO_GTERM_BLOCKING_EXCHANGES).
122 CALL DUMMY_IN_STEPPING( myTime, myIter, myThid )
123 cph I've commented this line since it may conflict with MITgcm's adjoint
124 cph However, need to check whether that's still consistent
125 cph with the ecco-branch (it should).
126 cph CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
127 #endif
128
129 #ifdef ALLOW_AUTODIFF_TAMC
130 c**************************************
131 #include "checkpoint_lev1_directives.h"
132 c**************************************
133 #endif
134
135 C-- Call external forcing package
136 #ifdef ALLOW_BULK_FORCE
137 IF ( useBulkForce ) THEN
138 #ifdef ALLOW_DEBUG
139 IF ( debugLevel .GE. debLevB )
140 & CALL DEBUG_CALL('BULKF_FIELDS_LOAD',myThid)
141 #endif
142 CALL TIMER_START('BULKF_FORCING [FORWARD_STEP]',mythid)
143 C- load all forcing fields at current time
144 CALL BULKF_FIELDS_LOAD( myTime, myIter, myThid )
145 C- calculate qnet and empmr (and wind stress)
146 CALL BULKF_FORCING( myTime, myIter, myThid )
147 CALL TIMER_STOP ('BULKF_FORCING [FORWARD_STEP]',mythid)
148 ELSE
149 #endif /* ALLOW_BULK_FORCE */
150
151 # ifdef ALLOW_EXF
152 # ifdef ALLOW_DEBUG
153 IF ( debugLevel .GE. debLevB )
154 & CALL DEBUG_CALL('EXF_GETFORCING',myThid)
155 # endif
156 CALL TIMER_START('EXF_GETFORCING [FORWARD_STEP]',mythid)
157 CALL EXF_GETFORCING( mytime, myiter, mythid )
158 CALL TIMER_STOP ('EXF_GETFORCING [FORWARD_STEP]',mythid)
159 # else /* ALLOW_EXF undef */
160 cph The following IF-statement creates an additional dependency
161 cph for the forcing fields requiring additional storing.
162 cph Therefore, the IF-statement will be put between CPP-OPTIONS,
163 cph assuming that ALLOW_SEAICE has not yet been differentiated.
164 # if (defined (ALLOW_SEAICE) || defined (ALLOW_EBM))
165 IF ( .NOT. useSEAICE .AND. .NOT. useEBM ) THEN
166 # endif
167 #ifdef ALLOW_DEBUG
168 IF ( debugLevel .GE. debLevB )
169 & CALL DEBUG_CALL('EXTERNAL_FIELDS_LOAD',myThid)
170 #endif
171 CALL TIMER_START('EXTERNAL_FIELDS_LOAD[FORWARD_STEP]',mythid)
172 CALL EXTERNAL_FIELDS_LOAD( mytime, myiter, mythid )
173 CALL TIMER_STOP ('EXTERNAL_FIELDS_LOAD[FORWARD_STEP]',mythid)
174 # if (defined (ALLOW_SEAICE) || defined (ALLOW_EBM))
175 ENDIF
176 # endif
177 # endif /* ALLOW_EXF */
178 #ifdef ALLOW_BULK_FORCE
179 C-- end of if/else block useBulfforce --
180 ENDIF
181 #endif /* ALLOW_BULK_FORCE */
182
183 #ifdef ALLOW_AUTODIFF
184 c-- Add control vector for forcing and parameter fields
185 if ( myiter .EQ. nIter0 )
186 & CALL CTRL_MAP_FORCING (mythid)
187 #endif
188
189 #ifdef ALLOW_THSICE
190 IF ( useThSIce .AND. buoyancyRelation(1:7) .EQ. 'OCEANIC' ) THEN
191 #ifdef ALLOW_DEBUG
192 IF ( debugLevel .GE. debLevB )
193 & CALL DEBUG_CALL('THSICE_MAIN',myThid)
194 #endif
195 C-- Step forward Therm.Sea-Ice variables
196 C and modify forcing terms including effects from ice
197 CALL TIMER_START('THSICE_MAIN [FORWARD_STEP]', myThid)
198 CALL THSICE_MAIN( myTime, myIter, myThid )
199 CALL TIMER_STOP( 'THSICE_MAIN [FORWARD_STEP]', myThid)
200 ENDIF
201 #endif /* ALLOW_THSICE */
202
203 # ifdef ALLOW_SEAICE
204 C-- Call sea ice model to compute forcing/external data fields. In
205 C addition to computing prognostic sea-ice variables and diagnosing the
206 C forcing/external data fields that drive the ocean model, SEAICE_MODEL
207 C also sets theta to the freezing point under sea-ice. The implied
208 C surface heat flux is then stored in variable surfaceTendencyTice,
209 C which is needed by KPP package (kpp_calc.F and kpp_transport_t.F)
210 C to diagnose surface buoyancy fluxes and for the non-local transport
211 C term. Because this call precedes model thermodynamics, temperature
212 C under sea-ice may not be "exactly" at the freezing point by the time
213 C theta is dumped or time-averaged.
214 IF ( useSEAICE ) THEN
215 #ifdef ALLOW_DEBUG
216 IF ( debugLevel .GE. debLevB )
217 & CALL DEBUG_CALL('SEAICE_MODEL',myThid)
218 #endif
219 CALL TIMER_START('SEAICE_MODEL [FORWARD_STEP]',myThid)
220 CALL SEAICE_MODEL( myTime, myIter, myThid )
221 CALL TIMER_STOP ('SEAICE_MODEL [FORWARD_STEP]',myThid)
222 ENDIF
223 # endif /* ALLOW_SEAICE */
224
225 C-- Freeze water at the surface
226 #ifdef ALLOW_AUTODIFF_TAMC
227 CADJ STORE theta = comlev1, key = ikey_dynamics
228 #endif
229 IF ( allowFreezing .AND. .NOT. useSEAICE
230 & .AND. .NOT. useThSIce ) THEN
231 CALL FREEZE_SURFACE( myTime, myIter, myThid )
232 ENDIF
233
234 #ifdef ALLOW_AUTODIFF_TAMC
235 # ifdef ALLOW_PTRACERS
236 cph this replaces _bibj storing of ptracer within thermodynamics
237 CADJ STORE ptracer = comlev1, key = ikey_dynamics
238 # endif
239 #endif
240
241 #ifdef ALLOW_PTRACERS
242 # ifdef ALLOW_GCHEM
243 CALL GCHEM_FIELDS_LOAD( mytime, myiter, mythid )
244 # endif
245 #endif
246
247 #ifdef COMPONENT_MODULE
248 IF ( useCoupler .AND. cpl_earlyExpImpCall ) THEN
249 C Post coupling data that I export.
250 C Read in coupling data that I import.
251 CALL TIMER_START('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid)
252 CALL CPL_EXPORT_MY_DATA( myIter, myTime, myThid )
253 CALL CPL_IMPORT_EXTERNAL_DATA( myIter, myTime, myThid )
254 CALL TIMER_STOP ('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid)
255 ENDIF
256 # ifndef ALLOW_AIM
257 C jmc: don't know precisely where to put this call. leave it here for now.
258 IF ( useCoupler ) THEN
259 CALL OCN_APPLY_IMPORT( .TRUE., myTime, myIter, myThid )
260 ENDIF
261 # endif
262 #endif /* COMPONENT_MODULE */
263
264 #ifdef ALLOW_EBM
265 IF ( useEBM ) THEN
266 # ifdef ALLOW_DEBUG
267 IF ( debugLevel .GE. debLevB )
268 & CALL DEBUG_CALL('EBM',myThid)
269 # endif
270 CALL TIMER_START('EBM [FORWARD_STEP]',mythid)
271 CALL EBM_DRIVER ( myTime, myIter, myThid )
272 CALL TIMER_STOP ('EBM [FORWARD_STEP]',mythid)
273 ENDIF
274 #endif
275
276 C-- Step forward fields and calculate time tendency terms.
277
278 #ifdef ALLOW_DEBUG
279 IF ( debugLevel .GE. debLevB )
280 & CALL DEBUG_CALL('DO_ATMOSPHERIC_PHYS',myThid)
281 #endif
282 CALL TIMER_START('DO_ATMOSPHERIC_PHYS [FORWARD_STEP]',mythid)
283 CALL DO_ATMOSPHERIC_PHYS( myTime, myIter, myThid )
284 CALL TIMER_STOP ('DO_ATMOSPHERIC_PHYS [FORWARD_STEP]',mythid)
285
286 #ifdef ALLOW_DEBUG
287 IF ( debugLevel .GE. debLevB )
288 & CALL DEBUG_CALL('DO_OCEANIC_PHYS',myThid)
289 #endif
290 CALL TIMER_START('DO_OCEANIC_PHYS [FORWARD_STEP]',mythid)
291 CALL DO_OCEANIC_PHYS( myTime, myIter, myThid )
292 CALL TIMER_STOP ('DO_OCEANIC_PHYS [FORWARD_STEP]',mythid)
293
294 IF ( .NOT.staggerTimeStep ) THEN
295 #ifdef ALLOW_DEBUG
296 IF ( debugLevel .GE. debLevB )
297 & CALL DEBUG_CALL('THERMODYNAMICS',myThid)
298 #endif
299 CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',mythid)
300 CALL THERMODYNAMICS( myTime, myIter, myThid )
301 CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',mythid)
302 C-- if not staggerTimeStep: end
303 ENDIF
304
305 #ifdef COMPONENT_MODULE
306 IF ( useCoupler .AND. .NOT.cpl_earlyExpImpCall ) THEN
307 C Post coupling data that I export.
308 C Read in coupling data that I import.
309 myItP1 = myIter + 1
310 CALL TIMER_START('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid)
311 CALL CPL_EXPORT_MY_DATA( myItP1, myTime, myThid )
312 CALL CPL_IMPORT_EXTERNAL_DATA( myItP1, myTime, myThid )
313 CALL TIMER_STOP ('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid)
314 # ifndef ALLOW_AIM
315 IF ( useRealFreshWaterFlux ) THEN
316 CALL OCN_APPLY_IMPORT( .FALSE., myTime, myIter, myThid )
317 ENDIF
318 # endif
319 ENDIF
320 #endif /* COMPONENT_MODULE */
321
322 C-- Step forward fields and calculate time tendency terms.
323 #ifndef ALLOW_AUTODIFF_TAMC
324 IF ( momStepping ) THEN
325 #endif
326 #ifdef ALLOW_DEBUG
327 IF ( debugLevel .GE. debLevB )
328 & CALL DEBUG_CALL('DYNAMICS',myThid)
329 #endif
330 CALL TIMER_START('DYNAMICS [FORWARD_STEP]',mythid)
331 CALL DYNAMICS( myTime, myIter, myThid )
332 CALL TIMER_STOP ('DYNAMICS [FORWARD_STEP]',mythid)
333 #ifndef ALLOW_AUTODIFF_TAMC
334 ENDIF
335 #endif
336
337 #ifdef ALLOW_NONHYDROSTATIC
338 C-- Step forward W field in N-H algorithm
339 IF ( momStepping .AND. nonHydrostatic ) THEN
340 #ifdef ALLOW_DEBUG
341 IF ( debugLevel .GE. debLevB )
342 & CALL DEBUG_CALL('CALC_GW',myThid)
343 #endif
344 CALL TIMER_START('CALC_GW [FORWARD_STEP]',myThid)
345 CALL CALC_GW(myThid)
346 CALL TIMER_STOP ('CALC_GW [FORWARD_STEP]',myThid)
347 ENDIF
348 #endif
349
350 C-- Update time-counter
351 myIter = nIter0 + iLoop
352 myTime = startTime + deltaTClock * float(iLoop)
353
354 C-- Update geometric factors:
355 #ifdef NONLIN_FRSURF
356 C- update hfacC,W,S and recip_hFac according to etaH(n+1) :
357 IF ( nonlinFreeSurf.GT.0) THEN
358 IF ( select_rStar.GT.0 ) THEN
359 CALL TIMER_START('UPDATE_R_STAR [FORWARD_STEP]',myThid)
360 CALL UPDATE_R_STAR( myTime, myIter, myThid )
361 CALL TIMER_STOP ('UPDATE_R_STAR [FORWARD_STEP]',myThid)
362 ELSE
363 CALL TIMER_START('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
364 CALL UPDATE_SURF_DR( myTime, myIter, myThid )
365 CALL TIMER_STOP ('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
366 ENDIF
367 ENDIF
368 C- update also CG2D matrix (and preconditioner)
369 IF ( momStepping .AND. nonlinFreeSurf.GT.2 ) THEN
370 CALL TIMER_START('UPDATE_CG2D [FORWARD_STEP]',myThid)
371 CALL UPDATE_CG2D( myTime, myIter, myThid )
372 CALL TIMER_STOP ('UPDATE_CG2D [FORWARD_STEP]',myThid)
373 ENDIF
374 #endif
375
376 C-- Apply Filters to u*,v* before SOLVE_FOR_PRESSURE
377 #ifdef ALLOW_SHAP_FILT
378 IF (useSHAP_FILT .AND. shap_filt_uvStar) THEN
379 CALL TIMER_START('SHAP_FILT [FORWARD_STEP]',myThid)
380 IF (implicDiv2Dflow.LT.1.) THEN
381 C-- Explicit+Implicit part of the Barotropic Flow Divergence
382 C => Filtering of uVel,vVel is necessary
383 CALL SHAP_FILT_APPLY_UV( uVel,vVel,
384 & myTime, myIter, myThid )
385 ENDIF
386 CALL SHAP_FILT_APPLY_UV( gU,gV,myTime,myIter,myThid)
387 CALL TIMER_STOP ('SHAP_FILT [FORWARD_STEP]',myThid)
388 ENDIF
389 #endif
390 #ifdef ALLOW_ZONAL_FILT
391 IF (useZONAL_FILT .AND. zonal_filt_uvStar) THEN
392 CALL TIMER_START('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid)
393 IF (implicDiv2Dflow.LT.1.) THEN
394 C-- Explicit+Implicit part of the Barotropic Flow Divergence
395 C => Filtering of uVel,vVel is necessary
396 CALL ZONAL_FILT_APPLY_UV( uVel, vVel, myThid )
397 ENDIF
398 CALL ZONAL_FILT_APPLY_UV( gU, gV, myThid )
399 CALL TIMER_STOP ('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid)
400 ENDIF
401 #endif
402
403 C-- Solve elliptic equation(s).
404 C Two-dimensional only for conventional hydrostatic or
405 C three-dimensional for non-hydrostatic and/or IGW scheme.
406 IF ( momStepping ) THEN
407 CALL TIMER_START('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
408 CALL SOLVE_FOR_PRESSURE(myTime, myIter, myThid)
409 CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
410 ENDIF
411
412 C-- Correct divergence in flow field and cycle time-stepping momentum
413 c IF ( momStepping ) THEN
414 CALL TIMER_START('UV_CORRECTION_STEP [FORWARD_STEP]',myThid)
415 CALL MOMENTUM_CORRECTION_STEP(myTime, myIter, myThid)
416 CALL TIMER_STOP ('UV_CORRECTION_STEP [FORWARD_STEP]',myThid)
417 c ENDIF
418
419 #ifdef EXACT_CONSERV
420 IF (exactConserv) THEN
421 C-- Update etaH(n+1) :
422 CALL TIMER_START('UPDATE_ETAH [FORWARD_STEP]',mythid)
423 CALL UPDATE_ETAH( myTime, myIter, myThid )
424 CALL TIMER_STOP ('UPDATE_ETAH [FORWARD_STEP]',mythid)
425 ENDIF
426 #endif /* EXACT_CONSERV */
427
428 #ifdef NONLIN_FRSURF
429 IF ( select_rStar.NE.0 ) THEN
430 C-- r* : compute the future level thickness according to etaH(n+1)
431 CALL TIMER_START('CALC_R_STAR [FORWARD_STEP]',mythid)
432 CALL CALC_R_STAR(etaH, myTime, myIter, myThid )
433 CALL TIMER_STOP ('CALC_R_STAR [FORWARD_STEP]',mythid)
434 ELSEIF ( nonlinFreeSurf.GT.0) THEN
435 C-- compute the future surface level thickness according to etaH(n+1)
436 CALL TIMER_START('CALC_SURF_DR [FORWARD_STEP]',mythid)
437 CALL CALC_SURF_DR(etaH, myTime, myIter, myThid )
438 CALL TIMER_STOP ('CALC_SURF_DR [FORWARD_STEP]',mythid)
439 ENDIF
440 #endif /* NONLIN_FRSURF */
441
442 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
443 IF ( staggerTimeStep ) THEN
444
445 #ifdef ALLOW_DEBUG
446 IF ( debugLevel .GE. debLevB )
447 & CALL DEBUG_CALL('THERMODYNAMICS',myThid)
448 #endif
449 CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',mythid)
450 CALL THERMODYNAMICS( myTime, myIter, myThid )
451 CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',mythid)
452
453 C-- if staggerTimeStep: end
454 ENDIF
455 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
456
457 #ifdef ALLOW_AUTODIFF_TAMC
458 cph This is needed because convective_adjustment calls
459 cph find_rho which may use pressure()
460 CADJ STORE totphihyd = comlev1, key = ikey_dynamics
461 #endif
462 C-- Cycle time-stepping Tracers arrays (T,S,+pTracers)
463 CALL TIMER_START('TS_CORRECTION_STEP [FORWARD_STEP]',myThid)
464 CALL TRACERS_CORRECTION_STEP(myTime, myIter, myThid)
465 CALL TIMER_STOP ('TS_CORRECTION_STEP [FORWARD_STEP]',myThid)
466
467 C-- Do "blocking" sends and receives for tendency "overlap" terms
468 c CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
469 c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )
470 c CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
471
472 C-- Do "blocking" sends and receives for field "overlap" terms
473 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
474 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
475 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
476
477 cswdptr -- add for seperate timestepping of chemical/biological/forcing
478 cswdptr of ptracers ---
479 #ifdef ALLOW_GCHEM
480 ceh3 This is broken -- this ifdef should not be visible!
481 #ifdef PTRACERS_SEPARATE_FORCING
482 ceh3 needs an IF ( use GCHEM ) THEN
483 call GCHEM_FORCING_SEP( myTime,myIter,myThid )
484 #endif /* PTRACERS_SEPARATE_FORCING */
485 #endif /* ALLOW_GCHEM */
486 cswdptr -- end add ---
487
488 C AMM
489 #ifdef ALLOW_GRIDALT
490 if (useGRIDALT) then
491 CALL GRIDALT_UPDATE(myThid)
492 endif
493 #endif
494 C AMM
495
496 C AMM
497 #ifdef ALLOW_FIZHI
498 if( useFIZHI) then
499 CALL TIMER_START('FIZHI [FORWARD_STEP]',mythid)
500 CALL STEP_FIZHI_CORR ( myTime, myIter, myThid )
501 CALL TIMER_STOP('FIZHI [FORWARD_STEP]',mythid)
502 endif
503 #endif
504 C AMM
505
506 #ifdef ALLOW_FLT
507 C-- Calculate float trajectories
508 IF (useFLT) THEN
509 CALL TIMER_START('FLOATS [FORWARD_STEP]',myThid)
510 CALL FLT_MAIN(myIter,myTime, myThid)
511 CALL TIMER_STOP ('FLOATS [FORWARD_STEP]',myThid)
512 ENDIF
513 #endif
514
515 C-- State-variables statistics (time-aver, diagnostics ...)
516 CALL TIMER_START('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
517 CALL DO_STATEVARS_DIAGS( myTime, myIter, myThid )
518 CALL TIMER_STOP ('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
519
520 #ifdef ALLOW_MONITOR
521 C-- Check status of solution (statistics, cfl, etc...)
522 CALL TIMER_START('MONITOR [FORWARD_STEP]',myThid)
523 CALL MONITOR( myIter, myTime, myThid )
524 CALL TIMER_STOP ('MONITOR [FORWARD_STEP]',myThid)
525 #endif /* ALLOW_MONITOR */
526
527 C-- Do IO if needed.
528 CALL TIMER_START('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
529 CALL DO_THE_MODEL_IO( myTime, myIter, myThid )
530 CALL TIMER_STOP ('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
531
532 C-- Save state for restarts
533 CALL TIMER_START('WRITE_CHECKPOINT [FORWARD_STEP]',myThid)
534 CALL PACKAGES_WRITE_PICKUP(
535 I .FALSE., myTime, myIter, myThid )
536 CALL WRITE_CHECKPOINT(
537 I .FALSE., myTime, myIter, myThid )
538 CALL TIMER_STOP ('WRITE_CHECKPOINT [FORWARD_STEP]',myThid)
539
540 #ifdef ALLOW_DEBUG
541 IF ( debugLevel .GE. debLevB )
542 & CALL DEBUG_LEAVE('FORWARD_STEP',myThid)
543 #endif
544
545 RETURN
546 END

  ViewVC Help
Powered by ViewVC 1.1.22