1 |
C $Header: /u/gcmpack/MITgcm/model/src/the_main_loop.F,v 1.65 2005/11/09 19:54:23 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "PACKAGES_CONFIG.h" |
5 |
#include "CPP_OPTIONS.h" |
6 |
#ifdef ALLOW_OBCS |
7 |
# include "OBCS_OPTIONS.h" |
8 |
#endif |
9 |
#ifdef ALLOW_SEAICE |
10 |
# include "SEAICE_OPTIONS.h" |
11 |
#endif |
12 |
#ifdef ALLOW_GMREDI |
13 |
# include "GMREDI_OPTIONS.h" |
14 |
#endif |
15 |
|
16 |
CBOP |
17 |
C !ROUTINE: THE_MAIN_LOOP |
18 |
C !INTERFACE: |
19 |
SUBROUTINE THE_MAIN_LOOP( myTime, myIter, myThid ) |
20 |
|
21 |
C !DESCRIPTION: \bv |
22 |
C *================================================================* |
23 |
C | SUBROUTINE the_main_loop |
24 |
C | o Run the ocean model and evaluate the specified cost function. |
25 |
C *================================================================* |
26 |
C | |
27 |
C | THE_MAIN_LOOP is the toplevel routine for the Tangent Linear and |
28 |
C | Adjoint Model Compiler (TAMC). For this purpose the initialization |
29 |
C | of the model was split into two parts. Those parameters that do |
30 |
C | not depend on a specific model run are set in INITIALISE_FIXED, |
31 |
C | whereas those that do depend on the specific realization are |
32 |
C | initialized in INITIALISE_VARIA. |
33 |
C | This routine is to be used in conjuction with the MITgcmuv |
34 |
C | checkpoint 37. |
35 |
C *================================================================* |
36 |
C \ev |
37 |
|
38 |
C !USES: |
39 |
IMPLICIT NONE |
40 |
C == Global variables == |
41 |
#include "SIZE.h" |
42 |
#include "EEPARAMS.h" |
43 |
#include "PARAMS.h" |
44 |
|
45 |
c************************************** |
46 |
#ifdef ALLOW_AUTODIFF_TAMC |
47 |
|
48 |
c These includes are needed for |
49 |
c AD-checkpointing. |
50 |
c They provide the fields to be stored. |
51 |
|
52 |
# include "GRID.h" |
53 |
# include "DYNVARS.h" |
54 |
# include "FFIELDS.h" |
55 |
# include "EOS.h" |
56 |
# include "GAD.h" |
57 |
# ifdef ALLOW_CD_CODE |
58 |
# include "CD_CODE_VARS.h" |
59 |
# endif |
60 |
# ifdef ALLOW_PTRACERS |
61 |
# include "PTRACERS_SIZE.h" |
62 |
# include "PTRACERS.h" |
63 |
# endif |
64 |
# if (defined (EXACT_CONSERV) || defined (NONLIN_FRSURF)) |
65 |
# include "SURFACE.h" |
66 |
# endif |
67 |
# ifdef ALLOW_OBCS |
68 |
# include "OBCS.h" |
69 |
# endif |
70 |
# ifdef ALLOW_EXF |
71 |
# include "exf_fields.h" |
72 |
# include "exf_clim_fields.h" |
73 |
# ifdef ALLOW_BULKFORMULAE |
74 |
# include "exf_constants.h" |
75 |
# endif |
76 |
# endif /* ALLOW_EXF */ |
77 |
# ifdef ALLOW_SEAICE |
78 |
# include "SEAICE.h" |
79 |
# endif |
80 |
# ifdef ALLOW_EBM |
81 |
# include "EBM.h" |
82 |
# endif |
83 |
# ifdef ALLOW_DIVIDED_ADJOINT_MPI |
84 |
# include "mpif.h" |
85 |
# endif |
86 |
|
87 |
# include "tamc.h" |
88 |
# include "ctrl.h" |
89 |
# include "ctrl_dummy.h" |
90 |
# include "cost.h" |
91 |
|
92 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
93 |
c************************************** |
94 |
|
95 |
C !INPUT/OUTPUT PARAMETERS: |
96 |
C == Routine arguments == |
97 |
C note: under the multi-threaded model myiter and |
98 |
C mytime are local variables passed around as routine |
99 |
C arguments. Although this is fiddly it saves the need to |
100 |
C impose additional synchronisation points when they are |
101 |
C updated. |
102 |
C myIter - iteration counter for this thread |
103 |
C myTime - time counter for this thread |
104 |
C myThid - thread number for this instance of the routine. |
105 |
INTEGER myThid |
106 |
INTEGER myIter |
107 |
_RL myTime |
108 |
|
109 |
C !FUNCTIONS: |
110 |
C == Functions == |
111 |
#ifdef ALLOW_RUNCLOCK |
112 |
LOGICAL RUNCLOCK_CONTINUE |
113 |
LOGICAL RC_CONT |
114 |
#endif |
115 |
|
116 |
C !LOCAL VARIABLES: |
117 |
C == Local variables == |
118 |
integer iloop |
119 |
#ifdef ALLOW_AUTODIFF_TAMC |
120 |
integer ilev_1 |
121 |
integer ilev_2 |
122 |
integer ilev_3 |
123 |
integer ilev_4 |
124 |
integer max_lev2 |
125 |
integer max_lev3 |
126 |
integer max_lev4 |
127 |
#endif |
128 |
CEOP |
129 |
|
130 |
#ifdef ALLOW_DEBUG |
131 |
IF (debugMode) CALL DEBUG_ENTER('THE_MAIN_LOOP',myThid) |
132 |
#endif |
133 |
|
134 |
#ifdef ALLOW_AUTODIFF_TAMC |
135 |
c-- Initialize storage for the cost function evaluation. |
136 |
CADJ INIT dummytape = common, 1 |
137 |
c-- Initialize storage for the outermost loop. |
138 |
CADJ INIT tapelev_ini_bibj_k = USER |
139 |
CADJ INIT tapelev_init = USER |
140 |
c |
141 |
#if (defined (AUTODIFF_2_LEVEL_CHECKPOINT)) |
142 |
CADJ INIT tapelev2 = USER |
143 |
#elif (defined (AUTODIFF_4_LEVEL_CHECKPOINT)) |
144 |
CADJ INIT tapelev4 = USER |
145 |
#else |
146 |
CADJ INIT tapelev3 = USER |
147 |
#endif |
148 |
|
149 |
nIter0 = NINT( (startTime-baseTime)/deltaTClock ) |
150 |
ikey_dynamics = 1 |
151 |
|
152 |
CALL TIMER_START('ADJOINT SPIN-UP', mythid) |
153 |
#endif |
154 |
|
155 |
#ifdef ALLOW_AUTODIFF_TAMC |
156 |
# ifdef NONLIN_FRSURF |
157 |
CADJ STORE hFacC = tapelev_init, key = 1 |
158 |
# endif |
159 |
#endif |
160 |
|
161 |
#ifdef ALLOW_DEBUG |
162 |
IF (debugMode) CALL DEBUG_CALL('INITIALISE_VARIA',myThid) |
163 |
#endif |
164 |
C-- Set initial conditions (variable arrays) |
165 |
CALL TIMER_START('INITIALISE_VARIA [THE_MAIN_LOOP]', mythid) |
166 |
CALL INITIALISE_VARIA( mythid ) |
167 |
CALL TIMER_STOP ('INITIALISE_VARIA [THE_MAIN_LOOP]', mythid) |
168 |
|
169 |
#ifdef ALLOW_MONITOR |
170 |
#ifdef ALLOW_DEBUG |
171 |
IF (debugMode) CALL DEBUG_CALL('MONITOR',myThid) |
172 |
#endif |
173 |
C-- Check status of solution (statistics, cfl, etc...) |
174 |
CALL TIMER_START('MONITOR [THE_MAIN_LOOP]', mythid) |
175 |
CALL MONITOR( myIter, myTime, myThid ) |
176 |
CALL TIMER_STOP ('MONITOR [THE_MAIN_LOOP]', mythid) |
177 |
#endif /* ALLOW_MONITOR */ |
178 |
|
179 |
C-- Do IO if needed (Dump for start state). |
180 |
#ifdef ALLOW_DEBUG |
181 |
IF (debugMode) CALL DEBUG_CALL('DO_THE_MODEL_IO',myThid) |
182 |
#endif |
183 |
|
184 |
#ifdef ALLOW_OFFLINE |
185 |
CALL TIMER_START('OFFLINE_MODEL_IO [FORWARD_STEP]',myThid) |
186 |
CALL OFFLINE_MODEL_IO( myTime, myIter, myThid ) |
187 |
CALL TIMER_STOP ('OFFLINE_MODEL_IO [FORWARD_STEP]',myThid) |
188 |
#else |
189 |
CALL TIMER_START('DO_THE_MODEL_IO [THE_MAIN_LOOP]', mythid) |
190 |
CALL DO_THE_MODEL_IO( myTime, myIter, mythid ) |
191 |
CALL TIMER_STOP ('DO_THE_MODEL_IO [THE_MAIN_LOOP]', mythid) |
192 |
#endif |
193 |
|
194 |
|
195 |
#ifdef ALLOW_AUTODIFF_TAMC |
196 |
CALL TIMER_STOP ('ADJOINT SPIN-UP', mythid) |
197 |
_BARRIER |
198 |
#endif |
199 |
|
200 |
c-- Do the model integration. |
201 |
CALL TIMER_START('MAIN LOOP [THE_MAIN_LOOP]', mythid) |
202 |
|
203 |
c >>>>>>>>>>>>>>>>>>>>>>>>>>> LOOP <<<<<<<<<<<<<<<<<<<<<<<<<<<< |
204 |
c >>>>>>>>>>>>>>>>>>>>>>>>>>> STARTS <<<<<<<<<<<<<<<<<<<<<<<<<<<< |
205 |
|
206 |
#ifdef ALLOW_AUTODIFF_TAMC |
207 |
#ifdef ALLOW_TAMC_CHECKPOINTING |
208 |
|
209 |
max_lev4=nTimeSteps/(nchklev_1*nchklev_2*nchklev_3)+1 |
210 |
max_lev3=nTimeSteps/(nchklev_1*nchklev_2)+1 |
211 |
max_lev2=nTimeSteps/nchklev_1+1 |
212 |
|
213 |
c************************************** |
214 |
#ifdef ALLOW_DIVIDED_ADJOINT |
215 |
CADJ loop = divided |
216 |
#endif |
217 |
c************************************** |
218 |
|
219 |
#ifdef AUTODIFF_4_LEVEL_CHECKPOINT |
220 |
do ilev_4 = 1,nchklev_4 |
221 |
if(ilev_4.le.max_lev4) then |
222 |
c************************************** |
223 |
#include "checkpoint_lev4_directives.h" |
224 |
c************************************** |
225 |
c-- Initialise storage for the middle loop. |
226 |
CADJ INIT tapelev3 = USER |
227 |
#endif /* AUTODIFF_4_LEVEL_CHECKPOINT */ |
228 |
|
229 |
#ifndef AUTODIFF_2_LEVEL_CHECKPOINT |
230 |
do ilev_3 = 1,nchklev_3 |
231 |
if(ilev_3.le.max_lev3) then |
232 |
c************************************** |
233 |
#include "checkpoint_lev3_directives.h" |
234 |
c************************************** |
235 |
c-- Initialise storage for the middle loop. |
236 |
CADJ INIT tapelev2 = USER |
237 |
#endif /* AUTODIFF_2_LEVEL_CHECKPOINT */ |
238 |
|
239 |
do ilev_2 = 1,nchklev_2 |
240 |
if(ilev_2.le.max_lev2) then |
241 |
c************************************** |
242 |
#include "checkpoint_lev2_directives.h" |
243 |
c************************************** |
244 |
|
245 |
c************************************** |
246 |
#ifdef ALLOW_AUTODIFF_TAMC |
247 |
c-- Initialize storage for the innermost loop. |
248 |
c-- Always check common block sizes for the checkpointing! |
249 |
c-- |
250 |
CADJ INIT comlev1 = COMMON,nchklev_1 |
251 |
CADJ INIT comlev1_bibj = COMMON,nchklev_1*nsx*nsy*nthreads_chkpt |
252 |
CADJ INIT comlev1_bibj_k = COMMON,nchklev_1*nsx*nsy*nr*nthreads_chkpt |
253 |
c-- |
254 |
#ifdef ALLOW_KPP |
255 |
CADJ INIT comlev1_kpp = COMMON,nchklev_1*nsx*nsy |
256 |
CADJ INIT comlev1_kpp_k = COMMON,nchklev_1*nsx*nsy*nr |
257 |
#endif /* ALLOW_KPP */ |
258 |
c-- |
259 |
#ifdef ALLOW_GMREDI |
260 |
CADJ INIT comlev1_gmredi_k_gad |
261 |
CADJ & = COMMON,nchklev_1*nsx*nsy*nr*nthreads_chkpt*maxpass |
262 |
#endif /* ALLOW_GMREDI */ |
263 |
c-- |
264 |
#ifdef ALLOW_PTRACERS |
265 |
CADJ INIT comlev1_bibj_ptracers = COMMON, |
266 |
CADJ & nchklev_1*nsx*nsy*nthreads_chkpt*PTRACERS_num |
267 |
#endif /* ALLOW_PTRACERS */ |
268 |
c-- |
269 |
#ifndef DISABLE_MULTIDIM_ADVECTION |
270 |
CADJ INIT comlev1_bibj_k_gad |
271 |
CADJ & = COMMON,nchklev_1*nsx*nsy*nr*nthreads_chkpt*maxpass |
272 |
CADJ INIT comlev1_bibj_k_gad_pass |
273 |
CADJ & = COMMON,nchklev_1*nsx*nsy*nr*nthreads_chkpt*maxpass*maxcube |
274 |
#endif /* DISABLE_MULTIDIM_ADVECTION */ |
275 |
c-- |
276 |
#if (defined (ALLOW_EXF) && defined (ALLOW_BULKFORMULAE)) |
277 |
CADJ INIT comlev1_exf_1 |
278 |
CADJ & = COMMON,nchklev_1*snx*nsx*sny*nsy*nthreads_chkpt |
279 |
CADJ INIT comlev1_exf_2 |
280 |
CADJ & = COMMON,niter_bulk*nchklev_1*snx*nsx*sny*nsy*nthreads_chkpt |
281 |
#endif /* ALLOW_BULKFORMULAE */ |
282 |
c-- |
283 |
#ifdef ALLOW_SEAICE |
284 |
# ifdef SEAICE_ALLOW_DYNAMICS |
285 |
CADJ INIT comlev1_lsr = COMMON,nchklev_1*2 |
286 |
# endif |
287 |
#endif /* ALLOW_SEAICE */ |
288 |
c-- |
289 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
290 |
c************************************** |
291 |
|
292 |
do ilev_1 = 1,nchklev_1 |
293 |
|
294 |
c-- The if-statement below introduces a some flexibility in the |
295 |
c-- choice of the 3-tupel ( nchklev_1, nchklev_2, nchklev_3 ). |
296 |
|
297 |
iloop = (ilev_2 - 1)*nchklev_1 + ilev_1 |
298 |
#ifndef AUTODIFF_2_LEVEL_CHECKPOINT |
299 |
& + (ilev_3 - 1)*nchklev_2*nchklev_1 |
300 |
#endif |
301 |
#ifdef AUTODIFF_4_LEVEL_CHECKPOINT |
302 |
& + (ilev_4 - 1)*nchklev_3*nchklev_2*nchklev_1 |
303 |
#endif |
304 |
|
305 |
if ( iloop .le. nTimeSteps ) then |
306 |
|
307 |
#else /* ALLOW_TAMC_CHECKPOINTING undefined */ |
308 |
c-- Initialise storage for reference trajectory without TAMC check- |
309 |
c-- pointing. |
310 |
CADJ INIT history = USER |
311 |
CADJ INIT comlev1_bibj = COMMON,nchklev_0*nsx*nsy*nthreads_chkpt |
312 |
CADJ INIT comlev1_bibj_k = COMMON,nchklev_0*nsx*nsy*nr*nthreads_chkpt |
313 |
CADJ INIT comlev1_kpp = COMMON,nchklev_0*nsx*nsy |
314 |
|
315 |
c-- Check the choice of the checkpointing parameters in relation |
316 |
c-- to nTimeSteps: (nchklev_0 .ge. nTimeSteps) |
317 |
if (nchklev_0 .lt. nTimeSteps) then |
318 |
print* |
319 |
print*, ' the_main_loop: TAMC checkpointing parameter ', |
320 |
& 'nchklev_0 = ', nchklev_0 |
321 |
print*, ' not consistent with nTimeSteps = ', |
322 |
& nTimeSteps |
323 |
stop ' ... stopped in the_main_loop.' |
324 |
endif |
325 |
|
326 |
DO iloop = 1, nTimeSteps |
327 |
|
328 |
#endif /* ALLOW_TAMC_CHECKPOINTING */ |
329 |
|
330 |
#else /* ALLOW_AUTODIFF_TAMC undefined */ |
331 |
|
332 |
c-- Start the main loop of adjoint_Objfunc. Automatic differentiation |
333 |
c-- NOT enabled. |
334 |
DO iloop = 1, nTimeSteps |
335 |
|
336 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
337 |
|
338 |
c-- >>> Loop body start <<< |
339 |
|
340 |
#ifdef ALLOW_AUTODIFF_TAMC |
341 |
nIter0 = NINT( (startTime-baseTime)/deltaTClock ) |
342 |
ikey_dynamics = ilev_1 |
343 |
CALL AUTODIFF_INADMODE_UNSET( myThid ) |
344 |
#endif |
345 |
|
346 |
#ifdef ALLOW_AUTODIFF_TAMC |
347 |
CALL AUTODIFF_INADMODE_UNSET( myThid ) |
348 |
#endif |
349 |
|
350 |
#ifdef ALLOW_DEBUG |
351 |
IF (debugMode) CALL DEBUG_CALL('FORWARD_STEP',myThid) |
352 |
#endif |
353 |
CALL TIMER_START('FORWARD_STEP [THE_MAIN_LOOP]',mythid) |
354 |
CALL FORWARD_STEP( iloop, mytime, myiter, mythid ) |
355 |
CALL TIMER_STOP ('FORWARD_STEP [THE_MAIN_LOOP]',mythid) |
356 |
|
357 |
#ifdef ALLOW_AUTODIFF_TAMC |
358 |
CALL AUTODIFF_INADMODE_SET( myThid ) |
359 |
#endif |
360 |
|
361 |
#ifdef ALLOW_RUNCLOCK |
362 |
IF (useRunClock) THEN |
363 |
RC_CONT=RUNCLOCK_CONTINUE( myThid ) |
364 |
IF (.NOT.RC_CONT) RETURN |
365 |
ENDIF |
366 |
#endif /* ALLOW_RUNCLOCK */ |
367 |
|
368 |
c-- >>> Loop body end <<< |
369 |
|
370 |
#ifdef ALLOW_AUTODIFF_TAMC |
371 |
CALL AUTODIFF_INADMODE_SET( myThid ) |
372 |
#ifdef ALLOW_TAMC_CHECKPOINTING |
373 |
endif |
374 |
enddo |
375 |
endif |
376 |
enddo |
377 |
#ifndef AUTODIFF_2_LEVEL_CHECKPOINT |
378 |
endif |
379 |
enddo |
380 |
#endif |
381 |
#ifdef AUTODIFF_4_LEVEL_CHECKPOINT |
382 |
endif |
383 |
enddo |
384 |
#endif |
385 |
#else /* ndef ALLOW_TAMC_CHECKPOINTING */ |
386 |
enddo |
387 |
#endif /* ALLOW_TAMC_CHECKPOINTING */ |
388 |
|
389 |
#else /* ndef ALLOW_AUTODIFF_TAMC */ |
390 |
enddo |
391 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
392 |
|
393 |
#ifdef ALLOW_COST |
394 |
c-- Sum all cost function contributions. |
395 |
call TIMER_START('COST_FINAL [ADJOINT SPIN-DOWN]', mythid) |
396 |
call COST_FINAL ( mythid ) |
397 |
call TIMER_STOP ('COST_FINAL [ADJOINT SPIN-DOWN]', mythid) |
398 |
#endif |
399 |
|
400 |
_BARRIER |
401 |
CALL TIMER_STOP ('MAIN LOOP [THE_MAIN_LOOP]', mythid) |
402 |
|
403 |
#ifdef ALLOW_DEBUG |
404 |
IF (debugMode) CALL DEBUG_LEAVE('THE_MAIN_LOOP',myThid) |
405 |
#endif |
406 |
|
407 |
END |