/[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.146 - (show annotations) (download)
Thu Dec 14 22:31:18 2006 UTC (17 years, 5 months ago) by heimbach
Branch: MAIN
Changes since 1.145: +4 -1 lines
Updating seaice adjoint, step 1 (everything, except SEAICE_EVP).

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

  ViewVC Help
Powered by ViewVC 1.1.22