/[MITgcm]/MITgcm_contrib/plumes/forward_step.F
ViewVC logotype

Contents of /MITgcm_contrib/plumes/forward_step.F

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


Revision 1.1 - (show annotations) (download)
Thu May 13 22:21:45 2004 UTC (19 years, 11 months ago) by molod
Branch: MAIN
CVS Tags: HEAD
More developing....

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

  ViewVC Help
Powered by ViewVC 1.1.22