/[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.44 - (hide annotations) (download)
Sat Dec 21 00:33:54 2002 UTC (21 years, 5 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint47e_post, branch-exfmods-tag
Branch point for: branch-exfmods-curt
Changes since 1.43: +2 -2 lines
add a call to TIMER_STOP that was missing.

1 jmc 1.44 C $Header: /u/gcmpack/MITgcm/model/src/forward_step.F,v 1.43 2002/12/06 20:10:58 cheisey Exp $
2 adcroft 1.15 C $Name: $
3 adcroft 1.1
4     #include "CPP_OPTIONS.h"
5    
6 cnh 1.22 CBOP
7     C !ROUTINE: FORWARD_STEP
8     C !INTERFACE:
9 adcroft 1.13 SUBROUTINE FORWARD_STEP( iloop, myTime, myIter, myThid )
10 heimbach 1.12
11 cnh 1.22 C !DESCRIPTION: \bv
12     C *==================================================================
13     C | SUBROUTINE forward_step
14     C | o Run the ocean model and, optionally, evaluate a cost function.
15     C *==================================================================
16     C |
17     C | THE_MAIN_LOOP is the toplevel routine for the Tangent Linear and
18     C | Adjoint Model Compiler (TAMC). For this purpose the initialization
19     C | of the model was split into two parts. Those parameters that do
20     C | not depend on a specific model run are set in INITIALISE_FIXED,
21     C | whereas those that do depend on the specific realization are
22     C | initialized in INITIALISE_VARIA.
23     C |
24     C *==================================================================
25     C \ev
26 adcroft 1.1
27 cnh 1.22 C !USES:
28     IMPLICIT NONE
29     C == Global variables ==
30 adcroft 1.1 #include "SIZE.h"
31     #include "EEPARAMS.h"
32     #include "PARAMS.h"
33     #include "DYNVARS.h"
34 heimbach 1.12 #include "FFIELDS.h"
35    
36 adcroft 1.1 #ifdef ALLOW_NONHYDROSTATIC
37     #include "CG3D.h"
38     #endif
39    
40 jmc 1.27 #ifdef ALLOW_SHAP_FILT
41     #include "SHAP_FILT.h"
42     #endif
43     #ifdef ALLOW_ZONAL_FILT
44     #include "ZONAL_FILT.h"
45     #endif
46    
47 heimbach 1.12 #ifdef ALLOW_AUTODIFF_TAMC
48 heimbach 1.30 # include "tamc.h"
49     # include "ctrl.h"
50     # include "ctrl_dummy.h"
51     # include "cost.h"
52 heimbach 1.37 # include "EOS.h"
53 heimbach 1.30 # ifdef INCLUDE_EXTERNAL_FORCING_PACKAGE
54     # include "exf_fields.h"
55 cheisey 1.40 # ifdef ALLOW_BULK_FORCE
56 heimbach 1.30 # include "exf_constants.h"
57     # endif
58     # endif
59     # ifdef ALLOW_OBCS
60     # include "OBCS.h"
61     # endif
62 heimbach 1.12 #endif
63    
64 cnh 1.22 C !LOCAL VARIABLES:
65     C == Routine arguments ==
66 adcroft 1.13 C note: under the multi-threaded model myiter and
67     C mytime are local variables passed around as routine
68     C arguments. Although this is fiddly it saves the need to
69     C impose additional synchronisation points when they are
70     C updated.
71     C myiter - iteration counter for this thread
72     C mytime - time counter for this thread
73     C mythid - thread number for this instance of the routine.
74 heimbach 1.12 integer iloop
75     integer mythid
76     integer myiter
77     _RL mytime
78 jmc 1.21 INTEGER bi,bj
79 heimbach 1.12
80 cnh 1.22 CEOP
81 heimbach 1.12
82    
83     #ifdef ALLOW_AUTODIFF_TAMC
84 heimbach 1.14 C-- Reset the model iteration counter and the model time.
85 heimbach 1.12 myiter = nIter0 + (iloop-1)
86     mytime = startTime + float(iloop-1)*deltaTclock
87     #endif
88    
89 heimbach 1.32 #if (defined (ALLOW_AUTODIFF_TAMC) && defined (ALLOW_AUTODIFF_MONITOR))
90 heimbach 1.16 C Include call to a dummy routine. Its adjoint will be
91     C called at the proper place in the adjoint code.
92     C The adjoint routine will print out adjoint values
93     C if requested. The location of the call is important,
94     C it has to be after the adjoint of the exchanges
95     C (DO_GTERM_BLOCKING_EXCHANGES).
96     CALL DUMMY_IN_STEPPING( myTime, myIter, myThid )
97 heimbach 1.26 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
98 heimbach 1.16 #endif
99    
100 jmc 1.21 #ifdef EXACT_CONSERV
101     IF (exactConserv) THEN
102     C-- Update etaH(n+1) :
103 heimbach 1.36 CALL TIMER_START('UPDATE_ETAH [FORWARD_STEP]',mythid)
104 jmc 1.35 CALL UPDATE_ETAH( myTime, myIter, myThid )
105 heimbach 1.36 CALL TIMER_STOP ('UPDATE_ETAH [FORWARD_STEP]',mythid)
106 jmc 1.21 ENDIF
107     #endif /* EXACT_CONSERV */
108    
109 jmc 1.18 #ifdef NONLIN_FRSURF
110     C-- compute the future surface level thickness
111 jmc 1.21 C according to etaH(n+1)
112 jmc 1.18 IF ( nonlinFreeSurf.GT.0) THEN
113 heimbach 1.36 CALL TIMER_START('CALC_SURF_DR [FORWARD_STEP]',mythid)
114 jmc 1.21 CALL CALC_SURF_DR(etaH, myTime, myIter, myThid )
115 heimbach 1.36 CALL TIMER_STOP ('CALC_SURF_DR [FORWARD_STEP]',mythid)
116 jmc 1.18 ENDIF
117 jmc 1.21 #endif /* NONLIN_FRSURF */
118 jmc 1.18
119 cheisey 1.38
120    
121    
122 adcroft 1.13 C-- Load forcing/external data fields.
123 heimbach 1.28 #ifdef ALLOW_AUTODIFF_TAMC
124     c**************************************
125     #include "checkpoint_lev1_directives.h"
126     c**************************************
127 heimbach 1.23 #endif
128 cheisey 1.38
129    
130 heimbach 1.37 C-- Call external forcing package
131 cheisey 1.38 cswdblk -- add ---
132 cheisey 1.40 #ifdef ALLOW_BULK_FORCE
133 cheisey 1.38 CALL TIMER_START('BULKF_FIELDS_LOAD[THE_MAIN_LOOP]',mythid)
134     CALL BULKF_FIELDS_LOAD( mytime, myiter, mythid )
135     CALL TIMER_STOP ('BULKF_FIELDS_LOAD[THE_MAIN_LOOP]',mythid)
136     c calculate qnet and empmr (and wind stress)
137     DO bj=myByLo(myThid),myByHi(myThid)
138     DO bi=myBxLo(myThid),myBxHi(myThid)
139     CALL BULKF_FORCING( bi,bj, mytime, myiter, mythid )
140     ENDDO
141     ENDDO
142     c Update the tile edges.
143     _EXCH_XY_R8(Qnet, mythid)
144     _EXCH_XY_R8(EmPmR, mythid)
145 cheisey 1.43 CALL EXCH_UV_XY_RS(fu, fv, .TRUE., myThid)
146     C _EXCH_XY_R8(fu , mythid)
147     C _EXCH_XY_R8(fv , mythid)
148 cheisey 1.38 cswdblk -- end add ---
149 cheisey 1.42 #else
150 cheisey 1.39 #ifdef INCLUDE_EXTERNAL_FORCING_PACKAGE
151 cheisey 1.42 IF ( .not. useSEAICE ) THEN
152 heimbach 1.36 CALL TIMER_START('EXF_GETFORCING [FORWARD_STEP]',mythid)
153     CALL EXF_GETFORCING( mytime, myiter, mythid )
154     CALL TIMER_STOP ('EXF_GETFORCING [FORWARD_STEP]',mythid)
155     ENDIF
156 heimbach 1.12 #else
157 cheisey 1.42 IF ( .not. useSEAICE ) THEN
158 heimbach 1.36 CALL TIMER_START('EXTERNAL_FIELDS_LOAD[FORWARD_STEP]',mythid)
159     CALL EXTERNAL_FIELDS_LOAD( mytime, myiter, mythid )
160     CALL TIMER_STOP ('EXTERNAL_FIELDS_LOAD[FORWARD_STEP]',mythid)
161     ENDIF
162     #endif /* INCLUDE_EXTERNAL_FORCING_PACKAGE */
163 cheisey 1.42 #endif
164 cheisey 1.38
165    
166    
167 heimbach 1.36 #ifdef ALLOW_SEAICE
168 heimbach 1.37 C-- Call sea ice model to compute forcing/external data fields.
169 heimbach 1.36 IF ( useSEAICE ) THEN
170     CALL TIMER_START('SEAICE_MODEL [FORWARD_STEP]',myThid)
171 jmc 1.41 CALL SEAICE_MODEL( myTime, myIter, myThid )
172 heimbach 1.36 CALL TIMER_STOP ('SEAICE_MODEL [FORWARD_STEP]',myThid)
173     ENDIF
174     #endif ALLOW_SEAICE
175 adcroft 1.15
176     C-- Step forward fields and calculate time tendency terms.
177 heimbach 1.36 CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',mythid)
178 adcroft 1.15 CALL THERMODYNAMICS( myTime, myIter, myThid )
179 heimbach 1.36 CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',mythid)
180 jmc 1.35
181     C-- do exchanges (needed for DYNAMICS) when using stagger time-step :
182 heimbach 1.36 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
183 jmc 1.35 CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid )
184 heimbach 1.36 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
185 jmc 1.24
186     #ifdef ALLOW_SHAP_FILT
187 jmc 1.27 IF (useSHAP_FILT .AND.
188     & staggerTimeStep .AND. shap_filt_TrStagg ) THEN
189 heimbach 1.36 CALL TIMER_START('SHAP_FILT [FORWARD_STEP]',myThid)
190 jmc 1.24 CALL SHAP_FILT_APPLY_TS( gT, gS, myTime, myIter, myThid )
191 heimbach 1.36 CALL TIMER_STOP ('SHAP_FILT [FORWARD_STEP]',myThid)
192 jmc 1.24 ENDIF
193     #endif
194 jmc 1.27 #ifdef ALLOW_ZONAL_FILT
195     IF (useZONAL_FILT .AND.
196     & staggerTimeStep .AND. zonal_filt_TrStagg ) THEN
197 heimbach 1.36 CALL TIMER_START('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid)
198 jmc 1.27 CALL ZONAL_FILT_APPLY_TS( gT, gS, myThid )
199 heimbach 1.36 CALL TIMER_STOP ('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid)
200 jmc 1.27 ENDIF
201     #endif
202 heimbach 1.12
203 adcroft 1.13 C-- Step forward fields and calculate time tendency terms.
204 adcroft 1.19 IF ( momStepping ) THEN
205 heimbach 1.36 CALL TIMER_START('DYNAMICS [FORWARD_STEP]',mythid)
206 heimbach 1.12 CALL DYNAMICS( myTime, myIter, myThid )
207 heimbach 1.36 CALL TIMER_STOP ('DYNAMICS [FORWARD_STEP]',mythid)
208 adcroft 1.19 ENDIF
209 heimbach 1.12
210 adcroft 1.1 #ifdef ALLOW_NONHYDROSTATIC
211     C-- Step forward W field in N-H algorithm
212 adcroft 1.19 IF ( momStepping .AND. nonHydrostatic ) THEN
213 heimbach 1.36 CALL TIMER_START('CALC_GW [FORWARD_STEP]',myThid)
214 heimbach 1.12 CALL CALC_GW(myThid)
215 heimbach 1.36 CALL TIMER_STOP ('CALC_GW [FORWARD_STEP]',myThid)
216 heimbach 1.12 ENDIF
217 adcroft 1.1 #endif
218 jmc 1.18
219     #ifdef NONLIN_FRSURF
220 jmc 1.21 C-- update hfacC,W,S and recip_hFac according to etaH(n+1) :
221 jmc 1.18 IF ( nonlinFreeSurf.GT.0) THEN
222 heimbach 1.36 CALL TIMER_START('UPDATE_SURF_DR [FORWARD_STEP]',mythid)
223 jmc 1.18 CALL UPDATE_SURF_DR( myTime, myIter, myThid )
224 heimbach 1.36 CALL TIMER_STOP ('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
225 jmc 1.18 ENDIF
226     C- update also CG2D matrix (and preconditioner)
227 jmc 1.33 IF ( momStepping .AND. nonlinFreeSurf.GT.2 ) THEN
228 heimbach 1.36 CALL TIMER_START('UPDATE_CG2D [FORWARD_STEP]',mythid)
229 jmc 1.18 CALL UPDATE_CG2D( myTime, myIter, myThid )
230 jmc 1.44 CALL TIMER_STOP ('UPDATE_CG2D [FORWARD_STEP]',mythid)
231 adcroft 1.19 ENDIF
232 jmc 1.18 #endif
233 adcroft 1.1
234 jmc 1.27 C-- Apply Filters to u*,v* before SOLVE_FOR_PRESSURE
235     #ifdef ALLOW_SHAP_FILT
236     IF (useSHAP_FILT .AND. shap_filt_uvStar) THEN
237 heimbach 1.36 CALL TIMER_START('SHAP_FILT [FORWARD_STEP]',myThid)
238 jmc 1.27 CALL SHAP_FILT_APPLY_UV( gUnm1,gVnm1, myTime,myIter,myThid )
239     IF (implicDiv2Dflow.LT.1.) THEN
240     C-- Explicit+Implicit part of the Barotropic Flow Divergence
241     C => Filtering of uVel,vVel is necessary
242     CALL SHAP_FILT_APPLY_UV( uVel,vVel, myTime,myIter,myThid )
243     ENDIF
244 heimbach 1.36 CALL TIMER_STOP ('SHAP_FILT [FORWARD_STEP]',myThid)
245 jmc 1.27 ENDIF
246     #endif
247     #ifdef ALLOW_ZONAL_FILT
248     IF (useZONAL_FILT .AND. zonal_filt_uvStar) THEN
249 heimbach 1.36 CALL TIMER_START('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid)
250 jmc 1.27 CALL ZONAL_FILT_APPLY_UV( gUnm1, gVnm1, myThid )
251     IF (implicDiv2Dflow.LT.1.) THEN
252     C-- Explicit+Implicit part of the Barotropic Flow Divergence
253     C => Filtering of uVel,vVel is necessary
254     CALL ZONAL_FILT_APPLY_UV( uVel, vVel, myThid )
255     ENDIF
256 heimbach 1.36 CALL TIMER_STOP ('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid)
257 jmc 1.27 ENDIF
258     #endif
259 heimbach 1.12
260 adcroft 1.1 C-- Solve elliptic equation(s).
261     C Two-dimensional only for conventional hydrostatic or
262     C three-dimensional for non-hydrostatic and/or IGW scheme.
263 adcroft 1.19 IF ( momStepping ) THEN
264 heimbach 1.36 CALL TIMER_START('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
265 jmc 1.31 CALL SOLVE_FOR_PRESSURE(myTime, myIter, myThid)
266 heimbach 1.36 CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
267 adcroft 1.19 ENDIF
268 adcroft 1.1
269 heimbach 1.37 #ifdef ALLOW_AUTODIFF_TAMC
270     cph This is needed because convective_adjustment calls
271     cph find_rho which may use pressure()
272     CADJ STORE pressure = comlev1, key = ikey_dynamics
273     #endif
274 adcroft 1.5 C-- Correct divergence in flow field and cycle time-stepping
275 jmc 1.7 C arrays (for all fields) ; update time-counter
276 heimbach 1.12 myIter = nIter0 + iLoop
277     myTime = startTime + deltaTClock * float(iLoop)
278 heimbach 1.36 CALL TIMER_START('THE_CORRECTION_STEP [FORWARD_STEP]',myThid)
279 heimbach 1.12 CALL THE_CORRECTION_STEP(myTime, myIter, myThid)
280 heimbach 1.36 CALL TIMER_STOP ('THE_CORRECTION_STEP [FORWARD_STEP]',myThid)
281 adcroft 1.5
282 adcroft 1.1 C-- Do "blocking" sends and receives for tendency "overlap" terms
283 heimbach 1.36 c CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
284 jmc 1.7 c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )
285 heimbach 1.36 c CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
286 adcroft 1.5
287     C-- Do "blocking" sends and receives for field "overlap" terms
288 heimbach 1.36 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
289 adcroft 1.5 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
290 heimbach 1.36 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
291 adcroft 1.20
292     #ifdef ALLOW_FLT
293     C-- Calculate float trajectories
294     IF (useFLT) THEN
295 heimbach 1.36 CALL TIMER_START('FLOATS [FORWARD_STEP]',myThid)
296 adcroft 1.20 CALL FLT_MAIN(myIter,myTime, myThid)
297 heimbach 1.36 CALL TIMER_STOP ('FLOATS [FORWARD_STEP]',myThid)
298 adcroft 1.20 ENDIF
299     #endif
300 heimbach 1.12
301     #ifndef EXCLUDE_MONITOR
302     C-- Check status of solution (statistics, cfl, etc...)
303 heimbach 1.36 CALL TIMER_START('MONITOR [FORWARD_STEP]',mythid)
304 heimbach 1.12 CALL MONITOR( myIter, myTime, myThid )
305 heimbach 1.36 CALL TIMER_STOP ('MONITOR [FORWARD_STEP]',myThid)
306 heimbach 1.12 #endif /* EXCLUDE_MONITOR */
307 adcroft 1.1
308 jmc 1.7 C-- Do IO if needed.
309 heimbach 1.36 CALL TIMER_START('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
310 heimbach 1.12 CALL DO_THE_MODEL_IO( myTime, myIter, myThid )
311 heimbach 1.36 CALL TIMER_STOP ('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
312 adcroft 1.1
313     C-- Save state for restarts
314 jmc 1.7 C Note: (jmc: is it still the case after ckp35 ?)
315 adcroft 1.1 C =====
316     C Because of the ordering of the timestepping code and
317     C tendency term code at end of loop model arrays hold
318     C U,V,T,S at "time-level" N but gu, gv, gs, gt, guNM1,...
319     C at "time-level" N+1/2 (guNM1 at "time-level" N+1/2 is
320 jmc 1.10 C gu at "time-level" N-1/2) and etaN at "time-level" N+1/2.
321 adcroft 1.1 C where N = I+timeLevBase-1
322     C Thus a checkpoint contains U.0000000000, GU.0000000001 and
323 jmc 1.10 C etaN.0000000001 in the indexing scheme used for the model
324 adcroft 1.1 C "state" files. This example is referred to as a checkpoint
325     C at time level 1
326 heimbach 1.36 CALL TIMER_START('WRITE_CHECKPOINT [FORWARD_STEP]',myThid)
327 adcroft 1.1 CALL WRITE_CHECKPOINT(
328 heimbach 1.12 & .FALSE., myTime, myIter, myThid )
329 heimbach 1.36 CALL TIMER_STOP ('WRITE_CHECKPOINT [FORWARD_STEP]',myThid)
330 adcroft 1.1
331     END

  ViewVC Help
Powered by ViewVC 1.1.22