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

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

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


Revision 1.18 - (hide annotations) (download)
Mon Aug 27 18:50:41 2001 UTC (22 years, 8 months ago) by jmc
Branch: MAIN
Changes since 1.17: +20 -1 lines
modified to incorporate NonLin-FreeSurf

1 jmc 1.18 C $Header: /u/gcmpack/models/MITgcmUV/model/src/forward_step.F,v 1.17 2001/08/14 00:20:48 heimbach Exp $
2 adcroft 1.15 C $Name: $
3 adcroft 1.1
4     #include "CPP_OPTIONS.h"
5    
6 adcroft 1.13 SUBROUTINE FORWARD_STEP( iloop, myTime, myIter, myThid )
7     IMPLICIT NONE
8 heimbach 1.12
9     c ==================================================================
10     c SUBROUTINE forward_step
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 adcroft 1.1
22 adcroft 1.13 C == global variables ==
23 adcroft 1.1 #include "SIZE.h"
24     #include "EEPARAMS.h"
25     #include "PARAMS.h"
26     #include "DYNVARS.h"
27 heimbach 1.12 #include "FFIELDS.h"
28    
29 adcroft 1.1 #ifdef ALLOW_NONHYDROSTATIC
30     #include "CG3D.h"
31     #endif
32    
33 heimbach 1.12 #ifdef ALLOW_AUTODIFF_TAMC
34     #include "tamc.h"
35     #include "ctrl.h"
36     #include "ctrl_dummy.h"
37     #include "cost.h"
38     #endif
39    
40 adcroft 1.13 C == routine arguments ==
41     C note: under the multi-threaded model myiter and
42     C mytime are local variables passed around as routine
43     C arguments. Although this is fiddly it saves the need to
44     C impose additional synchronisation points when they are
45     C updated.
46     C myiter - iteration counter for this thread
47     C mytime - time counter for this thread
48     C mythid - thread number for this instance of the routine.
49 heimbach 1.12 integer iloop
50     integer mythid
51     integer myiter
52     _RL mytime
53    
54 adcroft 1.13 C == local variables ==
55 heimbach 1.12
56    
57     #ifdef ALLOW_AUTODIFF_TAMC
58 heimbach 1.14 C-- Reset the model iteration counter and the model time.
59 heimbach 1.12 myiter = nIter0 + (iloop-1)
60     mytime = startTime + float(iloop-1)*deltaTclock
61     #endif
62    
63 heimbach 1.16 #ifdef ALLOW_AUTODIFF_TAMC
64     C Include call to a dummy routine. Its adjoint will be
65     C called at the proper place in the adjoint code.
66     C The adjoint routine will print out adjoint values
67     C if requested. The location of the call is important,
68     C it has to be after the adjoint of the exchanges
69     C (DO_GTERM_BLOCKING_EXCHANGES).
70     CALL DUMMY_IN_STEPPING( myTime, myIter, myThid )
71     #endif
72    
73 jmc 1.18 #ifdef NONLIN_FRSURF
74     C-- compute the future surface level thickness
75     C according to the current Eta field
76     IF ( nonlinFreeSurf.GT.0) THEN
77     CALL CALC_SURF_DR(etaN, myTime, myIter, myThid )
78     ENDIF
79     #endif
80    
81 adcroft 1.13 C-- Load forcing/external data fields.
82 heimbach 1.12 #ifdef INCLUDE_EXTERNAL_FORCING_PACKAGE
83 adcroft 1.13 C NOTE, that although the exf package is part of the
84     C distribution, it is not currently maintained, i.e.
85     C exf is disabled by default in genmake.
86 heimbach 1.12 CALL EXF_GETFORCING( mytime, myiter, mythid )
87     #else
88     CALL TIMER_START('EXTERNAL_FIELDS_LOAD[THE_MAIN_LOOP]',mythid)
89     CALL EXTERNAL_FIELDS_LOAD( mytime, myiter, mythid )
90     CALL TIMER_STOP ('EXTERNAL_FIELDS_LOAD[THE_MAIN_LOOP]',mythid)
91     #endif
92 adcroft 1.15
93     C-- Step forward fields and calculate time tendency terms.
94     CALL TIMER_START('THERMODYNAMICS [THE_MAIN_LOOP]',mythid)
95     CALL THERMODYNAMICS( myTime, myIter, myThid )
96     CALL TIMER_STOP ('THERMODYNAMICS [THE_MAIN_LOOP]',mythid)
97 heimbach 1.12
98 adcroft 1.13 C-- Step forward fields and calculate time tendency terms.
99 heimbach 1.12 CALL TIMER_START('DYNAMICS [THE_MAIN_LOOP]',mythid)
100     CALL DYNAMICS( myTime, myIter, myThid )
101     CALL TIMER_STOP ('DYNAMICS [THE_MAIN_LOOP]',mythid)
102    
103     #ifndef ALLOW_AUTODIFF_TAMC
104 adcroft 1.1 #ifdef ALLOW_NONHYDROSTATIC
105     C-- Step forward W field in N-H algorithm
106 heimbach 1.12 IF ( nonHydrostatic ) THEN
107     CALL TIMER_START('CALC_GW [THE_MAIN_LOOP]',myThid)
108     CALL CALC_GW(myThid)
109     CALL TIMER_STOP ('CALC_GW [THE_MAIN_LOOP]',myThid)
110     ENDIF
111 adcroft 1.1 #endif
112 adcroft 1.13 #endif /* ALLOW_AUTODIFF_TAMC */
113 jmc 1.18
114     #ifdef NONLIN_FRSURF
115     C-- update hfacC,W,S and recip_hFac according to the current Eta field
116     IF ( nonlinFreeSurf.GT.0) THEN
117     CALL UPDATE_SURF_DR( myTime, myIter, myThid )
118     ENDIF
119     C- update also CG2D matrix (and preconditioner)
120     IF ( nonlinFreeSurf.GT.2) THEN
121     CALL UPDATE_CG2D( myTime, myIter, myThid )
122     ENDIF
123     #endif
124 adcroft 1.1
125 adcroft 1.13 C-- This block has been moved to the_correction_step(). We have
126     C left this code here to indicate where the filters used to be
127     C in the algorithm before JMC moved them to after the pressure
128     C solver.
129     C
130     C#ifdef ALLOW_SHAP_FILT
131     CC-- Step forward all tiles, filter and exchange.
132     C CALL TIMER_START('SHAP_FILT [THE_MAIN_LOOP]',myThid)
133     C CALL SHAP_FILT_APPLY(
134     C I gUnm1, gVnm1, gTnm1, gSnm1,
135     C I myTime, myIter, myThid )
136     C IF (implicDiv2Dflow.LT.1.) THEN
137     CC-- Explicit+Implicit part of the Barotropic Flow Divergence
138     CC => Filtering of uVel,vVel is necessary
139     C CALL SHAP_FILT_UV( uVel, vVel, myTime, myThid )
140     C ENDIF
141     C CALL TIMER_STOP ('SHAP_FILT [THE_MAIN_LOOP]',myThid)
142     C#endif
143     C
144     C#ifdef ALLOW_ZONAL_FILT
145     C IF (zonal_filt_lat.LT.90.) THEN
146     C CALL TIMER_START('ZONAL_FILT_APPLY [THE_MAIN_LOOP]',myThid)
147     C CALL ZONAL_FILT_APPLY(
148     C U gUnm1, gVnm1, gTnm1, gSnm1,
149     C I myThid )
150     C CALL TIMER_STOP ('ZONAL_FILT_APPLY [THE_MAIN_LOOP]',myThid)
151     C ENDIF
152     C#endif
153 heimbach 1.12
154 adcroft 1.1 C-- Solve elliptic equation(s).
155     C Two-dimensional only for conventional hydrostatic or
156     C three-dimensional for non-hydrostatic and/or IGW scheme.
157 heimbach 1.12 CALL TIMER_START('SOLVE_FOR_PRESSURE [THE_MAIN_LOOP]',myThid)
158 adcroft 1.1 CALL SOLVE_FOR_PRESSURE( myThid )
159 heimbach 1.12 CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [THE_MAIN_LOOP]',myThid)
160 adcroft 1.1
161 adcroft 1.5 C-- Correct divergence in flow field and cycle time-stepping
162 jmc 1.7 C arrays (for all fields) ; update time-counter
163 heimbach 1.12 myIter = nIter0 + iLoop
164     myTime = startTime + deltaTClock * float(iLoop)
165     CALL TIMER_START('THE_CORRECTION_STEP [THE_MAIN_LOOP]',myThid)
166     CALL THE_CORRECTION_STEP(myTime, myIter, myThid)
167     CALL TIMER_STOP ('THE_CORRECTION_STEP [THE_MAIN_LOOP]',myThid)
168 adcroft 1.5
169 adcroft 1.1 C-- Do "blocking" sends and receives for tendency "overlap" terms
170 heimbach 1.12 c CALL TIMER_START('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
171 jmc 1.7 c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )
172 heimbach 1.12 c CALL TIMER_STOP ('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
173 adcroft 1.5
174     C-- Do "blocking" sends and receives for field "overlap" terms
175 heimbach 1.12 CALL TIMER_START('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
176 adcroft 1.5 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
177 heimbach 1.12 CALL TIMER_STOP ('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
178    
179     #ifndef EXCLUDE_MONITOR
180     C-- Check status of solution (statistics, cfl, etc...)
181     CALL MONITOR( myIter, myTime, myThid )
182     #endif /* EXCLUDE_MONITOR */
183 adcroft 1.1
184 heimbach 1.16 #ifndef ALLOW_AUTODIFF_TAMC
185 jmc 1.7 C-- Do IO if needed.
186 heimbach 1.12 CALL TIMER_START('DO_THE_MODEL_IO [THE_MAIN_LOOP]',myThid)
187     CALL DO_THE_MODEL_IO( myTime, myIter, myThid )
188     CALL TIMER_STOP ('DO_THE_MODEL_IO [THE_MAIN_LOOP]',myThid)
189 adcroft 1.1
190     C-- Save state for restarts
191 jmc 1.7 C Note: (jmc: is it still the case after ckp35 ?)
192 adcroft 1.1 C =====
193     C Because of the ordering of the timestepping code and
194     C tendency term code at end of loop model arrays hold
195     C U,V,T,S at "time-level" N but gu, gv, gs, gt, guNM1,...
196     C at "time-level" N+1/2 (guNM1 at "time-level" N+1/2 is
197 jmc 1.10 C gu at "time-level" N-1/2) and etaN at "time-level" N+1/2.
198 adcroft 1.1 C where N = I+timeLevBase-1
199     C Thus a checkpoint contains U.0000000000, GU.0000000001 and
200 jmc 1.10 C etaN.0000000001 in the indexing scheme used for the model
201 adcroft 1.1 C "state" files. This example is referred to as a checkpoint
202     C at time level 1
203 heimbach 1.12 CALL TIMER_START('WRITE_CHECKPOINT [THE_MAIN_LOOP]',myThid)
204 adcroft 1.1 CALL WRITE_CHECKPOINT(
205 heimbach 1.12 & .FALSE., myTime, myIter, myThid )
206     CALL TIMER_STOP ('WRITE_CHECKPOINT [THE_MAIN_LOOP]',myThid)
207    
208     #endif /* ALLOW_AUTODIFF_TAMC */
209 adcroft 1.1
210     END

  ViewVC Help
Powered by ViewVC 1.1.22