1 |
C $Header: /u/gcmpack/MITgcm/model/src/forward_step.F,v 1.92 2004/05/14 21:08:28 heimbach Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "PACKAGES_CONFIG.h" |
5 |
#include "CPP_OPTIONS.h" |
6 |
|
7 |
cswdptr -- add -- |
8 |
#ifdef ALLOW_GCHEM |
9 |
# include "GCHEM_OPTIONS.h" |
10 |
#endif |
11 |
cswdptr -- end add --- |
12 |
|
13 |
CBOP |
14 |
C !ROUTINE: FORWARD_STEP |
15 |
C !INTERFACE: |
16 |
SUBROUTINE FORWARD_STEP( iloop, myTime, myIter, myThid ) |
17 |
|
18 |
C !DESCRIPTION: \bv |
19 |
C *================================================================== |
20 |
C | SUBROUTINE forward_step |
21 |
C | o Run the ocean model and, optionally, evaluate a cost function. |
22 |
C *================================================================== |
23 |
C | |
24 |
C | THE_MAIN_LOOP is the toplevel routine for the Tangent Linear and |
25 |
C | Adjoint Model Compiler (TAMC). For this purpose the initialization |
26 |
C | of the model was split into two parts. Those parameters that do |
27 |
C | not depend on a specific model run are set in INITIALISE_FIXED, |
28 |
C | whereas those that do depend on the specific realization are |
29 |
C | initialized in INITIALISE_VARIA. |
30 |
C | |
31 |
C *================================================================== |
32 |
C \ev |
33 |
|
34 |
C !USES: |
35 |
IMPLICIT NONE |
36 |
C == Global variables == |
37 |
#include "SIZE.h" |
38 |
#include "EEPARAMS.h" |
39 |
#include "PARAMS.h" |
40 |
#include "DYNVARS.h" |
41 |
|
42 |
#ifdef ALLOW_SHAP_FILT |
43 |
# include "SHAP_FILT.h" |
44 |
#endif |
45 |
#ifdef ALLOW_ZONAL_FILT |
46 |
# include "ZONAL_FILT.h" |
47 |
#endif |
48 |
#ifdef COMPONENT_MODULE |
49 |
# include "CPL_PARAMS.h" |
50 |
#endif |
51 |
|
52 |
#ifdef ALLOW_AUTODIFF_TAMC |
53 |
# include "FFIELDS.h" |
54 |
|
55 |
# ifdef ALLOW_NONHYDROSTATIC |
56 |
# include "CG3D.h" |
57 |
# endif |
58 |
|
59 |
# include "tamc.h" |
60 |
# include "ctrl.h" |
61 |
# include "ctrl_dummy.h" |
62 |
# include "cost.h" |
63 |
# include "EOS.h" |
64 |
# ifdef ALLOW_EXF |
65 |
# include "exf_fields.h" |
66 |
# include "exf_clim_fields.h" |
67 |
# ifdef ALLOW_BULKFORMULAE |
68 |
# include "exf_constants.h" |
69 |
# endif |
70 |
# endif |
71 |
# ifdef ALLOW_OBCS |
72 |
# include "OBCS.h" |
73 |
# endif |
74 |
# ifdef ALLOW_PTRACERS |
75 |
# include "PTRACERS.h" |
76 |
# endif |
77 |
# ifdef ALLOW_CD_CODE |
78 |
# include "CD_CODE_VARS.h" |
79 |
# endif |
80 |
# ifdef ALLOW_EBM |
81 |
# include "EBM.h" |
82 |
# endif |
83 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
84 |
|
85 |
C !LOCAL VARIABLES: |
86 |
C == Routine arguments == |
87 |
C note: under the multi-threaded model myiter and |
88 |
C mytime are local variables passed around as routine |
89 |
C arguments. Although this is fiddly it saves the need to |
90 |
C impose additional synchronisation points when they are |
91 |
C updated. |
92 |
C myIter - iteration counter for this thread |
93 |
C myTime - time counter for this thread |
94 |
C myThid - thread number for this instance of the routine. |
95 |
INTEGER iloop |
96 |
INTEGER myThid |
97 |
INTEGER myIter |
98 |
_RL myTime |
99 |
|
100 |
C == Local variables == |
101 |
integer i,L |
102 |
INTEGER myItP1 |
103 |
CEOP |
104 |
|
105 |
#ifdef ALLOW_DEBUG |
106 |
IF ( debugLevel .GE. debLevB ) |
107 |
& CALL DEBUG_ENTER('FORWARD_STEP',myThid) |
108 |
#endif |
109 |
|
110 |
#ifdef ALLOW_AUTODIFF_TAMC |
111 |
C-- Reset the model iteration counter and the model time. |
112 |
myiter = nIter0 + (iloop-1) |
113 |
mytime = startTime + float(iloop-1)*deltaTclock |
114 |
#endif |
115 |
|
116 |
#if (defined (ALLOW_AUTODIFF_TAMC) && defined (ALLOW_AUTODIFF_MONITOR)) |
117 |
C Include call to a dummy routine. Its adjoint will be |
118 |
C called at the proper place in the adjoint code. |
119 |
C The adjoint routine will print out adjoint values |
120 |
C if requested. The location of the call is important, |
121 |
C it has to be after the adjoint of the exchanges |
122 |
C (DO_GTERM_BLOCKING_EXCHANGES). |
123 |
CALL DUMMY_IN_STEPPING( myTime, myIter, myThid ) |
124 |
cph I've commented this line since it may conflict with MITgcm's adjoint |
125 |
cph However, need to check whether that's still consistent |
126 |
cph with the ecco-branch (it should). |
127 |
cph CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid ) |
128 |
#endif |
129 |
|
130 |
#ifdef EXACT_CONSERV |
131 |
IF (exactConserv) THEN |
132 |
C-- Update etaH(n+1) : |
133 |
CALL TIMER_START('UPDATE_ETAH [FORWARD_STEP]',mythid) |
134 |
CALL UPDATE_ETAH( myTime, myIter, myThid ) |
135 |
CALL TIMER_STOP ('UPDATE_ETAH [FORWARD_STEP]',mythid) |
136 |
ENDIF |
137 |
#endif /* EXACT_CONSERV */ |
138 |
|
139 |
#ifdef NONLIN_FRSURF |
140 |
IF ( select_rStar.NE.0 ) THEN |
141 |
C-- r* : compute the future level thickness according to etaH(n+1) |
142 |
CALL TIMER_START('CALC_R_STAR [FORWARD_STEP]',mythid) |
143 |
CALL CALC_R_STAR(etaH, myTime, myIter, myThid ) |
144 |
CALL TIMER_STOP ('CALC_R_STAR [FORWARD_STEP]',mythid) |
145 |
ELSEIF ( nonlinFreeSurf.GT.0) THEN |
146 |
C-- compute the future surface level thickness according to etaH(n+1) |
147 |
CALL TIMER_START('CALC_SURF_DR [FORWARD_STEP]',mythid) |
148 |
CALL CALC_SURF_DR(etaH, myTime, myIter, myThid ) |
149 |
CALL TIMER_STOP ('CALC_SURF_DR [FORWARD_STEP]',mythid) |
150 |
ENDIF |
151 |
#endif /* NONLIN_FRSURF */ |
152 |
|
153 |
#ifdef ALLOW_AUTODIFF_TAMC |
154 |
c************************************** |
155 |
#include "checkpoint_lev1_directives.h" |
156 |
c************************************** |
157 |
#endif |
158 |
|
159 |
C-- Call external forcing package |
160 |
#ifdef ALLOW_BULK_FORCE |
161 |
IF ( useBulkForce ) THEN |
162 |
#ifdef ALLOW_DEBUG |
163 |
IF ( debugLevel .GE. debLevB ) |
164 |
& CALL DEBUG_CALL('BULKF_FIELDS_LOAD',myThid) |
165 |
#endif |
166 |
CALL TIMER_START('BULKF_FORCING [FORWARD_STEP]',mythid) |
167 |
C- load all forcing fields at current time |
168 |
CALL BULKF_FIELDS_LOAD( myTime, myIter, myThid ) |
169 |
C- calculate qnet and empmr (and wind stress) |
170 |
CALL BULKF_FORCING( myTime, myIter, myThid ) |
171 |
CALL TIMER_STOP ('BULKF_FORCING [FORWARD_STEP]',mythid) |
172 |
ELSE |
173 |
#endif /* ALLOW_BULK_FORCE */ |
174 |
|
175 |
# ifdef ALLOW_EXF |
176 |
# ifdef ALLOW_DEBUG |
177 |
IF ( debugLevel .GE. debLevB ) |
178 |
& CALL DEBUG_CALL('EXF_GETFORCING',myThid) |
179 |
# endif |
180 |
CALL TIMER_START('EXF_GETFORCING [FORWARD_STEP]',mythid) |
181 |
CALL EXF_GETFORCING( mytime, myiter, mythid ) |
182 |
CALL TIMER_STOP ('EXF_GETFORCING [FORWARD_STEP]',mythid) |
183 |
# else /* ALLOW_EXF undef */ |
184 |
cph The following IF-statement creates an additional dependency |
185 |
cph for the forcing fields requiring additional storing. |
186 |
cph Therefore, the IF-statement will be put between CPP-OPTIONS, |
187 |
cph assuming that ALLOW_SEAICE has not yet been differentiated. |
188 |
# if (defined (ALLOW_SEAICE) || defined (ALLOW_EBM)) |
189 |
IF ( .NOT. useSEAICE .AND. .NOT. useEBM ) THEN |
190 |
# endif |
191 |
#ifdef ALLOW_DEBUG |
192 |
IF ( debugLevel .GE. debLevB ) |
193 |
& CALL DEBUG_CALL('EXTERNAL_FIELDS_LOAD',myThid) |
194 |
#endif |
195 |
CALL TIMER_START('EXTERNAL_FIELDS_LOAD[FORWARD_STEP]',mythid) |
196 |
CALL EXTERNAL_FIELDS_LOAD( mytime, myiter, mythid ) |
197 |
CALL TIMER_STOP ('EXTERNAL_FIELDS_LOAD[FORWARD_STEP]',mythid) |
198 |
# if (defined (ALLOW_SEAICE) || defined (ALLOW_EBM)) |
199 |
ENDIF |
200 |
# endif |
201 |
# endif /* ALLOW_EXF */ |
202 |
#ifdef ALLOW_BULK_FORCE |
203 |
C-- end of if/else block useBulfforce -- |
204 |
ENDIF |
205 |
#endif /* ALLOW_BULK_FORCE */ |
206 |
|
207 |
#ifdef ALLOW_AUTODIFF |
208 |
c-- Add control vector for forcing and parameter fields |
209 |
if ( myiter .EQ. nIter0 ) |
210 |
& CALL CTRL_MAP_FORCING (mythid) |
211 |
#endif |
212 |
|
213 |
#ifdef ALLOW_THSICE |
214 |
IF ( useThSIce .AND. buoyancyRelation(1:7) .EQ. 'OCEANIC' ) THEN |
215 |
#ifdef ALLOW_DEBUG |
216 |
IF ( debugLevel .GE. debLevB ) |
217 |
& CALL DEBUG_CALL('THSICE_MAIN',myThid) |
218 |
#endif |
219 |
C-- Step forward Therm.Sea-Ice variables |
220 |
C and modify forcing terms including effects from ice |
221 |
CALL TIMER_START('THSICE_MAIN [FORWARD_STEP]', myThid) |
222 |
CALL THSICE_MAIN( myTime, myIter, myThid ) |
223 |
CALL TIMER_STOP( 'THSICE_MAIN [FORWARD_STEP]', myThid) |
224 |
ENDIF |
225 |
#endif /* ALLOW_THSICE */ |
226 |
|
227 |
# ifdef ALLOW_SEAICE |
228 |
C-- Call sea ice model to compute forcing/external data fields. In |
229 |
C addition to computing prognostic sea-ice variables and diagnosing the |
230 |
C forcing/external data fields that drive the ocean model, SEAICE_MODEL |
231 |
C also sets theta to the freezing point under sea-ice. The implied |
232 |
C surface heat flux is then stored in variable surfaceTendencyTice, |
233 |
C which is needed by KPP package (kpp_calc.F and kpp_transport_t.F) |
234 |
C to diagnose surface buoyancy fluxes and for the non-local transport |
235 |
C term. Because this call precedes model thermodynamics, temperature |
236 |
C under sea-ice may not be "exactly" at the freezing point by the time |
237 |
C theta is dumped or time-averaged. |
238 |
IF ( useSEAICE ) THEN |
239 |
#ifdef ALLOW_DEBUG |
240 |
IF ( debugLevel .GE. debLevB ) |
241 |
& CALL DEBUG_CALL('SEAICE_MODEL',myThid) |
242 |
#endif |
243 |
CALL TIMER_START('SEAICE_MODEL [FORWARD_STEP]',myThid) |
244 |
CALL SEAICE_MODEL( myTime, myIter, myThid ) |
245 |
CALL TIMER_STOP ('SEAICE_MODEL [FORWARD_STEP]',myThid) |
246 |
ENDIF |
247 |
# endif /* ALLOW_SEAICE */ |
248 |
|
249 |
C-- Freeze water at the surface |
250 |
#ifdef ALLOW_AUTODIFF_TAMC |
251 |
CADJ STORE theta = comlev1, key = ikey_dynamics |
252 |
#endif |
253 |
IF ( allowFreezing .AND. .NOT. useSEAICE |
254 |
& .AND. .NOT. useThSIce ) THEN |
255 |
CALL FREEZE_SURFACE( myTime, myIter, myThid ) |
256 |
ENDIF |
257 |
|
258 |
#ifdef ALLOW_AUTODIFF_TAMC |
259 |
# ifdef ALLOW_PTRACERS |
260 |
cph this replaces _bibj storing of ptracer within thermodynamics |
261 |
CADJ STORE ptracer = comlev1, key = ikey_dynamics |
262 |
# endif |
263 |
#endif |
264 |
|
265 |
#ifdef ALLOW_PTRACERS |
266 |
# ifdef ALLOW_GCHEM |
267 |
CALL GCHEM_FIELDS_LOAD( mytime, myiter, mythid ) |
268 |
# endif |
269 |
#endif |
270 |
|
271 |
#ifdef COMPONENT_MODULE |
272 |
IF ( useCoupler .AND. cpl_earlyExpImpCall ) THEN |
273 |
C Post coupling data that I export. |
274 |
C Read in coupling data that I import. |
275 |
CALL TIMER_START('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid) |
276 |
CALL CPL_EXPORT_MY_DATA( myIter, myTime, myThid ) |
277 |
CALL CPL_IMPORT_EXTERNAL_DATA( myIter, myTime, myThid ) |
278 |
CALL TIMER_STOP ('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid) |
279 |
ENDIF |
280 |
# ifndef ALLOW_AIM |
281 |
C jmc: don't know precisely where to put this call. leave it here for now. |
282 |
IF ( useCoupler ) THEN |
283 |
CALL OCN_APPLY_IMPORT( .TRUE., myTime, myIter, myThid ) |
284 |
ENDIF |
285 |
# endif |
286 |
#endif /* COMPONENT_MODULE */ |
287 |
|
288 |
C-- Step forward fields and calculate time tendency terms. |
289 |
|
290 |
C AMM |
291 |
#ifdef ALLOW_FIZHI |
292 |
if( useFIZHI) then |
293 |
CALL TIMER_START('FIZHI [FORWARD_STEP]',mythid) |
294 |
CALL UPDATE_OCEAN_EXPORTS ( myTime, myIter, myThid ) |
295 |
CALL UPDATE_EARTH_EXPORTS ( myTime, myIter, myThid ) |
296 |
CALL UPDATE_CHEMISTRY_EXPORTS ( myTime, myIter, myThid ) |
297 |
CALL FIZHI_WRAPPER ( myTime, myIter, myThid ) |
298 |
CALL STEP_FIZHI_FG ( myTime, myIter, myThid, deltaTtracer ) |
299 |
CALL TIMER_STOP ('FIZHI [FORWARD_STEP]',mythid) |
300 |
endif |
301 |
#endif |
302 |
C AMM |
303 |
|
304 |
#ifdef ALLOW_EBM |
305 |
IF ( useEBM ) THEN |
306 |
# ifdef ALLOW_DEBUG |
307 |
IF ( debugLevel .GE. debLevB ) |
308 |
& CALL DEBUG_CALL('EBM',myThid) |
309 |
# endif |
310 |
CALL TIMER_START('EBM [FORWARD_STEP]',mythid) |
311 |
CALL EBM_DRIVER ( myTime, myIter, myThid ) |
312 |
CALL TIMER_STOP ('EBM [FORWARD_STEP]',mythid) |
313 |
ENDIF |
314 |
#endif |
315 |
|
316 |
#ifdef ALLOW_DEBUG |
317 |
IF ( debugLevel .GE. debLevB ) |
318 |
& CALL DEBUG_CALL('THERMODYNAMICS',myThid) |
319 |
#endif |
320 |
CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',mythid) |
321 |
CALL THERMODYNAMICS( myTime, myIter, myThid ) |
322 |
CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',mythid) |
323 |
|
324 |
C-- do exchanges (needed for DYNAMICS) when using stagger time-step : |
325 |
#ifdef ALLOW_DEBUG |
326 |
IF ( debugLevel .GE. debLevB ) |
327 |
& CALL DEBUG_CALL('DO_STAGGER_FIELDS_EXCH.',myThid) |
328 |
#endif |
329 |
CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
330 |
CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid ) |
331 |
CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
332 |
|
333 |
#ifdef ALLOW_SHAP_FILT |
334 |
IF (useSHAP_FILT .AND. |
335 |
& staggerTimeStep .AND. shap_filt_TrStagg ) THEN |
336 |
#ifdef ALLOW_DEBUG |
337 |
IF ( debugLevel .GE. debLevB ) |
338 |
& CALL DEBUG_CALL('SHAP_FILT_APPLY_TS',myThid) |
339 |
#endif |
340 |
CALL TIMER_START('SHAP_FILT [FORWARD_STEP]',myThid) |
341 |
CALL SHAP_FILT_APPLY_TS(gT,gS,myTime+deltaT,myIter+1,myThid) |
342 |
CALL TIMER_STOP ('SHAP_FILT [FORWARD_STEP]',myThid) |
343 |
ENDIF |
344 |
#endif |
345 |
|
346 |
#ifdef ALLOW_ZONAL_FILT |
347 |
IF (useZONAL_FILT .AND. |
348 |
& staggerTimeStep .AND. zonal_filt_TrStagg ) THEN |
349 |
#ifdef ALLOW_DEBUG |
350 |
IF ( debugLevel .GE. debLevB ) |
351 |
& CALL DEBUG_CALL('ZONAL_FILT_APPLY_TS',myThid) |
352 |
#endif |
353 |
CALL TIMER_START('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid) |
354 |
CALL ZONAL_FILT_APPLY_TS( gT, gS, myThid ) |
355 |
CALL TIMER_STOP ('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid) |
356 |
ENDIF |
357 |
#endif |
358 |
|
359 |
#ifdef COMPONENT_MODULE |
360 |
IF ( useCoupler .AND. .NOT.cpl_earlyExpImpCall ) THEN |
361 |
C Post coupling data that I export. |
362 |
C Read in coupling data that I import. |
363 |
myItP1 = myIter + 1 |
364 |
CALL TIMER_START('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid) |
365 |
CALL CPL_EXPORT_MY_DATA( myItP1, myTime, myThid ) |
366 |
CALL CPL_IMPORT_EXTERNAL_DATA( myItP1, myTime, myThid ) |
367 |
CALL TIMER_STOP ('CPL_EXPORT-IMPORT [FORWARD_STEP]',myThid) |
368 |
# ifndef ALLOW_AIM |
369 |
IF ( useRealFreshWaterFlux ) THEN |
370 |
CALL OCN_APPLY_IMPORT( .FALSE., myTime, myIter, myThid ) |
371 |
ENDIF |
372 |
# endif |
373 |
ENDIF |
374 |
#endif /* COMPONENT_MODULE */ |
375 |
|
376 |
C-- Step forward fields and calculate time tendency terms. |
377 |
#ifndef ALLOW_AUTODIFF_TAMC |
378 |
IF ( momStepping ) THEN |
379 |
#endif |
380 |
#ifdef ALLOW_DEBUG |
381 |
IF ( debugLevel .GE. debLevB ) |
382 |
& CALL DEBUG_CALL('DYNAMICS',myThid) |
383 |
#endif |
384 |
CALL TIMER_START('DYNAMICS [FORWARD_STEP]',mythid) |
385 |
CALL DYNAMICS( myTime, myIter, myThid ) |
386 |
CALL TIMER_STOP ('DYNAMICS [FORWARD_STEP]',mythid) |
387 |
#ifndef ALLOW_AUTODIFF_TAMC |
388 |
ENDIF |
389 |
#endif |
390 |
|
391 |
#ifdef ALLOW_NONHYDROSTATIC |
392 |
C-- Step forward W field in N-H algorithm |
393 |
IF ( momStepping .AND. nonHydrostatic ) THEN |
394 |
#ifdef ALLOW_DEBUG |
395 |
IF ( debugLevel .GE. debLevB ) |
396 |
& CALL DEBUG_CALL('CALC_GW',myThid) |
397 |
#endif |
398 |
CALL TIMER_START('CALC_GW [FORWARD_STEP]',myThid) |
399 |
CALL CALC_GW(myThid) |
400 |
CALL TIMER_STOP ('CALC_GW [FORWARD_STEP]',myThid) |
401 |
ENDIF |
402 |
#endif |
403 |
|
404 |
#ifdef NONLIN_FRSURF |
405 |
C-- update hfacC,W,S and recip_hFac according to etaH(n+1) : |
406 |
IF ( nonlinFreeSurf.GT.0) THEN |
407 |
IF ( select_rStar.GT.0 ) THEN |
408 |
CALL TIMER_START('UPDATE_R_STAR [FORWARD_STEP]',myThid) |
409 |
CALL UPDATE_R_STAR( myTime, myIter, myThid ) |
410 |
CALL TIMER_STOP ('UPDATE_R_STAR [FORWARD_STEP]',myThid) |
411 |
ELSE |
412 |
CALL TIMER_START('UPDATE_SURF_DR [FORWARD_STEP]',myThid) |
413 |
CALL UPDATE_SURF_DR( myTime, myIter, myThid ) |
414 |
CALL TIMER_STOP ('UPDATE_SURF_DR [FORWARD_STEP]',myThid) |
415 |
ENDIF |
416 |
ENDIF |
417 |
C- update also CG2D matrix (and preconditioner) |
418 |
IF ( momStepping .AND. nonlinFreeSurf.GT.2 ) THEN |
419 |
CALL TIMER_START('UPDATE_CG2D [FORWARD_STEP]',myThid) |
420 |
CALL UPDATE_CG2D( myTime, myIter, myThid ) |
421 |
CALL TIMER_STOP ('UPDATE_CG2D [FORWARD_STEP]',myThid) |
422 |
ENDIF |
423 |
#endif |
424 |
|
425 |
C-- Apply Filters to u*,v* before SOLVE_FOR_PRESSURE |
426 |
#ifdef ALLOW_SHAP_FILT |
427 |
IF (useSHAP_FILT .AND. shap_filt_uvStar) THEN |
428 |
CALL TIMER_START('SHAP_FILT [FORWARD_STEP]',myThid) |
429 |
IF (implicDiv2Dflow.LT.1.) THEN |
430 |
C-- Explicit+Implicit part of the Barotropic Flow Divergence |
431 |
C => Filtering of uVel,vVel is necessary |
432 |
CALL SHAP_FILT_APPLY_UV( uVel,vVel, |
433 |
& myTime+deltaT, myIter+1, myThid ) |
434 |
ENDIF |
435 |
CALL SHAP_FILT_APPLY_UV( gU,gV,myTime+deltaT,myIter+1,myThid) |
436 |
CALL TIMER_STOP ('SHAP_FILT [FORWARD_STEP]',myThid) |
437 |
ENDIF |
438 |
#endif |
439 |
#ifdef ALLOW_ZONAL_FILT |
440 |
IF (useZONAL_FILT .AND. zonal_filt_uvStar) THEN |
441 |
CALL TIMER_START('ZONAL_FILT_APPLY [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 ZONAL_FILT_APPLY_UV( uVel, vVel, myThid ) |
446 |
ENDIF |
447 |
CALL ZONAL_FILT_APPLY_UV( gU, gV, myThid ) |
448 |
CALL TIMER_STOP ('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid) |
449 |
ENDIF |
450 |
#endif |
451 |
|
452 |
C-- Solve elliptic equation(s). |
453 |
C Two-dimensional only for conventional hydrostatic or |
454 |
C three-dimensional for non-hydrostatic and/or IGW scheme. |
455 |
IF ( momStepping ) THEN |
456 |
CALL TIMER_START('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid) |
457 |
CALL SOLVE_FOR_PRESSURE(myTime, myIter, myThid) |
458 |
CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid) |
459 |
ENDIF |
460 |
|
461 |
#ifdef ALLOW_AUTODIFF_TAMC |
462 |
cph This is needed because convective_adjustment calls |
463 |
cph find_rho which may use pressure() |
464 |
CADJ STORE totphihyd = comlev1, key = ikey_dynamics |
465 |
#endif |
466 |
C-- Correct divergence in flow field and cycle time-stepping |
467 |
C arrays (for all fields) ; update time-counter |
468 |
myIter = nIter0 + iLoop |
469 |
myTime = startTime + deltaTClock * float(iLoop) |
470 |
CALL TIMER_START('THE_CORRECTION_STEP [FORWARD_STEP]',myThid) |
471 |
CALL THE_CORRECTION_STEP(myTime, myIter, myThid) |
472 |
CALL TIMER_STOP ('THE_CORRECTION_STEP [FORWARD_STEP]',myThid) |
473 |
|
474 |
C-- Do "blocking" sends and receives for tendency "overlap" terms |
475 |
c CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
476 |
c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid ) |
477 |
c CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
478 |
|
479 |
C-- Do "blocking" sends and receives for field "overlap" terms |
480 |
CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
481 |
CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid ) |
482 |
CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid) |
483 |
|
484 |
cswdptr -- add for seperate timestepping of chemical/biological/forcing |
485 |
cswdptr of ptracers --- |
486 |
#ifdef ALLOW_GCHEM |
487 |
ceh3 This is broken -- this ifdef should not be visible! |
488 |
#ifdef PTRACERS_SEPARATE_FORCING |
489 |
ceh3 needs an IF ( use GCHEM ) THEN |
490 |
call GCHEM_FORCING_SEP( myTime,myIter,myThid ) |
491 |
#endif /* PTRACERS_SEPARATE_FORCING */ |
492 |
#endif /* ALLOW_GCHEM */ |
493 |
cswdptr -- end add --- |
494 |
|
495 |
C AMM |
496 |
#ifdef ALLOW_GRIDALT |
497 |
if (useGRIDALT) then |
498 |
CALL GRIDALT_UPDATE(myThid) |
499 |
endif |
500 |
#endif |
501 |
C AMM |
502 |
|
503 |
C AMM |
504 |
#ifdef ALLOW_FIZHI |
505 |
if( useFIZHI) then |
506 |
CALL TIMER_START('FIZHI [FORWARD_STEP]',mythid) |
507 |
CALL STEP_FIZHI_CORR ( myTime, myIter, myThid ) |
508 |
CALL TIMER_STOP('FIZHI [FORWARD_STEP]',mythid) |
509 |
endif |
510 |
#endif |
511 |
C AMM |
512 |
|
513 |
#ifdef ALLOW_FLT |
514 |
C-- Calculate float trajectories |
515 |
IF (useFLT) THEN |
516 |
CALL TIMER_START('FLOATS [FORWARD_STEP]',myThid) |
517 |
CALL FLT_MAIN(myIter,myTime, myThid) |
518 |
CALL TIMER_STOP ('FLOATS [FORWARD_STEP]',myThid) |
519 |
ENDIF |
520 |
#endif |
521 |
|
522 |
#ifdef ALLOW_MONITOR |
523 |
C-- Check status of solution (statistics, cfl, etc...) |
524 |
CALL TIMER_START('MONITOR [FORWARD_STEP]',myThid) |
525 |
CALL MONITOR( myIter, myTime, myThid ) |
526 |
CALL TIMER_STOP ('MONITOR [FORWARD_STEP]',myThid) |
527 |
#endif /* ALLOW_MONITOR */ |
528 |
|
529 |
C AMM -- Diagnostics package |
530 |
#ifdef ALLOW_DIAGNOSTICS |
531 |
|
532 |
if(usediagnostics)then |
533 |
call diagnostics_fill_state(myThid) |
534 |
endif |
535 |
|
536 |
# ifdef ALLOW_PTRACERS |
537 |
if(useptracers)then |
538 |
call diagnostics_fill_ptracers(myThid) |
539 |
endif |
540 |
# endif |
541 |
|
542 |
#endif |
543 |
C AMM |
544 |
|
545 |
C-- Do IO if needed. |
546 |
CALL TIMER_START('DO_THE_MODEL_IO [FORWARD_STEP]',myThid) |
547 |
CALL DO_THE_MODEL_IO( myTime, myIter, myThid ) |
548 |
CALL TIMER_STOP ('DO_THE_MODEL_IO [FORWARD_STEP]',myThid) |
549 |
|
550 |
C-- Save state for restarts |
551 |
C Note: (jmc: is it still the case after ckp35 ?) |
552 |
C ===== |
553 |
C Because of the ordering of the timestepping code and |
554 |
C tendency term code at end of loop model arrays hold |
555 |
C U,V,T,S at "time-level" N but gu, gv, gs, gt, guNM1,... |
556 |
C at "time-level" N+1/2 (guNM1 at "time-level" N+1/2 is |
557 |
C gu at "time-level" N-1/2) and etaN at "time-level" N+1/2. |
558 |
C where N = I+timeLevBase-1 |
559 |
C Thus a checkpoint contains U.0000000000, GU.0000000001 and |
560 |
C etaN.0000000001 in the indexing scheme used for the model |
561 |
C "state" files. This example is referred to as a checkpoint |
562 |
C at time level 1 |
563 |
CALL TIMER_START('WRITE_CHECKPOINT [FORWARD_STEP]',myThid) |
564 |
CALL PACKAGES_WRITE_PICKUP( |
565 |
I .FALSE., myTime, myIter, myThid ) |
566 |
CALL WRITE_CHECKPOINT( |
567 |
I .FALSE., myTime, myIter, myThid ) |
568 |
CALL TIMER_STOP ('WRITE_CHECKPOINT [FORWARD_STEP]',myThid) |
569 |
|
570 |
#ifdef ALLOW_DEBUG |
571 |
IF ( debugLevel .GE. debLevB ) |
572 |
& CALL DEBUG_LEAVE('FORWARD_STEP',myThid) |
573 |
#endif |
574 |
|
575 |
RETURN |
576 |
END |