/[MITgcm]/MITgcm_contrib/shelfice_remeshing/CPL1/code/forward_step.F
ViewVC logotype

Contents of /MITgcm_contrib/shelfice_remeshing/CPL1/code/forward_step.F

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


Revision 1.2 - (show annotations) (download)
Mon Nov 2 17:04:34 2015 UTC (9 years, 8 months ago) by dgoldberg
Branch: MAIN
CVS Tags: HEAD
Changes since 1.1: +1 -5 lines
*** empty log message ***

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

  ViewVC Help
Powered by ViewVC 1.1.22