1 |
C $Header: /u/gcmpack/MITgcm/model/src/forward_step.F,v 1.165 2009/02/13 21:56:48 heimbach Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "PACKAGES_CONFIG.h" |
5 |
#include "CPP_OPTIONS.h" |
6 |
|
7 |
#ifdef ALLOW_GMREDI |
8 |
# include "GMREDI_OPTIONS.h" |
9 |
#endif |
10 |
#ifdef ALLOW_OBCS |
11 |
# include "OBCS_OPTIONS.h" |
12 |
#endif |
13 |
#ifdef ALLOW_SEAICE |
14 |
# include "SEAICE_OPTIONS.h" |
15 |
#endif |
16 |
|
17 |
CBOP |
18 |
C !ROUTINE: FORWARD_STEP |
19 |
C !INTERFACE: |
20 |
SUBROUTINE FORWARD_STEP( iloop, myTime, myIter, myThid ) |
21 |
|
22 |
C !DESCRIPTION: \bv |
23 |
C *================================================================== |
24 |
C | SUBROUTINE forward_step |
25 |
C | o Run the ocean model and, optionally, evaluate a cost function. |
26 |
C *================================================================== |
27 |
C | |
28 |
C | THE_MAIN_LOOP is the toplevel routine for the Tangent Linear and |
29 |
C | Adjoint Model Compiler (TAMC). For this purpose the initialization |
30 |
C | of the model was split into two parts. Those parameters that do |
31 |
C | not depend on a specific model run are set in INITIALISE_FIXED, |
32 |
C | whereas those that do depend on the specific realization are |
33 |
C | initialized in INITIALISE_VARIA. |
34 |
C | |
35 |
C *================================================================== |
36 |
C \ev |
37 |
|
38 |
C !USES: |
39 |
IMPLICIT NONE |
40 |
C == Global variables == |
41 |
#include "SIZE.h" |
42 |
#include "EEPARAMS.h" |
43 |
#include "PARAMS.h" |
44 |
#include "DYNVARS.h" |
45 |
|
46 |
#ifdef ALLOW_MNC |
47 |
#include "MNC_PARAMS.h" |
48 |
#endif |
49 |
|
50 |
#ifdef HAVE_SIGREG |
51 |
#include "SIGREG.h" |
52 |
#endif |
53 |
|
54 |
#ifdef ALLOW_SHAP_FILT |
55 |
# include "SHAP_FILT.h" |
56 |
#endif |
57 |
#ifdef ALLOW_ZONAL_FILT |
58 |
# include "ZONAL_FILT.h" |
59 |
#endif |
60 |
#ifdef COMPONENT_MODULE |
61 |
# include "CPL_PARAMS.h" |
62 |
#endif |
63 |
|
64 |
#ifdef ALLOW_LONGSTEP |
65 |
# include "LONGSTEP_PARAMS.h" |
66 |
# include "LONGSTEP.h" |
67 |
#endif |
68 |
|
69 |
#ifdef ALLOW_AUTODIFF_TAMC |
70 |
# include "FFIELDS.h" |
71 |
# include "SURFACE.h" |
72 |
|
73 |
# include "tamc.h" |
74 |
# include "ctrl.h" |
75 |
# include "ctrl_dummy.h" |
76 |
# include "cost.h" |
77 |
# include "EOS.h" |
78 |
# if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL) |
79 |
# include "GRID.h" |
80 |
# endif |
81 |
# ifdef ALLOW_EXF |
82 |
# include "EXF_FIELDS.h" |
83 |
# ifdef ALLOW_BULKFORMULAE |
84 |
# include "EXF_CONSTANTS.h" |
85 |
# endif |
86 |
# endif |
87 |
# ifdef ALLOW_PTRACERS |
88 |
# include "PTRACERS_SIZE.h" |
89 |
# include "PTRACERS_FIELDS.h" |
90 |
# endif |
91 |
# ifdef ALLOW_GCHEM |
92 |
# include "GCHEM_FIELDS.h" |
93 |
# endif |
94 |
# ifdef ALLOW_CFC |
95 |
# include "CFC.h" |
96 |
# endif |
97 |
# ifdef ALLOW_DIC |
98 |
# include "DIC_VARS.h" |
99 |
# include "DIC_LOAD.h" |
100 |
# include "DIC_ATMOS.h" |
101 |
# endif |
102 |
# ifdef ALLOW_OBCS |
103 |
# include "OBCS.h" |
104 |
# ifdef ALLOW_PTRACERS |
105 |
# include "OBCS_PTRACERS.h" |
106 |
# endif |
107 |
# endif |
108 |
# ifdef ALLOW_CD_CODE |
109 |
# include "CD_CODE_VARS.h" |
110 |
# endif |
111 |
# ifdef ALLOW_THSICE |
112 |
# include "THSICE_VARS.h" |
113 |
# endif |
114 |
# ifdef ALLOW_SEAICE |
115 |
# include "SEAICE.h" |
116 |
# endif |
117 |
# ifdef ALLOW_SALT_PLUME |
118 |
# include "SALT_PLUME.h" |
119 |
# endif |
120 |
# ifdef ALLOW_EBM |
121 |
# include "EBM.h" |
122 |
# endif |
123 |
# ifdef ALLOW_KPP |
124 |
# include "KPP.h" |
125 |
# endif |
126 |
# ifdef ALLOW_GMREDI |
127 |
# include "GMREDI.h" |
128 |
# endif |
129 |
# ifdef ALLOW_RBCS |
130 |
# include "RBCS.h" |
131 |
# endif |
132 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
133 |
|
134 |
#ifdef ALLOW_MNC |
135 |
EXTERNAL DIFFERENT_MULTIPLE |
136 |
LOGICAL DIFFERENT_MULTIPLE |
137 |
#endif |
138 |
|
139 |
|
140 |
C !INPUT/OUTPUT PARAMETERS: |
141 |
C == Routine arguments == |
142 |
C note: under the multi-threaded model myIter and |
143 |
C myTime are local variables passed around as routine |
144 |
C arguments. Although this is fiddly it saves the need to |
145 |
C impose additional synchronisation points when they are |
146 |
C updated. |
147 |
C myTime :: time counter for this thread |
148 |
C myIter :: iteration counter for this thread |
149 |
C myThid :: thread number for this instance of the routine. |
150 |
INTEGER iloop |
151 |
_RL myTime |
152 |
INTEGER myIter |
153 |
INTEGER myThid |
154 |
|
155 |
C !LOCAL VARIABLES: |
156 |
C == Local variables == |
157 |
C modelEnd :: true if reaching the end of the run |
158 |
C myTimeBeg :: time at beginning of time step (needed by longstep) |
159 |
C myIterBeg :: iteration number at beginning of time step |
160 |
C myTimeTD :: time at call to thermodynamics (needed by longstep) |
161 |
C myIterTD :: iteration number at call to thermodynamics |
162 |
LOGICAL modelEnd |
163 |
#ifdef COMPONENT_MODULE |
164 |
INTEGER myItP1 |
165 |
#endif |
166 |
INTEGER myIterBeg, myIterTD |
167 |
_RL myTimeBeg, myTimeTD |
168 |
CEOP |
169 |
|
170 |
#ifdef ALLOW_DEBUG |
171 |
IF ( debugLevel .GE. debLevB ) |
172 |
& CALL DEBUG_ENTER('FORWARD_STEP',myThid) |
173 |
#endif |
174 |
|
175 |
#ifdef ALLOW_AUTODIFF_TAMC |
176 |
CALL AUTODIFF_INADMODE_UNSET( myThid ) |
177 |
#endif |
178 |
|
179 |
#ifdef ALLOW_AUTODIFF_TAMC |
180 |
C-- Reset the model iteration counter and the model time. |
181 |
myIter = nIter0 + (iloop-1) |
182 |
myTime = startTime + float(iloop-1)*deltaTclock |
183 |
#endif |
184 |
|
185 |
C store this for longstep_average with staggerTimeStep |
186 |
C and longstep_thermodynamics without staggerTimeStep |
187 |
C which are called after myIter and myTime are incremented |
188 |
C but need iter/time at beginning of time step |
189 |
myIterBeg = myIter |
190 |
myTimeBeg = myTime |
191 |
|
192 |
#ifdef ALLOW_AUTODIFF_TAMC |
193 |
c************************************** |
194 |
#include "checkpoint_lev1_directives.h" |
195 |
#include "checkpoint_lev1_template.h" |
196 |
c************************************** |
197 |
#endif |
198 |
|
199 |
C-- Switch on/off diagnostics for snap-shot output: |
200 |
#ifdef ALLOW_DIAGNOSTICS |
201 |
IF ( useDiagnostics ) THEN |
202 |
CALL DIAGNOSTICS_SWITCH_ONOFF( myTime, myIter, myThid ) |
203 |
C-- State-variables diagnostics |
204 |
CALL TIMER_START('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid) |
205 |
CALL DO_STATEVARS_DIAGS( myTime, 0, myIter, myThid ) |
206 |
CALL TIMER_STOP ('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid) |
207 |
ENDIF |
208 |
#endif |
209 |
|
210 |
#ifdef ALLOW_PROFILES |
211 |
#ifdef ALLOW_DEBUG |
212 |
IF (debugMode) CALL DEBUG_CALL('',myThid) |
213 |
#endif |
214 |
c-- Accumulate in-situ time averages of theta, salt, and SSH. |
215 |
CALL TIMER_START('PROFILES_INLOOP [THE_MAIN_LOOP]', mythid) |
216 |
CALL PROFILES_INLOOP( mytime, mythid ) |
217 |
CALL TIMER_STOP ('PROFILES_INLOOP [THE_MAIN_LOOP]', mythid) |
218 |
#endif |
219 |
|
220 |
C-- Call driver to load external forcing fields from file |
221 |
#ifdef ALLOW_DEBUG |
222 |
IF ( debugLevel .GE. debLevB ) |
223 |
& CALL DEBUG_CALL('LOAD_FIELDS_DRIVER',myThid) |
224 |
#endif |
225 |
#ifdef ALLOW_AUTODIFF_TAMC |
226 |
cph Important STORE that avoids hidden recomp. of load_fields_driver |
227 |
CADJ STORE theta = comlev1, key = ikey_dynamics, |
228 |
CADJ & kind = isbyte |
229 |
CADJ STORE uvel, vvel = comlev1, key = ikey_dynamics, |
230 |
CADJ & kind = isbyte |
231 |
#endif |
232 |
CALL TIMER_START('LOAD_FIELDS_DRIVER [FORWARD_STEP]',myThid) |
233 |
CALL LOAD_FIELDS_DRIVER( myTime, myIter, myThid ) |
234 |
CALL TIMER_STOP ('LOAD_FIELDS_DRIVER [FORWARD_STEP]',myThid) |
235 |
|
236 |
C-- Call Bulk-Formulae forcing package |
237 |
#ifdef ALLOW_BULK_FORCE |
238 |
IF ( useBulkForce ) THEN |
239 |
#ifdef ALLOW_DEBUG |
240 |
IF ( debugLevel .GE. debLevB ) |
241 |
& CALL DEBUG_CALL('BULKF_FORCING',myThid) |
242 |
#endif |
243 |
CALL TIMER_START('BULKF_FORCING [FORWARD_STEP]',myThid) |
244 |
C- calculate qnet and empmr (and wind stress) |
245 |
CALL BULKF_FORCING( myTime, myIter, myThid ) |
246 |
CALL TIMER_STOP ('BULKF_FORCING [FORWARD_STEP]',myThid) |
247 |
ENDIF |
248 |
#endif /* ALLOW_BULK_FORCE */ |
249 |
|
250 |
C-- Call external chepaml forcing package |
251 |
#ifdef ALLOW_CHEAPAML |
252 |
IF ( useCheapAML ) THEN |
253 |
#ifdef ALLOW_DEBUG |
254 |
IF ( debugLevel .GE. debLevB ) |
255 |
& CALL DEBUG_CALL('CHEAPAML_FIELDS_LOAD',myThid) |
256 |
#endif |
257 |
CALL TIMER_START('CHEAPAML [FORWARD_STEP]',mythid) |
258 |
C- calculate qnet (and wind stress) |
259 |
CALL CHEAPAML( myTime, myIter,myThid ) |
260 |
CALL TIMER_STOP ('CHEAPAML [FORWARD_STEP]',mythid) |
261 |
ENDIF |
262 |
#endif /*ALLOW_CHEAPAML */ |
263 |
|
264 |
|
265 |
#ifdef ALLOW_AUTODIFF |
266 |
c-- Add control vector for forcing and parameter fields |
267 |
IF ( myIter .EQ. nIter0 ) |
268 |
& CALL CTRL_MAP_FORCING (myThid) |
269 |
#endif |
270 |
|
271 |
#if (defined (ALLOW_AUTODIFF_TAMC) && defined (ALLOW_AUTODIFF_MONITOR)) |
272 |
CALL DUMMY_IN_STEPPING( myTime, myIter, myThid ) |
273 |
#endif |
274 |
|
275 |
#ifdef COMPONENT_MODULE |
276 |
IF ( useCoupler .AND. cpl_earlyExpImpCall ) THEN |
277 |
C Post coupling data that I export. |
278 |
C Read in coupling data that I import. |
279 |
CALL TIMER_START('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid) |
280 |
CALL CPL_EXPORT_MY_DATA( myIter, myTime, myThid ) |
281 |
CALL CPL_IMPORT_EXTERNAL_DATA( myIter, myTime, myThid ) |
282 |
CALL TIMER_STOP ('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid) |
283 |
ENDIF |
284 |
#endif /* COMPONENT_MODULE */ |
285 |
|
286 |
#ifdef ALLOW_EBM |
287 |
IF ( useEBM ) THEN |
288 |
# ifdef ALLOW_DEBUG |
289 |
IF ( debugLevel .GE. debLevB ) |
290 |
& CALL DEBUG_CALL('EBM',myThid) |
291 |
# endif |
292 |
CALL TIMER_START('EBM [FORWARD_STEP]',myThid) |
293 |
CALL EBM_DRIVER ( myTime, myIter, myThid ) |
294 |
CALL TIMER_STOP ('EBM [FORWARD_STEP]',myThid) |
295 |
ENDIF |
296 |
#endif /* ALLOW_EBM */ |
297 |
|
298 |
C-- Step forward fields and calculate time tendency terms. |
299 |
|
300 |
#ifdef ALLOW_DEBUG |
301 |
IF ( debugLevel .GE. debLevB ) |
302 |
& CALL DEBUG_CALL('DO_ATMOSPHERIC_PHYS',myThid) |
303 |
#endif |
304 |
CALL TIMER_START('DO_ATMOSPHERIC_PHYS [FORWARD_STEP]',myThid) |
305 |
CALL DO_ATMOSPHERIC_PHYS( myTime, myIter, myThid ) |
306 |
CALL TIMER_STOP ('DO_ATMOSPHERIC_PHYS [FORWARD_STEP]',myThid) |
307 |
|
308 |
#ifdef ALLOW_AUTODIFF_TAMC |
309 |
CADJ STORE surfaceforcingtice = comlev1, key = ikey_dynamics, |
310 |
CADJ & kind = isbyte |
311 |
# ifdef ALLOW_KPP |
312 |
CADJ STORE uvel = comlev1, key = ikey_dynamics, |
313 |
CADJ & kind = isbyte |
314 |
CADJ STORE vvel = comlev1, key = ikey_dynamics, |
315 |
CADJ & kind = isbyte |
316 |
# endif |
317 |
# ifdef EXACT_CONSERV |
318 |
cphCADJ STORE empmr = comlev1, key = ikey_dynamics |
319 |
cphCADJ STORE pmepr = comlev1, key = ikey_dynamics |
320 |
# endif |
321 |
# ifdef ALLOW_PTRACERS |
322 |
CADJ STORE ptracer = comlev1, key = ikey_dynamics, |
323 |
CADJ & kind = isbyte |
324 |
# endif |
325 |
# if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL) |
326 |
cph-test |
327 |
CADJ STORE hFacC = comlev1, key = ikey_dynamics, |
328 |
CADJ & kind = isbyte |
329 |
# ifndef DISABLE_RSTAR_CODE |
330 |
CADJ STORE rstarexpc = comlev1, key = ikey_dynamics, |
331 |
CADJ & kind = isbyte |
332 |
# endif |
333 |
# endif |
334 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
335 |
|
336 |
#ifndef ALLOW_AUTODIFF_TAMC |
337 |
IF ( .NOT. useOffLine ) THEN |
338 |
#endif |
339 |
#ifdef ALLOW_DEBUG |
340 |
IF ( debugLevel .GE. debLevB ) |
341 |
& CALL DEBUG_CALL('DO_OCEANIC_PHYS',myThid) |
342 |
#endif |
343 |
CALL TIMER_START('DO_OCEANIC_PHYS [FORWARD_STEP]',myThid) |
344 |
CALL DO_OCEANIC_PHYS( myTime, myIter, myThid ) |
345 |
CALL TIMER_STOP ('DO_OCEANIC_PHYS [FORWARD_STEP]',myThid) |
346 |
#ifdef ALLOW_AUTODIFF_TAMC |
347 |
CADJ STORE EmPmR = comlev1, key = ikey_dynamics, |
348 |
CADJ & kind = isbyte |
349 |
# ifdef EXACT_CONSERV |
350 |
CADJ STORE pmepr = comlev1, key = ikey_dynamics, |
351 |
CADJ & kind = isbyte |
352 |
# endif |
353 |
#else |
354 |
ENDIF |
355 |
#endif |
356 |
|
357 |
#ifdef ALLOW_AUTODIFF_TAMC |
358 |
# ifdef NONLIN_FRSURF |
359 |
cph-test |
360 |
CADJ STORE hFac_surfC = comlev1, key = ikey_dynamics, |
361 |
CADJ & kind = isbyte |
362 |
CADJ STORE hfac_surfs = comlev1, key = ikey_dynamics, |
363 |
CADJ & kind = isbyte |
364 |
CADJ STORE hfac_surfw = comlev1, key = ikey_dynamics, |
365 |
CADJ & kind = isbyte |
366 |
# endif |
367 |
# if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL) |
368 |
CADJ STORE hFacC, hFacS, hFacW |
369 |
CADJ & = comlev1, key = ikey_dynamics, |
370 |
CADJ & kind = isbyte |
371 |
CADJ STORE recip_hFacC, recip_hFacS, recip_hFacW |
372 |
CADJ & = comlev1, key = ikey_dynamics, |
373 |
CADJ & kind = isbyte |
374 |
c |
375 |
CADJ STORE surfaceforcingu = comlev1, key = ikey_dynamics, |
376 |
CADJ & kind = isbyte |
377 |
CADJ STORE surfaceforcingv = comlev1, key = ikey_dynamics, |
378 |
CADJ & kind = isbyte |
379 |
# endif |
380 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
381 |
|
382 |
#ifdef ALLOW_GCHEM |
383 |
#ifdef ALLOW_AUTODIFF_TAMC |
384 |
CADJ STORE ptracer = comlev1, key = ikey_dynamics, |
385 |
CADJ & kind = isbyte |
386 |
CADJ STORE theta = comlev1, key = ikey_dynamics, |
387 |
CADJ & kind = isbyte |
388 |
CADJ STORE salt = comlev1, key = ikey_dynamics, |
389 |
CADJ & kind = isbyte |
390 |
#endif |
391 |
IF ( useGCHEM ) THEN |
392 |
#ifdef ALLOW_DEBUG |
393 |
IF ( debugLevel .GE. debLevB ) |
394 |
& CALL DEBUG_CALL('GCHEM_CALC_TENDENCY',myThid) |
395 |
#endif |
396 |
CALL TIMER_START('GCHEM_CALC_TENDENCY [FORWARD_STEP]',myThid) |
397 |
CALL GCHEM_CALC_TENDENCY( myTime, myIter, myThid ) |
398 |
CALL TIMER_STOP ('GCHEM_CALC_TENDENCY [FORWARD_STEP]',myThid) |
399 |
ENDIF |
400 |
#endif /* ALLOW_GCHEM */ |
401 |
|
402 |
#ifdef ALLOW_AUTODIFF_TAMC |
403 |
cph needed to be moved here from do_oceanic_physics |
404 |
cph to be visible down the road |
405 |
c |
406 |
CADJ STORE rhoInSitu = comlev1, key = ikey_dynamics, |
407 |
CADJ & kind = isbyte |
408 |
CADJ STORE surfaceForcingS = comlev1, key = ikey_dynamics, |
409 |
CADJ & kind = isbyte |
410 |
CADJ STORE surfaceForcingT = comlev1, key = ikey_dynamics, |
411 |
CADJ & kind = isbyte |
412 |
CADJ STORE surfaceForcingTice = comlev1, key = ikey_dynamics, |
413 |
CADJ & kind = isbyte |
414 |
CADJ STORE IVDConvCount = comlev1, key = ikey_dynamics, |
415 |
CADJ & kind = isbyte |
416 |
# ifdef ALLOW_PTRACERS |
417 |
CADJ STORE surfaceForcingPTr = comlev1, key = ikey_dynamics, |
418 |
CADJ & kind = isbyte |
419 |
# endif |
420 |
c |
421 |
# ifdef ALLOW_GMREDI |
422 |
CADJ STORE Kwx = comlev1, key = ikey_dynamics, |
423 |
CADJ & kind = isbyte |
424 |
CADJ STORE Kwy = comlev1, key = ikey_dynamics, |
425 |
CADJ & kind = isbyte |
426 |
CADJ STORE Kwz = comlev1, key = ikey_dynamics, |
427 |
CADJ & kind = isbyte |
428 |
# ifdef GM_BOLUS_ADVEC |
429 |
CADJ STORE GM_PsiX = comlev1, key = ikey_dynamics, |
430 |
CADJ & kind = isbyte |
431 |
CADJ STORE GM_PsiY = comlev1, key = ikey_dynamics, |
432 |
CADJ & kind = isbyte |
433 |
# endif |
434 |
# endif |
435 |
c |
436 |
# ifdef ALLOW_KPP |
437 |
CADJ STORE KPPghat = comlev1, key = ikey_dynamics, |
438 |
CADJ & kind = isbyte |
439 |
CADJ STORE KPPfrac = comlev1, key = ikey_dynamics, |
440 |
CADJ & kind = isbyte |
441 |
CADJ STORE KPPdiffKzS = comlev1, key = ikey_dynamics, |
442 |
CADJ & kind = isbyte |
443 |
CADJ STORE KPPdiffKzT = comlev1, key = ikey_dynamics, |
444 |
CADJ & kind = isbyte |
445 |
# endif |
446 |
c |
447 |
# if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL) |
448 |
CADJ STORE theta,salt,wvel = comlev1, key = ikey_dynamics, |
449 |
CADJ & kind = isbyte |
450 |
CADJ STORE etaH = comlev1, key = ikey_dynamics, |
451 |
CADJ & kind = isbyte |
452 |
# ifdef ALLOW_CD_CODE |
453 |
CADJ STORE etanm1 = comlev1, key = ikey_dynamics, |
454 |
CADJ & kind = isbyte |
455 |
# endif |
456 |
# ifndef DISABLE_RSTAR_CODE |
457 |
cph-test |
458 |
CADJ STORE rstarexpc = comlev1, key = ikey_dynamics, |
459 |
CADJ & kind = isbyte |
460 |
# endif |
461 |
# endif |
462 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
463 |
|
464 |
#ifdef ALLOW_LONGSTEP |
465 |
IF ( usePTRACERS ) THEN |
466 |
IF ( .NOT.LS_afterTS .AND. .NOT.LS_staggerTimeStep ) THEN |
467 |
C Average all variables before advection (but after do_oceanic_phys |
468 |
C where Qsw, KPP and GMRedi stuff is computed). |
469 |
C This is like diagnostics package and will reproduce offline |
470 |
C results. |
471 |
#ifdef ALLOW_DEBUG |
472 |
IF ( debugLevel .GE. debLevB ) |
473 |
& CALL DEBUG_CALL('LONGSTEP_AVERAGE',myThid) |
474 |
#endif |
475 |
CALL TIMER_START('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid) |
476 |
CALL LONGSTEP_AVERAGE( myTime, myIter, myThid ) |
477 |
CALL TIMER_STOP ('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid) |
478 |
ENDIF |
479 |
ENDIF |
480 |
#endif /* ALLOW_LONGSTEP */ |
481 |
|
482 |
IF ( .NOT.staggerTimeStep ) THEN |
483 |
#ifdef ALLOW_DEBUG |
484 |
IF ( debugLevel .GE. debLevB ) |
485 |
& CALL DEBUG_CALL('THERMODYNAMICS',myThid) |
486 |
#endif |
487 |
CADJ STORE salt = comlev1, key = ikey_dynamics, |
488 |
CADJ & kind = isbyte |
489 |
CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',myThid) |
490 |
CALL THERMODYNAMICS( myTime, myIter, myThid ) |
491 |
CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',myThid) |
492 |
C-- if not staggerTimeStep: end |
493 |
ENDIF |
494 |
|
495 |
#ifdef ALLOW_LONGSTEP |
496 |
IF ( usePTRACERS ) THEN |
497 |
IF ( LS_afterTS .AND. .NOT.LS_staggerTimeStep ) THEN |
498 |
C Average T and S after thermodynamics. This will reproduce the |
499 |
C effect of gchem's separate timestep, and thus online results |
500 |
C for LS_nIter=1 |
501 |
#ifdef ALLOW_DEBUG |
502 |
IF ( debugLevel .GE. debLevB ) |
503 |
& CALL DEBUG_CALL('LONGSTEP_AVERAGE',myThid) |
504 |
#endif |
505 |
CALL TIMER_START('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid) |
506 |
CALL LONGSTEP_AVERAGE( myTime, myIter, myThid ) |
507 |
CALL TIMER_STOP ('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid) |
508 |
ENDIF |
509 |
ENDIF |
510 |
#endif /* ALLOW_LONGSTEP */ |
511 |
|
512 |
c #ifdef ALLOW_NONHYDROSTATIC |
513 |
IF ( implicitIntGravWave ) THEN |
514 |
CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
515 |
CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid ) |
516 |
CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
517 |
ENDIF |
518 |
c #endif |
519 |
|
520 |
#ifdef COMPONENT_MODULE |
521 |
IF ( useCoupler .AND. .NOT.cpl_earlyExpImpCall ) THEN |
522 |
C Post coupling data that I export. |
523 |
C Read in coupling data that I import. |
524 |
myItP1 = myIter + 1 |
525 |
CALL TIMER_START('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid) |
526 |
CALL CPL_EXPORT_MY_DATA( myItP1, myTime, myThid ) |
527 |
CALL CPL_IMPORT_EXTERNAL_DATA( myItP1, myTime, myThid ) |
528 |
CALL TIMER_STOP ('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid) |
529 |
# ifdef ALLOW_OCN_COMPON_INTERF |
530 |
IF ( useRealFreshWaterFlux ) THEN |
531 |
CALL OCN_APPLY_IMPORT( .FALSE., myTime, myIter, myThid ) |
532 |
ENDIF |
533 |
# endif /* ALLOW_OCN_COMPON_INTERF */ |
534 |
ENDIF |
535 |
#endif /* COMPONENT_MODULE */ |
536 |
|
537 |
#ifdef ALLOW_AUTODIFF_TAMC |
538 |
# if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL) |
539 |
CADJ STORE hFacC = comlev1, key = ikey_dynamics, |
540 |
CADJ & kind = isbyte |
541 |
CADJ STORE hFacS = comlev1, key = ikey_dynamics, |
542 |
CADJ & kind = isbyte |
543 |
CADJ STORE hFacW = comlev1, key = ikey_dynamics, |
544 |
CADJ & kind = isbyte |
545 |
CADJ STORE recip_hFacC = comlev1, key = ikey_dynamics, |
546 |
CADJ & kind = isbyte |
547 |
CADJ STORE recip_hFacS = comlev1, key = ikey_dynamics, |
548 |
CADJ & kind = isbyte |
549 |
CADJ STORE recip_hFacW = comlev1, key = ikey_dynamics, |
550 |
CADJ & kind = isbyte |
551 |
CADJ STORE etaN = comlev1, key = ikey_dynamics, |
552 |
CADJ & kind = isbyte |
553 |
c |
554 |
# ifndef DISABLE_RSTAR_CODE |
555 |
CADJ STORE rstarFacC = comlev1, key = ikey_dynamics, |
556 |
CADJ & kind = isbyte |
557 |
CADJ STORE rstarFacS = comlev1, key = ikey_dynamics, |
558 |
CADJ & kind = isbyte |
559 |
CADJ STORE rstarFacW = comlev1, key = ikey_dynamics, |
560 |
CADJ & kind = isbyte |
561 |
c |
562 |
CADJ STORE h0facc,h0facs,h0facw |
563 |
CADJ & = comlev1, key = ikey_dynamics, |
564 |
CADJ & kind = isbyte |
565 |
CADJ STORE rstardhcdt,rstardhsdt,rstardhwdt |
566 |
CADJ & = comlev1, key = ikey_dynamics, |
567 |
CADJ & kind = isbyte |
568 |
CADJ STORE rstarexpc,rstarexps,rstarexpw |
569 |
CADJ & = comlev1, key = ikey_dynamics, |
570 |
CADJ & kind = isbyte |
571 |
# endif |
572 |
# endif |
573 |
#endif |
574 |
C-- Step forward fields and calculate time tendency terms. |
575 |
#ifndef ALLOW_AUTODIFF_TAMC |
576 |
IF ( momStepping ) THEN |
577 |
#endif |
578 |
#ifdef ALLOW_DEBUG |
579 |
IF ( debugLevel .GE. debLevB ) |
580 |
& CALL DEBUG_CALL('DYNAMICS',myThid) |
581 |
#endif |
582 |
CALL TIMER_START('DYNAMICS [FORWARD_STEP]',myThid) |
583 |
CALL DYNAMICS( myTime, myIter, myThid ) |
584 |
CALL TIMER_STOP ('DYNAMICS [FORWARD_STEP]',myThid) |
585 |
#ifndef ALLOW_AUTODIFF_TAMC |
586 |
ENDIF |
587 |
#endif |
588 |
|
589 |
#ifdef ALLOW_AUTODIFF_TAMC |
590 |
# if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL) |
591 |
CADJ STORE gU, gV = comlev1, key = ikey_dynamics, |
592 |
CADJ & kind = isbyte |
593 |
# endif |
594 |
#endif |
595 |
|
596 |
C-- Update time-counter |
597 |
myIter = nIter0 + iLoop |
598 |
myTime = startTime + deltaTClock * float(iLoop) |
599 |
|
600 |
#ifdef ALLOW_MNC |
601 |
C Update the default next iter for MNC |
602 |
IF ( useMNC ) THEN |
603 |
CALL MNC_CW_CITER_SETG( 1, 1, -1, myIter , myThid ) |
604 |
|
605 |
C TODO: Logic should be added here so that users can specify, on |
606 |
C a per-citer-group basis, when it is time to update the |
607 |
C "current" (and not just the "next") iteration |
608 |
|
609 |
C TODO: the following is just a temporary band-aid (mostly, for |
610 |
C Baylor) until someone writes a routine that better handles time |
611 |
C boundaries such as weeks, months, years, etc. |
612 |
IF ( mnc_filefreq .GT. 0 ) THEN |
613 |
IF (DIFFERENT_MULTIPLE(mnc_filefreq,myTime,deltaTClock)) |
614 |
& THEN |
615 |
CALL MNC_CW_CITER_SETG( 1, 1, myIter, -1 , myThid ) |
616 |
ENDIF |
617 |
ENDIF |
618 |
ENDIF |
619 |
#endif /* ALLOW_MNC */ |
620 |
|
621 |
C-- Update geometric factors: |
622 |
#ifdef NONLIN_FRSURF |
623 |
C- update hfacC,W,S and recip_hFac according to etaH(n+1) : |
624 |
IF ( nonlinFreeSurf.GT.0) THEN |
625 |
IF ( select_rStar.GT.0 ) THEN |
626 |
# ifndef DISABLE_RSTAR_CODE |
627 |
# ifdef ALLOW_AUTODIFF_TAMC |
628 |
CADJ STORE hFacC = comlev1, key = ikey_dynamics, |
629 |
CADJ & kind = isbyte |
630 |
CADJ STORE hFacS = comlev1, key = ikey_dynamics, |
631 |
CADJ & kind = isbyte |
632 |
CADJ STORE hFacW = comlev1, key = ikey_dynamics, |
633 |
CADJ & kind = isbyte |
634 |
CADJ STORE recip_hFacC = comlev1, key = ikey_dynamics, |
635 |
CADJ & kind = isbyte |
636 |
CADJ STORE recip_hFacS = comlev1, key = ikey_dynamics, |
637 |
CADJ & kind = isbyte |
638 |
CADJ STORE recip_hFacW = comlev1, key = ikey_dynamics, |
639 |
CADJ & kind = isbyte |
640 |
c |
641 |
CADJ STORE rstarFacC = comlev1, key = ikey_dynamics, |
642 |
CADJ & kind = isbyte |
643 |
CADJ STORE rstarFacS = comlev1, key = ikey_dynamics, |
644 |
CADJ & kind = isbyte |
645 |
CADJ STORE rstarFacW = comlev1, key = ikey_dynamics, |
646 |
CADJ & kind = isbyte |
647 |
c |
648 |
CADJ STORE h0facc,h0facs,h0facw = comlev1, key = ikey_dynamics, |
649 |
CADJ & kind = isbyte |
650 |
# endif |
651 |
CALL TIMER_START('UPDATE_R_STAR [FORWARD_STEP]',myThid) |
652 |
CALL UPDATE_R_STAR( myTime, myIter, myThid ) |
653 |
CALL TIMER_STOP ('UPDATE_R_STAR [FORWARD_STEP]',myThid) |
654 |
# ifdef ALLOW_AUTODIFF_TAMC |
655 |
CADJ STORE hFacC = comlev1, key = ikey_dynamics, |
656 |
CADJ & kind = isbyte |
657 |
CADJ STORE hFacS = comlev1, key = ikey_dynamics, |
658 |
CADJ & kind = isbyte |
659 |
CADJ STORE hFacW = comlev1, key = ikey_dynamics, |
660 |
CADJ & kind = isbyte |
661 |
CADJ STORE recip_hFacC = comlev1, key = ikey_dynamics, |
662 |
CADJ & kind = isbyte |
663 |
CADJ STORE recip_hFacS = comlev1, key = ikey_dynamics, |
664 |
CADJ & kind = isbyte |
665 |
CADJ STORE recip_hFacW = comlev1, key = ikey_dynamics, |
666 |
CADJ & kind = isbyte |
667 |
# endif |
668 |
# endif /* DISABLE_RSTAR_CODE */ |
669 |
ELSE |
670 |
#ifdef ALLOW_AUTODIFF_TAMC |
671 |
CADJ STORE hFac_surfC, hFac_surfS, hFac_surfW |
672 |
CADJ & = comlev1, key = ikey_dynamics, |
673 |
CADJ & kind = isbyte |
674 |
#endif |
675 |
CALL TIMER_START('UPDATE_SURF_DR [FORWARD_STEP]',myThid) |
676 |
CALL UPDATE_SURF_DR( myTime, myIter, myThid ) |
677 |
CALL TIMER_STOP ('UPDATE_SURF_DR [FORWARD_STEP]',myThid) |
678 |
ENDIF |
679 |
ENDIF |
680 |
# ifdef ALLOW_AUTODIFF_TAMC |
681 |
CADJ STORE hFacC = comlev1, key = ikey_dynamics, |
682 |
CADJ & kind = isbyte |
683 |
CADJ STORE hFacS = comlev1, key = ikey_dynamics, |
684 |
CADJ & kind = isbyte |
685 |
CADJ STORE hFacW = comlev1, key = ikey_dynamics, |
686 |
CADJ & kind = isbyte |
687 |
CADJ STORE recip_hFacC = comlev1, key = ikey_dynamics, |
688 |
CADJ & kind = isbyte |
689 |
CADJ STORE recip_hFacS = comlev1, key = ikey_dynamics, |
690 |
CADJ & kind = isbyte |
691 |
CADJ STORE recip_hFacW = comlev1, key = ikey_dynamics, |
692 |
CADJ & kind = isbyte |
693 |
# endif |
694 |
C- update also CG2D matrix (and preconditioner) |
695 |
IF ( momStepping .AND. nonlinFreeSurf.GT.2 ) THEN |
696 |
CALL TIMER_START('UPDATE_CG2D [FORWARD_STEP]',myThid) |
697 |
CALL UPDATE_CG2D( myTime, myIter, myThid ) |
698 |
CALL TIMER_STOP ('UPDATE_CG2D [FORWARD_STEP]',myThid) |
699 |
ENDIF |
700 |
#endif /* NONLIN_FRSURF */ |
701 |
|
702 |
C-- Apply Filters to u*,v* before SOLVE_FOR_PRESSURE |
703 |
#ifdef ALLOW_SHAP_FILT |
704 |
IF (useSHAP_FILT .AND. shap_filt_uvStar) THEN |
705 |
CALL TIMER_START('SHAP_FILT_UV [FORWARD_STEP]',myThid) |
706 |
IF (implicDiv2Dflow.LT.1.) THEN |
707 |
C-- Explicit+Implicit part of the Barotropic Flow Divergence |
708 |
C => Filtering of uVel,vVel is necessary |
709 |
CALL SHAP_FILT_APPLY_UV( uVel,vVel, |
710 |
& myTime, myIter, myThid ) |
711 |
ENDIF |
712 |
CALL SHAP_FILT_APPLY_UV( gU,gV,myTime,myIter,myThid) |
713 |
CALL TIMER_STOP ('SHAP_FILT_UV [FORWARD_STEP]',myThid) |
714 |
ENDIF |
715 |
#endif |
716 |
#ifdef ALLOW_ZONAL_FILT |
717 |
IF (useZONAL_FILT .AND. zonal_filt_uvStar) THEN |
718 |
CALL TIMER_START('ZONAL_FILT_UV [FORWARD_STEP]',myThid) |
719 |
IF (implicDiv2Dflow.LT.1.) THEN |
720 |
C-- Explicit+Implicit part of the Barotropic Flow Divergence |
721 |
C => Filtering of uVel,vVel is necessary |
722 |
CALL ZONAL_FILT_APPLY_UV( uVel, vVel, myThid ) |
723 |
ENDIF |
724 |
CALL ZONAL_FILT_APPLY_UV( gU, gV, myThid ) |
725 |
CALL TIMER_STOP ('ZONAL_FILT_UV [FORWARD_STEP]',myThid) |
726 |
ENDIF |
727 |
#endif |
728 |
|
729 |
C-- Solve elliptic equation(s). |
730 |
C Two-dimensional only for conventional hydrostatic or |
731 |
C three-dimensional for non-hydrostatic and/or IGW scheme. |
732 |
IF ( momStepping ) THEN |
733 |
#ifdef ALLOW_AUTODIFF_TAMC |
734 |
# if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL) |
735 |
CADJ STORE uvel, vvel |
736 |
CADJ & = comlev1, key = ikey_dynamics, |
737 |
CADJ & kind = isbyte |
738 |
CADJ STORE empmr,hfacs,hfacw |
739 |
CADJ & = comlev1, key = ikey_dynamics, |
740 |
CADJ & kind = isbyte |
741 |
# endif |
742 |
#endif |
743 |
CALL TIMER_START('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid) |
744 |
CALL SOLVE_FOR_PRESSURE(myTime, myIter, myThid) |
745 |
CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid) |
746 |
ENDIF |
747 |
|
748 |
C-- Correct divergence in flow field and cycle time-stepping momentum |
749 |
#ifndef ALLOW_AUTODIFF_TAMC |
750 |
IF ( momStepping ) THEN |
751 |
#endif |
752 |
#ifdef ALLOW_AUTODIFF_TAMC |
753 |
# if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL) |
754 |
# ifndef DISABLE_RSTAR_CODE |
755 |
cph-test |
756 |
cph not clear, why this one |
757 |
CADJ STORE h0facc = comlev1, key = ikey_dynamics, |
758 |
CADJ & kind = isbyte |
759 |
# endif |
760 |
# endif |
761 |
# ifdef ALLOW_DEPTH_CONTROL |
762 |
CADJ STORE etan, uvel,vvel |
763 |
CADJ & = comlev1, key = ikey_dynamics |
764 |
# endif |
765 |
#endif |
766 |
CALL TIMER_START('MOM_CORRECTION_STEP [FORWARD_STEP]',myThid) |
767 |
CALL MOMENTUM_CORRECTION_STEP(myTime, myIter, myThid) |
768 |
CALL TIMER_STOP ('MOM_CORRECTION_STEP [FORWARD_STEP]',myThid) |
769 |
#ifndef ALLOW_AUTODIFF_TAMC |
770 |
ENDIF |
771 |
#endif |
772 |
|
773 |
#ifdef EXACT_CONSERV |
774 |
IF (exactConserv) THEN |
775 |
C-- Update etaH(n+1) : |
776 |
CALL TIMER_START('UPDATE_ETAH [FORWARD_STEP]',myThid) |
777 |
CALL UPDATE_ETAH( myTime, myIter, myThid ) |
778 |
CALL TIMER_STOP ('UPDATE_ETAH [FORWARD_STEP]',myThid) |
779 |
ENDIF |
780 |
#endif /* EXACT_CONSERV */ |
781 |
|
782 |
#ifdef NONLIN_FRSURF |
783 |
IF ( select_rStar.NE.0 ) THEN |
784 |
# ifndef DISABLE_RSTAR_CODE |
785 |
# ifdef ALLOW_AUTODIFF_TAMC |
786 |
cph-test |
787 |
CADJ STORE rstarfacc,rstarfacs,rstarfacw = comlev1, key = ikey_dynamics, |
788 |
CADJ & kind = isbyte |
789 |
# endif |
790 |
C-- r* : compute the future level thickness according to etaH(n+1) |
791 |
CALL TIMER_START('CALC_R_STAR [FORWARD_STEP]',myThid) |
792 |
CALL CALC_R_STAR(etaH, myTime, myIter, myThid ) |
793 |
CALL TIMER_STOP ('CALC_R_STAR [FORWARD_STEP]',myThid) |
794 |
# endif /* DISABLE_RSTAR_CODE */ |
795 |
ELSEIF ( nonlinFreeSurf.GT.0) THEN |
796 |
C-- compute the future surface level thickness according to etaH(n+1) |
797 |
# ifdef ALLOW_AUTODIFF_TAMC |
798 |
CADJ STORE etaH = comlev1, key = ikey_dynamics, |
799 |
CADJ & kind = isbyte |
800 |
# endif |
801 |
CALL TIMER_START('CALC_SURF_DR [FORWARD_STEP]',myThid) |
802 |
CALL CALC_SURF_DR(etaH, myTime, myIter, myThid ) |
803 |
CALL TIMER_STOP ('CALC_SURF_DR [FORWARD_STEP]',myThid) |
804 |
ENDIF |
805 |
# ifdef ALLOW_AUTODIFF_TAMC |
806 |
CADJ STORE hFac_surfC = comlev1, key = ikey_dynamics, |
807 |
CADJ & kind = isbyte |
808 |
CADJ STORE salt,theta,vvel = comlev1, key = ikey_dynamics, |
809 |
CADJ & kind = isbyte |
810 |
# endif |
811 |
#endif /* NONLIN_FRSURF */ |
812 |
|
813 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
814 |
IF ( staggerTimeStep ) THEN |
815 |
C-- do exchanges of U,V (needed for multiDim) when using stagger time-step : |
816 |
#ifdef ALLOW_DEBUG |
817 |
IF ( debugLevel .GE. debLevB ) |
818 |
& CALL DEBUG_CALL('DO_STAGGER_FIELDS_EXCH.',myThid) |
819 |
#endif |
820 |
CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
821 |
CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid ) |
822 |
CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
823 |
|
824 |
#ifdef ALLOW_DIAGNOSTICS |
825 |
C-- State-variables diagnostics |
826 |
IF ( useDiagnostics ) THEN |
827 |
CALL TIMER_START('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid) |
828 |
CALL DO_STATEVARS_DIAGS( myTime, 1, myIter, myThid ) |
829 |
CALL TIMER_STOP ('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid) |
830 |
ENDIF |
831 |
#endif |
832 |
|
833 |
#ifdef ALLOW_DEBUG |
834 |
IF ( debugLevel .GE. debLevB ) |
835 |
& CALL DEBUG_CALL('THERMODYNAMICS',myThid) |
836 |
#endif |
837 |
CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',myThid) |
838 |
CALL THERMODYNAMICS( myTime, myIter, myThid ) |
839 |
CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',myThid) |
840 |
|
841 |
C-- if staggerTimeStep: end |
842 |
ENDIF |
843 |
C---+--------+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
844 |
|
845 |
#ifdef ALLOW_AUTODIFF_TAMC |
846 |
cph This is needed because convective_adjustment calls |
847 |
cph find_rho which may use pressure() |
848 |
CADJ STORE totphihyd = comlev1, key = ikey_dynamics, |
849 |
CADJ & kind = isbyte |
850 |
#endif |
851 |
C-- Cycle time-stepping Tracers arrays (T,S,+pTracers) |
852 |
CALL TIMER_START('TRC_CORRECTION_STEP [FORWARD_STEP]',myThid) |
853 |
CALL TRACERS_CORRECTION_STEP(myTime, myIter, myThid) |
854 |
CALL TIMER_STOP ('TRC_CORRECTION_STEP [FORWARD_STEP]',myThid) |
855 |
|
856 |
#ifdef ALLOW_LONGSTEP |
857 |
IF ( usePTRACERS ) THEN |
858 |
IF ( LS_staggerTimeStep ) THEN |
859 |
C Average everything at the end of the timestep. This will |
860 |
C reproduce online results with staggerTimeStep=.TRUE. |
861 |
#ifdef ALLOW_DEBUG |
862 |
IF ( debugLevel .GE. debLevB ) |
863 |
& CALL DEBUG_CALL('LONGSTEP_AVERAGE',myThid) |
864 |
#endif |
865 |
CALL TIMER_START('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid) |
866 |
C myIter has been update after dynamics, but the averaging window |
867 |
C should be determined by myIter at beginning of timestep |
868 |
CALL LONGSTEP_AVERAGE( myTimeBeg, myIterBeg, myThid ) |
869 |
CALL TIMER_STOP ('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid) |
870 |
ENDIF |
871 |
|
872 |
#ifdef ALLOW_DEBUG |
873 |
IF ( debugLevel .GE. debLevB ) |
874 |
& CALL DEBUG_CALL('LONGSTEP_THERMODYNAMICS',myThid) |
875 |
#endif |
876 |
CALL TIMER_START('LONGSTEP_THERMODYNAMICS [FORWARD_STEP]', |
877 |
& myThid) |
878 |
IF ( LS_staggerTimeStep ) THEN |
879 |
myTimeTD = myTime |
880 |
myIterTD = myIter |
881 |
ELSE |
882 |
C set iter for forcing/obcs as at beginning of timestep |
883 |
myTimeTD = myTimeBeg |
884 |
myIterTD = myIterBeg |
885 |
ENDIF |
886 |
CALL LONGSTEP_THERMODYNAMICS( myTimeTD, myIterTD, myThid ) |
887 |
CALL TIMER_STOP ('LONGSTEP_THERMODYNAMICS [FORWARD_STEP]', |
888 |
& myThid) |
889 |
|
890 |
C-- Cycle time-stepping Tracers arrays (pTracers) |
891 |
CALL TIMER_START('LS_CORRECTION_STEP [FORWARD_STEP]',myThid) |
892 |
CALL LONGSTEP_CORRECTION_STEP(myTime, myIter, myThid) |
893 |
CALL TIMER_STOP ('LS_CORRECTION_STEP [FORWARD_STEP]',myThid) |
894 |
ENDIF |
895 |
#endif /* ALLOW_LONGSTEP */ |
896 |
|
897 |
#ifdef ALLOW_GCHEM |
898 |
C Add separate timestepping of chemical/biological/forcing |
899 |
C of ptracers here in GCHEM_FORCING_SEP |
900 |
#ifdef ALLOW_AUTODIFF_TAMC |
901 |
CADJ STORE ptracer = comlev1, key = ikey_dynamics, |
902 |
CADJ & kind = isbyte |
903 |
CADJ STORE theta = comlev1, key = ikey_dynamics, |
904 |
CADJ & kind = isbyte |
905 |
CADJ STORE salt = comlev1, key = ikey_dynamics, |
906 |
CADJ & kind = isbyte |
907 |
#endif |
908 |
|
909 |
#ifdef ALLOW_LONGSTEP |
910 |
IF ( LS_doTimeStep ) THEN |
911 |
#else |
912 |
IF ( .TRUE. ) THEN |
913 |
#endif |
914 |
IF ( useGCHEM ) THEN |
915 |
#ifdef ALLOW_DEBUG |
916 |
IF ( debugLevel .GE. debLevB ) |
917 |
& CALL DEBUG_CALL('GCHEM_FORCING_SEP',myThid) |
918 |
#endif /* ALLOW_DEBUG */ |
919 |
CALL TIMER_START('GCHEM_FORCING_SEP [FORWARD_STEP]',myThid) |
920 |
CALL GCHEM_FORCING_SEP( myTime,myIter,myThid ) |
921 |
CALL TIMER_STOP ('GCHEM_FORCING_SEP [FORWARD_STEP]',myThid) |
922 |
ENDIF |
923 |
C endif LS_doTimeStep |
924 |
ENDIF |
925 |
#endif /* ALLOW_GCHEM */ |
926 |
|
927 |
C-- Do "blocking" sends and receives for tendency "overlap" terms |
928 |
c CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
929 |
c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid ) |
930 |
c CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
931 |
|
932 |
C-- Do "blocking" sends and receives for field "overlap" terms |
933 |
CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
934 |
CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid ) |
935 |
CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
936 |
|
937 |
#ifdef ALLOW_DIAGNOSTICS |
938 |
IF ( useDiagnostics ) THEN |
939 |
CALL TIMER_START('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid) |
940 |
CALL DO_STATEVARS_DIAGS( myTime, 2, myIter, myThid ) |
941 |
CALL TIMER_STOP ('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid) |
942 |
ENDIF |
943 |
#endif |
944 |
|
945 |
#ifdef ALLOW_GRIDALT |
946 |
IF (useGRIDALT) THEN |
947 |
CALL GRIDALT_UPDATE(myThid) |
948 |
ENDIF |
949 |
#endif |
950 |
|
951 |
#ifdef ALLOW_FIZHI |
952 |
IF (useFIZHI) THEN |
953 |
CALL TIMER_START('FIZHI [FORWARD_STEP]',myThid) |
954 |
CALL STEP_FIZHI_CORR ( myTime, myIter, myThid, dTtracerLev(1) ) |
955 |
CALL TIMER_STOP ('FIZHI [FORWARD_STEP]',myThid) |
956 |
ENDIF |
957 |
#endif |
958 |
|
959 |
#ifdef ALLOW_FLT |
960 |
C-- Calculate float trajectories |
961 |
IF (useFLT) THEN |
962 |
CALL TIMER_START('FLOATS [FORWARD_STEP]',myThid) |
963 |
CALL FLT_MAIN( myTime, myIter, myThid ) |
964 |
CALL TIMER_STOP ('FLOATS [FORWARD_STEP]',myThid) |
965 |
ENDIF |
966 |
#endif |
967 |
|
968 |
#ifdef ALLOW_TIMEAVE |
969 |
C-- State-variables time-averaging |
970 |
CALL TIMER_START('DO_STATEVARS_TAVE [FORWARD_STEP]',myThid) |
971 |
CALL DO_STATEVARS_TAVE( myTime, myIter, myThid ) |
972 |
CALL TIMER_STOP ('DO_STATEVARS_TAVE [FORWARD_STEP]',myThid) |
973 |
#endif |
974 |
|
975 |
#ifdef ALLOW_MONITOR |
976 |
IF ( .NOT.useOffLine ) THEN |
977 |
C-- Check status of solution (statistics, cfl, etc...) |
978 |
CALL TIMER_START('MONITOR [FORWARD_STEP]',myThid) |
979 |
CALL MONITOR( myIter, myTime, myThid ) |
980 |
CALL TIMER_STOP ('MONITOR [FORWARD_STEP]',myThid) |
981 |
ENDIF |
982 |
#endif /* ALLOW_MONITOR */ |
983 |
|
984 |
#ifdef ALLOW_COST |
985 |
C-- compare model with data and compute cost function |
986 |
C-- this is done after exchanges to allow interpolation |
987 |
CALL TIMER_START('COST_TILE [FORWARD_STEP]',myThid) |
988 |
CALL COST_TILE ( myTime, myIter, myThid ) |
989 |
CALL TIMER_STOP ('COST_TILE [FORWARD_STEP]',myThid) |
990 |
#endif |
991 |
|
992 |
C-- Do IO if needed. |
993 |
CALL TIMER_START('DO_THE_MODEL_IO [FORWARD_STEP]',myThid) |
994 |
CALL DO_THE_MODEL_IO( myTime, myIter, myThid ) |
995 |
CALL TIMER_STOP ('DO_THE_MODEL_IO [FORWARD_STEP]',myThid) |
996 |
|
997 |
modelEnd = myTime.EQ.endTime .OR. myIter.EQ.nEndIter |
998 |
#ifdef HAVE_SIGREG |
999 |
IF ( useSIGREG ) THEN |
1000 |
modelEnd = modelEnd .OR. ( i_got_signal.GT.0 ) |
1001 |
ENDIF |
1002 |
#endif /* HAVE_SIGREG */ |
1003 |
|
1004 |
C-- Save state for restarts |
1005 |
CALL TIMER_START('DO_WRITE_PICKUP [FORWARD_STEP]',myThid) |
1006 |
CALL DO_WRITE_PICKUP( |
1007 |
I modelEnd, myTime, myIter, myThid ) |
1008 |
CALL TIMER_STOP ('DO_WRITE_PICKUP [FORWARD_STEP]',myThid) |
1009 |
|
1010 |
#ifdef HAVE_SIGREG |
1011 |
IF ( useSIGREG ) THEN |
1012 |
IF ( modelEnd .AND. i_got_signal.GT.0 ) THEN |
1013 |
STOP 'Checkpoint completed -- killed by signal handler' |
1014 |
ENDIF |
1015 |
ENDIF |
1016 |
#endif /* HAVE_SIGREG */ |
1017 |
|
1018 |
#ifdef ALLOW_AUTODIFF_TAMC |
1019 |
CALL AUTODIFF_INADMODE_SET( myThid ) |
1020 |
#endif |
1021 |
|
1022 |
#ifdef ALLOW_SHOWFLOPS |
1023 |
CALL TIMER_START('SHOWFLOPS_INLOOP [THE_MAIN_LOOP]', mythid) |
1024 |
CALL SHOWFLOPS_INLOOP( iloop, mythid ) |
1025 |
CALL TIMER_STOP ('SHOWFLOPS_INLOOP [THE_MAIN_LOOP]', mythid) |
1026 |
#endif |
1027 |
|
1028 |
#ifdef ALLOW_DEBUG |
1029 |
IF ( debugLevel .GE. debLevB ) |
1030 |
& CALL DEBUG_LEAVE('FORWARD_STEP',myThid) |
1031 |
#endif |
1032 |
|
1033 |
RETURN |
1034 |
END |