/[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.19 - (show annotations) (download)
Tue Sep 4 14:44:54 2001 UTC (22 years, 8 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint40pre9
Changes since 1.18: +8 -2 lines
Added some missing conditionals so that fixed/flow passive tracer
experiments don't do any dyunamics.

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/forward_step.F,v 1.18 2001/08/27 18:50:41 jmc 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
54 C == local variables ==
55
56
57 #ifdef ALLOW_AUTODIFF_TAMC
58 C-- Reset the model iteration counter and the model time.
59 myiter = nIter0 + (iloop-1)
60 mytime = startTime + float(iloop-1)*deltaTclock
61 #endif
62
63 #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 #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 C-- Load forcing/external data fields.
82 #ifdef INCLUDE_EXTERNAL_FORCING_PACKAGE
83 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 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
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
98 C-- Step forward fields and calculate time tendency terms.
99 IF ( momStepping ) THEN
100 CALL TIMER_START('DYNAMICS [THE_MAIN_LOOP]',mythid)
101 CALL DYNAMICS( myTime, myIter, myThid )
102 CALL TIMER_STOP ('DYNAMICS [THE_MAIN_LOOP]',mythid)
103 ENDIF
104
105 #ifndef ALLOW_AUTODIFF_TAMC
106 #ifdef ALLOW_NONHYDROSTATIC
107 C-- Step forward W field in N-H algorithm
108 IF ( momStepping .AND. nonHydrostatic ) THEN
109 CALL TIMER_START('CALC_GW [THE_MAIN_LOOP]',myThid)
110 CALL CALC_GW(myThid)
111 CALL TIMER_STOP ('CALC_GW [THE_MAIN_LOOP]',myThid)
112 ENDIF
113 #endif
114 #endif /* ALLOW_AUTODIFF_TAMC */
115
116 #ifdef NONLIN_FRSURF
117 C-- update hfacC,W,S and recip_hFac according to the current Eta field
118 IF ( momStepping ) THEN
119 IF ( nonlinFreeSurf.GT.0) THEN
120 CALL UPDATE_SURF_DR( myTime, myIter, myThid )
121 ENDIF
122 C- update also CG2D matrix (and preconditioner)
123 IF ( nonlinFreeSurf.GT.2) THEN
124 CALL UPDATE_CG2D( myTime, myIter, myThid )
125 ENDIF
126 ENDIF
127 #endif
128
129 C-- This block has been moved to the_correction_step(). We have
130 C left this code here to indicate where the filters used to be
131 C in the algorithm before JMC moved them to after the pressure
132 C solver.
133 C
134 C#ifdef ALLOW_SHAP_FILT
135 CC-- Step forward all tiles, filter and exchange.
136 C CALL TIMER_START('SHAP_FILT [THE_MAIN_LOOP]',myThid)
137 C CALL SHAP_FILT_APPLY(
138 C I gUnm1, gVnm1, gTnm1, gSnm1,
139 C I myTime, myIter, myThid )
140 C IF (implicDiv2Dflow.LT.1.) THEN
141 CC-- Explicit+Implicit part of the Barotropic Flow Divergence
142 CC => Filtering of uVel,vVel is necessary
143 C CALL SHAP_FILT_UV( uVel, vVel, myTime, myThid )
144 C ENDIF
145 C CALL TIMER_STOP ('SHAP_FILT [THE_MAIN_LOOP]',myThid)
146 C#endif
147 C
148 C#ifdef ALLOW_ZONAL_FILT
149 C IF (zonal_filt_lat.LT.90.) THEN
150 C CALL TIMER_START('ZONAL_FILT_APPLY [THE_MAIN_LOOP]',myThid)
151 C CALL ZONAL_FILT_APPLY(
152 C U gUnm1, gVnm1, gTnm1, gSnm1,
153 C I myThid )
154 C CALL TIMER_STOP ('ZONAL_FILT_APPLY [THE_MAIN_LOOP]',myThid)
155 C ENDIF
156 C#endif
157
158 C-- Solve elliptic equation(s).
159 C Two-dimensional only for conventional hydrostatic or
160 C three-dimensional for non-hydrostatic and/or IGW scheme.
161 IF ( momStepping ) THEN
162 CALL TIMER_START('SOLVE_FOR_PRESSURE [THE_MAIN_LOOP]',myThid)
163 CALL SOLVE_FOR_PRESSURE( myThid )
164 CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [THE_MAIN_LOOP]',myThid)
165 ENDIF
166
167 C-- Correct divergence in flow field and cycle time-stepping
168 C arrays (for all fields) ; update time-counter
169 myIter = nIter0 + iLoop
170 myTime = startTime + deltaTClock * float(iLoop)
171 CALL TIMER_START('THE_CORRECTION_STEP [THE_MAIN_LOOP]',myThid)
172 CALL THE_CORRECTION_STEP(myTime, myIter, myThid)
173 CALL TIMER_STOP ('THE_CORRECTION_STEP [THE_MAIN_LOOP]',myThid)
174
175 C-- Do "blocking" sends and receives for tendency "overlap" terms
176 c CALL TIMER_START('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
177 c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )
178 c CALL TIMER_STOP ('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
179
180 C-- Do "blocking" sends and receives for field "overlap" terms
181 CALL TIMER_START('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
182 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
183 CALL TIMER_STOP ('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
184
185 #ifndef EXCLUDE_MONITOR
186 C-- Check status of solution (statistics, cfl, etc...)
187 CALL MONITOR( myIter, myTime, myThid )
188 #endif /* EXCLUDE_MONITOR */
189
190 #ifndef ALLOW_AUTODIFF_TAMC
191 C-- Do IO if needed.
192 CALL TIMER_START('DO_THE_MODEL_IO [THE_MAIN_LOOP]',myThid)
193 CALL DO_THE_MODEL_IO( myTime, myIter, myThid )
194 CALL TIMER_STOP ('DO_THE_MODEL_IO [THE_MAIN_LOOP]',myThid)
195
196 C-- Save state for restarts
197 C Note: (jmc: is it still the case after ckp35 ?)
198 C =====
199 C Because of the ordering of the timestepping code and
200 C tendency term code at end of loop model arrays hold
201 C U,V,T,S at "time-level" N but gu, gv, gs, gt, guNM1,...
202 C at "time-level" N+1/2 (guNM1 at "time-level" N+1/2 is
203 C gu at "time-level" N-1/2) and etaN at "time-level" N+1/2.
204 C where N = I+timeLevBase-1
205 C Thus a checkpoint contains U.0000000000, GU.0000000001 and
206 C etaN.0000000001 in the indexing scheme used for the model
207 C "state" files. This example is referred to as a checkpoint
208 C at time level 1
209 CALL TIMER_START('WRITE_CHECKPOINT [THE_MAIN_LOOP]',myThid)
210 CALL WRITE_CHECKPOINT(
211 & .FALSE., myTime, myIter, myThid )
212 CALL TIMER_STOP ('WRITE_CHECKPOINT [THE_MAIN_LOOP]',myThid)
213
214 #endif /* ALLOW_AUTODIFF_TAMC */
215
216 END

  ViewVC Help
Powered by ViewVC 1.1.22