1 |
C $Header: /u/gcmpack/models/MITgcmUV/model/src/the_main_loop.F,v 1.16 2001/08/14 00:20:49 heimbach Exp $ |
2 |
|
3 |
#include "CPP_OPTIONS.h" |
4 |
|
5 |
subroutine the_main_loop( mytime, myiter, mythid ) |
6 |
|
7 |
c ================================================================== |
8 |
c SUBROUTINE the_main_loop |
9 |
c ================================================================== |
10 |
c |
11 |
c o Run the ocean model and evaluate the specified cost function. |
12 |
c |
13 |
c *the_main_loop* is the toplevel routine for the Tangent Linear and |
14 |
c Adjoint Model Compiler (TAMC). For this purpose the initialization |
15 |
c of the model was split into two parts. Those parameters that do |
16 |
c not depend on a specific model run are set in *initialise_fixed*, |
17 |
c whereas those that do depend on the specific realization are |
18 |
c initialized in *initialise_varia*. |
19 |
c This routine is to be used in conjuction with the MITgcmuv |
20 |
c checkpoint 37. |
21 |
c |
22 |
c ================================================================== |
23 |
c SUBROUTINE the_main_loop |
24 |
c ================================================================== |
25 |
|
26 |
implicit none |
27 |
|
28 |
c == global variables == |
29 |
|
30 |
#include "SIZE.h" |
31 |
#include "EEPARAMS.h" |
32 |
#include "PARAMS.h" |
33 |
|
34 |
#ifdef ALLOW_AUTODIFF_TAMC |
35 |
# include "tamc.h" |
36 |
# include "ctrl.h" |
37 |
# include "ctrl_dummy.h" |
38 |
# include "cost.h" |
39 |
# include "DYNVARS.h" |
40 |
# include "FFIELDS.h" |
41 |
# ifdef ALLOW_PASSIVE_TRACER |
42 |
# include "TR1.h" |
43 |
# endif |
44 |
# ifdef ALLOW_NONHYDROSTATIC |
45 |
# include "CG3D.h" |
46 |
# endif |
47 |
#endif |
48 |
|
49 |
c == routine arguments == |
50 |
c note: under the multi-threaded model myiter and |
51 |
c mytime are local variables passed around as routine |
52 |
c arguments. Although this is fiddly it saves the need to |
53 |
c impose additional synchronisation points when they are |
54 |
c updated. |
55 |
c myiter - iteration counter for this thread |
56 |
c mytime - time counter for this thread |
57 |
c mythid - thread number for this instance of the routine. |
58 |
integer mythid |
59 |
integer myiter |
60 |
_RL mytime |
61 |
|
62 |
c == local variables == |
63 |
|
64 |
integer iloop |
65 |
|
66 |
#ifdef ALLOW_TAMC_CHECKPOINTING |
67 |
integer ilev_1 |
68 |
integer ilev_2 |
69 |
integer ilev_3 |
70 |
integer max_lev2 |
71 |
integer max_lev3 |
72 |
#endif |
73 |
|
74 |
c-- == end of interface == |
75 |
|
76 |
#ifdef ALLOW_AUTODIFF_TAMC |
77 |
c-- Initialize storage for the cost function evaluation. |
78 |
CADJ INIT dummytape = common, 1 |
79 |
c-- Initialize storage for the outermost loop. |
80 |
CADJ INIT tapelev3 = USER |
81 |
#ifdef ALLOW_TAMC_CHECKPOINTING |
82 |
nIter0 = INT( startTime/deltaTClock ) |
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 |
CADJ STORE gu = tapelev3, key = ilev_3 |
157 |
CADJ STORE gv = tapelev3, key = ilev_3 |
158 |
#ifdef INCLUDE_CD_CODE |
159 |
CADJ STORE uveld = tapelev3, key = ilev_3 |
160 |
CADJ STORE vveld = tapelev3, key = ilev_3 |
161 |
CADJ STORE unm1 = tapelev3, key = ilev_3 |
162 |
CADJ STORE vnm1 = tapelev3, key = ilev_3 |
163 |
CADJ STORE gucd = tapelev3, key = ilev_3 |
164 |
CADJ STORE gvcd = tapelev3, key = ilev_3 |
165 |
#endif |
166 |
#ifdef ALLOW_COST_TRACER |
167 |
CADJ STORE tr1 = tapelev3, key = ilev_3 |
168 |
CADJ STORE gtr1nm1 = tapelev3, key = ilev_3 |
169 |
#endif |
170 |
|
171 |
c-- Initialise storage for the middle loop. |
172 |
CADJ INIT tapelev2 = USER |
173 |
|
174 |
do ilev_2 = 1,nchklev_2 |
175 |
if(ilev_2.le.max_lev2) then |
176 |
CADJ STORE gsnm1 = tapelev2, key = ilev_2 |
177 |
CADJ STORE gtnm1 = tapelev2, key = ilev_2 |
178 |
CADJ STORE gunm1 = tapelev2, key = ilev_2 |
179 |
CADJ STORE gvnm1 = tapelev2, key = ilev_2 |
180 |
CADJ STORE theta = tapelev2, key = ilev_2 |
181 |
CADJ STORE salt = tapelev2, key = ilev_2 |
182 |
CADJ STORE uvel = tapelev2, key = ilev_2 |
183 |
CADJ STORE vvel = tapelev2, key = ilev_2 |
184 |
CADJ STORE wvel = tapelev2, key = ilev_2 |
185 |
CADJ STORE etan = tapelev2, key = ilev_2 |
186 |
CADJ STORE etanm1 = tapelev2, key = ilev_2 |
187 |
CADJ STORE gu = tapelev2, key = ilev_2 |
188 |
CADJ STORE gv = tapelev2, key = ilev_2 |
189 |
#ifdef INCLUDE_CD_CODE |
190 |
CADJ STORE uveld = tapelev2, key = ilev_2 |
191 |
CADJ STORE vveld = tapelev2, key = ilev_2 |
192 |
CADJ STORE unm1 = tapelev2, key = ilev_2 |
193 |
CADJ STORE vnm1 = tapelev2, key = ilev_2 |
194 |
CADJ STORE gucd = tapelev2, key = ilev_2 |
195 |
CADJ STORE gvcd = tapelev2, key = ilev_2 |
196 |
#endif |
197 |
#ifdef ALLOW_COST_TRACER |
198 |
CADJ STORE tr1 = tapelev2, key = ilev_2 |
199 |
CADJ STORE gtr1nm1 = tapelev2, key = ilev_2 |
200 |
#endif |
201 |
|
202 |
c-- Initialize storage for the innermost loop. |
203 |
c-- Always check common block sizes for the checkpointing! |
204 |
CADJ INIT comlev1 = COMMON,nchklev_1 |
205 |
CADJ INIT comlev1_bibj = COMMON,nchklev_1*nsx*nsy*nthreads_chkpt |
206 |
CADJ INIT comlev1_bibj_k = COMMON,nchklev_1*nsx*nsy*nr*nthreads_chkpt |
207 |
CADJ INIT comlev1_kpp = COMMON,nchklev_1*nsx*nsy |
208 |
|
209 |
C-- RG replace 2 by max of num_v_smooth_Ri |
210 |
CADJ INIT comlev1_kpp_sm = COMMON,nchklev_1*nsx*nsy*2 |
211 |
|
212 |
do ilev_1 = 1,nchklev_1 |
213 |
|
214 |
c-- The if-statement below introduces a some flexibility in the |
215 |
c-- choice of the 3-tupel ( nchklev_1, nchklev_2, nchklev_3 ). |
216 |
c-- |
217 |
c-- Requirement: nchklev_1*nchklev_2*nchklev_3 .ge. nTimeSteps . |
218 |
|
219 |
iloop = (ilev_3 - 1)*nchklev_2*nchklev_1 + |
220 |
& (ilev_2 - 1)*nchklev_1 + ilev_1 |
221 |
|
222 |
if ( iloop .le. nTimeSteps ) then |
223 |
|
224 |
#else /* ALLOW_TAMC_CHECKPOINTING undefined */ |
225 |
c-- Initialise storage for reference trajectory without TAMC check- |
226 |
c-- pointing. |
227 |
CADJ INIT history = USER |
228 |
CADJ INIT comlev1_bibj = COMMON,nchklev_0*nsx*nsy*nthreads_chkpt |
229 |
CADJ INIT comlev1_bibj_k = COMMON,nchklev_0*nsx*nsy*nr*nthreads_chkpt |
230 |
CADJ INIT comlev1_kpp = COMMON,nchklev_0*nsx*nsy |
231 |
|
232 |
C-- RG replace 2 by max of num_v_smooth_Ri |
233 |
CADJ INIT comlev1_kpp_sm = COMMON,nchklev_0*nsx*nsy*2 |
234 |
|
235 |
c-- Check the choice of the checkpointing parameters in relation |
236 |
c-- to nTimeSteps: (nchklev_0 .ge. nTimeSteps) |
237 |
if (nchklev_0 .lt. nTimeSteps) then |
238 |
print* |
239 |
print*, ' the_main_loop: TAMC checkpointing parameter ', |
240 |
& nchklev_0 = ', nchklev_0 |
241 |
print*, ' not consistent with nTimeSteps = ', |
242 |
& nTimeSteps |
243 |
stop ' ... stopped in the_main_loop.' |
244 |
endif |
245 |
|
246 |
DO iloop = 1, nTimeSteps |
247 |
|
248 |
#endif /* ALLOW_TAMC_CHECKPOINTING */ |
249 |
|
250 |
#else /* ALLOW_AUTODIFF_TAMC undefined */ |
251 |
|
252 |
c-- Start the main loop of adjoint_Objfunc. Automatic differentiation |
253 |
c-- NOT enabled. |
254 |
DO iloop = 1, nTimeSteps |
255 |
|
256 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
257 |
|
258 |
c-- >>> Loop body start <<< |
259 |
|
260 |
#ifdef ALLOW_TAMC_CHECKPOINTING |
261 |
nIter0 = INT( startTime/deltaTClock ) |
262 |
ikey_dynamics = ilev_1 |
263 |
#endif |
264 |
|
265 |
|
266 |
CALL TIMER_START('FORWARD_STEP [THE_MAIN_LOOP]',mythid) |
267 |
CALL FORWARD_STEP( iloop, mytime, myiter, mythid ) |
268 |
CALL TIMER_STOP ('FORWARD_STEP [THE_MAIN_LOOP]',mythid) |
269 |
|
270 |
#ifdef ALLOW_COST |
271 |
C-- compare model with data and compute cost function |
272 |
C-- this is done after exchanges to allow interpolation |
273 |
CALL TIMER_START('COST_TILE [THE_MAIN_LOOP]',myThid) |
274 |
CALL COST_TILE ( myThid ) |
275 |
CALL TIMER_STOP ('COST_TILE [THE_MAIN_LOOP]',myThid) |
276 |
#endif |
277 |
|
278 |
c-- >>> Loop body end <<< |
279 |
|
280 |
#ifdef ALLOW_AUTODIFF_TAMC |
281 |
#ifdef ALLOW_TAMC_CHECKPOINTING |
282 |
endif |
283 |
enddo |
284 |
endif |
285 |
enddo |
286 |
endif |
287 |
enddo |
288 |
#else |
289 |
enddo |
290 |
#endif |
291 |
|
292 |
#else |
293 |
enddo |
294 |
#endif |
295 |
|
296 |
#ifdef ALLOW_COST |
297 |
c-- Sum all cost function contributions. |
298 |
call TIMER_START('COST_FINAL [ADJOINT SPIN-DOWN]', mythid) |
299 |
call COST_FINAL ( mythid ) |
300 |
call TIMER_STOP ('COST_FINAL [ADJOINT SPIN-DOWN]', mythid) |
301 |
#endif |
302 |
|
303 |
_BARRIER |
304 |
CALL TIMER_STOP ('MAIN LOOP [THE_MAIN_LOOP]', mythid) |
305 |
|
306 |
END |