/[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.150 - (show annotations) (download)
Mon Jun 4 21:37:26 2007 UTC (16 years, 11 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint59d
Changes since 1.149: +11 -1 lines
Further cleanup of top-level routines.

1 C $Header: /u/gcmpack/MITgcm/model/src/forward_step.F,v 1.149 2007/04/18 19:54:06 heimbach Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 #ifdef ALLOW_GMREDI
8 # include "GMREDI_OPTIONS.h"
9 #endif
10 #ifdef ALLOW_OBCS
11 # include "OBCS_OPTIONS.h"
12 #endif
13
14 CBOP
15 C !ROUTINE: FORWARD_STEP
16 C !INTERFACE:
17 SUBROUTINE FORWARD_STEP( iloop, myTime, myIter, myThid )
18
19 C !DESCRIPTION: \bv
20 C *==================================================================
21 C | SUBROUTINE forward_step
22 C | o Run the ocean model and, optionally, evaluate a cost function.
23 C *==================================================================
24 C |
25 C | THE_MAIN_LOOP is the toplevel routine for the Tangent Linear and
26 C | Adjoint Model Compiler (TAMC). For this purpose the initialization
27 C | of the model was split into two parts. Those parameters that do
28 C | not depend on a specific model run are set in INITIALISE_FIXED,
29 C | whereas those that do depend on the specific realization are
30 C | initialized in INITIALISE_VARIA.
31 C |
32 C *==================================================================
33 C \ev
34
35 C !USES:
36 IMPLICIT NONE
37 C == Global variables ==
38 #include "SIZE.h"
39 #include "EEPARAMS.h"
40 #include "PARAMS.h"
41 #include "DYNVARS.h"
42
43 #ifdef ALLOW_MNC
44 #include "MNC_PARAMS.h"
45 EXTERNAL DIFFERENT_MULTIPLE
46 LOGICAL DIFFERENT_MULTIPLE
47 #endif
48
49 #ifdef HAVE_SIGREG
50 #include "SIGREG.h"
51 #endif
52
53 #ifdef ALLOW_SHAP_FILT
54 # include "SHAP_FILT.h"
55 #endif
56 #ifdef ALLOW_ZONAL_FILT
57 # include "ZONAL_FILT.h"
58 #endif
59 #ifdef COMPONENT_MODULE
60 # include "CPL_PARAMS.h"
61 #endif
62
63 #ifdef ALLOW_AUTODIFF_TAMC
64 # include "FFIELDS.h"
65
66 # include "tamc.h"
67 # include "ctrl.h"
68 # include "ctrl_dummy.h"
69 # include "cost.h"
70 # include "EOS.h"
71 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
72 # include "GRID.h"
73 # endif
74 # ifdef ALLOW_EXF
75 # include "EXF_FIELDS.h"
76 # ifdef ALLOW_BULKFORMULAE
77 # include "EXF_CONSTANTS.h"
78 # endif
79 # endif
80 # ifdef ALLOW_PTRACERS
81 # include "PTRACERS_SIZE.h"
82 # include "PTRACERS.h"
83 # endif
84 # ifdef ALLOW_OBCS
85 # include "OBCS.h"
86 # ifdef ALLOW_PTRACERS
87 # include "OBCS_PTRACERS.h"
88 # endif
89 # endif
90 # ifdef ALLOW_CD_CODE
91 # include "CD_CODE_VARS.h"
92 # endif
93 # ifdef ALLOW_THSICE
94 # include "THSICE_VARS.h"
95 # endif
96 # ifdef ALLOW_EBM
97 # include "EBM.h"
98 # endif
99 # ifdef EXACT_CONSERV
100 # include "SURFACE.h"
101 # endif
102 # ifdef ALLOW_KPP
103 # include "KPP.h"
104 # endif
105 # ifdef ALLOW_GMREDI
106 # include "GMREDI.h"
107 # endif
108 # ifdef ALLOW_RBCS
109 # include "RBCS.h"
110 # endif
111 #endif /* ALLOW_AUTODIFF_TAMC */
112
113 C !LOCAL VARIABLES:
114 C == Routine arguments ==
115 C note: under the multi-threaded model myIter and
116 C myTime are local variables passed around as routine
117 C arguments. Although this is fiddly it saves the need to
118 C impose additional synchronisation points when they are
119 C updated.
120 C myIter - iteration counter for this thread
121 C myTime - time counter for this thread
122 C myThid - thread number for this instance of the routine.
123 INTEGER iloop
124 INTEGER myThid
125 INTEGER myIter
126 _RL myTime
127
128 C == Local variables ==
129 #ifdef COMPONENT_MODULE
130 INTEGER myItP1
131 #endif
132 CEOP
133
134 #ifdef ALLOW_DEBUG
135 IF ( debugLevel .GE. debLevB )
136 & CALL DEBUG_ENTER('FORWARD_STEP',myThid)
137 #endif
138
139 #ifdef ALLOW_AUTODIFF_TAMC
140 C-- Reset the model iteration counter and the model time.
141 myIter = nIter0 + (iloop-1)
142 myTime = startTime + float(iloop-1)*deltaTclock
143 #endif
144
145 #ifdef ALLOW_AUTODIFF_TAMC
146 c**************************************
147 #include "checkpoint_lev1_directives.h"
148 #include "checkpoint_lev1_template.h"
149 c**************************************
150 #endif
151
152 C-- Switch on/off diagnostics for snap-shot output:
153 #ifdef ALLOW_DIAGNOSTICS
154 IF ( useDiagnostics ) THEN
155 CALL DIAGNOSTICS_SWITCH_ONOFF( myTime, myIter, myThid )
156 C-- State-variables diagnostics
157 CALL TIMER_START('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
158 CALL DO_STATEVARS_DIAGS( myTime, 0, myIter, myThid )
159 CALL TIMER_STOP ('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
160 ENDIF
161 #endif
162
163 #ifdef ALLOW_PROFILES
164 #ifdef ALLOW_DEBUG
165 IF (debugMode) CALL DEBUG_CALL('',myThid)
166 #endif
167 c-- Accumulate in-situ time averages of theta, salt, and SSH.
168 CALL TIMER_START('PROFILES_INLOOP [THE_MAIN_LOOP]', mythid)
169 CALL PROFILES_INLOOP( mytime, mythid )
170 CALL TIMER_STOP ('PROFILES_INLOOP [THE_MAIN_LOOP]', mythid)
171 #endif
172
173 C-- Call driver to load external forcing fields from file
174 #ifdef ALLOW_DEBUG
175 IF ( debugLevel .GE. debLevB )
176 & CALL DEBUG_CALL('LOAD_FIELDS_DRIVER',myThid)
177 #endif
178 CALL TIMER_START('LOAD_FIELDS_DRIVER [FORWARD_STEP]',myThid)
179 CALL LOAD_FIELDS_DRIVER( myTime, myIter, myThid )
180 CALL TIMER_STOP ('LOAD_FIELDS_DRIVER [FORWARD_STEP]',myThid)
181
182 C-- Call Bulk-Formulae forcing package
183 #ifdef ALLOW_BULK_FORCE
184 IF ( useBulkForce ) THEN
185 #ifdef ALLOW_DEBUG
186 IF ( debugLevel .GE. debLevB )
187 & CALL DEBUG_CALL('BULKF_FORCING',myThid)
188 #endif
189 CALL TIMER_START('BULKF_FORCING [FORWARD_STEP]',myThid)
190 C- calculate qnet and empmr (and wind stress)
191 CALL BULKF_FORCING( myTime, myIter, myThid )
192 CALL TIMER_STOP ('BULKF_FORCING [FORWARD_STEP]',myThid)
193 ENDIF
194 #endif /* ALLOW_BULK_FORCE */
195
196 #ifdef ALLOW_AUTODIFF
197 c-- Add control vector for forcing and parameter fields
198 IF ( myIter .EQ. nIter0 )
199 & CALL CTRL_MAP_FORCING (myThid)
200 #endif
201
202 #if (defined (ALLOW_AUTODIFF_TAMC) && defined (ALLOW_AUTODIFF_MONITOR))
203 CALL DUMMY_IN_STEPPING( myTime, myIter, myThid )
204 #endif
205
206 #ifdef COMPONENT_MODULE
207 IF ( useCoupler .AND. cpl_earlyExpImpCall ) THEN
208 C Post coupling data that I export.
209 C Read in coupling data that I import.
210 CALL TIMER_START('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid)
211 CALL CPL_EXPORT_MY_DATA( myIter, myTime, myThid )
212 CALL CPL_IMPORT_EXTERNAL_DATA( myIter, myTime, myThid )
213 CALL TIMER_STOP ('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid)
214 ENDIF
215 #endif /* COMPONENT_MODULE */
216
217 #ifdef ALLOW_EBM
218 IF ( useEBM ) THEN
219 # ifdef ALLOW_DEBUG
220 IF ( debugLevel .GE. debLevB )
221 & CALL DEBUG_CALL('EBM',myThid)
222 # endif
223 CALL TIMER_START('EBM [FORWARD_STEP]',myThid)
224 CALL EBM_DRIVER ( myTime, myIter, myThid )
225 CALL TIMER_STOP ('EBM [FORWARD_STEP]',myThid)
226 ENDIF
227 #endif /* ALLOW_EBM */
228
229 C-- Step forward fields and calculate time tendency terms.
230
231 #ifdef ALLOW_DEBUG
232 IF ( debugLevel .GE. debLevB )
233 & CALL DEBUG_CALL('DO_ATMOSPHERIC_PHYS',myThid)
234 #endif
235 CALL TIMER_START('DO_ATMOSPHERIC_PHYS [FORWARD_STEP]',myThid)
236 CALL DO_ATMOSPHERIC_PHYS( myTime, myIter, myThid )
237 CALL TIMER_STOP ('DO_ATMOSPHERIC_PHYS [FORWARD_STEP]',myThid)
238
239 #ifdef ALLOW_AUTODIFF_TAMC
240 CADJ STORE surfaceforcingtice = comlev1, key = ikey_dynamics
241 # ifdef ALLOW_KPP
242 CADJ STORE uvel = comlev1, key = ikey_dynamics
243 CADJ STORE vvel = comlev1, key = ikey_dynamics
244 # endif
245 # ifdef EXACT_CONSERV
246 cphCADJ STORE empmr = comlev1, key = ikey_dynamics
247 cphCADJ STORE pmepr = comlev1, key = ikey_dynamics
248 # endif
249 # ifdef ALLOW_PTRACERS
250 CADJ STORE ptracer = comlev1, key = ikey_dynamics
251 # endif
252 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
253 cph-test
254 CADJ STORE hFacC = comlev1, key = ikey_dynamics
255 # ifndef DISABLE_RSTAR_CODE
256 CADJ STORE rstarexpc = comlev1, key = ikey_dynamics
257 # endif
258 # endif
259 #endif /* ALLOW_AUTODIFF_TAMC */
260
261 #ifndef ALLOW_AUTODIFF_TAMC
262 IF ( .NOT. useOffLine ) THEN
263 #endif
264 #ifdef ALLOW_DEBUG
265 IF ( debugLevel .GE. debLevB )
266 & CALL DEBUG_CALL('DO_OCEANIC_PHYS',myThid)
267 #endif
268 CALL TIMER_START('DO_OCEANIC_PHYS [FORWARD_STEP]',myThid)
269 CALL DO_OCEANIC_PHYS( myTime, myIter, myThid )
270 CALL TIMER_STOP ('DO_OCEANIC_PHYS [FORWARD_STEP]',myThid)
271 #ifdef ALLOW_AUTODIFF_TAMC
272 CADJ STORE EmPmR = comlev1, key = ikey_dynamics
273 # ifdef EXACT_CONSERV
274 CADJ STORE pmepr = comlev1, key = ikey_dynamics
275 # endif
276 #else
277 ENDIF
278 #endif
279
280 #ifdef ALLOW_AUTODIFF_TAMC
281 # ifdef NONLIN_FRSURF
282 cph-test
283 CADJ STORE hFac_surfC = comlev1, key = ikey_dynamics
284 CADJ STORE hfac_surfs = comlev1, key = ikey_dynamics
285 CADJ STORE hfac_surfw = comlev1, key = ikey_dynamics
286 # endif
287 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
288 CADJ STORE hFacC, hFacS, hFacW
289 CADJ & = comlev1, key = ikey_dynamics
290 CADJ STORE recip_hFacC, recip_hFacS, recip_hFacW
291 CADJ & = comlev1, key = ikey_dynamics
292 c
293 CADJ STORE surfaceforcingu = comlev1, key = ikey_dynamics
294 CADJ STORE surfaceforcingv = comlev1, key = ikey_dynamics
295 # endif
296 #endif /* ALLOW_AUTODIFF_TAMC */
297
298 #ifdef ALLOW_GCHEM
299 C GCHEM package is an interface for any bio-geochemical or
300 C ecosystem model you would like to include.
301 C If GCHEM_SEPARATE_FORCING is not defined, you are
302 C responsible for computing tendency terms for passive
303 C tracers and storing them on a 3DxNumPtracers-array called
304 C gchemTendency in GCHEM_CALC_TENDENCY. This tendency is then added
305 C to gPtr in ptracers_forcing later-on.
306 C If GCHEM_SEPARATE_FORCING is defined, you are reponsible for
307 C UPDATING ptracers directly in GCHEM_FORCING_SEP. This amounts
308 C to a completely separate time step that you have to implement
309 C yourself (Eulerian seems to be fine in most cases).
310 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
311 C CAVEAT: Up to now, when GCHEM is turned on the field ptracerForcingSurf,
312 C which is needed for KPP is not set properly. ptracerForcingSurf must
313 C be treated differently depending on whether GCHEM_SEPARATE_FORCING
314 C is define or not. TBD.
315 CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
316 IF ( useGCHEM ) THEN
317 #ifdef ALLOW_DEBUG
318 IF ( debugLevel .GE. debLevB )
319 & CALL DEBUG_CALL('GCHEM_CALC_TENDENCY',myThid)
320 #endif
321 CALL TIMER_START('GCHEM_CALC_TENDENCY [FORWARD_STEP]',myThid)
322 CALL GCHEM_CALC_TENDENCY( myTime, myIter, myThid )
323 CALL TIMER_STOP ('GCHEM_CALC_TENDENCY [FORWARD_STEP]',myThid)
324 ENDIF
325 #endif /* ALLOW_GCHEM */
326
327 #ifdef ALLOW_AUTODIFF_TAMC
328 cph needed to be moved here from do_oceanic_physics
329 cph to be visible down the road
330 c
331 CADJ STORE surfaceForcingS = comlev1, key = ikey_dynamics
332 CADJ STORE surfaceForcingT = comlev1, key = ikey_dynamics
333 CADJ STORE surfaceForcingTice = comlev1, key = ikey_dynamics
334 ctest(
335 CADJ STORE IVDConvCount = comlev1, key = ikey_dynamics
336 ctest)
337 # ifdef ALLOW_PTRACERS
338 CADJ STORE surfaceForcingPtr = comlev1, key = ikey_dynamics
339 # endif
340 c
341 # ifdef ALLOW_GMREDI
342 CADJ STORE Kwx = comlev1, key = ikey_dynamics
343 CADJ STORE Kwy = comlev1, key = ikey_dynamics
344 CADJ STORE Kwz = comlev1, key = ikey_dynamics
345 # ifdef GM_BOLUS_ADVEC
346 CADJ STORE GM_PsiX = comlev1, key = ikey_dynamics
347 CADJ STORE GM_PsiY = comlev1, key = ikey_dynamics
348 # endif
349 # endif
350 c
351 # ifdef ALLOW_KPP
352 CADJ STORE KPPghat = comlev1, key = ikey_dynamics
353 CADJ STORE KPPfrac = comlev1, key = ikey_dynamics
354 CADJ STORE KPPdiffKzS = comlev1, key = ikey_dynamics
355 CADJ STORE KPPdiffKzT = comlev1, key = ikey_dynamics
356 # endif
357 c
358 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
359 CADJ STORE etaH = comlev1, key = ikey_dynamics
360 # ifdef ALLOW_CD_CODE
361 CADJ STORE etanm1 = comlev1, key = ikey_dynamics
362 # endif
363 # ifndef DISABLE_RSTAR_CODE
364 cph-test
365 CADJ STORE rstarexpc = comlev1, key = ikey_dynamics
366 # endif
367 # endif
368 #endif /* ALLOW_AUTODIFF_TAMC */
369
370 IF ( .NOT.staggerTimeStep ) THEN
371 #ifdef ALLOW_DEBUG
372 IF ( debugLevel .GE. debLevB )
373 & CALL DEBUG_CALL('THERMODYNAMICS',myThid)
374 #endif
375 CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',myThid)
376 CALL THERMODYNAMICS( myTime, myIter, myThid )
377 CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',myThid)
378 C-- if not staggerTimeStep: end
379 ENDIF
380 c #ifdef ALLOW_NONHYDROSTATIC
381 IF ( implicitIntGravWave ) THEN
382 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
383 CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid )
384 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
385 ENDIF
386 c #endif
387
388 #ifdef COMPONENT_MODULE
389 IF ( useCoupler .AND. .NOT.cpl_earlyExpImpCall ) THEN
390 C Post coupling data that I export.
391 C Read in coupling data that I import.
392 myItP1 = myIter + 1
393 CALL TIMER_START('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid)
394 CALL CPL_EXPORT_MY_DATA( myItP1, myTime, myThid )
395 CALL CPL_IMPORT_EXTERNAL_DATA( myItP1, myTime, myThid )
396 CALL TIMER_STOP ('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid)
397 # ifdef ALLOW_OCN_COMPON_INTERF
398 IF ( useRealFreshWaterFlux ) THEN
399 CALL OCN_APPLY_IMPORT( .FALSE., myTime, myIter, myThid )
400 ENDIF
401 # endif /* ALLOW_OCN_COMPON_INTERF */
402 ENDIF
403 #endif /* COMPONENT_MODULE */
404
405 #ifdef ALLOW_AUTODIFF_TAMC
406 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
407 CADJ STORE hFacC = comlev1, key = ikey_dynamics
408 CADJ STORE hFacS = comlev1, key = ikey_dynamics
409 CADJ STORE hFacW = comlev1, key = ikey_dynamics
410 CADJ STORE recip_hFacC = comlev1, key = ikey_dynamics
411 CADJ STORE recip_hFacS = comlev1, key = ikey_dynamics
412 CADJ STORE recip_hFacW = comlev1, key = ikey_dynamics
413 CADJ STORE etaN = comlev1, key = ikey_dynamics
414 c
415 # ifndef DISABLE_RSTAR_CODE
416 CADJ STORE rstarFacC = comlev1, key = ikey_dynamics
417 CADJ STORE rstarFacS = comlev1, key = ikey_dynamics
418 CADJ STORE rstarFacW = comlev1, key = ikey_dynamics
419 c
420 CADJ STORE h0facc,h0facs,h0facw
421 CADJ & = comlev1, key = ikey_dynamics
422 CADJ STORE rstardhcdt,rstardhsdt,rstardhwdt
423 CADJ & = comlev1, key = ikey_dynamics
424 CADJ STORE rstarexpc,rstarexps,rstarexpw
425 CADJ & = comlev1, key = ikey_dynamics
426 # endif
427 # endif
428 #endif
429 C-- Step forward fields and calculate time tendency terms.
430 #ifndef ALLOW_AUTODIFF_TAMC
431 IF ( momStepping ) THEN
432 #endif
433 #ifdef ALLOW_DEBUG
434 IF ( debugLevel .GE. debLevB )
435 & CALL DEBUG_CALL('DYNAMICS',myThid)
436 #endif
437 CALL TIMER_START('DYNAMICS [FORWARD_STEP]',myThid)
438 CALL DYNAMICS( myTime, myIter, myThid )
439 CALL TIMER_STOP ('DYNAMICS [FORWARD_STEP]',myThid)
440 #ifndef ALLOW_AUTODIFF_TAMC
441 ENDIF
442 #endif
443
444 #ifdef ALLOW_AUTODIFF_TAMC
445 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
446 cph-test
447 CADJ STORE gU, gV = comlev1, key = ikey_dynamics
448 # endif
449 #endif
450
451 C-- Update time-counter
452 myIter = nIter0 + iLoop
453 myTime = startTime + deltaTClock * float(iLoop)
454
455 #ifdef ALLOW_MNC
456 C Update the default next iter for MNC
457 IF ( useMNC ) THEN
458 CALL MNC_CW_CITER_SETG( 1, 1, -1, myIter , myThid )
459
460 C TODO: Logic should be added here so that users can specify, on
461 C a per-citer-group basis, when it is time to update the
462 C "current" (and not just the "next") iteration
463
464 C TODO: the following is just a temporary band-aid (mostly, for
465 C Baylor) until someone writes a routine that better handles time
466 C boundaries such as weeks, months, years, etc.
467 IF ( mnc_filefreq .GT. 0 ) THEN
468 IF (DIFFERENT_MULTIPLE(mnc_filefreq,myTime,deltaTClock))
469 & THEN
470 CALL MNC_CW_CITER_SETG( 1, 1, myIter, -1 , myThid )
471 ENDIF
472 ENDIF
473 ENDIF
474 #endif /* ALLOW_MNC */
475
476 C-- Update geometric factors:
477 #ifdef NONLIN_FRSURF
478 C- update hfacC,W,S and recip_hFac according to etaH(n+1) :
479 IF ( nonlinFreeSurf.GT.0) THEN
480 IF ( select_rStar.GT.0 ) THEN
481 # ifndef DISABLE_RSTAR_CODE
482 # ifdef ALLOW_AUTODIFF_TAMC
483 cph-test
484 CADJ STORE hFacC = comlev1, key = ikey_dynamics
485 CADJ STORE hFacS = comlev1, key = ikey_dynamics
486 CADJ STORE hFacW = comlev1, key = ikey_dynamics
487 CADJ STORE recip_hFacC = comlev1, key = ikey_dynamics
488 CADJ STORE recip_hFacS = comlev1, key = ikey_dynamics
489 CADJ STORE recip_hFacW = comlev1, key = ikey_dynamics
490 c
491 CADJ STORE rstarFacC = comlev1, key = ikey_dynamics
492 CADJ STORE rstarFacS = comlev1, key = ikey_dynamics
493 CADJ STORE rstarFacW = comlev1, key = ikey_dynamics
494 c
495 CADJ STORE h0facc,h0facs,h0facw = comlev1, key = ikey_dynamics
496 # endif
497 CALL TIMER_START('UPDATE_R_STAR [FORWARD_STEP]',myThid)
498 CALL UPDATE_R_STAR( myTime, myIter, myThid )
499 CALL TIMER_STOP ('UPDATE_R_STAR [FORWARD_STEP]',myThid)
500 # ifdef ALLOW_AUTODIFF_TAMC
501 cph-test
502 CADJ STORE hFacC = comlev1, key = ikey_dynamics
503 CADJ STORE hFacS = comlev1, key = ikey_dynamics
504 CADJ STORE hFacW = comlev1, key = ikey_dynamics
505 CADJ STORE recip_hFacC = comlev1, key = ikey_dynamics
506 CADJ STORE recip_hFacS = comlev1, key = ikey_dynamics
507 CADJ STORE recip_hFacW = comlev1, key = ikey_dynamics
508 # endif
509 # endif /* DISABLE_RSTAR_CODE */
510 ELSE
511 #ifdef ALLOW_AUTODIFF_TAMC
512 CADJ STORE hFac_surfC, hFac_surfS, hFac_surfW
513 CADJ & = comlev1, key = ikey_dynamics
514 #endif
515 CALL TIMER_START('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
516 CALL UPDATE_SURF_DR( myTime, myIter, myThid )
517 CALL TIMER_STOP ('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
518 ENDIF
519 ENDIF
520 # ifdef ALLOW_AUTODIFF_TAMC
521 cph-test
522 CADJ STORE hFacC = comlev1, key = ikey_dynamics
523 CADJ STORE hFacS = comlev1, key = ikey_dynamics
524 CADJ STORE hFacW = comlev1, key = ikey_dynamics
525 CADJ STORE recip_hFacC = comlev1, key = ikey_dynamics
526 CADJ STORE recip_hFacS = comlev1, key = ikey_dynamics
527 CADJ STORE recip_hFacW = comlev1, key = ikey_dynamics
528 # endif
529 C- update also CG2D matrix (and preconditioner)
530 IF ( momStepping .AND. nonlinFreeSurf.GT.2 ) THEN
531 CALL TIMER_START('UPDATE_CG2D [FORWARD_STEP]',myThid)
532 CALL UPDATE_CG2D( myTime, myIter, myThid )
533 CALL TIMER_STOP ('UPDATE_CG2D [FORWARD_STEP]',myThid)
534 ENDIF
535 #endif /* NONLIN_FRSURF */
536
537 C-- Apply Filters to u*,v* before SOLVE_FOR_PRESSURE
538 #ifdef ALLOW_SHAP_FILT
539 IF (useSHAP_FILT .AND. shap_filt_uvStar) THEN
540 CALL TIMER_START('SHAP_FILT_UV [FORWARD_STEP]',myThid)
541 IF (implicDiv2Dflow.LT.1.) THEN
542 C-- Explicit+Implicit part of the Barotropic Flow Divergence
543 C => Filtering of uVel,vVel is necessary
544 CALL SHAP_FILT_APPLY_UV( uVel,vVel,
545 & myTime, myIter, myThid )
546 ENDIF
547 CALL SHAP_FILT_APPLY_UV( gU,gV,myTime,myIter,myThid)
548 CALL TIMER_STOP ('SHAP_FILT_UV [FORWARD_STEP]',myThid)
549 ENDIF
550 #endif
551 #ifdef ALLOW_ZONAL_FILT
552 IF (useZONAL_FILT .AND. zonal_filt_uvStar) THEN
553 CALL TIMER_START('ZONAL_FILT_UV [FORWARD_STEP]',myThid)
554 IF (implicDiv2Dflow.LT.1.) THEN
555 C-- Explicit+Implicit part of the Barotropic Flow Divergence
556 C => Filtering of uVel,vVel is necessary
557 CALL ZONAL_FILT_APPLY_UV( uVel, vVel, myThid )
558 ENDIF
559 CALL ZONAL_FILT_APPLY_UV( gU, gV, myThid )
560 CALL TIMER_STOP ('ZONAL_FILT_UV [FORWARD_STEP]',myThid)
561 ENDIF
562 #endif
563
564 C-- Solve elliptic equation(s).
565 C Two-dimensional only for conventional hydrostatic or
566 C three-dimensional for non-hydrostatic and/or IGW scheme.
567 IF ( momStepping ) THEN
568 #ifdef ALLOW_AUTODIFF_TAMC
569 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
570 CADJ STORE uvel, vvel
571 CADJ & = comlev1, key = ikey_dynamics
572 CADJ STORE empmr,hfacs,hfacw
573 CADJ & = comlev1, key = ikey_dynamics
574 # endif
575 #endif
576 CALL TIMER_START('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
577 CALL SOLVE_FOR_PRESSURE(myTime, myIter, myThid)
578 CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
579 ENDIF
580
581 C-- Correct divergence in flow field and cycle time-stepping momentum
582 #ifndef ALLOW_AUTODIFF_TAMC
583 IF ( momStepping ) THEN
584 #endif
585 #ifdef ALLOW_AUTODIFF_TAMC
586 # if (defined NONLIN_FRSURF) || (defined ALLOW_DEPTH_CONTROL)
587 # ifndef DISABLE_RSTAR_CODE
588 cph-test
589 cph not clear, why this one
590 CADJ STORE h0facc = comlev1, key = ikey_dynamics
591 # endif
592 # endif
593 # ifdef ALLOW_DEPTH_CONTROL
594 CADJ STORE etan, uvel,vvel
595 CADJ & = comlev1, key = ikey_dynamics
596 # endif
597 #endif
598 CALL TIMER_START('MOM_CORRECTION_STEP [FORWARD_STEP]',myThid)
599 CALL MOMENTUM_CORRECTION_STEP(myTime, myIter, myThid)
600 CALL TIMER_STOP ('MOM_CORRECTION_STEP [FORWARD_STEP]',myThid)
601 #ifndef ALLOW_AUTODIFF_TAMC
602 ENDIF
603 #endif
604
605 #ifdef EXACT_CONSERV
606 IF (exactConserv) THEN
607 C-- Update etaH(n+1) :
608 CALL TIMER_START('UPDATE_ETAH [FORWARD_STEP]',myThid)
609 CALL UPDATE_ETAH( myTime, myIter, myThid )
610 CALL TIMER_STOP ('UPDATE_ETAH [FORWARD_STEP]',myThid)
611 ENDIF
612 #endif /* EXACT_CONSERV */
613
614 #ifdef NONLIN_FRSURF
615 IF ( select_rStar.NE.0 ) THEN
616 # ifndef DISABLE_RSTAR_CODE
617 # ifdef ALLOW_AUTODIFF_TAMC
618 cph-test
619 CADJ STORE rstarfacc,rstarfacs,rstarfacw = comlev1, key = ikey_dynamics
620 # endif
621 C-- r* : compute the future level thickness according to etaH(n+1)
622 CALL TIMER_START('CALC_R_STAR [FORWARD_STEP]',myThid)
623 CALL CALC_R_STAR(etaH, myTime, myIter, myThid )
624 CALL TIMER_STOP ('CALC_R_STAR [FORWARD_STEP]',myThid)
625 # endif /* DISABLE_RSTAR_CODE */
626 ELSEIF ( nonlinFreeSurf.GT.0) THEN
627 C-- compute the future surface level thickness according to etaH(n+1)
628 # ifdef ALLOW_AUTODIFF_TAMC
629 CADJ STORE etaH = comlev1, key = ikey_dynamics
630 # endif
631 CALL TIMER_START('CALC_SURF_DR [FORWARD_STEP]',myThid)
632 CALL CALC_SURF_DR(etaH, myTime, myIter, myThid )
633 CALL TIMER_STOP ('CALC_SURF_DR [FORWARD_STEP]',myThid)
634 ENDIF
635 # ifdef ALLOW_AUTODIFF_TAMC
636 cph-test
637 CADJ STORE hFac_surfC = comlev1, key = ikey_dynamics
638 # endif
639 #endif /* NONLIN_FRSURF */
640
641 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
642 IF ( staggerTimeStep ) THEN
643 C-- do exchanges of U,V (needed for multiDim) when using stagger time-step :
644 #ifdef ALLOW_DEBUG
645 IF ( debugLevel .GE. debLevB )
646 & CALL DEBUG_CALL('DO_STAGGER_FIELDS_EXCH.',myThid)
647 #endif
648 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
649 CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid )
650 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
651
652 #ifdef ALLOW_DIAGNOSTICS
653 C-- State-variables diagnostics
654 IF ( useDiagnostics ) THEN
655 CALL TIMER_START('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
656 CALL DO_STATEVARS_DIAGS( myTime, 1, myIter, myThid )
657 CALL TIMER_STOP ('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
658 ENDIF
659 #endif
660
661 #ifdef ALLOW_DEBUG
662 IF ( debugLevel .GE. debLevB )
663 & CALL DEBUG_CALL('THERMODYNAMICS',myThid)
664 #endif
665 CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',myThid)
666 CALL THERMODYNAMICS( myTime, myIter, myThid )
667 CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',myThid)
668
669 C-- if staggerTimeStep: end
670 ENDIF
671 C---+--------+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
672
673 #ifdef ALLOW_AUTODIFF_TAMC
674 cph This is needed because convective_adjustment calls
675 cph find_rho which may use pressure()
676 CADJ STORE totphihyd = comlev1, key = ikey_dynamics
677 #endif
678 C-- Cycle time-stepping Tracers arrays (T,S,+pTracers)
679 CALL TIMER_START('TRC_CORRECTION_STEP [FORWARD_STEP]',myThid)
680 CALL TRACERS_CORRECTION_STEP(myTime, myIter, myThid)
681 CALL TIMER_STOP ('TRC_CORRECTION_STEP [FORWARD_STEP]',myThid)
682
683 #ifdef ALLOW_GCHEM
684 C Add separate timestepping of chemical/biological/forcing
685 C of ptracers here in GCHEM_FORCING_SEP
686 IF ( useGCHEM ) THEN
687 #ifdef ALLOW_DEBUG
688 IF ( debugLevel .GE. debLevB )
689 & CALL DEBUG_CALL('GCHEM_FORCING_SEP',myThid)
690 #endif /* ALLOW_DEBUG */
691 CALL TIMER_START('GCHEM_FORCING_SEP [FORWARD_STEP]',myThid)
692 CALL GCHEM_FORCING_SEP( myTime,myIter,myThid )
693 CALL TIMER_STOP ('GCHEM_FORCING_SEP [FORWARD_STEP]',myThid)
694 ENDIF
695 #endif /* ALLOW_GCHEM */
696
697 C-- Do "blocking" sends and receives for tendency "overlap" terms
698 c CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
699 c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )
700 c CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
701
702 C-- Do "blocking" sends and receives for field "overlap" terms
703 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
704 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
705 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
706
707 #ifdef ALLOW_DIAGNOSTICS
708 IF ( useDiagnostics ) THEN
709 CALL TIMER_START('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
710 CALL DO_STATEVARS_DIAGS( myTime, 2, myIter, myThid )
711 CALL TIMER_STOP ('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid)
712 ENDIF
713 #endif
714
715 #ifdef ALLOW_GRIDALT
716 IF (useGRIDALT) THEN
717 CALL GRIDALT_UPDATE(myThid)
718 ENDIF
719 #endif
720
721 #ifdef ALLOW_FIZHI
722 IF (useFIZHI) THEN
723 CALL TIMER_START('FIZHI [FORWARD_STEP]',myThid)
724 CALL STEP_FIZHI_CORR ( myTime, myIter, myThid, dTtracerLev(1) )
725 CALL TIMER_STOP ('FIZHI [FORWARD_STEP]',myThid)
726 ENDIF
727 #endif
728
729 #ifdef ALLOW_FLT
730 C-- Calculate float trajectories
731 IF (useFLT) THEN
732 CALL TIMER_START('FLOATS [FORWARD_STEP]',myThid)
733 CALL FLT_MAIN(myIter,myTime, myThid)
734 CALL TIMER_STOP ('FLOATS [FORWARD_STEP]',myThid)
735 ENDIF
736 #endif
737
738 #ifdef ALLOW_TIMEAVE
739 C-- State-variables time-averaging
740 CALL TIMER_START('DO_STATEVARS_TAVE [FORWARD_STEP]',myThid)
741 CALL DO_STATEVARS_TAVE( myTime, myIter, myThid )
742 CALL TIMER_STOP ('DO_STATEVARS_TAVE [FORWARD_STEP]',myThid)
743 #endif
744
745 #ifdef ALLOW_MONITOR
746 IF ( .NOT.useOffLine ) THEN
747 C-- Check status of solution (statistics, cfl, etc...)
748 CALL TIMER_START('MONITOR [FORWARD_STEP]',myThid)
749 CALL MONITOR( myIter, myTime, myThid )
750 CALL TIMER_STOP ('MONITOR [FORWARD_STEP]',myThid)
751 ENDIF
752 #endif /* ALLOW_MONITOR */
753
754 #ifdef ALLOW_COST
755 C-- compare model with data and compute cost function
756 C-- this is done after exchanges to allow interpolation
757 CALL TIMER_START('COST_TILE [FORWARD_STEP]',myThid)
758 CALL COST_TILE ( myTime, myIter, myThid )
759 CALL TIMER_STOP ('COST_TILE [FORWARD_STEP]',myThid)
760 #endif
761
762 C-- Do IO if needed.
763 CALL TIMER_START('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
764 CALL DO_THE_MODEL_IO( myTime, myIter, myThid )
765 CALL TIMER_STOP ('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
766
767 #ifdef HAVE_SIGREG
768 IF ( useSIGREG ) THEN
769 IF ( i_got_signal .GT. 0 ) THEN
770 CALL PACKAGES_WRITE_PICKUP(
771 I .TRUE., myTime, myIter, myThid )
772 CALL WRITE_PICKUP(
773 I .TRUE., myTime, myIter, myThid )
774 STOP 'Checkpoint completed -- killed by signal handler'
775 ENDIF
776 ENDIF
777 #endif /* HAVE_SIGREG */
778
779 C-- Save state for restarts
780 CALL TIMER_START('DO_WRITE_PICKUP [FORWARD_STEP]',myThid)
781 CALL DO_WRITE_PICKUP(
782 I .FALSE., myTime, myIter, myThid )
783 CALL TIMER_STOP ('DO_WRITE_PICKUP [FORWARD_STEP]',myThid)
784
785 #ifdef ALLOW_DEBUG
786 IF ( debugLevel .GE. debLevB )
787 & CALL DEBUG_LEAVE('FORWARD_STEP',myThid)
788 #endif
789
790 RETURN
791 END

  ViewVC Help
Powered by ViewVC 1.1.22