/[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.80 - (show annotations) (download)
Thu Jan 29 14:45:47 2004 UTC (20 years, 4 months ago) by molod
Branch: MAIN
Changes since 1.79: +19 -9 lines
Fizhi changes in forward step and diagnostics package_init

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

  ViewVC Help
Powered by ViewVC 1.1.22