1 |
C $Header: /u/gcmpack/MITgcm/model/src/forward_step.F,v 1.88 2004/03/29 22:07:40 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 |
# include "exf_clim_fields.h" |
63 |
# ifdef ALLOW_BULKFORMULAE |
64 |
# include "exf_constants.h" |
65 |
# endif |
66 |
# endif |
67 |
# ifdef ALLOW_OBCS |
68 |
# include "OBCS.h" |
69 |
# endif |
70 |
# ifdef ALLOW_PTRACERS |
71 |
# include "PTRACERS.h" |
72 |
# endif |
73 |
# ifdef ALLOW_CD_CODE |
74 |
# include "CD_CODE_VARS.h" |
75 |
# endif |
76 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
77 |
|
78 |
C !LOCAL VARIABLES: |
79 |
C == Routine arguments == |
80 |
C note: under the multi-threaded model myiter and |
81 |
C mytime are local variables passed around as routine |
82 |
C arguments. Although this is fiddly it saves the need to |
83 |
C impose additional synchronisation points when they are |
84 |
C updated. |
85 |
C myiter - iteration counter for this thread |
86 |
C mytime - time counter for this thread |
87 |
C mythid - thread number for this instance of the routine. |
88 |
integer iloop |
89 |
integer mythid |
90 |
integer myiter |
91 |
_RL mytime |
92 |
|
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 |
|
284 |
C AMM |
285 |
#ifdef ALLOW_FIZHI |
286 |
if( useFIZHI) then |
287 |
CALL TIMER_START('UPDATE EXPORTS [FORWARD_STEP]',mythid) |
288 |
CALL UPDATE_OCEAN_EXPORTS ( myTime, myIter, myThid ) |
289 |
CALL UPDATE_EARTH_EXPORTS ( myTime, myIter, myThid ) |
290 |
CALL UPDATE_CHEMISTRY_EXPORTS ( myTime, myIter, myThid ) |
291 |
CALL TIMER_STOP ('UPDATE EXPORTS [FORWARD_STEP]',mythid) |
292 |
CALL TIMER_START('FIZHI_WRAPPER [FORWARD_STEP]',mythid) |
293 |
CALL FIZHI_WRAPPER ( myTime, myIter, myThid ) |
294 |
CALL TIMER_STOP ('FIZHI_WRAPPER [FORWARD_STEP]',mythid) |
295 |
CALL TIMER_START('STEP_FIZHI_FG [FORWARD_STEP]',mythid) |
296 |
CALL STEP_FIZHI_FG ( myTime, myIter, myThid, deltaTtracer ) |
297 |
CALL TIMER_STOP ('STEP_FIZHI_FG [FORWARD_STEP]',mythid) |
298 |
endif |
299 |
#endif |
300 |
C AMM |
301 |
|
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 |
|
310 |
C-- do exchanges (needed for DYNAMICS) when using stagger time-step : |
311 |
#ifdef ALLOW_DEBUG |
312 |
IF ( debugLevel .GE. debLevB ) |
313 |
& CALL DEBUG_CALL('DO_STAGGER_FIELDS_EXCH.',myThid) |
314 |
#endif |
315 |
CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
316 |
CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid ) |
317 |
CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
318 |
|
319 |
#ifdef ALLOW_SHAP_FILT |
320 |
IF (useSHAP_FILT .AND. |
321 |
& staggerTimeStep .AND. shap_filt_TrStagg ) THEN |
322 |
#ifdef ALLOW_DEBUG |
323 |
IF ( debugLevel .GE. debLevB ) |
324 |
& CALL DEBUG_CALL('SHAP_FILT_APPLY_TS',myThid) |
325 |
#endif |
326 |
CALL TIMER_START('SHAP_FILT [FORWARD_STEP]',myThid) |
327 |
CALL SHAP_FILT_APPLY_TS(gT,gS,myTime+deltaT,myIter+1,myThid) |
328 |
CALL TIMER_STOP ('SHAP_FILT [FORWARD_STEP]',myThid) |
329 |
ENDIF |
330 |
#endif |
331 |
|
332 |
#ifdef ALLOW_ZONAL_FILT |
333 |
IF (useZONAL_FILT .AND. |
334 |
& staggerTimeStep .AND. zonal_filt_TrStagg ) THEN |
335 |
#ifdef ALLOW_DEBUG |
336 |
IF ( debugLevel .GE. debLevB ) |
337 |
& CALL DEBUG_CALL('ZONAL_FILT_APPLY_TS',myThid) |
338 |
#endif |
339 |
CALL TIMER_START('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid) |
340 |
CALL ZONAL_FILT_APPLY_TS( gT, gS, myThid ) |
341 |
CALL TIMER_STOP ('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid) |
342 |
ENDIF |
343 |
#endif |
344 |
|
345 |
C-- Step forward fields and calculate time tendency terms. |
346 |
#ifndef ALLOW_AUTODIFF_TAMC |
347 |
IF ( momStepping ) THEN |
348 |
#endif |
349 |
#ifdef ALLOW_DEBUG |
350 |
IF ( debugLevel .GE. debLevB ) |
351 |
& CALL DEBUG_CALL('DYNAMICS',myThid) |
352 |
#endif |
353 |
CALL TIMER_START('DYNAMICS [FORWARD_STEP]',mythid) |
354 |
CALL DYNAMICS( myTime, myIter, myThid ) |
355 |
CALL TIMER_STOP ('DYNAMICS [FORWARD_STEP]',mythid) |
356 |
#ifndef ALLOW_AUTODIFF_TAMC |
357 |
ENDIF |
358 |
#endif |
359 |
|
360 |
#ifdef ALLOW_NONHYDROSTATIC |
361 |
C-- Step forward W field in N-H algorithm |
362 |
IF ( momStepping .AND. nonHydrostatic ) THEN |
363 |
#ifdef ALLOW_DEBUG |
364 |
IF ( debugLevel .GE. debLevB ) |
365 |
& CALL DEBUG_CALL('CALC_GW',myThid) |
366 |
#endif |
367 |
CALL TIMER_START('CALC_GW [FORWARD_STEP]',myThid) |
368 |
CALL CALC_GW(myThid) |
369 |
CALL TIMER_STOP ('CALC_GW [FORWARD_STEP]',myThid) |
370 |
ENDIF |
371 |
#endif |
372 |
|
373 |
#ifdef NONLIN_FRSURF |
374 |
C-- update hfacC,W,S and recip_hFac according to etaH(n+1) : |
375 |
IF ( nonlinFreeSurf.GT.0) THEN |
376 |
IF ( select_rStar.GT.0 ) THEN |
377 |
CALL TIMER_START('UPDATE_R_STAR [FORWARD_STEP]',myThid) |
378 |
CALL UPDATE_R_STAR( myTime, myIter, myThid ) |
379 |
CALL TIMER_STOP ('UPDATE_R_STAR [FORWARD_STEP]',myThid) |
380 |
ELSE |
381 |
CALL TIMER_START('UPDATE_SURF_DR [FORWARD_STEP]',myThid) |
382 |
CALL UPDATE_SURF_DR( myTime, myIter, myThid ) |
383 |
CALL TIMER_STOP ('UPDATE_SURF_DR [FORWARD_STEP]',myThid) |
384 |
ENDIF |
385 |
ENDIF |
386 |
C- update also CG2D matrix (and preconditioner) |
387 |
IF ( momStepping .AND. nonlinFreeSurf.GT.2 ) THEN |
388 |
CALL TIMER_START('UPDATE_CG2D [FORWARD_STEP]',myThid) |
389 |
CALL UPDATE_CG2D( myTime, myIter, myThid ) |
390 |
CALL TIMER_STOP ('UPDATE_CG2D [FORWARD_STEP]',myThid) |
391 |
ENDIF |
392 |
#endif |
393 |
|
394 |
C-- Apply Filters to u*,v* before SOLVE_FOR_PRESSURE |
395 |
#ifdef ALLOW_SHAP_FILT |
396 |
IF (useSHAP_FILT .AND. shap_filt_uvStar) THEN |
397 |
CALL TIMER_START('SHAP_FILT [FORWARD_STEP]',myThid) |
398 |
IF (implicDiv2Dflow.LT.1.) THEN |
399 |
C-- Explicit+Implicit part of the Barotropic Flow Divergence |
400 |
C => Filtering of uVel,vVel is necessary |
401 |
CALL SHAP_FILT_APPLY_UV( uVel,vVel, |
402 |
& myTime+deltaT, myIter+1, myThid ) |
403 |
ENDIF |
404 |
CALL SHAP_FILT_APPLY_UV( gU,gV,myTime+deltaT,myIter+1,myThid) |
405 |
CALL TIMER_STOP ('SHAP_FILT [FORWARD_STEP]',myThid) |
406 |
ENDIF |
407 |
#endif |
408 |
#ifdef ALLOW_ZONAL_FILT |
409 |
IF (useZONAL_FILT .AND. zonal_filt_uvStar) THEN |
410 |
CALL TIMER_START('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid) |
411 |
IF (implicDiv2Dflow.LT.1.) THEN |
412 |
C-- Explicit+Implicit part of the Barotropic Flow Divergence |
413 |
C => Filtering of uVel,vVel is necessary |
414 |
CALL ZONAL_FILT_APPLY_UV( uVel, vVel, myThid ) |
415 |
ENDIF |
416 |
CALL ZONAL_FILT_APPLY_UV( gU, gV, myThid ) |
417 |
CALL TIMER_STOP ('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid) |
418 |
ENDIF |
419 |
#endif |
420 |
|
421 |
C-- Solve elliptic equation(s). |
422 |
C Two-dimensional only for conventional hydrostatic or |
423 |
C three-dimensional for non-hydrostatic and/or IGW scheme. |
424 |
IF ( momStepping ) THEN |
425 |
CALL TIMER_START('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid) |
426 |
CALL SOLVE_FOR_PRESSURE(myTime, myIter, myThid) |
427 |
CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid) |
428 |
ENDIF |
429 |
|
430 |
#ifdef ALLOW_AUTODIFF_TAMC |
431 |
cph This is needed because convective_adjustment calls |
432 |
cph find_rho which may use pressure() |
433 |
CADJ STORE totphihyd = comlev1, key = ikey_dynamics |
434 |
#endif |
435 |
C-- Correct divergence in flow field and cycle time-stepping |
436 |
C arrays (for all fields) ; update time-counter |
437 |
myIter = nIter0 + iLoop |
438 |
myTime = startTime + deltaTClock * float(iLoop) |
439 |
CALL TIMER_START('THE_CORRECTION_STEP [FORWARD_STEP]',myThid) |
440 |
CALL THE_CORRECTION_STEP(myTime, myIter, myThid) |
441 |
CALL TIMER_STOP ('THE_CORRECTION_STEP [FORWARD_STEP]',myThid) |
442 |
|
443 |
C-- Do "blocking" sends and receives for tendency "overlap" terms |
444 |
c CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
445 |
c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid ) |
446 |
c CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
447 |
|
448 |
C-- Do "blocking" sends and receives for field "overlap" terms |
449 |
CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
450 |
CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid ) |
451 |
CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
452 |
|
453 |
cswdptr -- add for seperate timestepping of chemical/biological/forcing |
454 |
cswdptr of ptracers --- |
455 |
#ifdef ALLOW_GCHEM |
456 |
ceh3 This is broken -- this ifdef should not be visible! |
457 |
#ifdef PTRACERS_SEPARATE_FORCING |
458 |
ceh3 needs an IF ( use GCHEM ) THEN |
459 |
call GCHEM_FORCING_SEP( myTime,myIter,myThid ) |
460 |
#endif /* PTRACERS_SEPARATE_FORCING */ |
461 |
#endif /* ALLOW_GCHEM */ |
462 |
cswdptr -- end add --- |
463 |
|
464 |
C AMM |
465 |
#ifdef ALLOW_GRIDALT |
466 |
if (useGRIDALT) then |
467 |
CALL TIMER_START('GRIDALT_UPDATE [FORWARD_STEP]',mythid) |
468 |
CALL GRIDALT_UPDATE(myThid) |
469 |
CALL TIMER_STOP ('GRIDALT_UPDATE [FORWARD_STEP]',mythid) |
470 |
endif |
471 |
#endif |
472 |
C AMM |
473 |
|
474 |
C AMM |
475 |
#ifdef ALLOW_FIZHI |
476 |
if( useFIZHI) then |
477 |
CALL TIMER_START('STEP_FIZHI_CORR [FORWARD_STEP]',mythid) |
478 |
CALL STEP_FIZHI_CORR ( myTime, myIter, myThid ) |
479 |
CALL TIMER_STOP ('STEP_FIZHI_CORR [FORWARD_STEP]',mythid) |
480 |
endif |
481 |
#endif |
482 |
C AMM |
483 |
|
484 |
#ifdef ALLOW_FLT |
485 |
C-- Calculate float trajectories |
486 |
IF (useFLT) THEN |
487 |
CALL TIMER_START('FLOATS [FORWARD_STEP]',myThid) |
488 |
CALL FLT_MAIN(myIter,myTime, myThid) |
489 |
CALL TIMER_STOP ('FLOATS [FORWARD_STEP]',myThid) |
490 |
ENDIF |
491 |
#endif |
492 |
|
493 |
#ifdef ALLOW_MONITOR |
494 |
C-- Check status of solution (statistics, cfl, etc...) |
495 |
CALL TIMER_START('MONITOR [FORWARD_STEP]',myThid) |
496 |
CALL MONITOR( myIter, myTime, myThid ) |
497 |
CALL TIMER_STOP ('MONITOR [FORWARD_STEP]',myThid) |
498 |
#endif /* ALLOW_MONITOR */ |
499 |
|
500 |
C AMM -- Diagnostics package |
501 |
#ifdef ALLOW_DIAGNOSTICS |
502 |
|
503 |
if(usediagnostics)then |
504 |
call diagnostics_fill_state(myThid) |
505 |
endif |
506 |
|
507 |
# ifdef ALLOW_PTRACERS |
508 |
if(useptracers)then |
509 |
call diagnostics_fill_ptracers(myThid) |
510 |
endif |
511 |
# endif |
512 |
|
513 |
#endif |
514 |
C AMM |
515 |
|
516 |
C-- Do IO if needed. |
517 |
CALL TIMER_START('DO_THE_MODEL_IO [FORWARD_STEP]',myThid) |
518 |
CALL DO_THE_MODEL_IO( myTime, myIter, myThid ) |
519 |
CALL TIMER_STOP ('DO_THE_MODEL_IO [FORWARD_STEP]',myThid) |
520 |
|
521 |
C-- Save state for restarts |
522 |
C Note: (jmc: is it still the case after ckp35 ?) |
523 |
C ===== |
524 |
C Because of the ordering of the timestepping code and |
525 |
C tendency term code at end of loop model arrays hold |
526 |
C U,V,T,S at "time-level" N but gu, gv, gs, gt, guNM1,... |
527 |
C at "time-level" N+1/2 (guNM1 at "time-level" N+1/2 is |
528 |
C gu at "time-level" N-1/2) and etaN at "time-level" N+1/2. |
529 |
C where N = I+timeLevBase-1 |
530 |
C Thus a checkpoint contains U.0000000000, GU.0000000001 and |
531 |
C etaN.0000000001 in the indexing scheme used for the model |
532 |
C "state" files. This example is referred to as a checkpoint |
533 |
C at time level 1 |
534 |
CALL TIMER_START('WRITE_CHECKPOINT [FORWARD_STEP]',myThid) |
535 |
CALL PACKAGES_WRITE_PICKUP( |
536 |
I .FALSE., myTime, myIter, myThid ) |
537 |
CALL WRITE_CHECKPOINT( |
538 |
I .FALSE., myTime, myIter, myThid ) |
539 |
CALL TIMER_STOP ('WRITE_CHECKPOINT [FORWARD_STEP]',myThid) |
540 |
|
541 |
#ifdef ALLOW_DEBUG |
542 |
IF ( debugLevel .GE. debLevB ) |
543 |
& CALL DEBUG_LEAVE('FORWARD_STEP',myThid) |
544 |
#endif |
545 |
|
546 |
RETURN |
547 |
END |