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