/[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.158 - (show annotations) (download)
Fri May 2 00:13:55 2008 UTC (16 years, 1 month ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint60, checkpoint59q, checkpoint59r
Changes since 1.157: +3 -8 lines
Clean up slightly.

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

  ViewVC Help
Powered by ViewVC 1.1.22