/[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.139 - (show annotations) (download)
Sun Apr 9 17:35:30 2006 UTC (18 years, 2 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint58d_post
Changes since 1.138: +4 -1 lines
Starting thsice adjoint

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

  ViewVC Help
Powered by ViewVC 1.1.22