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

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

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


Revision 1.7 - (hide annotations) (download)
Mon May 14 21:46:18 2001 UTC (23 years ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint39
Changes since 1.6: +13 -8 lines
Modifications/fixes to support TAMC differentiability
(mostly missing or wrong directives).

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

  ViewVC Help
Powered by ViewVC 1.1.22