1 |
C $Header: /u/gcmpack/MITgcm/model/src/forward_step.F,v 1.121 2005/09/14 18:32:07 edhill Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "PACKAGES_CONFIG.h" |
5 |
#include "CPP_OPTIONS.h" |
6 |
|
7 |
#ifdef ALLOW_OFFLINE |
8 |
# include "OFFLINE_OPTIONS.h" |
9 |
#endif |
10 |
#ifdef ALLOW_GMREDI |
11 |
# include "GMREDI_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 ALLOW_SHAP_FILT |
50 |
# include "SHAP_FILT.h" |
51 |
#endif |
52 |
#ifdef ALLOW_ZONAL_FILT |
53 |
# include "ZONAL_FILT.h" |
54 |
#endif |
55 |
#ifdef COMPONENT_MODULE |
56 |
# include "CPL_PARAMS.h" |
57 |
#endif |
58 |
|
59 |
#ifdef ALLOW_AUTODIFF_TAMC |
60 |
# include "FFIELDS.h" |
61 |
|
62 |
# ifdef ALLOW_NONHYDROSTATIC |
63 |
# include "CG3D.h" |
64 |
# endif |
65 |
|
66 |
# include "tamc.h" |
67 |
# include "ctrl.h" |
68 |
# include "ctrl_dummy.h" |
69 |
# include "cost.h" |
70 |
# include "EOS.h" |
71 |
# ifdef ALLOW_EXF |
72 |
# include "exf_fields.h" |
73 |
# include "exf_clim_fields.h" |
74 |
# ifdef ALLOW_BULKFORMULAE |
75 |
# include "exf_constants.h" |
76 |
# endif |
77 |
# endif |
78 |
# ifdef ALLOW_OBCS |
79 |
# include "OBCS.h" |
80 |
# endif |
81 |
# ifdef ALLOW_PTRACERS |
82 |
# include "PTRACERS_SIZE.h" |
83 |
# include "PTRACERS.h" |
84 |
# endif |
85 |
# ifdef ALLOW_CD_CODE |
86 |
# include "CD_CODE_VARS.h" |
87 |
# endif |
88 |
# ifdef ALLOW_EBM |
89 |
# include "EBM.h" |
90 |
# endif |
91 |
# ifdef EXACT_CONSERV |
92 |
# include "SURFACE.h" |
93 |
# endif |
94 |
# ifdef ALLOW_KPP |
95 |
# include "KPP.h" |
96 |
# endif |
97 |
# ifdef ALLOW_GMREDI |
98 |
# include "GMREDI.h" |
99 |
# endif |
100 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
101 |
|
102 |
C !LOCAL VARIABLES: |
103 |
C == Routine arguments == |
104 |
C note: under the multi-threaded model myiter and |
105 |
C mytime are local variables passed around as routine |
106 |
C arguments. Although this is fiddly it saves the need to |
107 |
C impose additional synchronisation points when they are |
108 |
C updated. |
109 |
C myIter - iteration counter for this thread |
110 |
C myTime - time counter for this thread |
111 |
C myThid - thread number for this instance of the routine. |
112 |
INTEGER iloop |
113 |
INTEGER myThid |
114 |
INTEGER myIter |
115 |
_RL myTime |
116 |
|
117 |
C == Local variables == |
118 |
INTEGER myItP1 |
119 |
CEOP |
120 |
|
121 |
#ifdef ALLOW_DEBUG |
122 |
IF ( debugLevel .GE. debLevB ) |
123 |
& CALL DEBUG_ENTER('FORWARD_STEP',myThid) |
124 |
#endif |
125 |
|
126 |
#ifdef ALLOW_AUTODIFF_TAMC |
127 |
C-- Reset the model iteration counter and the model time. |
128 |
myIter = nIter0 + (iloop-1) |
129 |
myTime = startTime + float(iloop-1)*deltaTclock |
130 |
#endif |
131 |
|
132 |
#ifdef ALLOW_AUTODIFF_TAMC |
133 |
c************************************** |
134 |
#include "checkpoint_lev1_directives.h" |
135 |
c************************************** |
136 |
#endif |
137 |
|
138 |
C-- Switch on/off diagnostics for snap-shot output: |
139 |
#ifdef ALLOW_DIAGNOSTICS |
140 |
IF ( useDiagnostics ) THEN |
141 |
CALL DIAGNOSTICS_SWITCH_ONOFF( myTime, myIter, myThid ) |
142 |
ENDIF |
143 |
#endif |
144 |
|
145 |
C-- State-variables diagnostics |
146 |
IF ( usediagnostics ) THEN |
147 |
CALL TIMER_START('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid) |
148 |
CALL DO_STATEVARS_DIAGS( myTime, 0, myIter, myThid ) |
149 |
CALL TIMER_STOP ('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid) |
150 |
ENDIF |
151 |
|
152 |
C-- Call external forcing package |
153 |
#ifdef ALLOW_BULK_FORCE |
154 |
IF ( useBulkForce ) THEN |
155 |
#ifdef ALLOW_DEBUG |
156 |
IF ( debugLevel .GE. debLevB ) |
157 |
& CALL DEBUG_CALL('BULKF_FIELDS_LOAD',myThid) |
158 |
#endif |
159 |
CALL TIMER_START('BULKF_FORCING [FORWARD_STEP]',mythid) |
160 |
C- load all forcing fields at current time |
161 |
CALL BULKF_FIELDS_LOAD( myTime, myIter, myThid ) |
162 |
C- calculate qnet and empmr (and wind stress) |
163 |
CALL BULKF_FORCING( myTime, myIter, myThid ) |
164 |
CALL TIMER_STOP ('BULKF_FORCING [FORWARD_STEP]',mythid) |
165 |
ELSE |
166 |
#endif /* ALLOW_BULK_FORCE */ |
167 |
|
168 |
# ifdef ALLOW_EXF |
169 |
# ifdef ALLOW_DEBUG |
170 |
IF ( debugLevel .GE. debLevB ) |
171 |
& CALL DEBUG_CALL('EXF_GETFORCING',myThid) |
172 |
# endif |
173 |
CALL TIMER_START('EXF_GETFORCING [FORWARD_STEP]',mythid) |
174 |
CALL EXF_GETFORCING( mytime, myiter, mythid ) |
175 |
CALL TIMER_STOP ('EXF_GETFORCING [FORWARD_STEP]',mythid) |
176 |
# else /* ALLOW_EXF undef */ |
177 |
cph The following IF-statement creates an additional dependency |
178 |
cph for the forcing fields requiring additional storing. |
179 |
cph Therefore, the IF-statement will be put between CPP-OPTIONS, |
180 |
cph assuming that ALLOW_SEAICE has not yet been differentiated. |
181 |
# if (defined (ALLOW_SEAICE) || defined (ALLOW_EBM)) |
182 |
IF ( .NOT. useSEAICE .AND. .NOT. useEBM ) THEN |
183 |
# endif |
184 |
#ifdef ALLOW_DEBUG |
185 |
IF ( debugLevel .GE. debLevB ) |
186 |
& CALL DEBUG_CALL('EXTERNAL_FIELDS_LOAD',myThid) |
187 |
#endif |
188 |
CALL TIMER_START('EXTERNAL_FIELDS_LOAD[FORWARD_STEP]',mythid) |
189 |
CALL EXTERNAL_FIELDS_LOAD( mytime, myiter, mythid ) |
190 |
CALL TIMER_STOP ('EXTERNAL_FIELDS_LOAD[FORWARD_STEP]',mythid) |
191 |
# if (defined (ALLOW_SEAICE) || defined (ALLOW_EBM)) |
192 |
ENDIF |
193 |
# endif |
194 |
# endif /* ALLOW_EXF */ |
195 |
#ifdef ALLOW_BULK_FORCE |
196 |
C-- end of if/else block useBulfforce -- |
197 |
ENDIF |
198 |
#endif /* ALLOW_BULK_FORCE */ |
199 |
|
200 |
#ifdef ALLOW_AUTODIFF |
201 |
c-- Add control vector for forcing and parameter fields |
202 |
if ( myiter .EQ. nIter0 ) |
203 |
& CALL CTRL_MAP_FORCING (mythid) |
204 |
#endif |
205 |
|
206 |
#if (defined (ALLOW_AUTODIFF_TAMC) && defined (ALLOW_AUTODIFF_MONITOR)) |
207 |
C Include call to a dummy routine. Its adjoint will be |
208 |
C called at the proper place in the adjoint code. |
209 |
C The adjoint routine will print out adjoint values |
210 |
C if requested. The location of the call is important, |
211 |
C it has to be after the adjoint of the exchanges |
212 |
C (DO_GTERM_BLOCKING_EXCHANGES). |
213 |
CALL DUMMY_IN_STEPPING( myTime, myIter, myThid ) |
214 |
cph I've commented this line since it may conflict with MITgcm's adjoint |
215 |
cph However, need to check whether that's still consistent |
216 |
cph with the ecco-branch (it should). |
217 |
cph CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid ) |
218 |
#endif |
219 |
|
220 |
# ifdef ALLOW_SEAICE |
221 |
C-- Call sea ice model to compute forcing/external data fields. In |
222 |
C addition to computing prognostic sea-ice variables and diagnosing the |
223 |
C forcing/external data fields that drive the ocean model, SEAICE_MODEL |
224 |
C also sets theta to the freezing point under sea-ice. The implied |
225 |
C surface heat flux is then stored in variable surfaceTendencyTice, |
226 |
C which is needed by KPP package (kpp_calc.F and kpp_transport_t.F) |
227 |
C to diagnose surface buoyancy fluxes and for the non-local transport |
228 |
C term. Because this call precedes model thermodynamics, temperature |
229 |
C under sea-ice may not be "exactly" at the freezing point by the time |
230 |
C theta is dumped or time-averaged. |
231 |
IF ( useSEAICE ) THEN |
232 |
#ifdef ALLOW_DEBUG |
233 |
IF ( debugLevel .GE. debLevB ) |
234 |
& CALL DEBUG_CALL('SEAICE_MODEL',myThid) |
235 |
#endif |
236 |
CALL TIMER_START('SEAICE_MODEL [FORWARD_STEP]',myThid) |
237 |
CALL SEAICE_MODEL( myTime, myIter, myThid ) |
238 |
CALL TIMER_STOP ('SEAICE_MODEL [FORWARD_STEP]',myThid) |
239 |
ENDIF |
240 |
# endif /* ALLOW_SEAICE */ |
241 |
|
242 |
#ifdef ALLOW_AUTODIFF_TAMC |
243 |
# ifdef ALLOW_PTRACERS |
244 |
cph this replaces _bibj storing of ptracer within thermodynamics |
245 |
CADJ STORE ptracer = comlev1, key = ikey_dynamics |
246 |
# endif |
247 |
#endif |
248 |
|
249 |
#ifdef ALLOW_OFFLINE |
250 |
call OFFLINE_FIELDS_LOAD( myTime, myIter, myThid ) |
251 |
#endif |
252 |
|
253 |
#ifdef ALLOW_PTRACERS |
254 |
# ifdef ALLOW_GCHEM |
255 |
IF ( useGCHEM ) THEN |
256 |
#ifdef ALLOW_DEBUG |
257 |
IF ( debugLevel .GE. debLevB ) |
258 |
& CALL DEBUG_CALL('GCHEM_FIELDS_LOAD',myThid) |
259 |
#endif /* ALLOW_DEBUG */ |
260 |
CALL GCHEM_FIELDS_LOAD( mytime, myiter, mythid ) |
261 |
ENDIF |
262 |
# endif |
263 |
#endif |
264 |
|
265 |
#ifdef COMPONENT_MODULE |
266 |
IF ( useCoupler .AND. cpl_earlyExpImpCall ) THEN |
267 |
C Post coupling data that I export. |
268 |
C Read in coupling data that I import. |
269 |
CALL TIMER_START('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid) |
270 |
CALL CPL_EXPORT_MY_DATA( myIter, myTime, myThid ) |
271 |
CALL CPL_IMPORT_EXTERNAL_DATA( myIter, myTime, myThid ) |
272 |
CALL TIMER_STOP ('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid) |
273 |
ENDIF |
274 |
#endif /* COMPONENT_MODULE */ |
275 |
|
276 |
#ifdef ALLOW_EBM |
277 |
IF ( useEBM ) THEN |
278 |
# ifdef ALLOW_DEBUG |
279 |
IF ( debugLevel .GE. debLevB ) |
280 |
& CALL DEBUG_CALL('EBM',myThid) |
281 |
# endif |
282 |
CALL TIMER_START('EBM [FORWARD_STEP]',mythid) |
283 |
CALL EBM_DRIVER ( myTime, myIter, myThid ) |
284 |
CALL TIMER_STOP ('EBM [FORWARD_STEP]',mythid) |
285 |
ENDIF |
286 |
#endif |
287 |
|
288 |
C-- Step forward fields and calculate time tendency terms. |
289 |
|
290 |
#ifdef ALLOW_DEBUG |
291 |
IF ( debugLevel .GE. debLevB ) |
292 |
& CALL DEBUG_CALL('DO_ATMOSPHERIC_PHYS',myThid) |
293 |
#endif |
294 |
CALL TIMER_START('DO_ATMOSPHERIC_PHYS [FORWARD_STEP]',mythid) |
295 |
CALL DO_ATMOSPHERIC_PHYS( myTime, myIter, myThid ) |
296 |
CALL TIMER_STOP ('DO_ATMOSPHERIC_PHYS [FORWARD_STEP]',mythid) |
297 |
|
298 |
#ifdef ALLOW_AUTODIFF_TAMC |
299 |
CADJ STORE theta = comlev1, key = ikey_dynamics |
300 |
CADJ STORE salt = comlev1, key = ikey_dynamics |
301 |
CADJ STORE totphihyd = comlev1, key = ikey_dynamics |
302 |
CADJ STORE surfaceforcingtice = comlev1, key = ikey_dynamics |
303 |
# ifdef ALLOW_KPP |
304 |
CADJ STORE uvel = comlev1, key = ikey_dynamics |
305 |
CADJ STORE vvel = comlev1, key = ikey_dynamics |
306 |
# endif |
307 |
# ifdef EXACT_CONSERV |
308 |
CADJ STORE empmr = comlev1, key = ikey_dynamics |
309 |
CADJ STORE pmepr = comlev1, key = ikey_dynamics |
310 |
# endif |
311 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
312 |
|
313 |
#ifndef ALLOW_OFFLINE |
314 |
#ifdef ALLOW_DEBUG |
315 |
IF ( debugLevel .GE. debLevB ) |
316 |
& CALL DEBUG_CALL('DO_OCEANIC_PHYS',myThid) |
317 |
#endif |
318 |
CALL TIMER_START('DO_OCEANIC_PHYS [FORWARD_STEP]',mythid) |
319 |
CALL DO_OCEANIC_PHYS( myTime, myIter, myThid ) |
320 |
CALL TIMER_STOP ('DO_OCEANIC_PHYS [FORWARD_STEP]',mythid) |
321 |
#endif |
322 |
|
323 |
#ifdef ALLOW_GCHEM |
324 |
C GCHEM package is an interface for any bio-geochemical or |
325 |
C ecosystem model you would like to include. |
326 |
C If GCHEM_SEPARATE_FORCING is not defined, you are |
327 |
C responsible for computing tendency terms for passive |
328 |
C tracers and storing them on a 3DxNumPtracers-array called |
329 |
C gchemTendency in GCHEM_CALC_TENDENCY. This tendency is then added |
330 |
C to gPtr in ptracers_forcing later-on. |
331 |
C If GCHEM_SEPARATE_FORCING is defined, you are reponsible for |
332 |
C UPDATING ptracers directly in GCHEM_FORCING_SEP. This amounts |
333 |
C to a completely separate time step that you have to implement |
334 |
C yourself (Eulerian seems to be fine in most cases). |
335 |
CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC |
336 |
C CAVEAT: Up to now, when GCHEM is turned on the field ptracerForcingSurf, |
337 |
C which is needed for KPP is not set properly. ptracerForcingSurf must |
338 |
C be treated differently depending on whether GCHEM_SEPARATE_FORCING |
339 |
C is define or not. TBD. |
340 |
CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC |
341 |
IF ( useGCHEM ) THEN |
342 |
#ifdef ALLOW_DEBUG |
343 |
IF ( debugLevel .GE. debLevB ) |
344 |
& CALL DEBUG_CALL('GCHEM_CALC_TENDENCY',myThid) |
345 |
#endif |
346 |
CALL TIMER_START('GCHEM_CALC_TENDENCY [FORWARD_STEP]',myThid) |
347 |
CALL GCHEM_CALC_TENDENCY( myTime, myIter, myThid ) |
348 |
CALL TIMER_STOP ('GCHEM_CALC_TENDENCY [FORWARD_STEP]',myThid) |
349 |
ENDIF |
350 |
#endif /* ALLOW_GCHEM */ |
351 |
|
352 |
#ifdef ALLOW_AUTODIFF_TAMC |
353 |
cph needed to be moved here from do_oceanic_physics |
354 |
cph to be visible down the road |
355 |
c |
356 |
CADJ STORE surfaceForcingS = comlev1, key = ikey_dynamics |
357 |
CADJ STORE surfaceForcingT = comlev1, key = ikey_dynamics |
358 |
CADJ STORE surfaceForcingTice = comlev1, key = ikey_dynamics |
359 |
ctest( |
360 |
CADJ STORE IVDConvCount = comlev1, key = ikey_dynamics |
361 |
ctest) |
362 |
# ifdef ALLOW_PTRACERS |
363 |
CADJ STORE surfaceForcingPtr = comlev1, key = ikey_dynamics |
364 |
# endif |
365 |
c |
366 |
# ifdef ALLOW_GMREDI |
367 |
CADJ STORE Kwx = comlev1, key = ikey_dynamics |
368 |
CADJ STORE Kwy = comlev1, key = ikey_dynamics |
369 |
CADJ STORE Kwz = comlev1, key = ikey_dynamics |
370 |
# ifdef GM_BOLUS_ADVEC |
371 |
CADJ STORE GM_PsiX = comlev1, key = ikey_dynamics |
372 |
CADJ STORE GM_PsiY = comlev1, key = ikey_dynamics |
373 |
# endif |
374 |
# endif |
375 |
c |
376 |
# ifdef ALLOW_KPP |
377 |
CADJ STORE KPPghat = comlev1, key = ikey_dynamics |
378 |
CADJ STORE KPPfrac = comlev1, key = ikey_dynamics |
379 |
CADJ STORE KPPdiffKzS = comlev1, key = ikey_dynamics |
380 |
CADJ STORE KPPdiffKzT = comlev1, key = ikey_dynamics |
381 |
# endif |
382 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
383 |
|
384 |
IF ( .NOT.staggerTimeStep ) THEN |
385 |
#ifdef ALLOW_DEBUG |
386 |
IF ( debugLevel .GE. debLevB ) |
387 |
& CALL DEBUG_CALL('THERMODYNAMICS',myThid) |
388 |
#endif |
389 |
CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',mythid) |
390 |
CALL THERMODYNAMICS( myTime, myIter, myThid ) |
391 |
CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',mythid) |
392 |
C-- if not staggerTimeStep: end |
393 |
ENDIF |
394 |
|
395 |
#ifdef COMPONENT_MODULE |
396 |
IF ( useCoupler .AND. .NOT.cpl_earlyExpImpCall ) THEN |
397 |
C Post coupling data that I export. |
398 |
C Read in coupling data that I import. |
399 |
myItP1 = myIter + 1 |
400 |
CALL TIMER_START('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid) |
401 |
CALL CPL_EXPORT_MY_DATA( myItP1, myTime, myThid ) |
402 |
CALL CPL_IMPORT_EXTERNAL_DATA( myItP1, myTime, myThid ) |
403 |
CALL TIMER_STOP ('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid) |
404 |
# ifndef ALLOW_AIM |
405 |
IF ( useRealFreshWaterFlux ) THEN |
406 |
CALL OCN_APPLY_IMPORT( .FALSE., myTime, myIter, myThid ) |
407 |
ENDIF |
408 |
# endif |
409 |
ENDIF |
410 |
#endif /* COMPONENT_MODULE */ |
411 |
|
412 |
C-- Step forward fields and calculate time tendency terms. |
413 |
#ifndef ALLOW_OFFLINE |
414 |
#ifndef ALLOW_AUTODIFF_TAMC |
415 |
IF ( momStepping ) THEN |
416 |
#endif |
417 |
#ifdef ALLOW_DEBUG |
418 |
IF ( debugLevel .GE. debLevB ) |
419 |
& CALL DEBUG_CALL('DYNAMICS',myThid) |
420 |
#endif |
421 |
CALL TIMER_START('DYNAMICS [FORWARD_STEP]',mythid) |
422 |
CALL DYNAMICS( myTime, myIter, myThid ) |
423 |
CALL TIMER_STOP ('DYNAMICS [FORWARD_STEP]',mythid) |
424 |
#ifndef ALLOW_AUTODIFF_TAMC |
425 |
ENDIF |
426 |
#endif |
427 |
#endif |
428 |
|
429 |
C-- Update time-counter |
430 |
myIter = nIter0 + iLoop |
431 |
myTime = startTime + deltaTClock * float(iLoop) |
432 |
|
433 |
#ifdef ALLOW_MNC |
434 |
C Update the default next iter for MNC |
435 |
IF ( useMNC ) THEN |
436 |
CALL MNC_CW_CITER_SETG( 1, 1, -1, myIter , myThid ) |
437 |
|
438 |
C TODO: Logic should be added here so that users can specify, on |
439 |
C a per-citer-group basis, when it is time to update the |
440 |
C "current" (and not just the "next") iteration |
441 |
|
442 |
C TODO: the following is just a temporary band-aid (mostly, for |
443 |
C Baylor) until someone writes a routine that better handles time |
444 |
C boundaries such as weeks, months, years, etc. |
445 |
IF ( mnc_filefreq .GT. 0 ) THEN |
446 |
IF (DIFFERENT_MULTIPLE(mnc_filefreq,myTime,deltaTClock)) |
447 |
& THEN |
448 |
CALL MNC_CW_CITER_SETG( 1, 1, myIter, -1 , myThid ) |
449 |
ENDIF |
450 |
ENDIF |
451 |
ENDIF |
452 |
#endif |
453 |
|
454 |
C-- Update geometric factors: |
455 |
#ifdef NONLIN_FRSURF |
456 |
C- update hfacC,W,S and recip_hFac according to etaH(n+1) : |
457 |
IF ( nonlinFreeSurf.GT.0) THEN |
458 |
IF ( select_rStar.GT.0 ) THEN |
459 |
CALL TIMER_START('UPDATE_R_STAR [FORWARD_STEP]',myThid) |
460 |
CALL UPDATE_R_STAR( myTime, myIter, myThid ) |
461 |
CALL TIMER_STOP ('UPDATE_R_STAR [FORWARD_STEP]',myThid) |
462 |
ELSE |
463 |
CALL TIMER_START('UPDATE_SURF_DR [FORWARD_STEP]',myThid) |
464 |
CALL UPDATE_SURF_DR( myTime, myIter, myThid ) |
465 |
CALL TIMER_STOP ('UPDATE_SURF_DR [FORWARD_STEP]',myThid) |
466 |
ENDIF |
467 |
ENDIF |
468 |
C- update also CG2D matrix (and preconditioner) |
469 |
IF ( momStepping .AND. nonlinFreeSurf.GT.2 ) THEN |
470 |
CALL TIMER_START('UPDATE_CG2D [FORWARD_STEP]',myThid) |
471 |
CALL UPDATE_CG2D( myTime, myIter, myThid ) |
472 |
CALL TIMER_STOP ('UPDATE_CG2D [FORWARD_STEP]',myThid) |
473 |
ENDIF |
474 |
#endif |
475 |
|
476 |
C-- Apply Filters to u*,v* before SOLVE_FOR_PRESSURE |
477 |
#ifdef ALLOW_SHAP_FILT |
478 |
IF (useSHAP_FILT .AND. shap_filt_uvStar) THEN |
479 |
CALL TIMER_START('SHAP_FILT_UV [FORWARD_STEP]',myThid) |
480 |
IF (implicDiv2Dflow.LT.1.) THEN |
481 |
C-- Explicit+Implicit part of the Barotropic Flow Divergence |
482 |
C => Filtering of uVel,vVel is necessary |
483 |
CALL SHAP_FILT_APPLY_UV( uVel,vVel, |
484 |
& myTime, myIter, myThid ) |
485 |
ENDIF |
486 |
CALL SHAP_FILT_APPLY_UV( gU,gV,myTime,myIter,myThid) |
487 |
CALL TIMER_STOP ('SHAP_FILT_UV [FORWARD_STEP]',myThid) |
488 |
ENDIF |
489 |
#endif |
490 |
#ifdef ALLOW_ZONAL_FILT |
491 |
IF (useZONAL_FILT .AND. zonal_filt_uvStar) THEN |
492 |
CALL TIMER_START('ZONAL_FILT_UV [FORWARD_STEP]',myThid) |
493 |
IF (implicDiv2Dflow.LT.1.) THEN |
494 |
C-- Explicit+Implicit part of the Barotropic Flow Divergence |
495 |
C => Filtering of uVel,vVel is necessary |
496 |
CALL ZONAL_FILT_APPLY_UV( uVel, vVel, myThid ) |
497 |
ENDIF |
498 |
CALL ZONAL_FILT_APPLY_UV( gU, gV, myThid ) |
499 |
CALL TIMER_STOP ('ZONAL_FILT_UV [FORWARD_STEP]',myThid) |
500 |
ENDIF |
501 |
#endif |
502 |
|
503 |
C-- Solve elliptic equation(s). |
504 |
C Two-dimensional only for conventional hydrostatic or |
505 |
C three-dimensional for non-hydrostatic and/or IGW scheme. |
506 |
#ifndef ALLOW_OFFLINE |
507 |
IF ( momStepping ) THEN |
508 |
CALL TIMER_START('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid) |
509 |
CALL SOLVE_FOR_PRESSURE(myTime, myIter, myThid) |
510 |
CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid) |
511 |
ENDIF |
512 |
#endif |
513 |
|
514 |
C-- Correct divergence in flow field and cycle time-stepping momentum |
515 |
c IF ( momStepping ) THEN |
516 |
#ifndef ALLOW_OFFLINE |
517 |
CALL TIMER_START('MOM_CORRECTION_STEP [FORWARD_STEP]',myThid) |
518 |
CALL MOMENTUM_CORRECTION_STEP(myTime, myIter, myThid) |
519 |
CALL TIMER_STOP ('MOM_CORRECTION_STEP [FORWARD_STEP]',myThid) |
520 |
#endif |
521 |
c ENDIF |
522 |
|
523 |
#ifdef EXACT_CONSERV |
524 |
IF (exactConserv) THEN |
525 |
C-- Update etaH(n+1) : |
526 |
CALL TIMER_START('UPDATE_ETAH [FORWARD_STEP]',mythid) |
527 |
CALL UPDATE_ETAH( myTime, myIter, myThid ) |
528 |
CALL TIMER_STOP ('UPDATE_ETAH [FORWARD_STEP]',mythid) |
529 |
ENDIF |
530 |
#endif /* EXACT_CONSERV */ |
531 |
|
532 |
#ifdef NONLIN_FRSURF |
533 |
IF ( select_rStar.NE.0 ) THEN |
534 |
C-- r* : compute the future level thickness according to etaH(n+1) |
535 |
CALL TIMER_START('CALC_R_STAR [FORWARD_STEP]',mythid) |
536 |
CALL CALC_R_STAR(etaH, myTime, myIter, myThid ) |
537 |
CALL TIMER_STOP ('CALC_R_STAR [FORWARD_STEP]',mythid) |
538 |
ELSEIF ( nonlinFreeSurf.GT.0) THEN |
539 |
C-- compute the future surface level thickness according to etaH(n+1) |
540 |
CALL TIMER_START('CALC_SURF_DR [FORWARD_STEP]',mythid) |
541 |
CALL CALC_SURF_DR(etaH, myTime, myIter, myThid ) |
542 |
CALL TIMER_STOP ('CALC_SURF_DR [FORWARD_STEP]',mythid) |
543 |
ENDIF |
544 |
#endif /* NONLIN_FRSURF */ |
545 |
|
546 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
547 |
IF ( staggerTimeStep ) THEN |
548 |
C-- do exchanges of U,V (needed for multiDim) when using stagger time-step : |
549 |
#ifdef ALLOW_DEBUG |
550 |
IF ( debugLevel .GE. debLevB ) |
551 |
& CALL DEBUG_CALL('DO_STAGGER_FIELDS_EXCH.',myThid) |
552 |
#endif |
553 |
CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
554 |
CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid ) |
555 |
CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
556 |
|
557 |
C-- State-variables diagnostics |
558 |
IF ( usediagnostics ) THEN |
559 |
CALL TIMER_START('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid) |
560 |
CALL DO_STATEVARS_DIAGS( myTime, 1, myIter, myThid ) |
561 |
CALL TIMER_STOP ('DO_STATEVARS_DIAGS [FORWARD_STEP]',myThid) |
562 |
ENDIF |
563 |
|
564 |
#ifdef ALLOW_DEBUG |
565 |
IF ( debugLevel .GE. debLevB ) |
566 |
& CALL DEBUG_CALL('THERMODYNAMICS',myThid) |
567 |
#endif |
568 |
CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',mythid) |
569 |
CALL THERMODYNAMICS( myTime, myIter, myThid ) |
570 |
CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',mythid) |
571 |
|
572 |
C-- if staggerTimeStep: end |
573 |
ENDIF |
574 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
575 |
|
576 |
#ifdef ALLOW_AUTODIFF_TAMC |
577 |
cph This is needed because convective_adjustment calls |
578 |
cph find_rho which may use pressure() |
579 |
CADJ STORE totphihyd = comlev1, key = ikey_dynamics |
580 |
#endif |
581 |
C-- Cycle time-stepping Tracers arrays (T,S,+pTracers) |
582 |
CALL TIMER_START('TRC_CORRECTION_STEP [FORWARD_STEP]',myThid) |
583 |
CALL TRACERS_CORRECTION_STEP(myTime, myIter, myThid) |
584 |
CALL TIMER_STOP ('TRC_CORRECTION_STEP [FORWARD_STEP]',myThid) |
585 |
|
586 |
#ifdef ALLOW_GCHEM |
587 |
C Add separate timestepping of chemical/biological/forcing |
588 |
C of ptracers here in GCHEM_FORCING_SEP |
589 |
IF ( useGCHEM ) THEN |
590 |
#ifdef ALLOW_DEBUG |
591 |
IF ( debugLevel .GE. debLevB ) |
592 |
& CALL DEBUG_CALL('GCHEM_FORCING_SEP',myThid) |
593 |
#endif /* ALLOW_DEBUG */ |
594 |
CALL TIMER_START('GCHEM_FORCING_SEP [FORWARD_STEP]',myThid) |
595 |
CALL GCHEM_FORCING_SEP( myTime,myIter,myThid ) |
596 |
CALL TIMER_STOP ('GCHEM_FORCING_SEP [FORWARD_STEP]',myThid) |
597 |
ENDIF |
598 |
#endif /* ALLOW_GCHEM */ |
599 |
|
600 |
C-- Do "blocking" sends and receives for tendency "overlap" terms |
601 |
c CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
602 |
c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid ) |
603 |
c CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
604 |
|
605 |
C-- Do "blocking" sends and receives for field "overlap" terms |
606 |
CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
607 |
CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid ) |
608 |
CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
609 |
|
610 |
|
611 |
C AMM |
612 |
#ifdef ALLOW_GRIDALT |
613 |
if (useGRIDALT) then |
614 |
CALL GRIDALT_UPDATE(myThid) |
615 |
endif |
616 |
#endif |
617 |
C AMM |
618 |
|
619 |
C AMM |
620 |
#ifdef ALLOW_FIZHI |
621 |
if( useFIZHI) then |
622 |
CALL TIMER_START('FIZHI [FORWARD_STEP]',mythid) |
623 |
CALL STEP_FIZHI_CORR ( myTime, myIter, myThid ) |
624 |
CALL TIMER_STOP('FIZHI [FORWARD_STEP]',mythid) |
625 |
endif |
626 |
#endif |
627 |
C AMM |
628 |
|
629 |
#ifdef ALLOW_FLT |
630 |
C-- Calculate float trajectories |
631 |
IF (useFLT) THEN |
632 |
CALL TIMER_START('FLOATS [FORWARD_STEP]',myThid) |
633 |
CALL FLT_MAIN(myIter,myTime, myThid) |
634 |
CALL TIMER_STOP ('FLOATS [FORWARD_STEP]',myThid) |
635 |
ENDIF |
636 |
#endif |
637 |
|
638 |
C-- State-variables time-averaging |
639 |
CALL TIMER_START('DO_STATEVARS_TAVE [FORWARD_STEP]',myThid) |
640 |
CALL DO_STATEVARS_TAVE( myTime, myIter, myThid ) |
641 |
CALL TIMER_STOP ('DO_STATEVARS_TAVE [FORWARD_STEP]',myThid) |
642 |
|
643 |
#ifndef ALLOW_OFFLINE |
644 |
#ifdef ALLOW_MONITOR |
645 |
C-- Check status of solution (statistics, cfl, etc...) |
646 |
CALL TIMER_START('MONITOR [FORWARD_STEP]',myThid) |
647 |
CALL MONITOR( myIter, myTime, myThid ) |
648 |
CALL TIMER_STOP ('MONITOR [FORWARD_STEP]',myThid) |
649 |
#endif /* ALLOW_MONITOR */ |
650 |
#endif |
651 |
|
652 |
#ifdef ALLOW_COST |
653 |
C-- compare model with data and compute cost function |
654 |
C-- this is done after exchanges to allow interpolation |
655 |
CALL TIMER_START('COST_TILE [FORWARD_STEP]',myThid) |
656 |
CALL COST_TILE ( mytime, myiter, myThid ) |
657 |
CALL TIMER_STOP ('COST_TILE [FORWARD_STEP]',myThid) |
658 |
#endif |
659 |
|
660 |
C-- Do IO if needed. |
661 |
#ifdef ALLOW_OFFLINE |
662 |
CALL TIMER_START('OFFLINE_MODEL_IO [FORWARD_STEP]',myThid) |
663 |
CALL OFFLINE_MODEL_IO( myTime, myIter, myThid ) |
664 |
CALL TIMER_STOP ('OFFLINE_MODEL_IO [FORWARD_STEP]',myThid) |
665 |
#else |
666 |
CALL TIMER_START('DO_THE_MODEL_IO [FORWARD_STEP]',myThid) |
667 |
CALL DO_THE_MODEL_IO( myTime, myIter, myThid ) |
668 |
CALL TIMER_STOP ('DO_THE_MODEL_IO [FORWARD_STEP]',myThid) |
669 |
#endif |
670 |
|
671 |
C-- Save state for restarts |
672 |
CALL TIMER_START('WRITE_CHECKPOINT [FORWARD_STEP]',myThid) |
673 |
CALL PACKAGES_WRITE_PICKUP( |
674 |
I .FALSE., myTime, myIter, myThid ) |
675 |
#ifndef ALLOW_OFFLINE |
676 |
CALL WRITE_CHECKPOINT( |
677 |
I .FALSE., myTime, myIter, myThid ) |
678 |
#endif |
679 |
CALL TIMER_STOP ('WRITE_CHECKPOINT [FORWARD_STEP]',myThid) |
680 |
|
681 |
#ifdef ALLOW_DEBUG |
682 |
IF ( debugLevel .GE. debLevB ) |
683 |
& CALL DEBUG_LEAVE('FORWARD_STEP',myThid) |
684 |
#endif |
685 |
|
686 |
RETURN |
687 |
END |