/[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.26 - (show annotations) (download)
Tue Nov 20 23:27:29 2001 UTC (22 years, 6 months ago) by heimbach
Branch: MAIN
Changes since 1.25: +2 -1 lines
o Bugfix in adcommon.h: commen blocks were adjusted to latest
  common block structure in DYNVARS.h
o placed a do_field_blocking_exchanges after dummy_in_stepping
  to ensure that addummy_in_stepping is preceded by exchanges.

1 C $Header: /u/gcmpack/models/MITgcmUV/model/src/forward_step.F,v 1.25 2001/11/20 21:13:09 heimbach Exp $
2 C $Name: $
3
4 #include "CPP_OPTIONS.h"
5
6 CBOP
7 C !ROUTINE: FORWARD_STEP
8 C !INTERFACE:
9 SUBROUTINE FORWARD_STEP( iloop, myTime, myIter, myThid )
10
11 C !DESCRIPTION: \bv
12 C *==================================================================
13 C | SUBROUTINE forward_step
14 C | o Run the ocean model and, optionally, evaluate a cost function.
15 C *==================================================================
16 C |
17 C | THE_MAIN_LOOP is the toplevel routine for the Tangent Linear and
18 C | Adjoint Model Compiler (TAMC). For this purpose the initialization
19 C | of the model was split into two parts. Those parameters that do
20 C | not depend on a specific model run are set in INITIALISE_FIXED,
21 C | whereas those that do depend on the specific realization are
22 C | initialized in INITIALISE_VARIA.
23 C |
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 #include "DYNVARS.h"
34 #include "FFIELDS.h"
35
36 #ifdef ALLOW_NONHYDROSTATIC
37 #include "CG3D.h"
38 #endif
39
40 #ifdef ALLOW_AUTODIFF_TAMC
41 #include "tamc.h"
42 #include "ctrl.h"
43 #include "ctrl_dummy.h"
44 #include "cost.h"
45 #endif
46
47 C !LOCAL VARIABLES:
48 C == Routine arguments ==
49 C note: under the multi-threaded model myiter and
50 C mytime are local variables passed around as routine
51 C arguments. Although this is fiddly it saves the need to
52 C impose additional synchronisation points when they are
53 C updated.
54 C myiter - iteration counter for this thread
55 C mytime - time counter for this thread
56 C mythid - thread number for this instance of the routine.
57 integer iloop
58 integer mythid
59 integer myiter
60 _RL mytime
61 INTEGER bi,bj
62
63 CEOP
64
65
66 #ifdef ALLOW_AUTODIFF_TAMC
67 C-- Reset the model iteration counter and the model time.
68 myiter = nIter0 + (iloop-1)
69 mytime = startTime + float(iloop-1)*deltaTclock
70 #endif
71
72 #ifdef ALLOW_AUTODIFF_TAMC
73 C Include call to a dummy routine. Its adjoint will be
74 C called at the proper place in the adjoint code.
75 C The adjoint routine will print out adjoint values
76 C if requested. The location of the call is important,
77 C it has to be after the adjoint of the exchanges
78 C (DO_GTERM_BLOCKING_EXCHANGES).
79 CALL DUMMY_IN_STEPPING( myTime, myIter, myThid )
80 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
81 #endif
82
83 #ifdef EXACT_CONSERV
84 IF (exactConserv) THEN
85 C-- Update etaH(n+1) :
86 DO bj=myByLo(myThid),myByHi(myThid)
87 DO bi=myBxLo(myThid),myBxHi(myThid)
88 CALL CALC_EXACT_ETA( .FALSE., bi,bj, uVel,vVel,
89 I startTime, nIter0, myThid )
90 ENDDO
91 ENDDO
92 IF (implicDiv2Dflow .NE. 1. _d 0 )
93 & _EXCH_XY_R8(etaH, myThid )
94 ENDIF
95 #endif /* EXACT_CONSERV */
96
97 #ifdef NONLIN_FRSURF
98 C-- compute the future surface level thickness
99 C according to etaH(n+1)
100 IF ( nonlinFreeSurf.GT.0) THEN
101 CALL CALC_SURF_DR(etaH, myTime, myIter, myThid )
102 ENDIF
103 #endif /* NONLIN_FRSURF */
104
105 C-- Load forcing/external data fields.
106 #ifdef INCLUDE_EXTERNAL_FORCING_PACKAGE
107 C NOTE, that although the exf package is part of the
108 C distribution, it is not currently maintained, i.e.
109 C exf is disabled by default in genmake.
110 #ifdef ALLOW_BULKFORMULAE
111 CADJ STORE theta = comlev1, key = ikey_dynamics
112 #endif
113 CALL EXF_GETFORCING( mytime, myiter, mythid )
114 #else
115 CALL TIMER_START('EXTERNAL_FIELDS_LOAD[THE_MAIN_LOOP]',mythid)
116 CALL EXTERNAL_FIELDS_LOAD( mytime, myiter, mythid )
117 CALL TIMER_STOP ('EXTERNAL_FIELDS_LOAD[THE_MAIN_LOOP]',mythid)
118 #endif
119
120 C-- Step forward fields and calculate time tendency terms.
121 CALL TIMER_START('THERMODYNAMICS [THE_MAIN_LOOP]',mythid)
122 CALL THERMODYNAMICS( myTime, myIter, myThid )
123 CALL TIMER_STOP ('THERMODYNAMICS [THE_MAIN_LOOP]',mythid)
124
125 #ifdef ALLOW_SHAP_FILT
126 IF (staggerTimeStep .AND. useSHAP_FILT) THEN
127 c CALL TIMER_START('SHAP_FILT [THE_MAIN_LOOP]',myThid)
128 CALL SHAP_FILT_APPLY_TS( gT, gS, myTime, myIter, myThid )
129 c CALL TIMER_STOP ('SHAP_FILT [THE_MAIN_LOOP]',myThid)
130 ENDIF
131 #endif
132
133 C-- Step forward fields and calculate time tendency terms.
134 IF ( momStepping ) THEN
135 CALL TIMER_START('DYNAMICS [THE_MAIN_LOOP]',mythid)
136 CALL DYNAMICS( myTime, myIter, myThid )
137 CALL TIMER_STOP ('DYNAMICS [THE_MAIN_LOOP]',mythid)
138 ENDIF
139
140 #ifdef ALLOW_NONHYDROSTATIC
141 C-- Step forward W field in N-H algorithm
142 IF ( momStepping .AND. nonHydrostatic ) THEN
143 CALL TIMER_START('CALC_GW [THE_MAIN_LOOP]',myThid)
144 CALL CALC_GW(myThid)
145 CALL TIMER_STOP ('CALC_GW [THE_MAIN_LOOP]',myThid)
146 ENDIF
147 #endif
148
149 #ifdef NONLIN_FRSURF
150 C-- update hfacC,W,S and recip_hFac according to etaH(n+1) :
151 IF ( momStepping ) THEN
152 IF ( nonlinFreeSurf.GT.0) THEN
153 CALL UPDATE_SURF_DR( myTime, myIter, myThid )
154 ENDIF
155 C- update also CG2D matrix (and preconditioner)
156 IF ( nonlinFreeSurf.GT.2) THEN
157 CALL UPDATE_CG2D( myTime, myIter, myThid )
158 ENDIF
159 ENDIF
160 #endif
161
162 C-- This block has been moved to the_correction_step(). We have
163 C left this code here to indicate where the filters used to be
164 C in the algorithm before JMC moved them to after the pressure
165 C solver.
166 C
167 C#ifdef ALLOW_SHAP_FILT
168 CC-- Step forward all tiles, filter and exchange.
169 C CALL TIMER_START('SHAP_FILT [THE_MAIN_LOOP]',myThid)
170 C CALL SHAP_FILT_APPLY(
171 C I gUnm1, gVnm1, gTnm1, gSnm1,
172 C I myTime, myIter, myThid )
173 C IF (implicDiv2Dflow.LT.1.) THEN
174 CC-- Explicit+Implicit part of the Barotropic Flow Divergence
175 CC => Filtering of uVel,vVel is necessary
176 C CALL SHAP_FILT_UV( uVel, vVel, myTime, myThid )
177 C ENDIF
178 C CALL TIMER_STOP ('SHAP_FILT [THE_MAIN_LOOP]',myThid)
179 C#endif
180 C
181 C#ifdef ALLOW_ZONAL_FILT
182 C IF (zonal_filt_lat.LT.90.) THEN
183 C CALL TIMER_START('ZONAL_FILT_APPLY [THE_MAIN_LOOP]',myThid)
184 C CALL ZONAL_FILT_APPLY(
185 C U gUnm1, gVnm1, gTnm1, gSnm1,
186 C I myThid )
187 C CALL TIMER_STOP ('ZONAL_FILT_APPLY [THE_MAIN_LOOP]',myThid)
188 C ENDIF
189 C#endif
190
191 C-- Solve elliptic equation(s).
192 C Two-dimensional only for conventional hydrostatic or
193 C three-dimensional for non-hydrostatic and/or IGW scheme.
194 IF ( momStepping ) THEN
195 CALL TIMER_START('SOLVE_FOR_PRESSURE [THE_MAIN_LOOP]',myThid)
196 CALL SOLVE_FOR_PRESSURE( myThid )
197 CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [THE_MAIN_LOOP]',myThid)
198 ENDIF
199
200 C-- Correct divergence in flow field and cycle time-stepping
201 C arrays (for all fields) ; update time-counter
202 myIter = nIter0 + iLoop
203 myTime = startTime + deltaTClock * float(iLoop)
204 CALL TIMER_START('THE_CORRECTION_STEP [THE_MAIN_LOOP]',myThid)
205 CALL THE_CORRECTION_STEP(myTime, myIter, myThid)
206 CALL TIMER_STOP ('THE_CORRECTION_STEP [THE_MAIN_LOOP]',myThid)
207
208 C-- Do "blocking" sends and receives for tendency "overlap" terms
209 c CALL TIMER_START('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
210 c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )
211 c CALL TIMER_STOP ('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
212
213 C-- Do "blocking" sends and receives for field "overlap" terms
214 CALL TIMER_START('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
215 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
216 CALL TIMER_STOP ('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid)
217
218 #ifdef ALLOW_FLT
219 C-- Calculate float trajectories
220 IF (useFLT) THEN
221 CALL TIMER_START('FLOATS [THE_MAIN_LOOP]',myThid)
222 CALL FLT_MAIN(myIter,myTime, myThid)
223 CALL TIMER_STOP ('FLOATS [THE_MAIN_LOOP]',myThid)
224 ENDIF
225 #endif
226
227 #ifndef EXCLUDE_MONITOR
228 C-- Check status of solution (statistics, cfl, etc...)
229 CALL MONITOR( myIter, myTime, myThid )
230 #endif /* EXCLUDE_MONITOR */
231
232 #ifndef ALLOW_AUTODIFF_TAMC
233 C-- Do IO if needed.
234 CALL TIMER_START('DO_THE_MODEL_IO [THE_MAIN_LOOP]',myThid)
235 CALL DO_THE_MODEL_IO( myTime, myIter, myThid )
236 CALL TIMER_STOP ('DO_THE_MODEL_IO [THE_MAIN_LOOP]',myThid)
237
238 C-- Save state for restarts
239 C Note: (jmc: is it still the case after ckp35 ?)
240 C =====
241 C Because of the ordering of the timestepping code and
242 C tendency term code at end of loop model arrays hold
243 C U,V,T,S at "time-level" N but gu, gv, gs, gt, guNM1,...
244 C at "time-level" N+1/2 (guNM1 at "time-level" N+1/2 is
245 C gu at "time-level" N-1/2) and etaN at "time-level" N+1/2.
246 C where N = I+timeLevBase-1
247 C Thus a checkpoint contains U.0000000000, GU.0000000001 and
248 C etaN.0000000001 in the indexing scheme used for the model
249 C "state" files. This example is referred to as a checkpoint
250 C at time level 1
251 CALL TIMER_START('WRITE_CHECKPOINT [THE_MAIN_LOOP]',myThid)
252 CALL WRITE_CHECKPOINT(
253 & .FALSE., myTime, myIter, myThid )
254 CALL TIMER_STOP ('WRITE_CHECKPOINT [THE_MAIN_LOOP]',myThid)
255
256 #endif /* ALLOW_AUTODIFF_TAMC */
257
258 END

  ViewVC Help
Powered by ViewVC 1.1.22