/[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.96 - (show annotations) (download)
Tue Jul 6 01:05:53 2004 UTC (19 years, 11 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint54a_post
Changes since 1.95: +72 -73 lines
re-write staggerTimeStep: step forward momentum 1rst and then T,S

1 C $Header: /u/gcmpack/MITgcm/model/src/forward_step.F,v 1.95 2004/06/14 22:42:00 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 C-- Step forward fields and calculate time tendency terms.
265
266 C AMM
267 #ifdef ALLOW_FIZHI
268 if( useFIZHI) then
269 CALL TIMER_START('FIZHI [FORWARD_STEP]',mythid)
270 CALL FIZHI_UPDATE_TIME ( myIter, deltaTclock )
271 CALL UPDATE_OCEAN_EXPORTS ( myTime, myIter, myThid )
272 CALL UPDATE_EARTH_EXPORTS ( myTime, myIter, myThid )
273 CALL UPDATE_CHEMISTRY_EXPORTS ( myTime, myIter, myThid )
274 CALL FIZHI_WRAPPER ( myTime, myIter, myThid )
275 CALL STEP_FIZHI_FG ( myTime, myIter, myThid, deltaTtracer )
276 CALL TIMER_STOP ('FIZHI [FORWARD_STEP]',mythid)
277 endif
278 #endif
279 C AMM
280
281 #ifdef ALLOW_EBM
282 IF ( useEBM ) THEN
283 # ifdef ALLOW_DEBUG
284 IF ( debugLevel .GE. debLevB )
285 & CALL DEBUG_CALL('EBM',myThid)
286 # endif
287 CALL TIMER_START('EBM [FORWARD_STEP]',mythid)
288 CALL EBM_DRIVER ( myTime, myIter, myThid )
289 CALL TIMER_STOP ('EBM [FORWARD_STEP]',mythid)
290 ENDIF
291 #endif
292
293 #ifdef ALLOW_DEBUG
294 IF ( debugLevel .GE. debLevB )
295 & CALL DEBUG_CALL('DO_OCEANIC_PHYS',myThid)
296 #endif
297 CALL TIMER_START('DO_OCEANIC_PHYS [FORWARD_STEP]',mythid)
298 CALL DO_OCEANIC_PHYS( myTime, myIter, myThid )
299 CALL TIMER_STOP ('DO_OCEANIC_PHYS [FORWARD_STEP]',mythid)
300
301 IF ( .NOT.staggerTimeStep ) THEN
302 #ifdef ALLOW_DEBUG
303 IF ( debugLevel .GE. debLevB )
304 & CALL DEBUG_CALL('THERMODYNAMICS',myThid)
305 #endif
306 CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',mythid)
307 CALL THERMODYNAMICS( myTime, myIter, myThid )
308 CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',mythid)
309 C-- if not staggerTimeStep: end
310 ENDIF
311
312 #ifdef COMPONENT_MODULE
313 IF ( useCoupler .AND. .NOT.cpl_earlyExpImpCall ) THEN
314 C Post coupling data that I export.
315 C Read in coupling data that I import.
316 myItP1 = myIter + 1
317 CALL TIMER_START('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid)
318 CALL CPL_EXPORT_MY_DATA( myItP1, myTime, myThid )
319 CALL CPL_IMPORT_EXTERNAL_DATA( myItP1, myTime, myThid )
320 CALL TIMER_STOP ('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid)
321 # ifndef ALLOW_AIM
322 IF ( useRealFreshWaterFlux ) THEN
323 CALL OCN_APPLY_IMPORT( .FALSE., myTime, myIter, myThid )
324 ENDIF
325 # endif
326 ENDIF
327 #endif /* COMPONENT_MODULE */
328
329 C-- Step forward fields and calculate time tendency terms.
330 #ifndef ALLOW_AUTODIFF_TAMC
331 IF ( momStepping ) THEN
332 #endif
333 #ifdef ALLOW_DEBUG
334 IF ( debugLevel .GE. debLevB )
335 & CALL DEBUG_CALL('DYNAMICS',myThid)
336 #endif
337 CALL TIMER_START('DYNAMICS [FORWARD_STEP]',mythid)
338 CALL DYNAMICS( myTime, myIter, myThid )
339 CALL TIMER_STOP ('DYNAMICS [FORWARD_STEP]',mythid)
340 #ifndef ALLOW_AUTODIFF_TAMC
341 ENDIF
342 #endif
343
344 #ifdef ALLOW_NONHYDROSTATIC
345 C-- Step forward W field in N-H algorithm
346 IF ( momStepping .AND. nonHydrostatic ) THEN
347 #ifdef ALLOW_DEBUG
348 IF ( debugLevel .GE. debLevB )
349 & CALL DEBUG_CALL('CALC_GW',myThid)
350 #endif
351 CALL TIMER_START('CALC_GW [FORWARD_STEP]',myThid)
352 CALL CALC_GW(myThid)
353 CALL TIMER_STOP ('CALC_GW [FORWARD_STEP]',myThid)
354 ENDIF
355 #endif
356
357 C-- Update time-counter
358 myIter = nIter0 + iLoop
359 myTime = startTime + deltaTClock * float(iLoop)
360
361 C-- Update geometric factors:
362 #ifdef NONLIN_FRSURF
363 C- update hfacC,W,S and recip_hFac according to etaH(n+1) :
364 IF ( nonlinFreeSurf.GT.0) THEN
365 IF ( select_rStar.GT.0 ) THEN
366 CALL TIMER_START('UPDATE_R_STAR [FORWARD_STEP]',myThid)
367 CALL UPDATE_R_STAR( myTime, myIter, myThid )
368 CALL TIMER_STOP ('UPDATE_R_STAR [FORWARD_STEP]',myThid)
369 ELSE
370 CALL TIMER_START('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
371 CALL UPDATE_SURF_DR( myTime, myIter, myThid )
372 CALL TIMER_STOP ('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
373 ENDIF
374 ENDIF
375 C- update also CG2D matrix (and preconditioner)
376 IF ( momStepping .AND. nonlinFreeSurf.GT.2 ) THEN
377 CALL TIMER_START('UPDATE_CG2D [FORWARD_STEP]',myThid)
378 CALL UPDATE_CG2D( myTime, myIter, myThid )
379 CALL TIMER_STOP ('UPDATE_CG2D [FORWARD_STEP]',myThid)
380 ENDIF
381 #endif
382
383 C-- Apply Filters to u*,v* before SOLVE_FOR_PRESSURE
384 #ifdef ALLOW_SHAP_FILT
385 IF (useSHAP_FILT .AND. shap_filt_uvStar) THEN
386 CALL TIMER_START('SHAP_FILT [FORWARD_STEP]',myThid)
387 IF (implicDiv2Dflow.LT.1.) THEN
388 C-- Explicit+Implicit part of the Barotropic Flow Divergence
389 C => Filtering of uVel,vVel is necessary
390 CALL SHAP_FILT_APPLY_UV( uVel,vVel,
391 & myTime, myIter, myThid )
392 ENDIF
393 CALL SHAP_FILT_APPLY_UV( gU,gV,myTime,myIter,myThid)
394 CALL TIMER_STOP ('SHAP_FILT [FORWARD_STEP]',myThid)
395 ENDIF
396 #endif
397 #ifdef ALLOW_ZONAL_FILT
398 IF (useZONAL_FILT .AND. zonal_filt_uvStar) THEN
399 CALL TIMER_START('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid)
400 IF (implicDiv2Dflow.LT.1.) THEN
401 C-- Explicit+Implicit part of the Barotropic Flow Divergence
402 C => Filtering of uVel,vVel is necessary
403 CALL ZONAL_FILT_APPLY_UV( uVel, vVel, myThid )
404 ENDIF
405 CALL ZONAL_FILT_APPLY_UV( gU, gV, myThid )
406 CALL TIMER_STOP ('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid)
407 ENDIF
408 #endif
409
410 C-- Solve elliptic equation(s).
411 C Two-dimensional only for conventional hydrostatic or
412 C three-dimensional for non-hydrostatic and/or IGW scheme.
413 IF ( momStepping ) THEN
414 CALL TIMER_START('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
415 CALL SOLVE_FOR_PRESSURE(myTime, myIter, myThid)
416 CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
417 ENDIF
418
419 C-- Correct divergence in flow field and cycle time-stepping momentum
420 c IF ( momStepping ) THEN
421 CALL TIMER_START('UV_CORRECTION_STEP [FORWARD_STEP]',myThid)
422 CALL MOMENTUM_CORRECTION_STEP(myTime, myIter, myThid)
423 CALL TIMER_STOP ('UV_CORRECTION_STEP [FORWARD_STEP]',myThid)
424 c ENDIF
425
426 #ifdef EXACT_CONSERV
427 IF (exactConserv) THEN
428 C-- Update etaH(n+1) :
429 CALL TIMER_START('UPDATE_ETAH [FORWARD_STEP]',mythid)
430 CALL UPDATE_ETAH( myTime, myIter, myThid )
431 CALL TIMER_STOP ('UPDATE_ETAH [FORWARD_STEP]',mythid)
432 ENDIF
433 #endif /* EXACT_CONSERV */
434
435 #ifdef NONLIN_FRSURF
436 IF ( select_rStar.NE.0 ) THEN
437 C-- r* : compute the future level thickness according to etaH(n+1)
438 CALL TIMER_START('CALC_R_STAR [FORWARD_STEP]',mythid)
439 CALL CALC_R_STAR(etaH, myTime, myIter, myThid )
440 CALL TIMER_STOP ('CALC_R_STAR [FORWARD_STEP]',mythid)
441 ELSEIF ( nonlinFreeSurf.GT.0) THEN
442 C-- compute the future surface level thickness according to etaH(n+1)
443 CALL TIMER_START('CALC_SURF_DR [FORWARD_STEP]',mythid)
444 CALL CALC_SURF_DR(etaH, myTime, myIter, myThid )
445 CALL TIMER_STOP ('CALC_SURF_DR [FORWARD_STEP]',mythid)
446 ENDIF
447 #endif /* NONLIN_FRSURF */
448
449 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
450 IF ( staggerTimeStep ) THEN
451
452 #ifdef ALLOW_DEBUG
453 IF ( debugLevel .GE. debLevB )
454 & CALL DEBUG_CALL('THERMODYNAMICS',myThid)
455 #endif
456 CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',mythid)
457 CALL THERMODYNAMICS( myTime, myIter, myThid )
458 CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',mythid)
459
460 C-- if staggerTimeStep: end
461 ENDIF
462 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
463
464 #ifdef ALLOW_AUTODIFF_TAMC
465 cph This is needed because convective_adjustment calls
466 cph find_rho which may use pressure()
467 CADJ STORE totphihyd = comlev1, key = ikey_dynamics
468 #endif
469 C-- Cycle time-stepping Tracers arrays (T,S,+pTracers)
470 CALL TIMER_START('TS_CORRECTION_STEP [FORWARD_STEP]',myThid)
471 CALL TRACERS_CORRECTION_STEP(myTime, myIter, myThid)
472 CALL TIMER_STOP ('TS_CORRECTION_STEP [FORWARD_STEP]',myThid)
473
474 C-- Do "blocking" sends and receives for tendency "overlap" terms
475 c CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
476 c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )
477 c CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
478
479 C-- Do "blocking" sends and receives for field "overlap" terms
480 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
481 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
482 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
483
484 cswdptr -- add for seperate timestepping of chemical/biological/forcing
485 cswdptr of ptracers ---
486 #ifdef ALLOW_GCHEM
487 ceh3 This is broken -- this ifdef should not be visible!
488 #ifdef PTRACERS_SEPARATE_FORCING
489 ceh3 needs an IF ( use GCHEM ) THEN
490 call GCHEM_FORCING_SEP( myTime,myIter,myThid )
491 #endif /* PTRACERS_SEPARATE_FORCING */
492 #endif /* ALLOW_GCHEM */
493 cswdptr -- end add ---
494
495 C AMM
496 #ifdef ALLOW_GRIDALT
497 if (useGRIDALT) then
498 CALL GRIDALT_UPDATE(myThid)
499 endif
500 #endif
501 C AMM
502
503 C AMM
504 #ifdef ALLOW_FIZHI
505 if( useFIZHI) then
506 CALL TIMER_START('FIZHI [FORWARD_STEP]',mythid)
507 CALL STEP_FIZHI_CORR ( myTime, myIter, myThid )
508 CALL TIMER_STOP('FIZHI [FORWARD_STEP]',mythid)
509 endif
510 #endif
511 C AMM
512
513 #ifdef ALLOW_FLT
514 C-- Calculate float trajectories
515 IF (useFLT) THEN
516 CALL TIMER_START('FLOATS [FORWARD_STEP]',myThid)
517 CALL FLT_MAIN(myIter,myTime, myThid)
518 CALL TIMER_STOP ('FLOATS [FORWARD_STEP]',myThid)
519 ENDIF
520 #endif
521
522 #ifdef ALLOW_MONITOR
523 C-- Check status of solution (statistics, cfl, etc...)
524 CALL TIMER_START('MONITOR [FORWARD_STEP]',myThid)
525 CALL MONITOR( myIter, myTime, myThid )
526 CALL TIMER_STOP ('MONITOR [FORWARD_STEP]',myThid)
527 #endif /* ALLOW_MONITOR */
528
529 C AMM -- Diagnostics package
530 #ifdef ALLOW_DIAGNOSTICS
531
532 if(usediagnostics)then
533 call diagnostics_fill_state(myThid)
534 endif
535
536 # ifdef ALLOW_PTRACERS
537 if(useptracers)then
538 CALL PTRACERS_FILL_DIAGNOSTICS(myThid)
539 endif
540 # endif
541
542 #endif
543 C AMM
544
545 C-- Do IO if needed.
546 CALL TIMER_START('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
547 CALL DO_THE_MODEL_IO( myTime, myIter, myThid )
548 CALL TIMER_STOP ('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
549
550 C-- Save state for restarts
551 C Note: (jmc: is it still the case after ckp35 ?)
552 C =====
553 C Because of the ordering of the timestepping code and
554 C tendency term code at end of loop model arrays hold
555 C U,V,T,S at "time-level" N but gu, gv, gs, gt, guNM1,...
556 C at "time-level" N+1/2 (guNM1 at "time-level" N+1/2 is
557 C gu at "time-level" N-1/2) and etaN at "time-level" N+1/2.
558 C where N = I+timeLevBase-1
559 C Thus a checkpoint contains U.0000000000, GU.0000000001 and
560 C etaN.0000000001 in the indexing scheme used for the model
561 C "state" files. This example is referred to as a checkpoint
562 C at time level 1
563 CALL TIMER_START('WRITE_CHECKPOINT [FORWARD_STEP]',myThid)
564 CALL PACKAGES_WRITE_PICKUP(
565 I .FALSE., myTime, myIter, myThid )
566 CALL WRITE_CHECKPOINT(
567 I .FALSE., myTime, myIter, myThid )
568 CALL TIMER_STOP ('WRITE_CHECKPOINT [FORWARD_STEP]',myThid)
569
570 #ifdef ALLOW_DEBUG
571 IF ( debugLevel .GE. debLevB )
572 & CALL DEBUG_LEAVE('FORWARD_STEP',myThid)
573 #endif
574
575 RETURN
576 END

  ViewVC Help
Powered by ViewVC 1.1.22