/[MITgcm]/MITgcm_contrib/sannino/OASIS_3.0_Coupler/code/forward_step.F
ViewVC logotype

Contents of /MITgcm_contrib/sannino/OASIS_3.0_Coupler/code/forward_step.F

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


Revision 1.1 - (show annotations) (download)
Thu Jul 20 21:08:16 2006 UTC (19 years ago) by sannino
Branch: MAIN
CVS Tags: HEAD
o Adding OASIS package
o Adding grid refinement package

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

  ViewVC Help
Powered by ViewVC 1.1.22