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

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

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


Revision 1.35 - (show annotations) (download)
Mon Oct 7 16:24:45 2002 UTC (21 years, 7 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint46l_post, checkpoint46l_pre, checkpoint46j_post, checkpoint46k_post, checkpoint46m_post
Changes since 1.34: +7 -9 lines
* split calc_exact_eta in 2 S/R : integr_continuity & update_etaH
* move wVel computation at the end of the time step, in S/R integr_continuity
* create specific S/R to exchange T,S before DYNAMICS (for stagger time step)
* update timeave pkg for wVel diagnostic ; put convertEmP2rUnit in PARAMS.h

1 C $Header: /u/gcmpack/MITgcm/model/src/forward_step.F,v 1.34 2002/07/13 04:59:42 heimbach Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: FORWARD_STEP
8 C !INTERFACE:
9 SUBROUTINE FORWARD_STEP( iloop, myTime, myIter, myThid )
10
11 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
27 C !USES:
28 IMPLICIT NONE
29 C == Global variables ==
30 #include "SIZE.h"
31 #include "EEPARAMS.h"
32 #include "PARAMS.h"
33 #include "DYNVARS.h"
34 #include "FFIELDS.h"
35
36 #ifdef ALLOW_NONHYDROSTATIC
37 #include "CG3D.h"
38 #endif
39
40 #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 #ifdef ALLOW_AUTODIFF_TAMC
48 # include "tamc.h"
49 # include "ctrl.h"
50 # include "ctrl_dummy.h"
51 # include "cost.h"
52 # ifdef INCLUDE_EXTERNAL_FORCING_PACKAGE
53 # include "exf_fields.h"
54 # ifdef ALLOW_BULKFORMULAE
55 # include "exf_constants.h"
56 # endif
57 # endif
58 # ifdef ALLOW_OBCS
59 # include "OBCS.h"
60 # endif
61 #endif
62
63 C !LOCAL VARIABLES:
64 C == Routine arguments ==
65 C note: under the multi-threaded model myiter and
66 C mytime are local variables passed around as routine
67 C arguments. Although this is fiddly it saves the need to
68 C impose additional synchronisation points when they are
69 C updated.
70 C myiter - iteration counter for this thread
71 C mytime - time counter for this thread
72 C mythid - thread number for this instance of the routine.
73 integer iloop
74 integer mythid
75 integer myiter
76 _RL mytime
77 INTEGER bi,bj
78
79 CEOP
80
81
82 #ifdef ALLOW_AUTODIFF_TAMC
83 C-- Reset the model iteration counter and the model time.
84 myiter = nIter0 + (iloop-1)
85 mytime = startTime + float(iloop-1)*deltaTclock
86 #endif
87
88 #if (defined (ALLOW_AUTODIFF_TAMC) && defined (ALLOW_AUTODIFF_MONITOR))
89 C Include call to a dummy routine. Its adjoint will be
90 C called at the proper place in the adjoint code.
91 C The adjoint routine will print out adjoint values
92 C if requested. The location of the call is important,
93 C it has to be after the adjoint of the exchanges
94 C (DO_GTERM_BLOCKING_EXCHANGES).
95 CALL DUMMY_IN_STEPPING( myTime, myIter, myThid )
96 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
97 #endif
98
99 #ifdef EXACT_CONSERV
100 IF (exactConserv) THEN
101 C-- Update etaH(n+1) :
102 CALL UPDATE_ETAH( myTime, myIter, myThid )
103 ENDIF
104 #endif /* EXACT_CONSERV */
105
106 #ifdef NONLIN_FRSURF
107 C-- compute the future surface level thickness
108 C according to etaH(n+1)
109 IF ( nonlinFreeSurf.GT.0) THEN
110 CALL CALC_SURF_DR(etaH, myTime, myIter, myThid )
111 ENDIF
112 #endif /* NONLIN_FRSURF */
113
114 C-- Load forcing/external data fields.
115 #ifdef INCLUDE_EXTERNAL_FORCING_PACKAGE
116 C NOTE, that although the exf package is part of the
117 C distribution, it is not currently maintained, i.e.
118 C exf is disabled by default in genmake.
119 #ifdef ALLOW_AUTODIFF_TAMC
120 c**************************************
121 #include "checkpoint_lev1_directives.h"
122 c**************************************
123 #endif
124 CALL EXF_GETFORCING( mytime, myiter, mythid )
125 #else
126 CALL TIMER_START('EXTERNAL_FIELDS_LOAD[THE_MAIN_LOOP]',mythid)
127 CALL EXTERNAL_FIELDS_LOAD( mytime, myiter, mythid )
128 CALL TIMER_STOP ('EXTERNAL_FIELDS_LOAD[THE_MAIN_LOOP]',mythid)
129 #endif
130
131 C-- Step forward fields and calculate time tendency terms.
132 CALL TIMER_START('THERMODYNAMICS [THE_MAIN_LOOP]',mythid)
133 CALL THERMODYNAMICS( myTime, myIter, myThid )
134 CALL TIMER_STOP ('THERMODYNAMICS [THE_MAIN_LOOP]',mythid)
135
136 C-- do exchanges (needed for DYNAMICS) when using stagger time-step :
137 CALL TIMER_START('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
138 CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid )
139 CALL TIMER_STOP ('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
140
141 #ifdef ALLOW_SHAP_FILT
142 IF (useSHAP_FILT .AND.
143 & staggerTimeStep .AND. shap_filt_TrStagg ) THEN
144 CALL TIMER_START('SHAP_FILT [THE_MAIN_LOOP]',myThid)
145 CALL SHAP_FILT_APPLY_TS( gT, gS, myTime, myIter, myThid )
146 CALL TIMER_STOP ('SHAP_FILT [THE_MAIN_LOOP]',myThid)
147 ENDIF
148 #endif
149 #ifdef ALLOW_ZONAL_FILT
150 IF (useZONAL_FILT .AND.
151 & staggerTimeStep .AND. zonal_filt_TrStagg ) THEN
152 CALL TIMER_START('ZONAL_FILT_APPLY [THE_MAIN_LOOP]',myThid)
153 CALL ZONAL_FILT_APPLY_TS( gT, gS, myThid )
154 CALL TIMER_STOP ('ZONAL_FILT_APPLY [THE_MAIN_LOOP]',myThid)
155 ENDIF
156 #endif
157
158 C-- Step forward fields and calculate time tendency terms.
159 IF ( momStepping ) THEN
160 CALL TIMER_START('DYNAMICS [THE_MAIN_LOOP]',mythid)
161 CALL DYNAMICS( myTime, myIter, myThid )
162 CALL TIMER_STOP ('DYNAMICS [THE_MAIN_LOOP]',mythid)
163 ENDIF
164
165 #ifdef ALLOW_NONHYDROSTATIC
166 C-- Step forward W field in N-H algorithm
167 IF ( momStepping .AND. nonHydrostatic ) THEN
168 CALL TIMER_START('CALC_GW [THE_MAIN_LOOP]',myThid)
169 CALL CALC_GW(myThid)
170 CALL TIMER_STOP ('CALC_GW [THE_MAIN_LOOP]',myThid)
171 ENDIF
172 #endif
173
174 #ifdef NONLIN_FRSURF
175 C-- update hfacC,W,S and recip_hFac according to etaH(n+1) :
176 IF ( nonlinFreeSurf.GT.0) THEN
177 CALL UPDATE_SURF_DR( myTime, myIter, myThid )
178 ENDIF
179 C- update also CG2D matrix (and preconditioner)
180 IF ( momStepping .AND. nonlinFreeSurf.GT.2 ) THEN
181 CALL UPDATE_CG2D( myTime, myIter, myThid )
182 ENDIF
183 #endif
184
185 C-- Apply Filters to u*,v* before SOLVE_FOR_PRESSURE
186 #ifdef ALLOW_SHAP_FILT
187 IF (useSHAP_FILT .AND. shap_filt_uvStar) THEN
188 CALL TIMER_START('SHAP_FILT [THE_MAIN_LOOP]',myThid)
189 CALL SHAP_FILT_APPLY_UV( gUnm1,gVnm1, myTime,myIter,myThid )
190 IF (implicDiv2Dflow.LT.1.) THEN
191 C-- Explicit+Implicit part of the Barotropic Flow Divergence
192 C => Filtering of uVel,vVel is necessary
193 CALL SHAP_FILT_APPLY_UV( uVel,vVel, myTime,myIter,myThid )
194 ENDIF
195 CALL TIMER_STOP ('SHAP_FILT [THE_MAIN_LOOP]',myThid)
196 ENDIF
197 #endif
198 #ifdef ALLOW_ZONAL_FILT
199 IF (useZONAL_FILT .AND. zonal_filt_uvStar) THEN
200 CALL TIMER_START('ZONAL_FILT_APPLY [THE_MAIN_LOOP]',myThid)
201 CALL ZONAL_FILT_APPLY_UV( gUnm1, gVnm1, myThid )
202 IF (implicDiv2Dflow.LT.1.) THEN
203 C-- Explicit+Implicit part of the Barotropic Flow Divergence
204 C => Filtering of uVel,vVel is necessary
205 CALL ZONAL_FILT_APPLY_UV( uVel, vVel, myThid )
206 ENDIF
207 CALL TIMER_STOP ('ZONAL_FILT_APPLY [THE_MAIN_LOOP]',myThid)
208 ENDIF
209 #endif
210
211 C-- Solve elliptic equation(s).
212 C Two-dimensional only for conventional hydrostatic or
213 C three-dimensional for non-hydrostatic and/or IGW scheme.
214 IF ( momStepping ) THEN
215 CALL TIMER_START('SOLVE_FOR_PRESSURE [THE_MAIN_LOOP]',myThid)
216 CALL SOLVE_FOR_PRESSURE(myTime, myIter, myThid)
217 CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [THE_MAIN_LOOP]',myThid)
218 ENDIF
219
220 C-- Correct divergence in flow field and cycle time-stepping
221 C arrays (for all fields) ; update time-counter
222 myIter = nIter0 + iLoop
223 myTime = startTime + deltaTClock * float(iLoop)
224 CALL TIMER_START('THE_CORRECTION_STEP [THE_MAIN_LOOP]',myThid)
225 CALL THE_CORRECTION_STEP(myTime, myIter, myThid)
226 CALL TIMER_STOP ('THE_CORRECTION_STEP [THE_MAIN_LOOP]',myThid)
227
228 C-- Do "blocking" sends and receives for tendency "overlap" terms
229 c CALL TIMER_START('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
230 c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )
231 c CALL TIMER_STOP ('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
232
233 C-- Do "blocking" sends and receives for field "overlap" terms
234 CALL TIMER_START('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
235 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
236 CALL TIMER_STOP ('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
237
238 #ifdef ALLOW_FLT
239 C-- Calculate float trajectories
240 IF (useFLT) THEN
241 CALL TIMER_START('FLOATS [THE_MAIN_LOOP]',myThid)
242 CALL FLT_MAIN(myIter,myTime, myThid)
243 CALL TIMER_STOP ('FLOATS [THE_MAIN_LOOP]',myThid)
244 ENDIF
245 #endif
246
247 #ifndef EXCLUDE_MONITOR
248 C-- Check status of solution (statistics, cfl, etc...)
249 CALL MONITOR( myIter, myTime, myThid )
250 #endif /* EXCLUDE_MONITOR */
251
252 C-- Do IO if needed.
253 CALL TIMER_START('DO_THE_MODEL_IO [THE_MAIN_LOOP]',myThid)
254 CALL DO_THE_MODEL_IO( myTime, myIter, myThid )
255 CALL TIMER_STOP ('DO_THE_MODEL_IO [THE_MAIN_LOOP]',myThid)
256
257 C-- Save state for restarts
258 C Note: (jmc: is it still the case after ckp35 ?)
259 C =====
260 C Because of the ordering of the timestepping code and
261 C tendency term code at end of loop model arrays hold
262 C U,V,T,S at "time-level" N but gu, gv, gs, gt, guNM1,...
263 C at "time-level" N+1/2 (guNM1 at "time-level" N+1/2 is
264 C gu at "time-level" N-1/2) and etaN at "time-level" N+1/2.
265 C where N = I+timeLevBase-1
266 C Thus a checkpoint contains U.0000000000, GU.0000000001 and
267 C etaN.0000000001 in the indexing scheme used for the model
268 C "state" files. This example is referred to as a checkpoint
269 C at time level 1
270 CALL TIMER_START('WRITE_CHECKPOINT [THE_MAIN_LOOP]',myThid)
271 CALL WRITE_CHECKPOINT(
272 & .FALSE., myTime, myIter, myThid )
273 CALL TIMER_STOP ('WRITE_CHECKPOINT [THE_MAIN_LOOP]',myThid)
274
275 END

  ViewVC Help
Powered by ViewVC 1.1.22