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

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

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


Revision 1.21 - (show annotations) (download)
Wed Sep 19 13:58:08 2001 UTC (22 years, 8 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint40
Changes since 1.20: +20 -5 lines
"Volume exact-Conservation" modified for
non-linear free-surface + Crank-Nickelson

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

  ViewVC Help
Powered by ViewVC 1.1.22