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