/[MITgcm]/MITgcm/model/src/main_do_loop.F
ViewVC logotype

Annotation of /MITgcm/model/src/main_do_loop.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.6 - (hide annotations) (download)
Thu Mar 21 18:11:58 2013 UTC (11 years, 2 months ago) by jahn
Branch: MAIN
Changes since 1.5: +9 -1 lines
include ptracers headers also with OpenAD

1 jahn 1.6 C $Header: /u/gcmpack/MITgcm/model/src/main_do_loop.F,v 1.5 2013/02/23 15:57:57 jmc Exp $
2 heimbach 1.1 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_GENERIC_ADVDIFF
13     # include "GAD_OPTIONS.h"
14     #endif
15     #ifdef ALLOW_GMREDI
16     # include "GMREDI_OPTIONS.h"
17     #endif
18     #ifdef ALLOW_STREAMICE
19     # include "STREAMICE_OPTIONS.h"
20     #endif
21     #ifdef ALLOW_GGL90
22     # include "GGL90_OPTIONS.h"
23     #endif
24     #ifdef ALLOW_EXF
25     # include "EXF_OPTIONS.h"
26     #endif
27     #ifdef ALLOW_CTRL
28     # include "CTRL_OPTIONS.h"
29     #endif
30    
31     CBOP
32     C !ROUTINE: MAIN_DO_LOOP
33     C !INTERFACE:
34 jmc 1.2 SUBROUTINE MAIN_DO_LOOP( myTime, myIter, myThid )
35 heimbach 1.1
36     C !DESCRIPTION: \bv
37     C *================================================================*
38     C | SUBROUTINE the_loop_body
39     C | o Run the ocean model and evaluate the specified cost function.
40     C *================================================================*
41     C |
42     C | MAIN_DO_LOOP is the toplevel routine for the Tangent Linear and
43     C | Adjoint Model Compiler (TAMC).
44     C | For this purpose the initialization
45     C | of the model was split into two parts. Those parameters that do
46     C | not depend on a specific model run are set in INITIALISE_FIXED,
47     C | whereas those that do depend on the specific realization are
48     C | initialized in INITIALISE_VARIA.
49     C | This routine is to be used in conjuction with the MITgcmuv
50     C | checkpoint 37.
51     C *================================================================*
52     C \ev
53    
54     C !USES:
55     IMPLICIT NONE
56     C == Global variables ==
57     #include "SIZE.h"
58     #include "EEPARAMS.h"
59     #include "PARAMS.h"
60    
61     c**************************************
62     #ifdef ALLOW_AUTODIFF
63     # ifndef ALLOW_AUTODIFF_OPENAD
64    
65     c These includes are needed for
66     c AD-checkpointing.
67     c They provide the fields to be stored.
68    
69     # include "AUTODIFF_MYFIELDS.h"
70     # include "GRID.h"
71     # include "DYNVARS.h"
72     # include "SURFACE.h"
73     # include "FFIELDS.h"
74     # include "EOS.h"
75     # include "AUTODIFF.h"
76    
77     # ifdef ALLOW_GENERIC_ADVDIFF
78     # include "GAD.h"
79     # include "GAD_SOM_VARS.h"
80     # endif
81     # ifdef ALLOW_MOM_FLUXFORM
82     # include "MOM_FLUXFORM.h"
83     # endif
84     # ifdef ALLOW_CD_CODE
85     # include "CD_CODE_VARS.h"
86     # endif
87     # ifdef ALLOW_PTRACERS
88     # include "PTRACERS_SIZE.h"
89     # include "PTRACERS_FIELDS.h"
90     # include "PTRACERS_START.h"
91     # endif
92     # ifdef ALLOW_GCHEM
93     # include "GCHEM_FIELDS.h"
94     # endif
95     # ifdef ALLOW_CFC
96     # include "CFC.h"
97     # endif
98     # ifdef ALLOW_DIC
99     # include "DIC_VARS.h"
100     # include "DIC_LOAD.h"
101     # include "DIC_ATMOS.h"
102     # include "DIC_CTRL.h"
103     # include "DIC_COST.h"
104     # endif
105     # ifdef ALLOW_OBCS
106     # include "OBCS_PARAMS.h"
107     # include "OBCS_FIELDS.h"
108     # include "OBCS_SEAICE.h"
109     # ifdef ALLOW_PTRACERS
110     # include "OBCS_PTRACERS.h"
111     # endif
112     # endif
113     # ifdef ALLOW_EXF
114     # include "EXF_FIELDS.h"
115     # ifdef ALLOW_BULKFORMULAE
116     # include "EXF_CONSTANTS.h"
117     # endif
118     # endif /* ALLOW_EXF */
119     # ifdef ALLOW_SEAICE
120     # include "SEAICE_SIZE.h"
121     # include "SEAICE.h"
122     # include "SEAICE_PARAMS.h"
123     # include "SEAICE_COST.h"
124     # include "SEAICE_TRACER.h"
125     # endif
126     # ifdef ALLOW_SALT_PLUME
127     # include "SALT_PLUME.h"
128     # endif
129     # ifdef ALLOW_THSICE
130     # include "THSICE_SIZE.h"
131     # include "THSICE_VARS.h"
132     # endif
133     # ifdef ALLOW_SHELFICE
134     # include "SHELFICE.h"
135     # include "SHELFICE_COST.h"
136     # endif
137     # ifdef ALLOW_STREAMICE
138     # include "STREAMICE.h"
139     # include "STREAMICE_ADV.h"
140     # include "STREAMICE_BDRY.h"
141     # include "STREAMICE_CG.h"
142     # endif
143     # ifdef ALLOW_EBM
144     # include "EBM.h"
145     # endif
146     # ifdef ALLOW_RBCS
147     # include "RBCS_SIZE.h"
148     # include "RBCS_FIELDS.h"
149     # endif
150     # ifdef ALLOW_OFFLINE
151     # include "OFFLINE.h"
152     # endif
153     # ifdef ALLOW_CG2D_NSA
154     # include "CG2D.h"
155     # endif
156     # ifdef ALLOW_DIVIDED_ADJOINT_MPI
157     # include "mpif.h"
158     # endif
159    
160     # include "tamc.h"
161    
162     # ifdef ALLOW_GGL90
163     # include "GGL90.h"
164     # endif
165     # ifdef ALLOW_PROFILES
166     # include "profiles.h"
167     # endif
168    
169     # ifdef ALLOW_ECCO_EVOLUTION
170     # ifdef ALLOW_ECCO
171     # include "ecco_cost.h"
172     # endif
173     # endif
174    
175 jahn 1.6 # else /* undef ALLOW_AUTODIFF_OPENAD */
176    
177     # ifdef ALLOW_PTRACERS
178     # include "PTRACERS_SIZE.h"
179     # include "PTRACERS_FIELDS.h"
180     # include "PTRACERS_START.h"
181     # endif
182    
183 heimbach 1.1 # endif /* undef ALLOW_AUTODIFF_OPENAD */
184    
185     # ifdef ALLOW_CTRL
186     # include "CTRL_SIZE.h"
187     # include "ctrl.h"
188     # include "ctrl_dummy.h"
189     # include "CTRL_GENARR.h"
190     # endif
191     # ifdef ALLOW_COST
192     # include "cost.h"
193     # endif
194    
195     #endif /* ALLOW_AUTODIFF */
196     c**************************************
197    
198     C !INPUT/OUTPUT PARAMETERS:
199     C == Routine arguments ==
200     C note: under the multi-threaded model myiter and
201     C mytime are local variables passed around as routine
202     C arguments. Although this is fiddly it saves the need to
203     C impose additional synchronisation points when they are
204     C updated.
205     C myTime :: time counter for this thread
206     C myIter :: iteration counter for this thread
207     C myThid :: thread number for this instance of the routine.
208     _RL myTime
209     INTEGER myIter
210     INTEGER myThid
211    
212     C !FUNCTIONS:
213     C == Functions ==
214    
215     C !LOCAL VARIABLES:
216     C == Local variables ==
217 jmc 1.4 INTEGER iloop
218 heimbach 1.1 #ifdef ALLOW_AUTODIFF_TAMC
219     #ifdef STORE_LOADEDREC_TEST
220 jmc 1.4 INTEGER bi,bj
221 heimbach 1.1 #endif /* STORE_LOADEDREC_TEST */
222     #endif
223    
224     CEOP
225    
226     #ifdef ALLOW_DEBUG
227     IF (debugMode) CALL DEBUG_ENTER('MAIN_DO_LOOP',myThid)
228     #endif
229    
230     c >>>>>>>>>>>>>>>>>>>>>>>>>>> LOOP <<<<<<<<<<<<<<<<<<<<<<<<<<<<
231     c >>>>>>>>>>>>>>>>>>>>>>>>>>> STARTS <<<<<<<<<<<<<<<<<<<<<<<<<<<<
232    
233     c++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
234     #ifndef ALLOW_AUTODIFF_OPENAD
235     c++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
236     # ifdef ALLOW_AUTODIFF
237     # ifdef ALLOW_TAMC_CHECKPOINTING
238    
239     max_lev4=nTimeSteps/(nchklev_1*nchklev_2*nchklev_3)+1
240     max_lev3=nTimeSteps/(nchklev_1*nchklev_2)+1
241     max_lev2=nTimeSteps/nchklev_1+1
242    
243     c**************************************
244     # ifdef ALLOW_DIVIDED_ADJOINT
245     CADJ loop = divided
246     # endif
247     c**************************************
248    
249     # ifdef AUTODIFF_4_LEVEL_CHECKPOINT
250     do ilev_4 = 1,nchklev_4
251     if(ilev_4.le.max_lev4) then
252     c**************************************
253     #ifdef ALLOW_AUTODIFF_WHTAPEIO
254     CALL AUTODIFF_WHTAPEIO_SYNC( 4 , 0, mythid )
255     #endif
256     CALL AUTODIFF_STORE( myThid )
257     #include "checkpoint_lev4_directives.h"
258     CALL AUTODIFF_RESTORE( myThid )
259     #ifdef ALLOW_AUTODIFF_WHTAPEIO
260     CALL AUTODIFF_WHTAPEIO_SYNC( 4 , 1, mythid )
261     #endif
262     c**************************************
263     c-- Initialise storage for the middle loop.
264     CADJ INIT tapelev3 = USER
265     # endif /* AUTODIFF_4_LEVEL_CHECKPOINT */
266    
267     # ifndef AUTODIFF_2_LEVEL_CHECKPOINT
268     do ilev_3 = 1,nchklev_3
269     if(ilev_3.le.max_lev3) then
270     c**************************************
271     #ifdef ALLOW_AUTODIFF_WHTAPEIO
272     CALL AUTODIFF_WHTAPEIO_SYNC( 3 , 0, mythid )
273     #endif
274     CALL AUTODIFF_STORE( myThid )
275     #include "checkpoint_lev3_directives.h"
276     CALL AUTODIFF_RESTORE( myThid )
277     #ifdef ALLOW_AUTODIFF_WHTAPEIO
278     CALL AUTODIFF_WHTAPEIO_SYNC( 3 , 1, mythid )
279     #endif
280     c**************************************
281     c-- Initialise storage for the middle loop.
282     CADJ INIT tapelev2 = USER
283     # endif /* AUTODIFF_2_LEVEL_CHECKPOINT */
284    
285     do ilev_2 = 1,nchklev_2
286     if(ilev_2.le.max_lev2) then
287     c**************************************
288     #ifdef ALLOW_AUTODIFF_WHTAPEIO
289     CALL AUTODIFF_WHTAPEIO_SYNC( 2 , 0, mythid )
290     #endif
291     CALL AUTODIFF_STORE( myThid )
292     #include "checkpoint_lev2_directives.h"
293     CALL AUTODIFF_RESTORE( myThid )
294     #ifdef ALLOW_AUTODIFF_WHTAPEIO
295     CALL AUTODIFF_WHTAPEIO_SYNC( 2 , 1, mythid )
296     #endif
297     c**************************************
298    
299     # endif /* ALLOW_TAMC_CHECKPOINTING */
300    
301     c**************************************
302     c--
303     c-- Initialize storage for the innermost loop.
304     c-- Always check common block sizes for the checkpointing!
305     c--
306     CADJ INIT comlev1 = COMMON,nchklev_1
307     CADJ INIT comlev1_bibj = COMMON,nchklev_1*nsx*nsy*nthreads_chkpt
308     CADJ INIT comlev1_bibj_k = COMMON,nchklev_1*nsx*nsy*nr*nthreads_chkpt
309     c--
310     # ifdef ALLOW_KPP
311     CADJ INIT comlev1_kpp = COMMON,nchklev_1*nsx*nsy
312     CADJ INIT comlev1_kpp_k = COMMON,nchklev_1*nsx*nsy*nr
313     # endif /* ALLOW_KPP */
314     c--
315     # ifdef ALLOW_GMREDI
316     CADJ INIT comlev1_gmredi_k_gad
317     CADJ & = COMMON,nchklev_1*nsx*nsy*nr*nthreads_chkpt*maxpass
318     # endif /* ALLOW_GMREDI */
319     c--
320     # ifdef ALLOW_PTRACERS
321     CADJ INIT comlev1_bibj_ptracers = COMMON,
322     CADJ & nchklev_1*nsx*nsy*nthreads_chkpt*PTRACERS_num
323     CADJ INIT comlev1_bibj_k_ptracers = COMMON,
324     CADJ & nchklev_1*nsx*nsy*nthreads_chkpt*PTRACERS_num*nr
325     # endif /* ALLOW_PTRACERS */
326     c--
327     # ifndef DISABLE_MULTIDIM_ADVECTION
328     CADJ INIT comlev1_bibj_gad = COMMON,
329     CADJ & nchklev_1*nsx*nsy*nthreads_chkpt*maxpass
330     CADJ INIT comlev1_bibj_k_gad = COMMON,
331     CADJ & nchklev_1*nsx*nsy*nr*nthreads_chkpt*maxpass
332     CADJ INIT comlev1_bibj_k_gad_pass = COMMON,
333     CADJ & nchklev_1*nsx*nsy*nr*nthreads_chkpt*maxpass*maxpass
334     # endif /* DISABLE_MULTIDIM_ADVECTION */
335     c--
336     # ifdef ALLOW_MOM_COMMON
337     # ifndef AUTODIFF_DISABLE_LEITH
338     CADJ INIT comlev1_mom_ijk_loop
339     CADJ & = COMMON,nchklev_1*
340     CADJ & (snx+2*olx)*nsx*(sny+2*oly)*nsy*nr*nthreads_chkpt
341     # endif /* AUTODIFF_DISABLE_LEITH */
342     # endif /* ALLOW_MOM_COMMON */
343     c--
344     # if (defined (ALLOW_EXF) && defined (ALLOW_BULKFORMULAE))
345     CADJ INIT comlev1_exf_1
346     CADJ & = COMMON,nchklev_1*snx*nsx*sny*nsy*nthreads_chkpt
347     CADJ INIT comlev1_exf_2
348     CADJ & = COMMON,niter_bulk*nchklev_1*snx*nsx*sny*nsy*nthreads_chkpt
349     # endif /* ALLOW_BULKFORMULAE */
350     c--
351     # ifdef ALLOW_SEAICE
352     # ifdef SEAICE_ALLOW_DYNAMICS
353     CADJ INIT comlev1_dynsol = COMMON,nchklev_1*MPSEUDOTIMESTEPS
354     # ifdef SEAICE_LSR_ADJOINT_ITER
355     CADJ INIT comlev1_dyniter =
356     CADJ & COMMON,nchklev_1*MPSEUDOTIMESTEPS*SOLV_MAX_FIXED
357     # endif
358     # endif
359     # ifdef SEAICE_ALLOW_EVP
360     CADJ INIT comlev1_evp = COMMON,nEVPstepMax*nchklev_1
361     # endif
362     # ifdef SEAICE_MULTICATEGORY
363     CADJ INIT comlev1_multdim
364     CADJ & = COMMON,nchklev_1*nsx*nsy*nthreads_chkpt*multdim
365     # endif
366     # ifndef DISABLE_MULTIDIM_ADVECTION
367     CADJ INIT comlev1_bibj_k_gadice = COMMON,
368     CADJ & nchklev_1*nsx*nsy*nthreads_chkpt*maxpass
369     CADJ INIT comlev1_bibj_k_gadice_pass = COMMON,
370     CADJ & nchklev_1*nsx*nsy*nthreads_chkpt*maxpass*maxpass
371     # endif /* DISABLE_MULTIDIM_ADVECTION */
372     # endif /* ALLOW_SEAICE */
373     c--
374     # ifdef ALLOW_THSICE
375     CADJ INIT comlev1_thsice_1
376     CADJ & = COMMON,nchklev_1*snx*nsx*sny*nsy*nthreads_chkpt
377     CADJ INIT comlev1_thsice_2
378     CADJ & = COMMON,nchklev_1*snx*nsx*sny*nsy*nlyr*nthreads_chkpt
379     CADJ INIT comlev1_thsice_3
380     CADJ & = COMMON,nchklev_1*snx*nsx*sny*nsy*MaxTsf*nthreads_chkpt
381     CADJ INIT comlev1_thsice_4
382     CADJ & = COMMON,nchklev_1*nsx*nsy*maxpass*nthreads_chkpt
383     # endif /* ALLOW_THSICE */
384     c--
385     # ifdef ALLOW_STREAMICE
386     CADJ INIT comlev1_stream_nl = COMMON,nchklev_1*streamice_max_nl
387     CADJ INIT comlev1_stream_front = COMMON,nchklev_1*4
388     CADJ INIT comlev1_stream_ij
389     CADJ & = COMMON,nchklev_1*4*(snx+2)*nsx*(sny+2)*nsy
390     CADJ INIT comlev1_stream_hybrid
391     CADJ & = COMMON,nchklev_1*snx*nsx*sny*nsy*nr*nthreads_chkpt
392     # endif
393     c--
394     # ifdef ALLOW_CG2D_NSA
395     CADJ INIT comlev1_cg2d
396     CADJ & = COMMON,nchklev_1*nthreads_chkpt
397     CADJ INIT comlev1_cg2d_iter
398     CADJ & = COMMON,nchklev_1*nthreads_chkpt*numItersMax
399     # endif
400     c--
401     c**************************************
402    
403     #ifdef STORE_LOADEDREC_TEST
404     DO bj = myByLo(myThid), myByHi(myThid)
405     DO bi = myBxLo(myThid), myBxHi(myThid)
406     loadedRec(bi,bj) = 0
407     ENDDO
408     ENDDO
409     #endif /* STORE_LOADEDREC_TEST */
410    
411     #ifdef ALLOW_TAMC_CHECKPOINTING
412    
413     do ilev_1 = 1,nchklev_1
414    
415     c-- The if-statement below introduces a some flexibility in the
416     c-- choice of the 3-tupel ( nchklev_1, nchklev_2, nchklev_3 ).
417    
418     iloop = (ilev_2 - 1)*nchklev_1 + ilev_1
419     # ifndef AUTODIFF_2_LEVEL_CHECKPOINT
420     & + (ilev_3 - 1)*nchklev_2*nchklev_1
421     # endif
422     # ifdef AUTODIFF_4_LEVEL_CHECKPOINT
423     & + (ilev_4 - 1)*nchklev_3*nchklev_2*nchklev_1
424     # endif
425    
426     if ( iloop .le. nTimeSteps ) then
427    
428     # else /* ALLOW_TAMC_CHECKPOINTING undefined */
429    
430     DO iloop = 1, nTimeSteps
431    
432     # endif /* ALLOW_TAMC_CHECKPOINTING */
433     # endif /* ALLOW_AUTODIFF */
434    
435     #endif /* undef ALLOW_AUTODIFF_OPENAD */
436    
437     #ifdef ALLOW_AUTODIFF_OPENAD
438 utke 1.3 DO iloop = 1, nTimeSteps
439 heimbach 1.1 #endif
440    
441     #ifndef ALLOW_AUTODIFF
442    
443     c-- Start the main loop of adjoint_Objfunc. Automatic differentiation
444     c-- NOT enabled.
445     DO iloop = 1, nTimeSteps
446    
447     #endif /* ALLOW_AUTODIFF */
448    
449     c-- >>> Loop body start <<<
450    
451     #ifdef ALLOW_AUTODIFF_TAMC
452     nIter0 = NINT( (startTime-baseTime)/deltaTClock )
453 utke 1.3 # ifndef ALLOW_AUTODIFF_OPENAD
454 heimbach 1.1 ikey_dynamics = ilev_1
455 utke 1.3 # endif
456 heimbach 1.1 #endif
457    
458     #ifdef ALLOW_ECCO
459     #ifdef ALLOW_ECCO_EVOLUTION
460     #ifdef ALLOW_DEBUG
461     IF (debugMode) CALL DEBUG_CALL('cost_averagesfields',myThid)
462     #endif
463     c-- Accumulate time averages of temperature, salinity
464     #ifdef ALLOW_AUTODIFF
465     C-- Reset the model iteration counter and the model time.
466     myIter = nIter0 + (iloop-1)
467 jmc 1.2 myTime = startTime + float(iloop-1)*deltaTClock
468 heimbach 1.1 #endif
469     CALL TIMER_START('COST_AVERAGESFIELDS [MAIN_DO_LOOP]',mythid)
470     CALL COST_AVERAGESFIELDS( mytime, mythid )
471     CALL TIMER_STOP ('COST_AVERAGESFIELDS [MAIN_DO_LOOP]',mythid)
472     #endif /* ALLOW_ECCO_EVOLUTION */
473     #endif /* ALLOW_ECCO */
474    
475     #ifdef ALLOW_PROFILES
476     IF (usePROFILES) THEN
477     #ifdef ALLOW_DEBUG
478     IF (debugMode) CALL DEBUG_CALL('profiles_inloop',myThid)
479     #endif
480     c-- Accumulate in-situ time averages of theta, salt, and SSH.
481     #ifdef ALLOW_AUTODIFF
482     C-- Reset the model iteration counter and the model time.
483     myIter = nIter0 + (iloop-1)
484 jmc 1.2 myTime = startTime + float(iloop-1)*deltaTClock
485 heimbach 1.1 #endif
486     CALL TIMER_START('PROFILES_INLOOP [MAIN_DO_LOOP]', mythid)
487     CALL PROFILES_INLOOP( mytime, mythid )
488     CALL TIMER_STOP ('PROFILES_INLOOP [MAIN_DO_LOOP]', mythid)
489     ENDIF
490     #endif
491    
492     #ifdef ALLOW_DEBUG
493     IF (debugMode) CALL DEBUG_CALL('FORWARD_STEP',myThid)
494     #endif
495    
496     #ifdef ALLOW_ATM2D
497     CALL TIMER_START('FORWARD_STEP_ATM2D [MAIN_DO_LOOP]',mythid)
498     CALL FORWARD_STEP_ATM2D( iloop, mytime, myiter, mythid )
499     CALL TIMER_STOP ('FORWARD_STEP_ATM2D [MAIN_DO_LOOP]',mythid)
500     #else
501     CALL TIMER_START('FORWARD_STEP [MAIN_DO_LOOP]',mythid)
502     CALL FORWARD_STEP( iloop, mytime, myiter, mythid )
503     CALL TIMER_STOP ('FORWARD_STEP [MAIN_DO_LOOP]',mythid)
504     #endif
505    
506     c-- >>> Loop body end <<<
507     #ifdef ALLOW_AUTODIFF
508     # ifndef ALLOW_AUTODIFF_OPENAD
509    
510     # ifdef ALLOW_TAMC_CHECKPOINTING
511     endif
512     enddo
513     endif
514     enddo
515     # ifndef AUTODIFF_2_LEVEL_CHECKPOINT
516     endif
517     enddo
518     # endif
519     # ifdef AUTODIFF_4_LEVEL_CHECKPOINT
520     endif
521     enddo
522     # endif
523     # else /* ndef ALLOW_TAMC_CHECKPOINTING */
524 jmc 1.4 ENDDO
525 heimbach 1.1 # endif /* ALLOW_TAMC_CHECKPOINTING */
526 jmc 1.5 # else /* ndef ALLOW_AUTODIFF_OPENAD */
527 jmc 1.4 ENDDO
528 heimbach 1.1 # endif /* ALLOW_AUTODIFF_OPENAD */
529     #else /* ALLOW_AUTODIFF */
530 jmc 1.4 ENDDO
531 heimbach 1.1 #endif /* ALLOW_AUTODIFF */
532    
533     #ifdef ALLOW_DEBUG
534     IF (debugMode) CALL DEBUG_LEAVE('MAIN_DO_LOOP',myThid)
535     #endif
536    
537     RETURN
538     END

  ViewVC Help
Powered by ViewVC 1.1.22