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