1 |
C $Header: /u/gcmpack/MITgcm/model/src/the_main_loop.F,v 1.23 2001/11/20 21:14:41 heimbach Exp $ |
2 |
|
3 |
#include "CPP_OPTIONS.h" |
4 |
|
5 |
CBOP |
6 |
C !ROUTINE: THE_MAIN_LOOP |
7 |
C !INTERFACE: |
8 |
SUBROUTINE THE_MAIN_LOOP( mytime, myiter, mythid ) |
9 |
|
10 |
C !DESCRIPTION: \bv |
11 |
C *================================================================* |
12 |
C | SUBROUTINE the_main_loop |
13 |
C | o Run the ocean model and evaluate the specified cost function. |
14 |
C *================================================================* |
15 |
C | |
16 |
C | THE_MAIN_LOOP is the toplevel routine for the Tangent Linear and |
17 |
C | Adjoint Model Compiler (TAMC). For this purpose the initialization |
18 |
C | of the model was split into two parts. Those parameters that do |
19 |
C | not depend on a specific model run are set in INITIALISE_FIXED, |
20 |
C | whereas those that do depend on the specific realization are |
21 |
C | initialized in INITIALISE_VARIA. |
22 |
C | This routine is to be used in conjuction with the MITgcmuv |
23 |
C | checkpoint 37. |
24 |
C *================================================================* |
25 |
C \ev |
26 |
|
27 |
C !USES: |
28 |
IMPLICIT NONE |
29 |
C == Global variables == |
30 |
#include "SIZE.h" |
31 |
#include "EEPARAMS.h" |
32 |
#include "PARAMS.h" |
33 |
#ifdef ALLOW_AUTODIFF_TAMC |
34 |
# include "tamc.h" |
35 |
# include "ctrl.h" |
36 |
# include "ctrl_dummy.h" |
37 |
# include "cost.h" |
38 |
# include "DYNVARS.h" |
39 |
# include "FFIELDS.h" |
40 |
# include "GAD.h" |
41 |
# ifdef ALLOW_PASSIVE_TRACER |
42 |
# include "TR1.h" |
43 |
# endif |
44 |
# ifdef ALLOW_NONHYDROSTATIC |
45 |
# include "CG3D.h" |
46 |
# endif |
47 |
# ifdef EXACT_CONSERV |
48 |
# include "SURFACE.h" |
49 |
# endif |
50 |
#endif |
51 |
|
52 |
C !INPUT/OUTPUT PARAMETERS: |
53 |
C == Routine arguments == |
54 |
C note: under the multi-threaded model myiter and |
55 |
C mytime are local variables passed around as routine |
56 |
C arguments. Although this is fiddly it saves the need to |
57 |
C impose additional synchronisation points when they are |
58 |
C updated. |
59 |
C myIter - iteration counter for this thread |
60 |
C myTime - time counter for this thread |
61 |
C myThid - thread number for this instance of the routine. |
62 |
INTEGER myThid |
63 |
INTEGER myIter |
64 |
_RL myTime |
65 |
|
66 |
C !LOCAL VARIABLES: |
67 |
C == Local variables == |
68 |
integer iloop |
69 |
#ifdef ALLOW_TAMC_CHECKPOINTING |
70 |
integer ilev_1 |
71 |
integer ilev_2 |
72 |
integer ilev_3 |
73 |
integer max_lev2 |
74 |
integer max_lev3 |
75 |
#endif |
76 |
CEOP |
77 |
|
78 |
#ifdef ALLOW_AUTODIFF_TAMC |
79 |
c-- Initialize storage for the cost function evaluation. |
80 |
CADJ INIT dummytape = common, 1 |
81 |
c-- Initialize storage for the outermost loop. |
82 |
CADJ INIT tapelev3 = USER |
83 |
#ifdef ALLOW_TAMC_CHECKPOINTING |
84 |
nIter0 = INT( startTime/deltaTClock ) |
85 |
ikey_dynamics = 1 |
86 |
#endif |
87 |
CALL TIMER_START('ADJOINT SPIN-UP', mythid) |
88 |
#endif |
89 |
|
90 |
C-- Set initial conditions (variable arrays) |
91 |
CALL TIMER_START('INITIALISE_VARIA [THE_MAIN_LOOP]', mythid) |
92 |
CALL INITIALISE_VARIA( mythid ) |
93 |
CALL TIMER_STOP ('INITIALISE_VARIA [THE_MAIN_LOOP]', mythid) |
94 |
|
95 |
#ifndef ALLOW_AUTODIFF_TAMC |
96 |
c-- Dump for start state. |
97 |
CALL TIMER_START('WRITE_STATE [THE_MAIN_LOOP]', mythid) |
98 |
CALL WRITE_STATE( mytime, myiter, mythid ) |
99 |
CALL TIMER_STOP ('WRITE_STATE [THE_MAIN_LOOP]', mythid) |
100 |
#endif |
101 |
|
102 |
#ifndef EXCLUDE_MONITOR |
103 |
C-- Check status of solution (statistics, cfl, etc...) |
104 |
CALL TIMER_START('MONITOR [THE_MAIN_LOOP]', mythid) |
105 |
CALL MONITOR( myIter, myTime, myThid ) |
106 |
CALL TIMER_STOP ('MONITOR [THE_MAIN_LOOP]', mythid) |
107 |
#endif /* EXCLUDE_MONITOR */ |
108 |
|
109 |
#ifdef ALLOW_ADJOINT_RUN |
110 |
c-- Add control vector for forcing and parameter fields |
111 |
CALL CTRL_MAP_FORCING (mythid) |
112 |
#endif |
113 |
|
114 |
#ifdef ALLOW_AUTODIFF_TAMC |
115 |
CALL TIMER_STOP ('ADJOINT SPIN-UP', mythid) |
116 |
_BARRIER |
117 |
#endif |
118 |
|
119 |
c-- Do the model integration. |
120 |
CALL TIMER_START('MAIN LOOP [THE_MAIN_LOOP]', mythid) |
121 |
|
122 |
c >>>>>>>>>>>>>>>>>>>>>>>>>>> LOOP <<<<<<<<<<<<<<<<<<<<<<<<<<<< |
123 |
c >>>>>>>>>>>>>>>>>>>>>>>>>>> STARTS <<<<<<<<<<<<<<<<<<<<<<<<<<<< |
124 |
|
125 |
#ifdef ALLOW_AUTODIFF_TAMC |
126 |
#ifdef ALLOW_TAMC_CHECKPOINTING |
127 |
c-- Implement a three level checkpointing. For a two level |
128 |
c-- checkpointing delete the middle loop; for n levels (n > 3) |
129 |
c-- insert more loops. |
130 |
|
131 |
c-- Check the choice of the checkpointing parameters in relation |
132 |
c-- to nTimeSteps: (nchklev_1*nchklev_2*nchklev_3 .ge. nTimeSteps) |
133 |
if (nchklev_1*nchklev_2*nchklev_3 .lt. nTimeSteps) then |
134 |
print* |
135 |
print*, ' the_main_loop: TAMC checkpointing parameters' |
136 |
print*, ' nchklev_1*nchklev_2*nchklev_3 = ', |
137 |
& nchklev_1*nchklev_2*nchklev_3 |
138 |
print*, ' are not consistent with nTimeSteps = ', |
139 |
& nTimeSteps |
140 |
stop ' ... stopped in the_main_loop.' |
141 |
endif |
142 |
max_lev3=nTimeSteps/(nchklev_1*nchklev_2)+1 |
143 |
max_lev2=nTimeSteps/nchklev_1+1 |
144 |
|
145 |
do ilev_3 = 1,nchklev_3 |
146 |
if(ilev_3.le.max_lev3) then |
147 |
c************************************** |
148 |
#include "checkpoint_lev3_directives.h" |
149 |
c************************************** |
150 |
|
151 |
c-- Initialise storage for the middle loop. |
152 |
CADJ INIT tapelev2 = USER |
153 |
|
154 |
do ilev_2 = 1,nchklev_2 |
155 |
if(ilev_2.le.max_lev2) then |
156 |
c************************************** |
157 |
#include "checkpoint_lev2_directives.h" |
158 |
c************************************** |
159 |
|
160 |
c-- Initialize storage for the innermost loop. |
161 |
c-- Always check common block sizes for the checkpointing! |
162 |
CADJ INIT comlev1 = COMMON,nchklev_1 |
163 |
CADJ INIT comlev1_bibj = COMMON,nchklev_1*nsx*nsy*nthreads_chkpt |
164 |
CADJ INIT comlev1_bibj_k = COMMON,nchklev_1*nsx*nsy*nr*nthreads_chkpt |
165 |
CADJ INIT comlev1_kpp = COMMON,nchklev_1*nsx*nsy |
166 |
#ifndef DISABLE_MULTIDIM_ADVECTION |
167 |
CADJ INIT comlev1_bibj_pass |
168 |
CADJ & = COMMON,nchklev_1*nsx*nsy*nr*nthreads_chkpt*maxpass |
169 |
#endif /* DISABLE_MULTIDIM_ADVECTION */ |
170 |
#ifdef ALLOW_BULKFORMULAE |
171 |
CADJ INIT comlev1_exf_1 |
172 |
CADJ & = COMMON,nchklev_1*snx*nsx*sny*nsy*nthreads_chkpt |
173 |
CADJ INIT comlev1_exf_2 |
174 |
CADJ & = COMMON,niter_bulk*nchklev_1*snx*nsx*sny*nsy*nthreads_chkpt |
175 |
#endif |
176 |
|
177 |
do ilev_1 = 1,nchklev_1 |
178 |
|
179 |
c-- The if-statement below introduces a some flexibility in the |
180 |
c-- choice of the 3-tupel ( nchklev_1, nchklev_2, nchklev_3 ). |
181 |
c-- |
182 |
c-- Requirement: nchklev_1*nchklev_2*nchklev_3 .ge. nTimeSteps . |
183 |
|
184 |
iloop = (ilev_3 - 1)*nchklev_2*nchklev_1 + |
185 |
& (ilev_2 - 1)*nchklev_1 + ilev_1 |
186 |
|
187 |
if ( iloop .le. nTimeSteps ) then |
188 |
|
189 |
#else /* ALLOW_TAMC_CHECKPOINTING undefined */ |
190 |
c-- Initialise storage for reference trajectory without TAMC check- |
191 |
c-- pointing. |
192 |
CADJ INIT history = USER |
193 |
CADJ INIT comlev1_bibj = COMMON,nchklev_0*nsx*nsy*nthreads_chkpt |
194 |
CADJ INIT comlev1_bibj_k = COMMON,nchklev_0*nsx*nsy*nr*nthreads_chkpt |
195 |
CADJ INIT comlev1_kpp = COMMON,nchklev_0*nsx*nsy |
196 |
|
197 |
C-- RG replace 2 by max of num_v_smooth_Ri |
198 |
CADJ INIT comlev1_kpp_sm = COMMON,nchklev_0*nsx*nsy*2 |
199 |
|
200 |
c-- Check the choice of the checkpointing parameters in relation |
201 |
c-- to nTimeSteps: (nchklev_0 .ge. nTimeSteps) |
202 |
if (nchklev_0 .lt. nTimeSteps) then |
203 |
print* |
204 |
print*, ' the_main_loop: TAMC checkpointing parameter ', |
205 |
& 'nchklev_0 = ', nchklev_0 |
206 |
print*, ' not consistent with nTimeSteps = ', |
207 |
& nTimeSteps |
208 |
stop ' ... stopped in the_main_loop.' |
209 |
endif |
210 |
|
211 |
DO iloop = 1, nTimeSteps |
212 |
|
213 |
#endif /* ALLOW_TAMC_CHECKPOINTING */ |
214 |
|
215 |
#else /* ALLOW_AUTODIFF_TAMC undefined */ |
216 |
|
217 |
c-- Start the main loop of adjoint_Objfunc. Automatic differentiation |
218 |
c-- NOT enabled. |
219 |
DO iloop = 1, nTimeSteps |
220 |
|
221 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
222 |
|
223 |
c-- >>> Loop body start <<< |
224 |
|
225 |
#ifdef ALLOW_TAMC_CHECKPOINTING |
226 |
nIter0 = INT( startTime/deltaTClock ) |
227 |
ikey_dynamics = ilev_1 |
228 |
#endif |
229 |
|
230 |
|
231 |
CALL TIMER_START('FORWARD_STEP [THE_MAIN_LOOP]',mythid) |
232 |
CALL FORWARD_STEP( iloop, mytime, myiter, mythid ) |
233 |
CALL TIMER_STOP ('FORWARD_STEP [THE_MAIN_LOOP]',mythid) |
234 |
|
235 |
#ifdef ALLOW_COST |
236 |
C-- compare model with data and compute cost function |
237 |
C-- this is done after exchanges to allow interpolation |
238 |
CALL TIMER_START('COST_TILE [THE_MAIN_LOOP]',myThid) |
239 |
CALL COST_TILE ( myThid ) |
240 |
CALL TIMER_STOP ('COST_TILE [THE_MAIN_LOOP]',myThid) |
241 |
#endif |
242 |
|
243 |
c-- >>> Loop body end <<< |
244 |
|
245 |
#ifdef ALLOW_AUTODIFF_TAMC |
246 |
#ifdef ALLOW_TAMC_CHECKPOINTING |
247 |
endif |
248 |
enddo |
249 |
endif |
250 |
enddo |
251 |
endif |
252 |
enddo |
253 |
#else |
254 |
enddo |
255 |
#endif |
256 |
|
257 |
#else |
258 |
enddo |
259 |
#endif |
260 |
|
261 |
#ifdef ALLOW_COST |
262 |
c-- Sum all cost function contributions. |
263 |
call TIMER_START('COST_FINAL [ADJOINT SPIN-DOWN]', mythid) |
264 |
call COST_FINAL ( mythid ) |
265 |
call TIMER_STOP ('COST_FINAL [ADJOINT SPIN-DOWN]', mythid) |
266 |
#endif |
267 |
|
268 |
_BARRIER |
269 |
CALL TIMER_STOP ('MAIN LOOP [THE_MAIN_LOOP]', mythid) |
270 |
|
271 |
END |