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

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

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


Revision 1.9 - (show annotations) (download)
Mon Jun 4 13:25:35 2001 UTC (23 years ago) by adcroft
Branch: MAIN
Changes since 1.8: +6 -1 lines
Added Kinetic energy monitoring.

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/the_main_loop.F,v 1.8 2001/05/30 19:33:18 adcroft Exp $
2
3 #include "CPP_OPTIONS.h"
4
5
6 subroutine the_main_loop( mytime, myiter, mythid )
7
8 c ==================================================================
9 c SUBROUTINE the_main_loop
10 c ==================================================================
11 c
12 c o Run the ocean model and evaluate the specified cost function.
13 c
14 c *the_main_loop* is the toplevel routine for the Tangent Linear and
15 c Adjoint Model Compiler (TAMC). For this purpose the initialization
16 c of the model was split into two parts. Those parameters that do
17 c not depend on a specific model run are set in *initialise_fixed*,
18 c whereas those that do depend on the specific realization are
19 c initialized in *initialise_varia*.
20 c This routine is to be used in conjuction with the MITgcmuv
21 c checkpoint 37.
22 c
23 c ==================================================================
24 c SUBROUTINE the_main_loop
25 c ==================================================================
26
27 implicit none
28
29 c == global variables ==
30
31 #include "SIZE.h"
32 #include "EEPARAMS.h"
33 #include "PARAMS.h"
34 #include "DYNVARS.h"
35 #include "FFIELDS.h"
36
37 #ifdef ALLOW_NONHYDROSTATIC
38 #include "CG3D.h"
39 #endif
40
41 #ifdef ALLOW_AUTODIFF_TAMC
42 #include "tamc.h"
43 #include "ctrl.h"
44 #include "ctrl_dummy.h"
45 #include "cost.h"
46 #endif
47
48 c == routine arguments ==
49 c note: under the multi-threaded model myiter and
50 c mytime are local variables passed around as routine
51 c arguments. Although this is fiddly it saves the need to
52 c impose additional synchronisation points when they are
53 c updated.
54 c myiter - iteration counter for this thread
55 c mytime - time counter for this thread
56 c mythid - thread number for this instance of the routine.
57 integer mythid
58 integer myiter
59 _RL mytime
60
61 c == local variables ==
62
63 integer iloop
64
65 #ifdef ALLOW_TAMC_CHECKPOINTING
66 integer ilev_1
67 integer ilev_2
68 integer ilev_3
69 integer max_lev2
70 integer max_lev3
71 #endif
72
73 c-- == end of interface ==
74
75 #ifdef ALLOW_AUTODIFF_TAMC
76 c-- Initialize storage for the cost function evaluation.
77 CADJ INIT dummytape = common, 1
78 c-- Initialize storage for the outermost loop.
79 CADJ INIT tapelev3 = USER
80 #ifdef ALLOW_TAMC_CHECKPOINTING
81 ikey_dynamics = 1
82 #endif
83 #endif
84
85 CALL TIMER_START('ADJOINT SPIN-UP', mythid)
86
87 C-- Set initial conditions (variable arrays)
88 CALL TIMER_START('INITIALISE_VARIA [ADJOINT SPIN-UP]', mythid)
89 CALL INITIALISE_VARIA( mythid )
90 CALL TIMER_STOP ('INITIALISE_VARIA [ADJOINT SPIN-UP]', mythid)
91
92 #ifndef ALLOW_AUTODIFF_TAMC
93 c-- Dump for start state.
94 CALL TIMER_START('I/O (WRITE) [ADJOINT SPIN-UP]', mythid)
95 CALL WRITE_STATE( mytime, myiter, mythid )
96 CALL TIMER_STOP ('I/O (WRITE) [ADJOINT SPIN-UP]', mythid)
97 #endif
98
99 #ifndef EXCLUDE_MONITOR
100 C-- Check status of solution (statistics, cfl, etc...)
101 CALL MONITOR( myIter, myTime, myThid )
102 #endif /* EXCLUDE_MONITOR */
103
104 #ifdef ALLOW_COST_TEST
105 c-- Add control vector for forcing and parameter fields
106 CALL CTRL_MAP_FORCING (mythid)
107 #endif
108
109 CALL TIMER_STOP ('ADJOINT SPIN-UP', mythid)
110 _BARRIER
111
112 c-- Do the model integration.
113 call timer_start('ADJOINT MAIN LOOP',mythid)
114
115 c >>>>>>>>>>>>>>>>>>>>>>>>>>> LOOP <<<<<<<<<<<<<<<<<<<<<<<<<<<<
116 c >>>>>>>>>>>>>>>>>>>>>>>>>>> STARTS <<<<<<<<<<<<<<<<<<<<<<<<<<<<
117
118 #ifdef ALLOW_AUTODIFF_TAMC
119 #ifdef ALLOW_TAMC_CHECKPOINTING
120 c-- Implement a three level checkpointing. For a two level
121 c-- checkpointing delete the middle loop; for n levels (n > 3)
122 c-- insert more loops.
123
124 c-- Check the choice of the checkpointing parameters in relation
125 c-- to nTimeSteps: (nchklev_1*nchklev_2*nchklev_3 .ge. nTimeSteps)
126 if (nchklev_1*nchklev_2*nchklev_3 .lt. nTimeSteps) then
127 print*
128 print*, ' the_main_loop: TAMC checkpointing parameters'
129 print*, ' nchklev_1*nchklev_2*nchklev_3 = ',
130 & nchklev_1*nchklev_2*nchklev_3
131 print*, ' are not consistent with nTimeSteps = ',
132 & nTimeSteps
133 stop ' ... stopped in the_main_loop.'
134 endif
135 max_lev3=nTimeSteps/(nchklev_1*nchklev_2)+1
136 max_lev2=nTimeSteps/nchklev_1+1
137
138 do ilev_3 = 1,nchklev_3
139 if(ilev_3.le.max_lev3) then
140 CADJ STORE gsnm1 = tapelev3, key = ilev_3
141 CADJ STORE gtnm1 = tapelev3, key = ilev_3
142 CADJ STORE gunm1 = tapelev3, key = ilev_3
143 CADJ STORE gvnm1 = tapelev3, key = ilev_3
144 CADJ STORE theta = tapelev3, key = ilev_3
145 CADJ STORE salt = tapelev3, key = ilev_3
146 CADJ STORE uvel = tapelev3, key = ilev_3
147 CADJ STORE vvel = tapelev3, key = ilev_3
148 CADJ STORE wvel = tapelev3, key = ilev_3
149 CADJ STORE etan = tapelev3, key = ilev_3
150 CADJ STORE etanm1 = tapelev3, key = ilev_3
151 #ifdef INCLUDE_CD_CODE
152 CADJ STORE uveld = tapelev3, key = ilev_3
153 CADJ STORE vveld = tapelev3, key = ilev_3
154 CADJ STORE unm1 = tapelev3, key = ilev_3
155 CADJ STORE vnm1 = tapelev3, key = ilev_3
156 #endif
157
158 c-- Initialise storage for the middle loop.
159 CADJ INIT tapelev2 = USER
160
161 do ilev_2 = 1,nchklev_2
162 if(ilev_2.le.max_lev2) then
163 CADJ STORE gsnm1 = tapelev2, key = ilev_2
164 CADJ STORE gtnm1 = tapelev2, key = ilev_2
165 CADJ STORE gunm1 = tapelev2, key = ilev_2
166 CADJ STORE gvnm1 = tapelev2, key = ilev_2
167 CADJ STORE theta = tapelev2, key = ilev_2
168 CADJ STORE salt = tapelev2, key = ilev_2
169 CADJ STORE uvel = tapelev2, key = ilev_2
170 CADJ STORE vvel = tapelev2, key = ilev_2
171 CADJ STORE wvel = tapelev2, key = ilev_2
172 CADJ STORE etan = tapelev2, key = ilev_2
173 CADJ STORE etanm1 = tapelev2, key = ilev_2
174 #ifdef INCLUDE_CD_CODE
175 CADJ STORE uveld = tapelev2, key = ilev_2
176 CADJ STORE vveld = tapelev2, key = ilev_2
177 CADJ STORE unm1 = tapelev2, key = ilev_2
178 CADJ STORE vnm1 = tapelev2, key = ilev_2
179 #endif
180
181 c-- Initialize storage for the innermost loop.
182 c-- Always check common block sizes for the checkpointing!
183 CADJ INIT comlev1 = COMMON,nchklev_1
184 CADJ INIT comlev1_bibj = COMMON,nchklev_1*nsx*nsy*nthreads_chkpt
185 CADJ INIT comlev1_bibj_k = COMMON,nchklev_1*nsx*nsy*nr*nthreads_chkpt
186 cphCADJ INIT comlev1_impl = COMMON,nchklev_1*nsx*nsy*6
187 CADJ INIT comlev1_kpp = COMMON,nchklev_1*nsx*nsy
188
189 C-- RG replace 2 by max of num_v_smooth_Ri
190 CADJ INIT comlev1_kpp_sm = COMMON,nchklev_1*nsx*nsy*2
191
192 do ilev_1 = 1,nchklev_1
193
194 c-- The if-statement below introduces a some flexibility in the
195 c-- choice of the 3-tupel ( nchklev_1, nchklev_2, nchklev_3 ).
196 c--
197 c-- Requirement: nchklev_1*nchklev_2*nchklev_3 .ge. nTimeSteps .
198
199 iloop = (ilev_3 - 1)*nchklev_2*nchklev_1 +
200 & (ilev_2 - 1)*nchklev_1 + ilev_1
201
202 if ( iloop .le. nTimeSteps ) then
203
204 #else /* ALLOW_TAMC_CHECKPOINTING undefined */
205 c-- Initialise storage for reference trajectory without TAMC check-
206 c-- pointing.
207 CADJ INIT history = USER
208 CADJ INIT comlev1_bibj = COMMON,nchklev_0*nsx*nsy*nthreads_chkpt
209 CADJ INIT comlev1_bibj_k = COMMON,nchklev_0*nsx*nsy*nr*nthreads_chkpt
210 CADJ INIT comlev1_impl = COMMON,nchklev_0*nsx*nsy*6
211 CADJ INIT comlev1_impl_k = COMMON,nchklev_0*nsx*nsy*(nr-2)*6
212 CADJ INIT comlev1_kpp = COMMON,nchklev_0*nsx*nsy
213
214 C-- RG replace 2 by max of num_v_smooth_Ri
215 CADJ INIT comlev1_kpp_sm = COMMON,nchklev_0*nsx*nsy*2
216
217 c-- Check the choice of the checkpointing parameters in relation
218 c-- to nTimeSteps: (nchklev_0 .ge. nTimeSteps)
219 if (nchklev_0 .lt. nTimeSteps) then
220 print*
221 print*, ' the_main_loop: TAMC checkpointing parameter ',
222 & nchklev_0 = ', nchklev_0
223 print*, ' not consistent with nTimeSteps = ',
224 & nTimeSteps
225 stop ' ... stopped in the_main_loop.'
226 endif
227
228 do iloop = 1, nTimeSteps
229
230 #endif /* ALLOW_TAMC_CHECKPOINTING */
231
232 #else /* ALLOW_AUTODIFF_TAMC undefined */
233
234 c-- Start the main loop of adjoint_Objfunc. Automatic differentiation
235 c-- NOT enabled.
236 do iloop = 1, nTimeSteps
237
238 #endif /* ALLOW_AUTODIFF_TAMC */
239
240 c-- >>> Loop body start <<<
241
242 #ifdef ALLOW_AUTODIFF_TAMC
243 c-- Set the model iteration counter and the model time.
244 myiter = nIter0 + (iloop-1)
245 mytime = startTime + float(iloop-1)*deltaTclock
246 #endif
247
248 c-- Load forcing/external data fields.
249 call timer_start('EXTERNAL_FIELDS_LOAD [ADJOINT]',mythid)
250 #ifdef INCLUDE_EXTERNAL_FORCING_PACKAGE
251 c NOTE, that although the exf package is part of the
252 c distribution, it is not currently maintained, i.e.
253 c exf is disabled by default in genmake.
254 call exf_getforcing( mytime, myiter, mythid )
255 #else
256 call external_fields_load( mytime, myiter, mythid )
257 #endif
258 call timer_stop ('EXTERNAL_FIELDS_LOAD [ADJOINT]',mythid)
259
260 #ifdef ALLOW_TAMC_CHECKPOINTING
261 ikey_dynamics = ilev_1
262 #endif
263 c-- Step forward fields and calculate time tendency terms.
264 call timer_start('DYNAMICS [ADJOINT]',mythid)
265 call dynamics( mytime, myiter, mythid )
266 call timer_stop ('DYNAMICS [ADJOINT]',mythid)
267
268 #ifndef ALLOW_AUTODIFF_TAMC
269
270 #ifdef ALLOW_NONHYDROSTATIC
271 C-- Step forward W field in N-H algorithm
272 IF ( nonHydrostatic ) THEN
273 CALL TIMER_START('CALC_GW [ADJOINT]',myThid)
274 CALL CALC_GW( myThid)
275 CALL TIMER_STOP ('CALC_GW [ADJOINT]',myThid)
276 ENDIF
277 #endif
278
279 #ifdef ALLOW_SHAP_FILT
280 C-- Step forward all tiles, filter and exchange.
281 CALL TIMER_START('SHAP_FILT [ADJOINT]',myThid)
282 CALL SHAP_FILT_APPLY(
283 I gUnm1, gVnm1, gTnm1, gSnm1,
284 I myTime, myIter, myThid )
285 IF (implicDiv2Dflow.LT.1.) THEN
286 C-- Explicit+Implicit part of the Barotropic Flow Divergence
287 C => Filtering of uVel,vVel is necessary
288 CALL SHAP_FILT_UV( uVel, vVel, myTime, myThid )
289 ENDIF
290 CALL TIMER_STOP ('SHAP_FILT [ADJOINT]',myThid)
291 #endif
292
293 #ifdef ALLOW_ZONAL_FILT
294 IF (zonal_filt_lat.LT.90.) THEN
295 CALL ZONAL_FILT_APPLY(
296 U gUnm1, gVnm1, gTnm1, gSnm1,
297 I myThid )
298 ENDIF
299 #endif
300
301 #endif /* ALLOW_AUTODIFF_TAMC */
302
303 C-- Solve elliptic equation(s).
304 C Two-dimensional only for conventional hydrostatic or
305 C three-dimensional for non-hydrostatic and/or IGW scheme.
306 CALL TIMER_START('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
307 CALL SOLVE_FOR_PRESSURE( myThid )
308 CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
309
310 #ifdef ALLOW_AUTODIFF_TAMC
311 c Include call to a dummy routine. Its adjoint will be
312 c called at the proper place in the adjoint code.
313 c The adjoint routine will print out adjoint values
314 c if requested. The location of the call is important,
315 c it has to be after the adjoint of the exchanges
316 c (DO_GTERM_BLOCKING_EXCHANGES).
317 call dummy_in_stepping( myTime, myIter, myThid )
318 #endif
319
320 C-- Correct divergence in flow field and cycle time-stepping
321 C arrays (for all fields) ; update time-counter
322 myIter = nIter0 + iLoop
323 myTime = startTime + deltaTClock * float(iLoop)
324 CALL THE_CORRECTION_STEP(myTime, myIter, myThid)
325
326 C-- Do "blocking" sends and receives for tendency "overlap" terms
327 c CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
328 c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )
329 c CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
330
331 C-- Do "blocking" sends and receives for field "overlap" terms
332 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
333 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
334 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
335
336 #ifndef EXCLUDE_MONITOR
337 C-- Check status of solution (statistics, cfl, etc...)
338 CALL MONITOR( myIter, myTime, myThid )
339 #endif /* EXCLUDE_MONITOR */
340
341 #ifndef ALLOW_AUTODIFF_TAMC
342 C-- Do IO if needed.
343 CALL TIMER_START('I/O (WRITE) [FORWARD_STEP]',myThid)
344 CALL DO_THE_MODEL_IO( myTime, myIter, myThid )
345 CALL TIMER_STOP ('I/O (WRITE) [FORWARD_STEP]',myThid)
346
347 C-- Save state for restarts
348 C Note: (jmc: is it still the case after ckp35 ?)
349 C =====
350 C Because of the ordering of the timestepping code and
351 C tendency term code at end of loop model arrays hold
352 C U,V,T,S at "time-level" N but gu, gv, gs, gt, guNM1,...
353 C at "time-level" N+1/2 (guNM1 at "time-level" N+1/2 is
354 C gu at "time-level" N-1/2) and etaN at "time-level" N+1/2.
355 C where N = I+timeLevBase-1
356 C Thus a checkpoint contains U.0000000000, GU.0000000001 and
357 C etaN.0000000001 in the indexing scheme used for the model
358 C "state" files. This example is referred to as a checkpoint
359 C at time level 1
360 CALL TIMER_START('I/O (WRITE) [FORWARD_STEP]',myThid)
361 CALL WRITE_CHECKPOINT(
362 & .FALSE., myTime, myIter, myThid )
363 CALL TIMER_STOP ('I/O (WRITE) [FORWARD_STEP]',myThid)
364
365 #endif /* ALLOW_AUTODIFF_TAMC */
366
367 c-- >>> Loop body end <<<
368
369 #ifdef ALLOW_AUTODIFF_TAMC
370 #ifdef ALLOW_TAMC_CHECKPOINTING
371 endif
372 enddo
373 endif
374 enddo
375 endif
376 enddo
377 #else
378 enddo
379 #endif
380
381 #else
382 enddo
383 #endif
384
385 _BARRIER
386 call timer_stop ('ADJOINT MAIN LOOP', mythid)
387
388 call timer_start('ADJOINT FINALIZE COST', mythid)
389
390 #ifdef ALLOW_COST_TEST
391 c-- Veronique's test case
392 call timer_start('cost_test [ADJOINT SPIN-DOWN]', mythid)
393 call cost_test( myThid )
394 call timer_stop ('cost_test [ADJOINT SPIN-DOWN]', mythid)
395
396 c-- Sum all cost function contributions.
397 call timer_start('COST_FINAL [ADJOINT SPIN-DOWN]', mythid)
398 call cost_Final( mythid )
399 call timer_stop ('COST_FINAL [ADJOINT SPIN-DOWN]', mythid)
400 #endif
401
402 call timer_stop ('ADJOINT FINALIZE COST', mythid)
403
404 end
405

  ViewVC Help
Powered by ViewVC 1.1.22