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