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