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