/[MITgcm]/MITgcm/model/src/forward_step.F
ViewVC logotype

Annotation of /MITgcm/model/src/forward_step.F

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


Revision 1.106 - (hide annotations) (download)
Tue Nov 16 05:42:11 2004 UTC (19 years, 6 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint56
Changes since 1.105: +15 -15 lines
More on dsvd vs. MITgcm interfacing
o handling of g_, ad, via admtlm_vector (mds...vector)
o use ctrl_pack/unpack for admtlm_vector I/O
o use optimcycle for dsvd iteration
o make sure norm is w.r.t. derived quantities

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

  ViewVC Help
Powered by ViewVC 1.1.22