/[MITgcm]/MITgcm/model/src/forward_step.F
ViewVC logotype

Contents of /MITgcm/model/src/forward_step.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.243 - (show annotations) (download)
Fri Dec 29 19:30:22 2017 UTC (6 years, 4 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint66o, checkpoint66n, HEAD
Changes since 1.242: +8 -1 lines
- include GCHEM_SIZE.h before GCHEM_FIELDS.h
- add missing (since Apr 2008) option files: GCHEM_OPTIONS.h & DIC_OPTIONS.h

1 C $Header: /u/gcmpack/MITgcm/model/src/forward_step.F,v 1.242 2017/03/29 15:44:00 mmazloff Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 #ifdef ALLOW_AUTODIFF
8 # include "AUTODIFF_OPTIONS.h"
9 #endif
10 #ifdef ALLOW_GENERIC_ADVDIFF
11 # include "GAD_OPTIONS.h"
12 #endif
13 #ifdef ALLOW_GGL90
14 # include "GGL90_OPTIONS.h"
15 #endif
16 #ifdef ALLOW_GMREDI
17 # include "GMREDI_OPTIONS.h"
18 #endif
19 #ifdef ALLOW_OBCS
20 # include "OBCS_OPTIONS.h"
21 #endif
22 #ifdef ALLOW_THSICE
23 # include "THSICE_OPTIONS.h"
24 #endif
25 #ifdef ALLOW_SEAICE
26 # include "SEAICE_OPTIONS.h"
27 #endif
28 #ifdef ALLOW_PTRACERS
29 # include "PTRACERS_OPTIONS.h"
30 #endif
31 #ifdef ALLOW_GCHEM
32 # include "GCHEM_OPTIONS.h"
33 #endif
34 #ifdef ALLOW_DIC
35 # include "DIC_OPTIONS.h"
36 #endif
37 #ifdef ALLOW_EXF
38 # include "EXF_OPTIONS.h"
39 #endif
40 #ifdef ALLOW_STREAMICE
41 # include "STREAMICE_OPTIONS.h"
42 #endif
43 #ifdef ALLOW_COST
44 # include "COST_OPTIONS.h"
45 #endif
46 #ifdef ALLOW_CTRL
47 # include "CTRL_OPTIONS.h"
48 #endif
49 #ifdef ALLOW_ECCO
50 # include "ECCO_OPTIONS.h"
51 #endif
52
53 #define ALLOW_MOM_STEPPING
54 #if ( defined (ALLOW_AUTODIFF) && defined (ALLOW_OFFLINE) )
55 # undef ALLOW_MOM_STEPPING
56 #endif
57
58 CBOP
59 C !ROUTINE: FORWARD_STEP
60 C !INTERFACE:
61 SUBROUTINE FORWARD_STEP( iloop, myTime, myIter, myThid )
62
63 C !DESCRIPTION: \bv
64 C *=================================================================
65 C | SUBROUTINE forward_step
66 C | o Step forward in time the model variables for one time-step
67 C *=================================================================
68 C | The algorithm...
69 C |
70 C | "Calculation of Gs"
71 C | ===================
72 C | This is where all the accelerations and tendencies (ie.
73 C | physics, parameterizations etc...) are calculated
74 C | rho = rho ( theta[n], salt[n] )
75 C | b = b(rho, theta)
76 C | K31 = K31 ( rho )
77 C | Gu[n] = Gu( u[n], v[n], wVel, b, ... )
78 C | Gv[n] = Gv( u[n], v[n], wVel, b, ... )
79 C | Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
80 C | Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
81 C |
82 C | "Time-stepping" or "Prediction"
83 C | ================================
84 C | The models variables are stepped forward with the appropriate
85 C | time-stepping scheme (currently we use Adams-Bashforth II)
86 C | - For momentum, the result is always *only* a "prediction"
87 C | in that the flow may be divergent and will be "corrected"
88 C | later with a surface pressure gradient.
89 C | - Normally for tracers the result is the new field at time
90 C | level [n+1} *BUT* in the case of implicit diffusion the result
91 C | is also *only* a prediction.
92 C | - We denote "predictors" with an asterisk (*).
93 C | U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
94 C | V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
95 C | theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
96 C | salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
97 C | With implicit diffusion:
98 C | theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
99 C | salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
100 C | (1 + dt * K * d_zz) theta[n+1] = theta*
101 C | (1 + dt * K * d_zz) salt[n+1] = salt*
102 C |
103 C | "Correction Step"
104 C | =================
105 C | Here we update the horizontal velocities with the surface
106 C | pressure such that the resulting flow is either consistent
107 C | with the free-surface evolution or the rigid-lid:
108 C | U[n] = U* + dt x d/dx P
109 C | V[n] = V* + dt x d/dy P
110 C | W[n] = W* + dt x d/dz P (NH mode)
111 C *=================================================================
112 C \ev
113
114 C !CALLING SEQUENCE:
115 C FORWARD_STEP
116 C |
117 C |-- AUTODIFF_INADMODE_UNSET
118 C |
119 C |-- RESET_NLFS_VARS
120 C |-- UPDATE_R_STAR
121 C |-- UPDATE_SURF_DR
122 C |
123 C |-- PTRACERS_SWITCH_ONOFF
124 C |
125 C |-- DIAGNOSTICS_SWITCH_ONOFF
126 C |-- DO_STATEVARS_DIAGS
127 C |
128 C |-- NEST_CHILD_SETMEMO
129 C |-- NEST_PARENT_IO_1
130 C |
131 C |-- LOAD_FIELDS_DRIVER
132 C |
133 C |-- BULKF_FORCING
134 C |
135 C |-- CHEAPAML
136 C |
137 C |-- CTRL_MAP_FORCING
138 C |-- DUMMY_IN_STEPPING
139 C |
140 C |-- CPL_EXPORT_IMPORT_DATA
141 C |
142 C |-- OASIS_PUT
143 C |-- OASIS_GET
144 C |
145 C |-- EBM_DRIVER
146 C |
147 C |-- DO_ATMOSPHERIC_PHYS
148 C |
149 C |-- DO_OCEANIC_PHYS
150 C |
151 C |-- STREAMICE_TIMESTEP
152 C |
153 C |-- GCHEM_CALC_TENDENCY
154 C |
155 C |-- LONGSTEP_AVERAGE
156 C |-- LONGSTEP_THERMODYNAMICS
157 C |
158 C |-- THERMODYNAMICS
159 C |
160 C |-- LONGSTEP_AVERAGE
161 C |-- LONGSTEP_THERMODYNAMICS
162 C |
163 C |-- DO_STAGGER_FIELDS_EXCHANGES
164 C |
165 C |-- DYNAMICS
166 C |
167 C |-- MNC_UPDATE_TIME
168 C |
169 C |-- OFFLINE_FIELDS_LOAD
170 C |
171 C |-- UPDATE_R_STAR
172 C |-- UPDATE_SIGMA
173 C |-- UPDATE_SURF_DR
174 C |-- UPDATE_CG2D
175 C |
176 C |-- SHAP_FILT_APPLY_UV
177 C |-- ZONAL_FILT_APPLY_UV
178 C |
179 C |-- SOLVE_FOR_PRESSURE
180 C |
181 C |-- MOMENTUM_CORRECTION_STEP
182 C |
183 C |-- INTEGR_CONTINUITY
184 C |
185 C |-- CALC_R_STAR
186 C |-- CALC_SURF_DR
187 C |
188 C |-- DO_STAGGER_FIELDS_EXCHANGES
189 C |
190 C |-- DO_STATEVARS_DIAGS
191 C |
192 C |-- THERMODYNAMICS
193 C |
194 C |-- TRACERS_CORRECTION_STEP
195 C |
196 C |-- LONGSTEP_AVERAGE
197 C |-- LONGSTEP_THERMODYNAMICS
198 C |
199 C |-- GCHEM_FORCING_SEP
200 C |
201 C |-- DO_FIELDS_BLOCKING_EXCHANGES
202 C |
203 C |-- DO_STATEVARS_DIAGS
204 C |
205 C |-- GRIDALT_UPDATE
206 C |-- STEP_FIZHI_CORR
207 C |
208 C |-- FLT_MAIN
209 C |
210 C |-- DO_STATEVARS_TAVE
211 C |
212 C |-- NEST_PARENT_IO_2
213 C |-- NEST_CHILD_TRANSP
214 C |
215 C |-- MONITOR
216 C |
217 C |-- COST_TILE
218 C |
219 C |-- DO_THE_MODEL_IO
220 C |
221 C |-- PTRACERS_RESET
222 C |
223 C |-- DO_WRITE_PICKUP
224 C |
225 C |-- AUTODIFF_INADMODE_SET
226 C |
227 C |-- SHOWFLOPS_INLOOP
228
229 C !USES:
230 IMPLICIT NONE
231 C == Global variables ==
232 #include "SIZE.h"
233 #include "EEPARAMS.h"
234 #include "PARAMS.h"
235 #include "DYNVARS.h"
236
237 #ifdef HAVE_SIGREG
238 #include "SIGREG.h"
239 #endif
240
241 #ifdef ALLOW_SHAP_FILT
242 # include "SHAP_FILT.h"
243 #endif
244 #ifdef ALLOW_ZONAL_FILT
245 # include "ZONAL_FILT.h"
246 #endif
247
248 #ifdef ALLOW_LONGSTEP
249 # include "LONGSTEP_PARAMS.h"
250 # include "LONGSTEP.h"
251 #endif
252
253 #ifdef ALLOW_AUTODIFF
254 # include "AUTODIFF_MYFIELDS.h"
255 # include "FFIELDS.h"
256 # include "SURFACE.h"
257
258 # include "tamc.h"
259 # ifdef ALLOW_CTRL
260 # include "CTRL_SIZE.h"
261 # include "ctrl.h"
262 # include "ctrl_dummy.h"
263 # include "CTRL_GENARR.h"
264 # include "CTRL_OBCS.h"
265 # endif
266 # ifdef ALLOW_COST
267 # include "cost.h"
268 # endif
269 # ifdef ALLOW_ECCO
270 # include "ecco_cost.h"
271 # endif
272 # include "EOS.h"
273 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
274 # include "GRID.h"
275 # endif
276 # ifdef ALLOW_GMREDI
277 # include "GMREDI.h"
278 # endif
279 # ifdef ALLOW_EXF
280 # include "EXF_FIELDS.h"
281 # include "EXF_PARAM.h"
282 # ifdef ALLOW_BULKFORMULAE
283 # include "EXF_CONSTANTS.h"
284 # endif
285 # endif
286 # ifdef ALLOW_CD_CODE
287 # include "CD_CODE_VARS.h"
288 # endif
289 # ifdef ALLOW_GENERIC_ADVDIFF
290 # include "GAD.h"
291 # include "GAD_SOM_VARS.h"
292 # endif
293 # ifdef ALLOW_GGL90
294 # include "GGL90.h"
295 # endif
296 # ifdef ALLOW_PTRACERS
297 # include "PTRACERS_SIZE.h"
298 # include "PTRACERS_FIELDS.h"
299 # endif
300 # ifdef ALLOW_GCHEM
301 # include "GCHEM_SIZE.h"
302 # include "GCHEM_FIELDS.h"
303 # endif
304 # ifdef ALLOW_CFC
305 # include "CFC.h"
306 # endif
307 # ifdef ALLOW_DIC
308 # include "DIC_VARS.h"
309 # include "DIC_LOAD.h"
310 # include "DIC_ATMOS.h"
311 # include "DIC_COST.h"
312 # endif
313 # ifdef ALLOW_BLING
314 # include "BLING_VARS.h"
315 # include "BLING_LOAD.h"
316 # endif
317 # ifdef ALLOW_OBCS
318 # include "OBCS_PARAMS.h"
319 # include "OBCS_FIELDS.h"
320 # include "OBCS_SEAICE.h"
321 # ifdef ALLOW_PTRACERS
322 # include "OBCS_PTRACERS.h"
323 # endif
324 # endif
325 # ifdef ALLOW_THSICE
326 # include "THSICE_PARAMS.h"
327 # include "THSICE_SIZE.h"
328 # include "THSICE_VARS.h"
329 # include "THSICE_COST.h"
330 # endif
331 # ifdef ALLOW_SEAICE
332 # include "SEAICE_SIZE.h"
333 # include "SEAICE.h"
334 # include "SEAICE_COST.h"
335 # endif
336 # ifdef ALLOW_SALT_PLUME
337 # include "SALT_PLUME.h"
338 # endif
339 # ifdef ALLOW_SHELFICE
340 # include "SHELFICE.h"
341 # include "SHELFICE_COST.h"
342 # endif
343 # ifdef ALLOW_STREAMICE
344 # include "STREAMICE.h"
345 # include "STREAMICE_ADV.h"
346 # include "STREAMICE_BDRY.h"
347 # include "STREAMICE_CG.h"
348 # endif
349 # ifdef ALLOW_EBM
350 # include "EBM.h"
351 # endif
352 # ifdef ALLOW_KPP
353 # include "KPP.h"
354 # endif
355 # ifdef ALLOW_RBCS
356 # include "RBCS_SIZE.h"
357 # include "RBCS_FIELDS.h"
358 # endif
359 # ifdef ALLOW_OFFLINE
360 # include "OFFLINE.h"
361 # endif
362 # ifdef ALLOW_CG2D_NSA
363 # include "CG2D.h"
364 # endif
365 #endif /* ALLOW_AUTODIFF */
366
367 C !INPUT/OUTPUT PARAMETERS:
368 C == Routine arguments ==
369 C note: under the multi-threaded model myIter and
370 C myTime are local variables passed around as routine
371 C arguments. Although this is fiddly it saves the need to
372 C impose additional synchronisation points when they are
373 C updated.
374 C myTime :: time counter for this thread
375 C myIter :: iteration counter for this thread
376 C myThid :: thread number for this instance of the routine.
377 INTEGER iloop
378 _RL myTime
379 INTEGER myIter
380 INTEGER myThid
381
382 C !LOCAL VARIABLES:
383 C == Local variables ==
384 C modelEnd :: true if reaching the end of the run
385 C myTimeBeg :: time at beginning of time step (needed by longstep)
386 C myIterBeg :: iteration number at beginning of time step
387 LOGICAL modelEnd
388 #ifdef ALLOW_LONGSTEP
389 INTEGER myIterBeg
390 _RL myTimeBeg
391 #endif /* ALLOW_LONGSTEP */
392 CEOP
393
394 #ifdef ALLOW_DEBUG
395 IF (debugMode) CALL DEBUG_ENTER('FORWARD_STEP',myThid)
396 #endif
397
398 #ifdef ALLOW_AUTODIFF
399 CALL AUTODIFF_INADMODE_UNSET( myThid )
400 #endif
401
402 #ifdef ALLOW_AUTODIFF
403 C-- Reset the model iteration counter and the model time.
404 myIter = nIter0 + (iloop-1)
405 myTime = startTime + deltaTClock*(iLoop-1)
406 #endif
407
408 #ifdef ALLOW_LONGSTEP
409 C store this for longstep_average with staggerTimeStep
410 C which is called after myIter and myTime are incremented
411 C but needs iter/time at beginning of time step
412 myIterBeg = myIter
413 myTimeBeg = myTime
414 #endif /* ALLOW_LONGSTEP */
415
416 #ifdef ALLOW_AUTODIFF_TAMC
417 c**************************************
418 #include "checkpoint_lev1_directives.h"
419 #include "checkpoint_lev1_template.h"
420 c**************************************
421 #endif
422
423 C-- Reset geometric factors (hFacC,W,S & recip_hFac) to their current values:
424 C added to simplify adjoint derivation - no effect in forward run
425 #ifdef NONLIN_FRSURF
426 #ifndef ALLOW_AUTODIFF
427 IF ( doResetHFactors ) THEN
428 #endif
429 CALL RESET_NLFS_VARS( myTime, myIter, myThid )
430 IF ( select_rStar.GT.0 ) THEN
431 # ifndef DISABLE_RSTAR_CODE
432 # ifdef ALLOW_AUTODIFF_TAMC
433 CADJ STORE rStarFacC, rStarFacS, rStarFacW =
434 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
435 # endif
436 CALL TIMER_START('UPDATE_R_STAR [FORWARD_STEP]',myThid)
437 CALL UPDATE_R_STAR( .FALSE., myTime, myIter, myThid )
438 CALL TIMER_STOP ('UPDATE_R_STAR [FORWARD_STEP]',myThid)
439 # endif /* DISABLE_RSTAR_CODE */
440 ELSE
441 #ifdef ALLOW_AUTODIFF_TAMC
442 CADJ STORE hFac_surfC, hFac_surfS, hFac_surfW
443 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
444 #endif
445 CALL TIMER_START('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
446 CALL UPDATE_SURF_DR( .FALSE., myTime, myIter, myThid )
447 CALL TIMER_STOP ('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
448 ENDIF
449 #ifdef ALLOW_AUTODIFF_TAMC
450 CADJ STORE hFacC, hFacS, hFacW =
451 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
452 CADJ STORE recip_hFacC, recip_hFacS, recip_hFacW =
453 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
454 #endif
455 #ifndef ALLOW_AUTODIFF
456 ENDIF
457 #endif
458 #endif /* NONLIN_FRSURF */
459
460 #ifdef ALLOW_PTRACERS
461 C-- Switch on/off individual tracer time-stepping
462 IF ( usePTRACERS ) THEN
463 CALL PTRACERS_SWITCH_ONOFF( myTime, myIter, myThid )
464 ENDIF
465 #endif /* ALLOW_PTRACERS */
466
467 C-- Switch on/off diagnostics for snap-shot output:
468 #ifdef ALLOW_DIAGNOSTICS
469 IF ( useDiagnostics ) THEN
470 CALL DIAGNOSTICS_SWITCH_ONOFF( myTime, myIter, myThid )
471 C-- State-variables diagnostics
472 CALL TIMER_START('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
473 CALL DO_STATEVARS_DIAGS( myTime, 0, myIter, myThid )
474 CALL TIMER_STOP ('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
475 ENDIF
476 #endif /* ALLOW_DIAGNOSTICS */
477
478 #ifdef ALLOW_NEST_CHILD
479 IF ( useNEST_CHILD) THEN
480 CALL NEST_CHILD_SETMEMO( myTime, myIter, myThid )
481 ENDIF
482 #endif /* ALLOW_NEST_CHILD */
483
484 #ifdef ALLOW_NEST_PARENT
485 IF ( useNEST_PARENT) THEN
486 CALL NEST_PARENT_IO_1( myTime, myIter, myThid )
487 ENDIF
488 #endif /* ALLOW_NEST_PARENT */
489
490 C-- Call driver to load external forcing fields from file
491 #ifdef ALLOW_DEBUG
492 IF (debugMode) CALL DEBUG_CALL('LOAD_FIELDS_DRIVER',myThid)
493 #endif
494 #ifdef ALLOW_AUTODIFF_TAMC
495 cph Important STORE that avoids hidden recomp. of load_fields_driver
496 CADJ STORE theta = comlev1, key = ikey_dynamics,
497 CADJ & kind = isbyte
498 CADJ STORE uVel, vVel = comlev1, key = ikey_dynamics,
499 CADJ & kind = isbyte
500 #endif
501 CALL TIMER_START('LOAD_FIELDS_DRIVER [FORWARD_STEP]',myThid)
502 CALL LOAD_FIELDS_DRIVER( myTime, myIter, myThid )
503 CALL TIMER_STOP ('LOAD_FIELDS_DRIVER [FORWARD_STEP]',myThid)
504
505 C-- Call Bulk-Formulae forcing package
506 #ifdef ALLOW_BULK_FORCE
507 IF ( useBulkForce ) THEN
508 #ifdef ALLOW_DEBUG
509 IF (debugMode) CALL DEBUG_CALL('BULKF_FORCING',myThid)
510 #endif
511 CALL TIMER_START('BULKF_FORCING [FORWARD_STEP]',myThid)
512 C- calculate qnet and empmr (and wind stress)
513 CALL BULKF_FORCING( myTime, myIter, myThid )
514 CALL TIMER_STOP ('BULKF_FORCING [FORWARD_STEP]',myThid)
515 ENDIF
516 #endif /* ALLOW_BULK_FORCE */
517
518 C-- Call external chepaml forcing package
519 #ifdef ALLOW_CHEAPAML
520 IF ( useCheapAML ) THEN
521 #ifdef ALLOW_DEBUG
522 IF (debugMode) CALL DEBUG_CALL('CHEAPAML',myThid)
523 #endif
524 CALL TIMER_START('CHEAPAML [FORWARD_STEP]',myThid)
525 C- calculate qnet (and wind stress)
526 CALL CHEAPAML( myTime, myIter,myThid )
527 CALL TIMER_STOP ('CHEAPAML [FORWARD_STEP]',myThid)
528 ENDIF
529 #endif /*ALLOW_CHEAPAML */
530
531 #ifdef ALLOW_CTRL
532 C-- Add control vector for forcing and parameter fields
533 IF ( useCTRL ) THEN
534 CALL TIMER_START('CTRL_MAP_FORCING [FORWARD_STEP]',myThid)
535 CALL CTRL_MAP_FORCING( myTime, myIter, myThid )
536 CALL TIMER_STOP ('CTRL_MAP_FORCING [FORWARD_STEP]',myThid)
537 ENDIF
538 #endif
539
540 #ifdef ALLOW_AUTODIFF_MONITOR
541 CALL DUMMY_IN_STEPPING( myTime, myIter, myThid )
542 #endif
543
544 #ifdef COMPONENT_MODULE
545 IF ( useCoupler ) THEN
546 C Post coupling data that I export.
547 C Read in coupling data that I import.
548 CALL TIMER_START('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid)
549 CALL CPL_EXPORT_IMPORT_DATA( myTime, myIter, myThid )
550 CALL TIMER_STOP ('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid)
551 ENDIF
552 #endif /* COMPONENT_MODULE */
553 #ifdef ALLOW_OASIS
554 IF ( useOASIS ) THEN
555 CALL TIMER_START('OASIS_PUT-GET [FORWARD_STEP]',myThid)
556 C Post coupling data that I export.
557 CALL OASIS_PUT( myTime, myIter, myThid )
558 C Read in coupling data that I import.
559 CALL OASIS_GET( myTime, myIter, myThid )
560 CALL TIMER_STOP ('OASIS_PUT-GET [FORWARD_STEP]',myThid)
561 ENDIF
562 #endif /* ALLOW_OASIS */
563
564 #ifdef ALLOW_EBM
565 IF ( useEBM ) THEN
566 # ifdef ALLOW_DEBUG
567 IF (debugMode) CALL DEBUG_CALL('EBM',myThid)
568 # endif
569 CALL TIMER_START('EBM [FORWARD_STEP]',myThid)
570 CALL EBM_DRIVER ( myTime, myIter, myThid )
571 CALL TIMER_STOP ('EBM [FORWARD_STEP]',myThid)
572 ENDIF
573 #endif /* ALLOW_EBM */
574
575 C-- Step forward fields and calculate time tendency terms.
576
577 #ifdef ALLOW_DEBUG
578 IF (debugMode) CALL DEBUG_CALL('DO_ATMOSPHERIC_PHYS',myThid)
579 #endif
580 CALL TIMER_START('DO_ATMOSPHERIC_PHYS [FORWARD_STEP]',myThid)
581 CALL DO_ATMOSPHERIC_PHYS( myTime, myIter, myThid )
582 CALL TIMER_STOP ('DO_ATMOSPHERIC_PHYS [FORWARD_STEP]',myThid)
583
584 #ifdef ALLOW_AUTODIFF_TAMC
585 CADJ STORE surfaceForcingTice = comlev1, key = ikey_dynamics,
586 CADJ & kind = isbyte
587 # ifdef ALLOW_KPP
588 CADJ STORE uVel, vVel = comlev1, key = ikey_dynamics,
589 CADJ & kind = isbyte
590 # endif /* ALLOW_KPP */
591 # ifdef EXACT_CONSERV
592 CADJ STORE EmPmR = comlev1, key=ikey_dynamics, kind=isbyte
593 CADJ STORE PmEpR = comlev1, key=ikey_dynamics, kind=isbyte
594 # endif
595 # ifdef ALLOW_OBCS
596 CADJ STORE salt = comlev1, key=ikey_dynamics, kind=isbyte
597 CADJ STORE totphihyd = comlev1, key=ikey_dynamics, kind=isbyte
598 # ifdef ALLOW_OBCS_STEVENS
599 CADJ STORE gsNm1 = comlev1, key=ikey_dynamics, kind=isbyte
600 CADJ STORE gtNm1 = comlev1, key=ikey_dynamics, kind=isbyte
601 # endif
602 # endif /* ALLOW_OBCS */
603 # ifdef ALLOW_PTRACERS
604 CADJ STORE pTracer = comlev1, key = ikey_dynamics,
605 CADJ & kind = isbyte
606 # endif /* ALLOW_PTRACERS */
607 # ifdef ALLOW_DEPTH_CONTROL
608 CADJ STORE hFacC = comlev1, key = ikey_dynamics, kind = isbyte
609 # endif
610 #endif /* ALLOW_AUTODIFF_TAMC */
611
612 #ifdef ALLOW_DEBUG
613 IF (debugMode) CALL DEBUG_CALL('DO_OCEANIC_PHYS',myThid)
614 #endif
615 CALL TIMER_START('DO_OCEANIC_PHYS [FORWARD_STEP]',myThid)
616 CALL DO_OCEANIC_PHYS( myTime, myIter, myThid )
617 CALL TIMER_STOP ('DO_OCEANIC_PHYS [FORWARD_STEP]',myThid)
618
619 #ifdef ALLOW_STREAMICE
620 IF (useStreamIce) THEN
621 CALL STREAMICE_TIMESTEP ( myThid, myIter,
622 & iLoop, myTime )
623 ENDIF
624 #endif
625
626 #ifdef ALLOW_AUTODIFF_TAMC
627 CADJ STORE EmPmR = comlev1, key = ikey_dynamics,
628 CADJ & kind = isbyte
629 # ifdef EXACT_CONSERV
630 CADJ STORE PmEpR = comlev1, key = ikey_dynamics,
631 CADJ & kind = isbyte
632 # endif
633 cph-test(
634 CADJ STORE surfaceForcingU = comlev1, key=ikey_dynamics, kind=isbyte
635 CADJ STORE surfaceForcingV = comlev1, key=ikey_dynamics, kind=isbyte
636 CADJ STORE qsw = comlev1, key = ikey_dynamics, kind = isbyte
637 # ifdef ATMOSPHERIC_LOADING
638 CADJ STORE phi0surf = comlev1, key = ikey_dynamics, kind = isbyte
639 # endif
640 cph-test)
641
642 # ifdef ALLOW_DEPTH_CONTROL
643 CADJ STORE hFacC, hFacS, hFacW
644 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
645 CADJ STORE recip_hFacC, recip_hFacS, recip_hFacW
646 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
647 CADJ STORE surfaceForcingU, surfaceForcingV =
648 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
649 # endif
650 #endif /* ALLOW_AUTODIFF_TAMC */
651
652 #ifdef ALLOW_GCHEM
653 IF ( useGCHEM ) THEN
654 #ifdef ALLOW_AUTODIFF_TAMC
655 CADJ STORE pTracer = comlev1, key = ikey_dynamics,
656 CADJ & kind = isbyte
657 CADJ STORE theta, salt = comlev1, key = ikey_dynamics,
658 CADJ & kind = isbyte
659 #endif
660 #ifdef ALLOW_DEBUG
661 IF (debugMode) CALL DEBUG_CALL('GCHEM_CALC_TENDENCY',myThid)
662 #endif
663 CALL TIMER_START('GCHEM_CALC_TENDENCY [FORWARD_STEP]',myThid)
664 CALL GCHEM_CALC_TENDENCY( myTime, myIter, myThid )
665 CALL TIMER_STOP ('GCHEM_CALC_TENDENCY [FORWARD_STEP]',myThid)
666 ENDIF
667 #endif /* ALLOW_GCHEM */
668
669 #ifdef ALLOW_AUTODIFF_TAMC
670 cph needed to be moved here from do_oceanic_physics
671 cph to be visible down the road
672
673 CADJ STORE rhoInSitu = comlev1, key = ikey_dynamics,
674 CADJ & kind = isbyte
675 CADJ STORE surfaceForcingS = comlev1, key = ikey_dynamics,
676 CADJ & kind = isbyte
677 CADJ STORE surfaceForcingT = comlev1, key = ikey_dynamics,
678 CADJ & kind = isbyte
679 CADJ STORE surfaceForcingTice = comlev1, key = ikey_dynamics,
680 CADJ & kind = isbyte
681 CADJ STORE IVDConvCount = comlev1, key = ikey_dynamics,
682 CADJ & kind = isbyte
683 # ifdef ALLOW_PTRACERS
684 CADJ STORE surfaceForcingPTr = comlev1, key = ikey_dynamics,
685 CADJ & kind = isbyte
686 # endif
687
688 # ifdef ALLOW_GMREDI
689 CADJ STORE Kwx = comlev1, key = ikey_dynamics,
690 CADJ & kind = isbyte
691 CADJ STORE Kwy = comlev1, key = ikey_dynamics,
692 CADJ & kind = isbyte
693 CADJ STORE Kwz = comlev1, key = ikey_dynamics,
694 CADJ & kind = isbyte
695 # ifdef GM_BOLUS_ADVEC
696 CADJ STORE GM_PsiX = comlev1, key = ikey_dynamics,
697 CADJ & kind = isbyte
698 CADJ STORE GM_PsiY = comlev1, key = ikey_dynamics,
699 CADJ & kind = isbyte
700 # endif
701 # endif
702
703 # ifdef ALLOW_KPP
704 CADJ STORE KPPghat = comlev1, key = ikey_dynamics,
705 CADJ & kind = isbyte
706 CADJ STORE KPPfrac = comlev1, key = ikey_dynamics,
707 CADJ & kind = isbyte
708 CADJ STORE KPPdiffKzS = comlev1, key = ikey_dynamics,
709 CADJ & kind = isbyte
710 CADJ STORE KPPdiffKzT = comlev1, key = ikey_dynamics,
711 CADJ & kind = isbyte
712 # endif
713
714 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
715 CADJ STORE theta,salt = comlev1, key = ikey_dynamics, kind = isbyte
716 CADJ STORE etaH = comlev1, key = ikey_dynamics, kind = isbyte
717 # ifdef ALLOW_CD_CODE
718 CADJ STORE etaNm1 = comlev1, key = ikey_dynamics, kind = isbyte
719 # endif
720 # ifndef DISABLE_RSTAR_CODE
721 CADJ STORE rStarExpC = comlev1, key = ikey_dynamics, kind = isbyte
722 # endif
723 # endif
724 #endif /* ALLOW_AUTODIFF_TAMC */
725
726 #ifdef ALLOW_LONGSTEP
727 IF ( usePTRACERS .AND. LS_whenToSample .EQ. 0 ) THEN
728 C Average all variables before advection (but after do_oceanic_phys
729 C where Qsw, KPP and GMRedi stuff is computed).
730 C This is like diagnostics package and will reproduce offline
731 C results.
732 #ifdef ALLOW_DEBUG
733 IF (debugMode) CALL DEBUG_CALL('LONGSTEP_AVERAGE',myThid)
734 #endif
735 CALL TIMER_START('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid)
736 CALL LONGSTEP_AVERAGE( myTime, myIter, myThid )
737 CALL TIMER_STOP ('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid)
738
739 #ifdef ALLOW_DEBUG
740 IF (debugMode)
741 & CALL DEBUG_CALL('LONGSTEP_THERMODYNAMICS',myThid)
742 #endif
743 CALL TIMER_START('LONGSTEP_THERMODYNAMICS [FORWARD_STEP]',
744 & myThid)
745 CALL LONGSTEP_THERMODYNAMICS( myTime, myIter, myThid )
746 CALL TIMER_STOP ('LONGSTEP_THERMODYNAMICS [FORWARD_STEP]',
747 & myThid)
748 ENDIF
749 #endif /* ALLOW_LONGSTEP */
750
751 IF ( .NOT.staggerTimeStep ) THEN
752 #ifdef ALLOW_AUTODIFF_TAMC
753 CADJ STORE wVel = comlev1, key = ikey_dynamics, kind = isbyte
754 #endif
755 #ifdef ALLOW_DEBUG
756 IF (debugMode) CALL DEBUG_CALL('THERMODYNAMICS',myThid)
757 #endif
758 CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',myThid)
759 CALL THERMODYNAMICS( myTime, myIter, myThid )
760 CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',myThid)
761 C-- if not staggerTimeStep: end
762 ENDIF
763
764 #ifdef ALLOW_LONGSTEP
765 IF ( usePTRACERS .AND. LS_whenToSample .EQ. 1 ) THEN
766 C Average T and S after thermodynamics, but U,V,W before dynamics.
767 C This will reproduce online results with staggerTimeStep=.FALSE.
768 C for LS_nIter=1
769 #ifdef ALLOW_DEBUG
770 IF (debugMode) CALL DEBUG_CALL('LONGSTEP_AVERAGE',myThid)
771 #endif
772 CALL TIMER_START('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid)
773 CALL LONGSTEP_AVERAGE( myTime, myIter, myThid )
774 CALL TIMER_STOP ('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid)
775
776 #ifdef ALLOW_DEBUG
777 IF (debugMode)
778 & CALL DEBUG_CALL('LONGSTEP_THERMODYNAMICS',myThid)
779 #endif
780 CALL TIMER_START('LONGSTEP_THERMODYNAMICS [FORWARD_STEP]',
781 & myThid)
782 CALL LONGSTEP_THERMODYNAMICS( myTime, myIter, myThid )
783 CALL TIMER_STOP ('LONGSTEP_THERMODYNAMICS [FORWARD_STEP]',
784 & myThid)
785 ENDIF
786 #endif /* ALLOW_LONGSTEP */
787
788 c #ifdef ALLOW_NONHYDROSTATIC
789 IF ( implicitIntGravWave ) THEN
790 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
791 CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid )
792 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
793 ENDIF
794 c #endif
795
796 #ifdef ALLOW_AUTODIFF_TAMC
797 CADJ STORE etaN = comlev1, key = ikey_dynamics, kind = isbyte
798 # ifdef ALLOW_DEPTH_CONTROL
799 CADJ STORE hFacC, hFacS, hFacW
800 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
801 CADJ STORE recip_hFacC, recip_hFacS, recip_hFacW
802 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
803 # endif /* ALLOW_DEPTH_CONTROL */
804 #endif /* ALLOW_AUTODIFF_TAMC */
805
806 C-- Step forward fields and calculate time tendency terms.
807 #ifdef ALLOW_MOM_STEPPING
808 #ifndef ALLOW_AUTODIFF
809 IF ( momStepping ) THEN
810 #endif
811 #ifdef ALLOW_DEBUG
812 IF (debugMode) CALL DEBUG_CALL('DYNAMICS',myThid)
813 #endif
814 CALL TIMER_START('DYNAMICS [FORWARD_STEP]',myThid)
815 CALL DYNAMICS( myTime, myIter, myThid )
816 CALL TIMER_STOP ('DYNAMICS [FORWARD_STEP]',myThid)
817 #ifndef ALLOW_AUTODIFF
818 ENDIF
819 #endif
820 #endif /* ALLOW_MOM_STEPPING */
821
822 #ifdef ALLOW_AUTODIFF_TAMC
823 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
824 CADJ STORE gU, gV = comlev1, key = ikey_dynamics,
825 CADJ & kind = isbyte
826 # endif
827 #endif
828
829 C-- Update time-counter
830 myIter = nIter0 + iLoop
831 myTime = startTime + deltaTClock*iLoop
832
833 #ifdef ALLOW_MNC
834 C Update MNC time information
835 IF ( useMNC ) THEN
836 CALL MNC_UPDATE_TIME( myTime, myIter, myThid )
837 ENDIF
838 #endif /* ALLOW_MNC */
839
840 #ifdef ALLOW_OFFLINE
841 C Load new Offline fields and update state-variable
842 IF ( useOffLine ) THEN
843 #ifdef ALLOW_DEBUG
844 IF (debugMode) CALL DEBUG_CALL('OFFLINE_FIELDS_LOAD',myThid)
845 #endif /* ALLOW_DEBUG */
846 CALL TIMER_START('OFFLINE_FLDS_LOAD [FORWARD_STEP]',myThid)
847 CALL OFFLINE_FIELDS_LOAD( myTime, myIter, myThid )
848 CALL TIMER_STOP ('OFFLINE_FLDS_LOAD [FORWARD_STEP]',myThid)
849 ENDIF
850 #endif /* ALLOW_OFFLINE */
851
852 C-- Update geometric factors:
853 #ifdef NONLIN_FRSURF
854 C- update hfacC,W,S and recip_hFac according to etaH(n+1) :
855 IF ( select_rStar.GT.0 ) THEN
856 # ifndef DISABLE_RSTAR_CODE
857 # ifdef ALLOW_AUTODIFF_TAMC
858 CADJ STORE rStarFacC, rStarFacS, rStarFacW =
859 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
860 # endif
861 CALL TIMER_START('UPDATE_R_STAR [FORWARD_STEP]',myThid)
862 CALL UPDATE_R_STAR( .TRUE., myTime, myIter, myThid )
863 CALL TIMER_STOP ('UPDATE_R_STAR [FORWARD_STEP]',myThid)
864 # endif /* DISABLE_RSTAR_CODE */
865 ELSEIF ( selectSigmaCoord.NE.0 ) THEN
866 # ifndef DISABLE_SIGMA_CODE
867 CALL UPDATE_SIGMA( etaH, myTime, myIter, myThid )
868 # endif /* DISABLE_RSTAR_CODE */
869 ELSE
870 # ifdef ALLOW_AUTODIFF_TAMC
871 CADJ STORE hFac_surfC, hFac_surfS, hFac_surfW
872 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
873 # endif
874 CALL TIMER_START('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
875 CALL UPDATE_SURF_DR( .TRUE., myTime, myIter, myThid )
876 CALL TIMER_STOP ('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
877 ENDIF
878 # ifdef ALLOW_AUTODIFF_TAMC
879 CADJ STORE hFacC, hFacS, hFacW =
880 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
881 CADJ STORE recip_hFacC, recip_hFacS, recip_hFacW =
882 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
883 # ifdef ALLOW_CG2D_NSA
884 CADJ STORE aW2d, aS2d, aC2d =
885 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
886 CADJ STORE pC, pS, pW =
887 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
888 # endif
889 # endif
890 #endif /* NONLIN_FRSURF */
891
892 #if ( defined NONLIN_FRSURF || defined ALLOW_SOLVE4_PS_AND_DRAG )
893 C- update CG2D matrix (and preconditioner)
894 IF ( momStepping .AND.
895 & ( nonlinFreeSurf.GT.2 .OR. selectImplicitDrag.EQ.2 ) ) THEN
896 CALL TIMER_START('UPDATE_CG2D [FORWARD_STEP]',myThid)
897 CALL UPDATE_CG2D( myTime, myIter, myThid )
898 CALL TIMER_STOP ('UPDATE_CG2D [FORWARD_STEP]',myThid)
899 ENDIF
900 #endif /* NONLIN_FRSURF or ALLOW_SOLVE4_PS_AND_DRAG */
901
902 C-- Apply Filters to u*,v* before SOLVE_FOR_PRESSURE
903 #ifdef ALLOW_SHAP_FILT
904 IF (useSHAP_FILT .AND. shap_filt_uvStar) THEN
905 CALL TIMER_START('SHAP_FILT_UV [FORWARD_STEP]',myThid)
906 IF (implicDiv2DFlow.LT.1.) THEN
907 C-- Explicit+Implicit part of the Barotropic Flow Divergence
908 C => Filtering of uVel,vVel is necessary
909 CALL SHAP_FILT_APPLY_UV( uVel,vVel,
910 & myTime, myIter, myThid )
911 ENDIF
912 CALL SHAP_FILT_APPLY_UV( gU,gV,myTime,myIter,myThid)
913 CALL TIMER_STOP ('SHAP_FILT_UV [FORWARD_STEP]',myThid)
914 ENDIF
915 #endif
916 #ifdef ALLOW_ZONAL_FILT
917 IF (useZONAL_FILT .AND. zonal_filt_uvStar) THEN
918 CALL TIMER_START('ZONAL_FILT_UV [FORWARD_STEP]',myThid)
919 IF (implicDiv2DFlow.LT.1.) THEN
920 C-- Explicit+Implicit part of the Barotropic Flow Divergence
921 C => Filtering of uVel,vVel is necessary
922 CALL ZONAL_FILT_APPLY_UV( uVel, vVel, myThid )
923 ENDIF
924 CALL ZONAL_FILT_APPLY_UV( gU, gV, myThid )
925 CALL TIMER_STOP ('ZONAL_FILT_UV [FORWARD_STEP]',myThid)
926 ENDIF
927 #endif
928
929 C-- Solve elliptic equation(s).
930 C Two-dimensional only for conventional hydrostatic or
931 C three-dimensional for non-hydrostatic and/or IGW scheme.
932 IF ( momStepping ) THEN
933 #ifdef ALLOW_AUTODIFF_TAMC
934 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
935 CADJ STORE uVel, vVel
936 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
937 CADJ STORE EmPmR,hFacS,hFacW
938 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
939 # endif
940 #endif
941 CALL TIMER_START('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
942 CALL SOLVE_FOR_PRESSURE(myTime, myIter, myThid)
943 CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
944 ENDIF
945
946 C-- Correct divergence in flow field and cycle time-stepping momentum
947 #ifdef ALLOW_MOM_STEPPING
948 #ifndef ALLOW_AUTODIFF
949 IF ( momStepping ) THEN
950 #endif
951 #ifdef ALLOW_AUTODIFF_TAMC
952 # ifdef ALLOW_DEPTH_CONTROL
953 CADJ STORE etaN, uVel,vVel
954 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
955 # endif /* ALLOW_DEPTH_CONTROL */
956 #endif /* ALLOW_AUTODIFF_TAMC */
957 CALL TIMER_START('MOM_CORRECTION_STEP [FORWARD_STEP]',myThid)
958 CALL MOMENTUM_CORRECTION_STEP(myTime, myIter, myThid)
959 CALL TIMER_STOP ('MOM_CORRECTION_STEP [FORWARD_STEP]',myThid)
960 #ifndef ALLOW_AUTODIFF
961 ENDIF
962 #endif
963 #endif /* ALLOW_MOM_STEPPING */
964 #ifdef ALLOW_AUTODIFF_TAMC
965 CADJ STORE uVel, vVel = comlev1, key = ikey_dynamics, kind = isbyte
966 #endif
967
968 IF ( calc_wVelocity ) THEN
969 C-- Integrate continuity vertically for vertical velocity
970 C (+ update "etaN" & "etaH", exact volume conservation):
971 CALL TIMER_START('INTEGR_CONTINUITY [FORWARD_STEP]',myThid)
972 CALL INTEGR_CONTINUITY( uVel, vVel, myTime, myIter, myThid)
973 CALL TIMER_STOP ('INTEGR_CONTINUITY [FORWARD_STEP]',myThid)
974 ENDIF
975
976 #ifdef NONLIN_FRSURF
977 IF ( select_rStar.NE.0 ) THEN
978 # ifndef DISABLE_RSTAR_CODE
979 # ifdef ALLOW_AUTODIFF_TAMC
980 CADJ STORE etaH
981 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
982 CADJ STORE rStarFacC,rStarFacS,rStarFacW
983 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
984 # endif
985 C-- r* : compute the future level thickness according to etaH(n+1)
986 CALL TIMER_START('CALC_R_STAR [FORWARD_STEP]',myThid)
987 CALL CALC_R_STAR(etaH, myTime, myIter, myThid )
988 CALL TIMER_STOP ('CALC_R_STAR [FORWARD_STEP]',myThid)
989 # endif /* DISABLE_RSTAR_CODE */
990 ELSEIF ( nonlinFreeSurf.GT.0 .AND. selectSigmaCoord.EQ.0 ) THEN
991 C-- compute the future surface level thickness according to etaH(n+1)
992 # ifdef ALLOW_AUTODIFF_TAMC
993 CADJ STORE etaH = comlev1, key = ikey_dynamics,
994 CADJ & kind = isbyte
995 # endif
996 CALL TIMER_START('CALC_SURF_DR [FORWARD_STEP]',myThid)
997 CALL CALC_SURF_DR(etaH, myTime, myIter, myThid )
998 CALL TIMER_STOP ('CALC_SURF_DR [FORWARD_STEP]',myThid)
999 ENDIF
1000 # ifdef ALLOW_AUTODIFF_TAMC
1001 CADJ STORE rStarExpC
1002 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
1003 CADJ STORE hFac_surfC = comlev1, key = ikey_dynamics,
1004 CADJ & kind = isbyte
1005 CADJ STORE salt,theta = comlev1, key = ikey_dynamics,
1006 CADJ & kind = isbyte
1007 # endif
1008 #endif /* NONLIN_FRSURF */
1009
1010 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
1011 IF ( staggerTimeStep ) THEN
1012 C-- do exchanges of U,V (needed for multiDim) when using stagger time-step :
1013 #ifdef ALLOW_DEBUG
1014 IF (debugMode)
1015 & CALL DEBUG_CALL('DO_STAGGER_FIELDS_EXCH.',myThid)
1016 #endif
1017 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
1018 CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid )
1019 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
1020
1021 #ifdef ALLOW_DIAGNOSTICS
1022 C-- State-variables diagnostics
1023 IF ( useDiagnostics ) THEN
1024 CALL TIMER_START('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
1025 CALL DO_STATEVARS_DIAGS( myTime, 1, myIter, myThid )
1026 CALL TIMER_STOP ('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
1027 ENDIF
1028 #endif
1029
1030 #ifdef ALLOW_AUTODIFF_TAMC
1031 CADJ STORE wVel = comlev1, key = ikey_dynamics, kind = isbyte
1032 #endif
1033 #ifdef ALLOW_DEBUG
1034 IF (debugMode) CALL DEBUG_CALL('THERMODYNAMICS',myThid)
1035 #endif
1036 CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',myThid)
1037 CALL THERMODYNAMICS( myTime, myIter, myThid )
1038 CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',myThid)
1039
1040 C-- if staggerTimeStep: end
1041 ENDIF
1042 C---+--------+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
1043
1044 #ifdef ALLOW_AUTODIFF_TAMC
1045 cph This is needed because convective_adjustment calls
1046 cph find_rho which may use pressure()
1047 CADJ STORE totPhiHyd = comlev1, key = ikey_dynamics,
1048 CADJ & kind = isbyte
1049 #endif
1050 C-- Apply adjustments to Tracers arrays (T,S,+pTracers)
1051 CALL TIMER_START('TRC_CORRECTION_STEP [FORWARD_STEP]',myThid)
1052 CALL TRACERS_CORRECTION_STEP(myTime, myIter, myThid)
1053 CALL TIMER_STOP ('TRC_CORRECTION_STEP [FORWARD_STEP]',myThid)
1054
1055 #ifdef ALLOW_LONGSTEP
1056 IF ( usePTRACERS ) THEN
1057 IF ( LS_whenToSample .EQ. 2 ) THEN
1058 C Average everything at the end of the timestep. This will
1059 C reproduce online results with staggerTimeStep=.TRUE.
1060 C when LS_nIter=1
1061 #ifdef ALLOW_DEBUG
1062 IF (debugMode) CALL DEBUG_CALL('LONGSTEP_AVERAGE',myThid)
1063 #endif
1064 CALL TIMER_START('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid)
1065 C myIter has been update after dynamics, but the averaging window
1066 C should be determined by myIter at beginning of timestep
1067 CALL LONGSTEP_AVERAGE( myTimeBeg, myIterBeg, myThid )
1068 CALL TIMER_STOP ('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid)
1069
1070 #ifdef ALLOW_DEBUG
1071 IF (debugMode)
1072 & CALL DEBUG_CALL('LONGSTEP_THERMODYNAMICS',myThid)
1073 #endif
1074 CALL TIMER_START('LONGSTEP_THERMODYNAMICS [FORWARD_STEP]',
1075 & myThid)
1076 CALL LONGSTEP_THERMODYNAMICS( myTime, myIter, myThid )
1077 CALL TIMER_STOP ('LONGSTEP_THERMODYNAMICS [FORWARD_STEP]',
1078 & myThid)
1079 C-- if LS_whenToSample.EQ.2: end
1080 ENDIF
1081
1082 C-- Apply adjustments to passive Tracers arrays (pTracers)
1083 c CALL TIMER_START('LS_CORRECTION_STEP [FORWARD_STEP]',myThid)
1084 c CALL LONGSTEP_CORRECTION_STEP(myTime, myIter, myThid)
1085 c CALL TIMER_STOP ('LS_CORRECTION_STEP [FORWARD_STEP]',myThid)
1086 C-- if usePTRACERS: end
1087 ENDIF
1088 #endif /* ALLOW_LONGSTEP */
1089
1090 #ifdef ALLOW_GCHEM
1091 C Add separate timestepping of chemical/biological/forcing
1092 C of ptracers here in GCHEM_FORCING_SEP
1093 #ifdef ALLOW_LONGSTEP
1094 IF ( useGCHEM .AND. LS_doTimeStep ) THEN
1095 #else
1096 IF ( useGCHEM ) THEN
1097 #endif
1098 #ifdef ALLOW_AUTODIFF_TAMC
1099 CADJ STORE pTracer = comlev1, key = ikey_dynamics,
1100 CADJ & kind = isbyte
1101 CADJ STORE theta, salt = comlev1, key = ikey_dynamics,
1102 CADJ & kind = isbyte
1103 #endif
1104 #ifdef ALLOW_DEBUG
1105 IF (debugMode) CALL DEBUG_CALL('GCHEM_FORCING_SEP',myThid)
1106 #endif /* ALLOW_DEBUG */
1107 CALL TIMER_START('GCHEM_FORCING_SEP [FORWARD_STEP]',myThid)
1108 CALL GCHEM_FORCING_SEP( myTime,myIter,myThid )
1109 CALL TIMER_STOP ('GCHEM_FORCING_SEP [FORWARD_STEP]',myThid)
1110 ENDIF
1111 #endif /* ALLOW_GCHEM */
1112
1113 C-- Do "blocking" sends and receives for tendency "overlap" terms
1114 c CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
1115 c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )
1116 c CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
1117
1118 C-- Do "blocking" sends and receives for field "overlap" terms
1119 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
1120 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
1121 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
1122
1123 #ifdef ALLOW_DIAGNOSTICS
1124 IF ( useDiagnostics ) THEN
1125 CALL TIMER_START('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
1126 CALL DO_STATEVARS_DIAGS( myTime, 2, myIter, myThid )
1127 CALL TIMER_STOP ('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
1128 ENDIF
1129 #endif
1130
1131 #ifdef ALLOW_GRIDALT
1132 IF (useGRIDALT) THEN
1133 CALL GRIDALT_UPDATE(myThid)
1134 ENDIF
1135 #endif
1136
1137 #ifdef ALLOW_FIZHI
1138 IF (useFIZHI) THEN
1139 CALL TIMER_START('FIZHI [FORWARD_STEP]',myThid)
1140 CALL STEP_FIZHI_CORR ( myTime, myIter, myThid, dTtracerLev(1) )
1141 CALL TIMER_STOP ('FIZHI [FORWARD_STEP]',myThid)
1142 ENDIF
1143 #endif
1144
1145 #ifdef ALLOW_FLT
1146 C-- Calculate float trajectories
1147 IF (useFLT) THEN
1148 CALL TIMER_START('FLOATS [FORWARD_STEP]',myThid)
1149 CALL FLT_MAIN( myTime, myIter, myThid )
1150 CALL TIMER_STOP ('FLOATS [FORWARD_STEP]',myThid)
1151 ENDIF
1152 #endif
1153
1154 #ifdef ALLOW_TIMEAVE
1155 C-- State-variables time-averaging
1156 CALL TIMER_START('DO_STATEVARS_TAVE [FORWARD_STEP]',myThid)
1157 CALL DO_STATEVARS_TAVE( myTime, myIter, myThid )
1158 CALL TIMER_STOP ('DO_STATEVARS_TAVE [FORWARD_STEP]',myThid)
1159 #endif
1160
1161 #ifdef ALLOW_NEST_PARENT
1162 IF ( useNEST_PARENT) THEN
1163 CALL NEST_PARENT_IO_2( myTime, myIter, myThid )
1164 ENDIF
1165 #endif /* ALLOW_NEST_PARENT */
1166
1167 #ifdef ALLOW_NEST_CHILD
1168 IF ( useNEST_CHILD) THEN
1169 CALL NEST_CHILD_TRANSP( myTime, myIter, myThid )
1170 ENDIF
1171 #endif /* ALLOW_NEST_CHILD */
1172
1173 #ifdef ALLOW_MONITOR
1174 IF ( monitorFreq.GT.0. .OR. adjMonitorFreq.GT.0. ) THEN
1175 C-- Check status of solution (statistics, cfl, etc...)
1176 CALL TIMER_START('MONITOR [FORWARD_STEP]',myThid)
1177 CALL MONITOR( myTime, myIter, myThid )
1178 CALL TIMER_STOP ('MONITOR [FORWARD_STEP]',myThid)
1179 ENDIF
1180 #endif /* ALLOW_MONITOR */
1181
1182 #ifdef ALLOW_COST
1183 C-- compare model with data and compute cost function
1184 C-- this is done after exchanges to allow interpolation
1185 CALL TIMER_START('COST_TILE [FORWARD_STEP]',myThid)
1186 CALL COST_TILE ( myTime, myIter, myThid )
1187 CALL TIMER_STOP ('COST_TILE [FORWARD_STEP]',myThid)
1188 #endif
1189
1190 #ifdef ALLOW_ECCO
1191 C-- Diagnoze variables for pkg/ecco averaging and cost function purposes
1192 IF ( useECCO ) CALL ECCO_PHYS( myThid )
1193 #endif
1194
1195 C-- Check if it has reached the end of simulation
1196 modelEnd = myTime.EQ.endTime .OR. myIter.EQ.nEndIter
1197 #ifdef HAVE_SIGREG
1198 IF ( useSIGREG ) THEN
1199 modelEnd = modelEnd .OR. ( i_got_signal.GT.0 )
1200 ENDIF
1201 #endif /* HAVE_SIGREG */
1202
1203 C-- Do IO if needed.
1204 CALL TIMER_START('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
1205 CALL DO_THE_MODEL_IO( modelEnd, myTime, myIter, myThid )
1206 CALL TIMER_STOP ('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
1207
1208 #ifdef ALLOW_PTRACERS
1209 C Reset the ptracers (but after the io is done)
1210 IF ( usePTRACERS ) THEN
1211 CALL TIMER_START('PTRACERS_RESET [FORWARD_STEP]',myThid)
1212 CALL PTRACERS_RESET( myTime, myIter, myThid )
1213 CALL TIMER_STOP ('PTRACERS_RESET [FORWARD_STEP]',myThid)
1214 ENDIF
1215 #endif /* ALLOW_PTRACERS */
1216
1217 C-- Save state for restarts
1218 CALL TIMER_START('DO_WRITE_PICKUP [FORWARD_STEP]',myThid)
1219 CALL DO_WRITE_PICKUP( modelEnd, myTime, myIter, myThid )
1220 CALL TIMER_STOP ('DO_WRITE_PICKUP [FORWARD_STEP]',myThid)
1221
1222 #ifdef HAVE_SIGREG
1223 IF ( useSIGREG ) THEN
1224 IF ( modelEnd .AND. i_got_signal.GT.0 ) THEN
1225 STOP 'Checkpoint completed -- killed by signal handler'
1226 ENDIF
1227 ENDIF
1228 #endif /* HAVE_SIGREG */
1229
1230 #ifdef ALLOW_AUTODIFF
1231 CALL AUTODIFF_INADMODE_SET( myThid )
1232 #endif
1233
1234 #ifdef ALLOW_SHOWFLOPS
1235 CALL TIMER_START('SHOWFLOPS_INLOOP [THE_MAIN_LOOP]', myThid)
1236 CALL SHOWFLOPS_INLOOP( iloop, myThid )
1237 CALL TIMER_STOP ('SHOWFLOPS_INLOOP [THE_MAIN_LOOP]', myThid)
1238 #endif
1239
1240 #ifdef ALLOW_DEBUG
1241 IF (debugMode) CALL DEBUG_LEAVE('FORWARD_STEP',myThid)
1242 #endif
1243
1244 RETURN
1245 END

  ViewVC Help
Powered by ViewVC 1.1.22