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

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

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


Revision 1.1 - (show annotations) (download)
Fri Apr 1 10:19:37 2016 UTC (8 years, 2 months ago) by dgoldberg
Branch: MAIN
CVS Tags: HEAD
Added rough code to dig ice shelf to make continuous ocean

1 C $Header: /u/gcmpack/MITgcm_contrib/verification_other/shelfice_remeshing/code/forward_step.F,v 1.2 2016/01/22 16:09:33 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 #define ALLOW_MOM_STEPPING
48 #if ( defined (ALLOW_AUTODIFF) && defined (ALLOW_OFFLINE) )
49 # undef ALLOW_MOM_STEPPING
50 #endif
51
52 CBOP
53 C !ROUTINE: FORWARD_STEP
54 C !INTERFACE:
55 SUBROUTINE FORWARD_STEP( iloop, myTime, myIter, myThid )
56
57 C !DESCRIPTION: \bv
58 C *=================================================================
59 C | SUBROUTINE forward_step
60 C | o Step forward in time the model variables for one time-step
61 C *=================================================================
62 C | The algorithm...
63 C |
64 C | "Calculation of Gs"
65 C | ===================
66 C | This is where all the accelerations and tendencies (ie.
67 C | physics, parameterizations etc...) are calculated
68 C | rho = rho ( theta[n], salt[n] )
69 C | b = b(rho, theta)
70 C | K31 = K31 ( rho )
71 C | Gu[n] = Gu( u[n], v[n], wVel, b, ... )
72 C | Gv[n] = Gv( u[n], v[n], wVel, b, ... )
73 C | Gt[n] = Gt( theta[n], u[n], v[n], wVel, K31, ... )
74 C | Gs[n] = Gs( salt[n], u[n], v[n], wVel, K31, ... )
75 C |
76 C | "Time-stepping" or "Prediction"
77 C | ================================
78 C | The models variables are stepped forward with the appropriate
79 C | time-stepping scheme (currently we use Adams-Bashforth II)
80 C | - For momentum, the result is always *only* a "prediction"
81 C | in that the flow may be divergent and will be "corrected"
82 C | later with a surface pressure gradient.
83 C | - Normally for tracers the result is the new field at time
84 C | level [n+1} *BUT* in the case of implicit diffusion the result
85 C | is also *only* a prediction.
86 C | - We denote "predictors" with an asterisk (*).
87 C | U* = U[n] + dt x ( 3/2 Gu[n] - 1/2 Gu[n-1] )
88 C | V* = V[n] + dt x ( 3/2 Gv[n] - 1/2 Gv[n-1] )
89 C | theta[n+1] = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
90 C | salt[n+1] = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
91 C | With implicit diffusion:
92 C | theta* = theta[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
93 C | salt* = salt[n] + dt x ( 3/2 Gt[n] - 1/2 atG[n-1] )
94 C | (1 + dt * K * d_zz) theta[n+1] = theta*
95 C | (1 + dt * K * d_zz) salt[n+1] = salt*
96 C |
97 C | "Correction Step"
98 C | =================
99 C | Here we update the horizontal velocities with the surface
100 C | pressure such that the resulting flow is either consistent
101 C | with the free-surface evolution or the rigid-lid:
102 C | U[n] = U* + dt x d/dx P
103 C | V[n] = V* + dt x d/dy P
104 C | W[n] = W* + dt x d/dz P (NH mode)
105 C *=================================================================
106 C \ev
107
108 C !CALLING SEQUENCE:
109 C FORWARD_STEP
110 C |
111 C |-- SHELFICE_REMESHING
112 C |
113 C |-- AUTODIFF_INADMODE_UNSET
114 C |
115 C |-- RESET_NLFS_VARS
116 C |-- UPDATE_R_STAR
117 C |-- UPDATE_SURF_DR
118 C |
119 C |-- PTRACERS_SWITCH_ONOFF
120 C |
121 C |-- DIAGNOSTICS_SWITCH_ONOFF
122 C |-- DO_STATEVARS_DIAGS
123 C |
124 C |-- NEST_CHILD_SETMEMO
125 C |-- NEST_PARENT_IO_1
126 C |
127 C |-- LOAD_FIELDS_DRIVER
128 C |
129 C |-- BULKF_FORCING
130 C |
131 C |-- CHEAPAML
132 C |
133 C |-- CTRL_MAP_FORCING
134 C |-- DUMMY_IN_STEPPING
135 C |
136 C |-- CPL_EXPORT_IMPORT_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
233 #ifdef HAVE_SIGREG
234 #include "SIGREG.h"
235 #endif
236
237 #ifdef ALLOW_SHAP_FILT
238 # include "SHAP_FILT.h"
239 #endif
240 #ifdef ALLOW_ZONAL_FILT
241 # include "ZONAL_FILT.h"
242 #endif
243
244 #ifdef ALLOW_LONGSTEP
245 # include "LONGSTEP_PARAMS.h"
246 # include "LONGSTEP.h"
247 #endif
248
249 #ifdef ALLOW_AUTODIFF
250 # include "AUTODIFF_MYFIELDS.h"
251 # include "FFIELDS.h"
252 # include "SURFACE.h"
253
254 # include "tamc.h"
255 # ifdef ALLOW_CTRL
256 # include "CTRL_SIZE.h"
257 # include "ctrl.h"
258 # include "ctrl_dummy.h"
259 # include "CTRL_GENARR.h"
260 # include "CTRL_OBCS.h"
261 # endif
262 # ifdef ALLOW_COST
263 # include "cost.h"
264 # endif
265 # ifdef ALLOW_ECCO
266 # include "ecco_cost.h"
267 # endif
268 # include "EOS.h"
269 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
270 # include "GRID.h"
271 # endif
272 # ifdef ALLOW_GMREDI
273 # include "GMREDI.h"
274 # endif
275 # ifdef ALLOW_EXF
276 # include "EXF_FIELDS.h"
277 # include "EXF_PARAM.h"
278 # ifdef ALLOW_BULKFORMULAE
279 # include "EXF_CONSTANTS.h"
280 # endif
281 # endif
282 # ifdef ALLOW_CD_CODE
283 # include "CD_CODE_VARS.h"
284 # endif
285 # ifdef ALLOW_GENERIC_ADVDIFF
286 # include "GAD.h"
287 # include "GAD_SOM_VARS.h"
288 # endif
289 # ifdef ALLOW_GGL90
290 # include "GGL90.h"
291 # endif
292 # ifdef ALLOW_PTRACERS
293 # include "PTRACERS_SIZE.h"
294 # include "PTRACERS_FIELDS.h"
295 # endif
296 # ifdef ALLOW_GCHEM
297 # include "GCHEM_FIELDS.h"
298 # endif
299 # ifdef ALLOW_CFC
300 # include "CFC.h"
301 # endif
302 # ifdef ALLOW_DIC
303 # include "DIC_VARS.h"
304 # include "DIC_LOAD.h"
305 # include "DIC_ATMOS.h"
306 # include "DIC_COST.h"
307 # endif
308 # ifdef ALLOW_OBCS
309 # include "OBCS_PARAMS.h"
310 # include "OBCS_FIELDS.h"
311 # include "OBCS_SEAICE.h"
312 # ifdef ALLOW_PTRACERS
313 # include "OBCS_PTRACERS.h"
314 # endif
315 # endif
316 # ifdef ALLOW_THSICE
317 # include "THSICE_PARAMS.h"
318 # include "THSICE_SIZE.h"
319 # include "THSICE_VARS.h"
320 # include "THSICE_COST.h"
321 # endif
322 # ifdef ALLOW_SEAICE
323 # include "SEAICE_SIZE.h"
324 # include "SEAICE.h"
325 # include "SEAICE_COST.h"
326 # endif
327 # ifdef ALLOW_SALT_PLUME
328 # include "SALT_PLUME.h"
329 # endif
330 # ifdef ALLOW_SHELFICE
331 # include "SHELFICE.h"
332 # include "SHELFICE_COST.h"
333 # endif
334 # ifdef ALLOW_STREAMICE
335 # include "STREAMICE.h"
336 # include "STREAMICE_ADV.h"
337 # include "STREAMICE_BDRY.h"
338 # include "STREAMICE_CG.h"
339 # endif
340 # ifdef ALLOW_EBM
341 # include "EBM.h"
342 # endif
343 # ifdef ALLOW_KPP
344 # include "KPP.h"
345 # endif
346 # ifdef ALLOW_RBCS
347 # include "RBCS_SIZE.h"
348 # include "RBCS_FIELDS.h"
349 # endif
350 # ifdef ALLOW_OFFLINE
351 # include "OFFLINE.h"
352 # endif
353 # ifdef ALLOW_CG2D_NSA
354 # include "CG2D.h"
355 # endif
356 #endif /* ALLOW_AUTODIFF */
357
358 C !INPUT/OUTPUT PARAMETERS:
359 C == Routine arguments ==
360 C note: under the multi-threaded model myIter and
361 C myTime are local variables passed around as routine
362 C arguments. Although this is fiddly it saves the need to
363 C impose additional synchronisation points when they are
364 C updated.
365 C myTime :: time counter for this thread
366 C myIter :: iteration counter for this thread
367 C myThid :: thread number for this instance of the routine.
368 INTEGER iloop
369 _RL myTime
370 INTEGER myIter
371 INTEGER myThid
372
373 C !LOCAL VARIABLES:
374 C == Local variables ==
375 C modelEnd :: true if reaching the end of the run
376 C myTimeBeg :: time at beginning of time step (needed by longstep)
377 C myIterBeg :: iteration number at beginning of time step
378 LOGICAL modelEnd
379 #ifdef ALLOW_LONGSTEP
380 INTEGER myIterBeg
381 _RL myTimeBeg
382 #endif /* ALLOW_LONGSTEP */
383 CEOP
384
385 #ifdef ALLOW_SHELFICE
386 C-- Remesh shelfice
387 CALL TIMER_START('SHELFICE_REMESHING [FORWARD_STEP]',myThid)
388 CALL SHELFICE_REMESHING (myTime, myIter, myThid )
389 CALL TIMER_STOP('SHELFICE_REMESHING [FORWARD_STEP]',myThid)
390 #endif /* ALLOW_SHELFICE */
391
392
393
394
395
396 #ifdef ALLOW_DEBUG
397 IF (debugMode) CALL DEBUG_ENTER('FORWARD_STEP',myThid)
398 #endif
399
400 #ifdef ALLOW_AUTODIFF
401 CALL AUTODIFF_INADMODE_UNSET( myThid )
402 #endif
403
404 #ifdef ALLOW_AUTODIFF
405 C-- Reset the model iteration counter and the model time.
406 myIter = nIter0 + (iloop-1)
407 myTime = startTime + deltaTClock*(iLoop-1)
408 #endif
409
410 #ifdef ALLOW_LONGSTEP
411 C store this for longstep_average with staggerTimeStep
412 C which is called after myIter and myTime are incremented
413 C but needs iter/time at beginning of time step
414 myIterBeg = myIter
415 myTimeBeg = myTime
416 #endif /* ALLOW_LONGSTEP */
417
418 #ifdef ALLOW_AUTODIFF_TAMC
419 c**************************************
420 #include "checkpoint_lev1_directives.h"
421 #include "checkpoint_lev1_template.h"
422 c**************************************
423 #endif
424
425 C-- Reset geometric factors (hFacC,W,S & recip_hFac) to their current values:
426 C added to simplify adjoint derivation - no effect in forward run
427 #ifdef NONLIN_FRSURF
428 #ifndef ALLOW_AUTODIFF
429 IF ( doResetHFactors ) THEN
430 #endif
431 CALL RESET_NLFS_VARS( myTime, myIter, myThid )
432 IF ( select_rStar.GT.0 ) THEN
433 # ifndef DISABLE_RSTAR_CODE
434 # ifdef ALLOW_AUTODIFF_TAMC
435 CADJ STORE rStarFacC, rStarFacS, rStarFacW =
436 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
437 # endif
438 CALL TIMER_START('UPDATE_R_STAR [FORWARD_STEP]',myThid)
439 CALL UPDATE_R_STAR( .FALSE., myTime, myIter, myThid )
440 CALL TIMER_STOP ('UPDATE_R_STAR [FORWARD_STEP]',myThid)
441 # endif /* DISABLE_RSTAR_CODE */
442 ELSE
443 #ifdef ALLOW_AUTODIFF_TAMC
444 CADJ STORE hFac_surfC, hFac_surfS, hFac_surfW
445 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
446 #endif
447 CALL TIMER_START('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
448 CALL UPDATE_SURF_DR( .FALSE., myTime, myIter, myThid )
449 CALL TIMER_STOP ('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
450 ENDIF
451
452 #ifdef ALLOW_AUTODIFF_TAMC
453 CADJ STORE hFacC, hFacS, hFacW =
454 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
455 CADJ STORE recip_hFacC, recip_hFacS, recip_hFacW =
456 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
457 #endif
458 #ifndef ALLOW_AUTODIFF
459 ENDIF
460 #endif
461 #endif /* NONLIN_FRSURF */
462
463 #ifdef ALLOW_PTRACERS
464 C-- Switch on/off individual tracer time-stepping
465 IF ( usePTRACERS ) THEN
466 CALL PTRACERS_SWITCH_ONOFF( myTime, myIter, myThid )
467 ENDIF
468 #endif /* ALLOW_PTRACERS */
469
470 C-- Switch on/off diagnostics for snap-shot output:
471 #ifdef ALLOW_DIAGNOSTICS
472 IF ( useDiagnostics ) THEN
473 CALL DIAGNOSTICS_SWITCH_ONOFF( myTime, myIter, myThid )
474 C-- State-variables diagnostics
475 CALL TIMER_START('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
476 CALL DO_STATEVARS_DIAGS( myTime, 0, myIter, myThid )
477 CALL TIMER_STOP ('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
478 ENDIF
479 #endif /* ALLOW_DIAGNOSTICS */
480
481 #ifdef ALLOW_NEST_CHILD
482 IF ( useNEST_CHILD) THEN
483 CALL NEST_CHILD_SETMEMO( myTime, myIter, myThid )
484 ENDIF
485 #endif /* ALLOW_NEST_CHILD */
486
487 #ifdef ALLOW_NEST_PARENT
488 IF ( useNEST_PARENT) THEN
489 CALL NEST_PARENT_IO_1( myTime, myIter, myThid )
490 ENDIF
491 #endif /* ALLOW_NEST_PARENT */
492
493 C-- Call driver to load external forcing fields from file
494 #ifdef ALLOW_DEBUG
495 IF (debugMode) CALL DEBUG_CALL('LOAD_FIELDS_DRIVER',myThid)
496 #endif
497 #ifdef ALLOW_AUTODIFF_TAMC
498 cph Important STORE that avoids hidden recomp. of load_fields_driver
499 CADJ STORE theta = comlev1, key = ikey_dynamics,
500 CADJ & kind = isbyte
501 CADJ STORE uVel, vVel = comlev1, key = ikey_dynamics,
502 CADJ & kind = isbyte
503 #endif
504 CALL TIMER_START('LOAD_FIELDS_DRIVER [FORWARD_STEP]',myThid)
505 CALL LOAD_FIELDS_DRIVER( myTime, myIter, myThid )
506 CALL TIMER_STOP ('LOAD_FIELDS_DRIVER [FORWARD_STEP]',myThid)
507
508 C-- Call Bulk-Formulae forcing package
509 #ifdef ALLOW_BULK_FORCE
510 IF ( useBulkForce ) THEN
511 #ifdef ALLOW_DEBUG
512 IF (debugMode) CALL DEBUG_CALL('BULKF_FORCING',myThid)
513 #endif
514 CALL TIMER_START('BULKF_FORCING [FORWARD_STEP]',myThid)
515 C- calculate qnet and empmr (and wind stress)
516 CALL BULKF_FORCING( myTime, myIter, myThid )
517 CALL TIMER_STOP ('BULKF_FORCING [FORWARD_STEP]',myThid)
518 ENDIF
519 #endif /* ALLOW_BULK_FORCE */
520
521 C-- Call external chepaml forcing package
522 #ifdef ALLOW_CHEAPAML
523 IF ( useCheapAML ) THEN
524 #ifdef ALLOW_DEBUG
525 IF (debugMode) CALL DEBUG_CALL('CHEAPAML',myThid)
526 #endif
527 CALL TIMER_START('CHEAPAML [FORWARD_STEP]',myThid)
528 C- calculate qnet (and wind stress)
529 CALL CHEAPAML( myTime, myIter,myThid )
530 CALL TIMER_STOP ('CHEAPAML [FORWARD_STEP]',myThid)
531 ENDIF
532 #endif /*ALLOW_CHEAPAML */
533
534 #ifdef ALLOW_CTRL
535 C-- Add control vector for forcing and parameter fields
536 IF ( ( useCTRL ).AND.( myIter .EQ. nIter0 ) )
537 & CALL CTRL_MAP_FORCING (myThid)
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
872
873 CADJ STORE hFac_surfC, hFac_surfS, hFac_surfW
874 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
875 # endif
876 CALL TIMER_START('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
877 CALL UPDATE_SURF_DR( .TRUE., myTime, myIter, myThid )
878 CALL TIMER_STOP ('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
879 ENDIF
880 # ifdef ALLOW_AUTODIFF_TAMC
881
882
883 CADJ STORE hFacC, hFacS, hFacW =
884 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
885 CADJ STORE recip_hFacC, recip_hFacS, recip_hFacW =
886 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
887 # ifdef ALLOW_CG2D_NSA
888 CADJ STORE aW2d, aS2d, aC2d =
889 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
890 CADJ STORE pC, pS, pW =
891 CADJ & comlev1, key = ikey_dynamics, kind = isbyte
892 # endif
893 # endif
894
895 C- update also CG2D matrix (and preconditioner)
896 IF ( momStepping .AND. nonlinFreeSurf.GT.2 ) THEN
897 CALL TIMER_START('UPDATE_CG2D [FORWARD_STEP]',myThid)
898 CALL UPDATE_CG2D( myTime, myIter, myThid )
899 CALL TIMER_STOP ('UPDATE_CG2D [FORWARD_STEP]',myThid)
900 ENDIF
901 #endif /* NONLIN_FRSURF */
902
903 C-- Apply Filters to u*,v* before SOLVE_FOR_PRESSURE
904 #ifdef ALLOW_SHAP_FILT
905 IF (useSHAP_FILT .AND. shap_filt_uvStar) THEN
906 CALL TIMER_START('SHAP_FILT_UV [FORWARD_STEP]',myThid)
907 IF (implicDiv2Dflow.LT.1.) THEN
908 C-- Explicit+Implicit part of the Barotropic Flow Divergence
909 C => Filtering of uVel,vVel is necessary
910 CALL SHAP_FILT_APPLY_UV( uVel,vVel,
911 & myTime, myIter, myThid )
912 ENDIF
913 CALL SHAP_FILT_APPLY_UV( gU,gV,myTime,myIter,myThid)
914 CALL TIMER_STOP ('SHAP_FILT_UV [FORWARD_STEP]',myThid)
915 ENDIF
916 #endif
917 #ifdef ALLOW_ZONAL_FILT
918 IF (useZONAL_FILT .AND. zonal_filt_uvStar) THEN
919 CALL TIMER_START('ZONAL_FILT_UV [FORWARD_STEP]',myThid)
920 IF (implicDiv2Dflow.LT.1.) THEN
921 C-- Explicit+Implicit part of the Barotropic Flow Divergence
922 C => Filtering of uVel,vVel is necessary
923 CALL ZONAL_FILT_APPLY_UV( uVel, vVel, myThid )
924 ENDIF
925 CALL ZONAL_FILT_APPLY_UV( gU, gV, myThid )
926 CALL TIMER_STOP ('ZONAL_FILT_UV [FORWARD_STEP]',myThid)
927 ENDIF
928 #endif
929
930 C-- Solve elliptic equation(s).
931 C Two-dimensional only for conventional hydrostatic or
932 C three-dimensional for non-hydrostatic and/or IGW scheme.
933 IF ( momStepping ) THEN
934 #ifdef ALLOW_AUTODIFF_TAMC
935 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
936 CADJ STORE uVel, vVel
937 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
938 CADJ STORE EmPmR,hFacS,hFacW
939 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
940 # endif
941 #endif
942 CALL TIMER_START('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
943 CALL SOLVE_FOR_PRESSURE(myTime, myIter, myThid)
944 CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
945 ENDIF
946
947 C-- Correct divergence in flow field and cycle time-stepping momentum
948 #ifdef ALLOW_MOM_STEPPING
949 #ifndef ALLOW_AUTODIFF
950 IF ( momStepping ) THEN
951 #endif
952 #ifdef ALLOW_AUTODIFF_TAMC
953 # ifdef ALLOW_DEPTH_CONTROL
954 CADJ STORE etaN, uVel,vVel
955 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
956 # endif /* ALLOW_DEPTH_CONTROL */
957 #endif /* ALLOW_AUTODIFF_TAMC */
958 CALL TIMER_START('MOM_CORRECTION_STEP [FORWARD_STEP]',myThid)
959 CALL MOMENTUM_CORRECTION_STEP(myTime, myIter, myThid)
960 CALL TIMER_STOP ('MOM_CORRECTION_STEP [FORWARD_STEP]',myThid)
961 #ifndef ALLOW_AUTODIFF
962 ENDIF
963 #endif
964 #endif /* ALLOW_MOM_STEPPING */
965 #ifdef ALLOW_AUTODIFF_TAMC
966 CADJ STORE uVel, vVel = comlev1, key = ikey_dynamics, kind = isbyte
967 #endif
968
969 IF ( calc_wVelocity ) THEN
970 C-- Integrate continuity vertically for vertical velocity
971 C (+ update "etaN" & "etaH", exact volume conservation):
972 CALL TIMER_START('INTEGR_CONTINUITY [FORWARD_STEP]',myThid)
973 CALL INTEGR_CONTINUITY( uVel, vVel, myTime, myIter, myThid)
974 CALL TIMER_STOP ('INTEGR_CONTINUITY [FORWARD_STEP]',myThid)
975 ENDIF
976
977 #ifdef NONLIN_FRSURF
978 IF ( select_rStar.NE.0 ) THEN
979 # ifndef DISABLE_RSTAR_CODE
980 # ifdef ALLOW_AUTODIFF_TAMC
981 CADJ STORE etaH
982 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
983 CADJ STORE rStarFacC,rStarFacS,rStarFacW
984 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
985 # endif
986 C-- r* : compute the future level thickness according to etaH(n+1)
987 CALL TIMER_START('CALC_R_STAR [FORWARD_STEP]',myThid)
988 CALL CALC_R_STAR(etaH, myTime, myIter, myThid )
989 CALL TIMER_STOP ('CALC_R_STAR [FORWARD_STEP]',myThid)
990 # endif /* DISABLE_RSTAR_CODE */
991 ELSEIF ( nonlinFreeSurf.GT.0 .AND. selectSigmaCoord.EQ.0 ) THEN
992 C-- compute the future surface level thickness according to etaH(n+1)
993 # ifdef ALLOW_AUTODIFF_TAMC
994 CADJ STORE etaH = comlev1, key = ikey_dynamics,
995 CADJ & kind = isbyte
996 # endif
997 CALL TIMER_START('CALC_SURF_DR [FORWARD_STEP]',myThid)
998 CALL CALC_SURF_DR(etaH, myTime, myIter, myThid )
999 CALL TIMER_STOP ('CALC_SURF_DR [FORWARD_STEP]',myThid)
1000 ENDIF
1001 # ifdef ALLOW_AUTODIFF_TAMC
1002 CADJ STORE rStarExpC
1003 CADJ & = comlev1, key = ikey_dynamics, kind = isbyte
1004 CADJ STORE hFac_surfC = comlev1, key = ikey_dynamics,
1005 CADJ & kind = isbyte
1006 CADJ STORE salt,theta = comlev1, key = ikey_dynamics,
1007 CADJ & kind = isbyte
1008 # endif
1009 #endif /* NONLIN_FRSURF */
1010
1011 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
1012 IF ( staggerTimeStep ) THEN
1013 C-- do exchanges of U,V (needed for multiDim) when using stagger time-step :
1014 #ifdef ALLOW_DEBUG
1015 IF (debugMode)
1016 & CALL DEBUG_CALL('DO_STAGGER_FIELDS_EXCH.',myThid)
1017 #endif
1018 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
1019 CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid )
1020 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
1021
1022 #ifdef ALLOW_DIAGNOSTICS
1023 C-- State-variables diagnostics
1024 IF ( useDiagnostics ) THEN
1025 CALL TIMER_START('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
1026 CALL DO_STATEVARS_DIAGS( myTime, 1, myIter, myThid )
1027 CALL TIMER_STOP ('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
1028 ENDIF
1029 #endif
1030
1031
1032
1033
1034
1035
1036 #ifdef ALLOW_AUTODIFF_TAMC
1037 CADJ STORE wVel = comlev1, key = ikey_dynamics, kind = isbyte
1038 #endif
1039 #ifdef ALLOW_DEBUG
1040 IF (debugMode) CALL DEBUG_CALL('THERMODYNAMICS',myThid)
1041 #endif
1042 CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',myThid)
1043 CALL THERMODYNAMICS( myTime, myIter, myThid )
1044 CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',myThid)
1045
1046 C-- if staggerTimeStep: end
1047 ENDIF
1048 C---+--------+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
1049
1050 #ifdef ALLOW_AUTODIFF_TAMC
1051 cph This is needed because convective_adjustment calls
1052 cph find_rho which may use pressure()
1053 CADJ STORE totPhiHyd = comlev1, key = ikey_dynamics,
1054 CADJ & kind = isbyte
1055 #endif
1056 C-- Apply adjustments to Tracers arrays (T,S,+pTracers)
1057 CALL TIMER_START('TRC_CORRECTION_STEP [FORWARD_STEP]',myThid)
1058 CALL TRACERS_CORRECTION_STEP(myTime, myIter, myThid)
1059 CALL TIMER_STOP ('TRC_CORRECTION_STEP [FORWARD_STEP]',myThid)
1060
1061 #ifdef ALLOW_LONGSTEP
1062 IF ( usePTRACERS ) THEN
1063 IF ( LS_whenToSample .EQ. 2 ) THEN
1064 C Average everything at the end of the timestep. This will
1065 C reproduce online results with staggerTimeStep=.TRUE.
1066 C when LS_nIter=1
1067 #ifdef ALLOW_DEBUG
1068 IF (debugMode) CALL DEBUG_CALL('LONGSTEP_AVERAGE',myThid)
1069 #endif
1070 CALL TIMER_START('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid)
1071 C myIter has been update after dynamics, but the averaging window
1072 C should be determined by myIter at beginning of timestep
1073 CALL LONGSTEP_AVERAGE( myTimeBeg, myIterBeg, myThid )
1074 CALL TIMER_STOP ('LONGSTEP_AVERAGE [FORWARD_STEP]',myThid)
1075
1076 #ifdef ALLOW_DEBUG
1077 IF (debugMode)
1078 & CALL DEBUG_CALL('LONGSTEP_THERMODYNAMICS',myThid)
1079 #endif
1080 CALL TIMER_START('LONGSTEP_THERMODYNAMICS [FORWARD_STEP]',
1081 & myThid)
1082 CALL LONGSTEP_THERMODYNAMICS( myTime, myIter, myThid )
1083 CALL TIMER_STOP ('LONGSTEP_THERMODYNAMICS [FORWARD_STEP]',
1084 & myThid)
1085 C-- if LS_whenToSample.EQ.2: end
1086 ENDIF
1087
1088 C-- Apply adjustments to passive Tracers arrays (pTracers)
1089 c CALL TIMER_START('LS_CORRECTION_STEP [FORWARD_STEP]',myThid)
1090 c CALL LONGSTEP_CORRECTION_STEP(myTime, myIter, myThid)
1091 c CALL TIMER_STOP ('LS_CORRECTION_STEP [FORWARD_STEP]',myThid)
1092 C-- if usePTRACERS: end
1093 ENDIF
1094 #endif /* ALLOW_LONGSTEP */
1095
1096 #ifdef ALLOW_GCHEM
1097 C Add separate timestepping of chemical/biological/forcing
1098 C of ptracers here in GCHEM_FORCING_SEP
1099 #ifdef ALLOW_LONGSTEP
1100 IF ( useGCHEM .AND. LS_doTimeStep ) THEN
1101 #else
1102 IF ( useGCHEM ) THEN
1103 #endif
1104 #ifdef ALLOW_AUTODIFF_TAMC
1105 CADJ STORE pTracer = comlev1, key = ikey_dynamics,
1106 CADJ & kind = isbyte
1107 CADJ STORE theta, salt = comlev1, key = ikey_dynamics,
1108 CADJ & kind = isbyte
1109 #endif
1110 #ifdef ALLOW_DEBUG
1111 IF (debugMode) CALL DEBUG_CALL('GCHEM_FORCING_SEP',myThid)
1112 #endif /* ALLOW_DEBUG */
1113 CALL TIMER_START('GCHEM_FORCING_SEP [FORWARD_STEP]',myThid)
1114 CALL GCHEM_FORCING_SEP( myTime,myIter,myThid )
1115 CALL TIMER_STOP ('GCHEM_FORCING_SEP [FORWARD_STEP]',myThid)
1116 ENDIF
1117 #endif /* ALLOW_GCHEM */
1118
1119 C-- Do "blocking" sends and receives for tendency "overlap" terms
1120 c CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
1121 c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )
1122 c CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
1123
1124 C-- Do "blocking" sends and receives for field "overlap" terms
1125 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
1126 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
1127 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
1128
1129 #ifdef ALLOW_DIAGNOSTICS
1130 IF ( useDiagnostics ) THEN
1131 CALL TIMER_START('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
1132 CALL DO_STATEVARS_DIAGS( myTime, 2, myIter, myThid )
1133 CALL TIMER_STOP ('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
1134 ENDIF
1135 #endif
1136
1137 #ifdef ALLOW_GRIDALT
1138 IF (useGRIDALT) THEN
1139 CALL GRIDALT_UPDATE(myThid)
1140 ENDIF
1141 #endif
1142
1143 #ifdef ALLOW_FIZHI
1144 IF (useFIZHI) THEN
1145 CALL TIMER_START('FIZHI [FORWARD_STEP]',myThid)
1146 CALL STEP_FIZHI_CORR ( myTime, myIter, myThid, dTtracerLev(1) )
1147 CALL TIMER_STOP ('FIZHI [FORWARD_STEP]',myThid)
1148 ENDIF
1149 #endif
1150
1151 #ifdef ALLOW_FLT
1152 C-- Calculate float trajectories
1153 IF (useFLT) THEN
1154 CALL TIMER_START('FLOATS [FORWARD_STEP]',myThid)
1155 CALL FLT_MAIN( myTime, myIter, myThid )
1156 CALL TIMER_STOP ('FLOATS [FORWARD_STEP]',myThid)
1157 ENDIF
1158 #endif
1159
1160 #ifdef ALLOW_TIMEAVE
1161 C-- State-variables time-averaging
1162 CALL TIMER_START('DO_STATEVARS_TAVE [FORWARD_STEP]',myThid)
1163 CALL DO_STATEVARS_TAVE( myTime, myIter, myThid )
1164 CALL TIMER_STOP ('DO_STATEVARS_TAVE [FORWARD_STEP]',myThid)
1165 #endif
1166
1167 #ifdef ALLOW_NEST_PARENT
1168 IF ( useNEST_PARENT) THEN
1169 CALL NEST_PARENT_IO_2( myTime, myIter, myThid )
1170 ENDIF
1171 #endif /* ALLOW_NEST_PARENT */
1172
1173 #ifdef ALLOW_NEST_CHILD
1174 IF ( useNEST_CHILD) THEN
1175 CALL NEST_CHILD_TRANSP( myTime, myIter, myThid )
1176 ENDIF
1177 #endif /* ALLOW_NEST_CHILD */
1178
1179 #ifdef ALLOW_MONITOR
1180 IF ( monitorFreq.GT.0. .OR. adjMonitorFreq.GT.0. ) THEN
1181 C-- Check status of solution (statistics, cfl, etc...)
1182 CALL TIMER_START('MONITOR [FORWARD_STEP]',myThid)
1183 CALL MONITOR( myTime, myIter, myThid )
1184 CALL TIMER_STOP ('MONITOR [FORWARD_STEP]',myThid)
1185 ENDIF
1186 #endif /* ALLOW_MONITOR */
1187
1188 #ifdef ALLOW_COST
1189 C-- compare model with data and compute cost function
1190 C-- this is done after exchanges to allow interpolation
1191 CALL TIMER_START('COST_TILE [FORWARD_STEP]',myThid)
1192 CALL COST_TILE ( myTime, myIter, myThid )
1193 CALL TIMER_STOP ('COST_TILE [FORWARD_STEP]',myThid)
1194 #endif
1195
1196 C-- Check if it has reached the end of simulation
1197 modelEnd = myTime.EQ.endTime .OR. myIter.EQ.nEndIter
1198 #ifdef HAVE_SIGREG
1199 IF ( useSIGREG ) THEN
1200 modelEnd = modelEnd .OR. ( i_got_signal.GT.0 )
1201 ENDIF
1202 #endif /* HAVE_SIGREG */
1203
1204 C-- Do IO if needed.
1205 CALL TIMER_START('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
1206 CALL DO_THE_MODEL_IO( modelEnd, myTime, myIter, myThid )
1207 CALL TIMER_STOP ('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
1208
1209 #ifdef ALLOW_PTRACERS
1210 C Reset the ptracers (but after the io is done)
1211 IF ( usePTRACERS ) THEN
1212 CALL TIMER_START('PTRACERS_RESET [FORWARD_STEP]',myThid)
1213 CALL PTRACERS_RESET( myTime, myIter, myThid )
1214 CALL TIMER_STOP ('PTRACERS_RESET [FORWARD_STEP]',myThid)
1215 ENDIF
1216 #endif /* ALLOW_PTRACERS */
1217
1218 C-- Save state for restarts
1219 CALL TIMER_START('DO_WRITE_PICKUP [FORWARD_STEP]',myThid)
1220 CALL DO_WRITE_PICKUP( modelEnd, myTime, myIter, myThid )
1221 CALL TIMER_STOP ('DO_WRITE_PICKUP [FORWARD_STEP]',myThid)
1222
1223 #ifdef HAVE_SIGREG
1224 IF ( useSIGREG ) THEN
1225 IF ( modelEnd .AND. i_got_signal.GT.0 ) THEN
1226 STOP 'Checkpoint completed -- killed by signal handler'
1227 ENDIF
1228 ENDIF
1229 #endif /* HAVE_SIGREG */
1230
1231 #ifdef ALLOW_AUTODIFF
1232 CALL AUTODIFF_INADMODE_SET( myThid )
1233 #endif
1234
1235 #ifdef ALLOW_SHOWFLOPS
1236 CALL TIMER_START('SHOWFLOPS_INLOOP [THE_MAIN_LOOP]', myThid)
1237 CALL SHOWFLOPS_INLOOP( iloop, myThid )
1238 CALL TIMER_STOP ('SHOWFLOPS_INLOOP [THE_MAIN_LOOP]', myThid)
1239 #endif
1240
1241 #ifdef ALLOW_DEBUG
1242 IF (debugMode) CALL DEBUG_LEAVE('FORWARD_STEP',myThid)
1243 #endif
1244
1245
1246
1247
1248
1249 RETURN
1250 END

  ViewVC Help
Powered by ViewVC 1.1.22