/[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.48 - (hide annotations) (download)
Sun Jan 26 21:07:25 2003 UTC (21 years, 4 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint48b_post
Changes since 1.47: +15 -5 lines
r* coordinate added in #ifdef NONLIN_FRSURF block.
  (modification to pressure gradient not yet implemented)

1 jmc 1.48 C $Header: /u/gcmpack/MITgcm/model/src/forward_step.F,v 1.47 2003/01/22 16:09:24 jmc Exp $
2 adcroft 1.15 C $Name: $
3 adcroft 1.1
4     #include "CPP_OPTIONS.h"
5    
6 cnh 1.22 CBOP
7     C !ROUTINE: FORWARD_STEP
8     C !INTERFACE:
9 adcroft 1.13 SUBROUTINE FORWARD_STEP( iloop, myTime, myIter, myThid )
10 heimbach 1.12
11 cnh 1.22 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 adcroft 1.1
27 cnh 1.22 C !USES:
28     IMPLICIT NONE
29     C == Global variables ==
30 adcroft 1.1 #include "SIZE.h"
31     #include "EEPARAMS.h"
32     #include "PARAMS.h"
33     #include "DYNVARS.h"
34 heimbach 1.12 #include "FFIELDS.h"
35    
36 adcroft 1.1 #ifdef ALLOW_NONHYDROSTATIC
37     #include "CG3D.h"
38     #endif
39    
40 jmc 1.27 #ifdef ALLOW_SHAP_FILT
41     #include "SHAP_FILT.h"
42     #endif
43     #ifdef ALLOW_ZONAL_FILT
44     #include "ZONAL_FILT.h"
45     #endif
46    
47 heimbach 1.12 #ifdef ALLOW_AUTODIFF_TAMC
48 heimbach 1.30 # include "tamc.h"
49     # include "ctrl.h"
50     # include "ctrl_dummy.h"
51     # include "cost.h"
52 heimbach 1.37 # include "EOS.h"
53 heimbach 1.30 # ifdef INCLUDE_EXTERNAL_FORCING_PACKAGE
54     # include "exf_fields.h"
55 dimitri 1.45 # if (defined (ALLOW_BULKFORMULAE) || defined (ALLOW_BULK_FORCE))
56 heimbach 1.30 # include "exf_constants.h"
57     # endif
58     # endif
59     # ifdef ALLOW_OBCS
60     # include "OBCS.h"
61     # endif
62 heimbach 1.12 #endif
63    
64 cnh 1.22 C !LOCAL VARIABLES:
65     C == Routine arguments ==
66 adcroft 1.13 C note: under the multi-threaded model myiter and
67     C mytime are local variables passed around as routine
68     C arguments. Although this is fiddly it saves the need to
69     C impose additional synchronisation points when they are
70     C updated.
71     C myiter - iteration counter for this thread
72     C mytime - time counter for this thread
73     C mythid - thread number for this instance of the routine.
74 heimbach 1.12 integer iloop
75     integer mythid
76     integer myiter
77     _RL mytime
78 dimitri 1.45 #ifdef ALLOW_BULK_FORCE
79 jmc 1.21 INTEGER bi,bj
80 dimitri 1.45 #endif ALLOW_BULK_FORCE
81 heimbach 1.12
82 cnh 1.22 CEOP
83 heimbach 1.12
84     #ifdef ALLOW_AUTODIFF_TAMC
85 dimitri 1.45 C-- Reset the model iteration counter and the model time.
86     myiter = nIter0 + (iloop-1)
87     mytime = startTime + float(iloop-1)*deltaTclock
88 heimbach 1.12 #endif
89    
90 heimbach 1.32 #if (defined (ALLOW_AUTODIFF_TAMC) && defined (ALLOW_AUTODIFF_MONITOR))
91 dimitri 1.45 C Include call to a dummy routine. Its adjoint will be
92     C called at the proper place in the adjoint code.
93     C The adjoint routine will print out adjoint values
94     C if requested. The location of the call is important,
95     C it has to be after the adjoint of the exchanges
96     C (DO_GTERM_BLOCKING_EXCHANGES).
97     CALL DUMMY_IN_STEPPING( myTime, myIter, myThid )
98 heimbach 1.46 cph I've commented this line since it may conflict with MITgcm's adjoint
99     cph However, need to check whether that's still consistent
100     cph with the ecco-branch (it should).
101     cph CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
102 heimbach 1.16 #endif
103    
104 jmc 1.21 #ifdef EXACT_CONSERV
105     IF (exactConserv) THEN
106     C-- Update etaH(n+1) :
107 heimbach 1.36 CALL TIMER_START('UPDATE_ETAH [FORWARD_STEP]',mythid)
108 jmc 1.35 CALL UPDATE_ETAH( myTime, myIter, myThid )
109 heimbach 1.36 CALL TIMER_STOP ('UPDATE_ETAH [FORWARD_STEP]',mythid)
110 jmc 1.21 ENDIF
111     #endif /* EXACT_CONSERV */
112    
113 jmc 1.18 #ifdef NONLIN_FRSURF
114 jmc 1.48 IF ( select_rStar.NE.0 ) THEN
115     C-- r* : compute the future level thickness according to etaH(n+1)
116     CALL TIMER_START('CALC_R_STAR [FORWARD_STEP]',mythid)
117     CALL CALC_R_STAR(etaH, myTime, myIter, myThid )
118     CALL TIMER_STOP ('CALC_R_STAR [FORWARD_STEP]',mythid)
119     ELSEIF ( nonlinFreeSurf.GT.0) THEN
120     C-- compute the future surface level thickness according to etaH(n+1)
121 heimbach 1.36 CALL TIMER_START('CALC_SURF_DR [FORWARD_STEP]',mythid)
122 jmc 1.21 CALL CALC_SURF_DR(etaH, myTime, myIter, myThid )
123 heimbach 1.36 CALL TIMER_STOP ('CALC_SURF_DR [FORWARD_STEP]',mythid)
124 jmc 1.48 ENDIF
125 jmc 1.21 #endif /* NONLIN_FRSURF */
126 jmc 1.18
127 dimitri 1.45 C-- Load forcing/external data fields.
128 heimbach 1.28 #ifdef ALLOW_AUTODIFF_TAMC
129     c**************************************
130     #include "checkpoint_lev1_directives.h"
131     c**************************************
132 heimbach 1.23 #endif
133 cheisey 1.38
134    
135 heimbach 1.37 C-- Call external forcing package
136 cheisey 1.38 cswdblk -- add ---
137 cheisey 1.40 #ifdef ALLOW_BULK_FORCE
138 cheisey 1.38 CALL TIMER_START('BULKF_FIELDS_LOAD[THE_MAIN_LOOP]',mythid)
139     CALL BULKF_FIELDS_LOAD( mytime, myiter, mythid )
140     CALL TIMER_STOP ('BULKF_FIELDS_LOAD[THE_MAIN_LOOP]',mythid)
141     c calculate qnet and empmr (and wind stress)
142     DO bj=myByLo(myThid),myByHi(myThid)
143     DO bi=myBxLo(myThid),myBxHi(myThid)
144     CALL BULKF_FORCING( bi,bj, mytime, myiter, mythid )
145     ENDDO
146     ENDDO
147     c Update the tile edges.
148     _EXCH_XY_R8(Qnet, mythid)
149     _EXCH_XY_R8(EmPmR, mythid)
150 cheisey 1.43 CALL EXCH_UV_XY_RS(fu, fv, .TRUE., myThid)
151     C _EXCH_XY_R8(fu , mythid)
152     C _EXCH_XY_R8(fv , mythid)
153 cheisey 1.38 cswdblk -- end add ---
154 cheisey 1.42 #else
155 cheisey 1.39 #ifdef INCLUDE_EXTERNAL_FORCING_PACKAGE
156 dimitri 1.45 C NOTE, that although the exf package is part of the
157     C distribution, it is not currently maintained, i.e.
158     C exf is disabled by default in genmake.
159     CALL TIMER_START('EXF_GETFORCING [FORWARD_STEP]',mythid)
160     CALL EXF_GETFORCING( mytime, myiter, mythid )
161     CALL TIMER_STOP ('EXF_GETFORCING [FORWARD_STEP]',mythid)
162 heimbach 1.12 #else
163 dimitri 1.45 IF ( .NOT. useSEAICE ) THEN
164     CALL TIMER_START('EXTERNAL_FIELDS_LOAD[FORWARD_STEP]',mythid)
165     CALL EXTERNAL_FIELDS_LOAD( mytime, myiter, mythid )
166     CALL TIMER_STOP ('EXTERNAL_FIELDS_LOAD[FORWARD_STEP]',mythid)
167     ENDIF
168 heimbach 1.36 #endif /* INCLUDE_EXTERNAL_FORCING_PACKAGE */
169 cheisey 1.38
170 heimbach 1.36 #ifdef ALLOW_SEAICE
171 dimitri 1.45 C-- Call sea ice model to compute forcing/external data fields. In
172     C addition to computing prognostic sea-ice variables and diagnosing the
173     C forcing/external data fields that drive the ocean model, SEAICE_MODEL
174     C also sets theta to the freezing point under sea-ice. The implied
175     C surface heat flux is then stored in variable surfaceTendencyTice,
176     C which is needed by KPP package (kpp_calc.F and kpp_transport_t.F)
177     C to diagnose surface buoyancy fluxes and for the non-local transport
178     C term. Because this call precedes model thermodynamics, temperature
179     C under sea-ice may not be "exactly" at the freezing point by the time
180     C theta is dumped or time-averaged.
181     IF ( useSEAICE ) THEN
182 heimbach 1.36 CALL TIMER_START('SEAICE_MODEL [FORWARD_STEP]',myThid)
183 jmc 1.41 CALL SEAICE_MODEL( myTime, myIter, myThid )
184 heimbach 1.36 CALL TIMER_STOP ('SEAICE_MODEL [FORWARD_STEP]',myThid)
185 dimitri 1.45 ENDIF
186 heimbach 1.36 #endif ALLOW_SEAICE
187 dimitri 1.45 #endif ALLOW_BULK_FORCE
188 adcroft 1.15
189     C-- Step forward fields and calculate time tendency terms.
190 heimbach 1.36 CALL TIMER_START('THERMODYNAMICS [FORWARD_STEP]',mythid)
191 adcroft 1.15 CALL THERMODYNAMICS( myTime, myIter, myThid )
192 heimbach 1.36 CALL TIMER_STOP ('THERMODYNAMICS [FORWARD_STEP]',mythid)
193 jmc 1.35
194     C-- do exchanges (needed for DYNAMICS) when using stagger time-step :
195 heimbach 1.36 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
196 jmc 1.35 CALL DO_STAGGER_FIELDS_EXCHANGES( myTime, myIter, myThid )
197 heimbach 1.36 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
198 jmc 1.24
199     #ifdef ALLOW_SHAP_FILT
200 dimitri 1.45 IF (useSHAP_FILT .AND.
201 jmc 1.27 & staggerTimeStep .AND. shap_filt_TrStagg ) THEN
202 heimbach 1.36 CALL TIMER_START('SHAP_FILT [FORWARD_STEP]',myThid)
203 jmc 1.24 CALL SHAP_FILT_APPLY_TS( gT, gS, myTime, myIter, myThid )
204 heimbach 1.36 CALL TIMER_STOP ('SHAP_FILT [FORWARD_STEP]',myThid)
205 dimitri 1.45 ENDIF
206 jmc 1.24 #endif
207 jmc 1.27 #ifdef ALLOW_ZONAL_FILT
208 dimitri 1.45 IF (useZONAL_FILT .AND.
209 jmc 1.27 & staggerTimeStep .AND. zonal_filt_TrStagg ) THEN
210 heimbach 1.36 CALL TIMER_START('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid)
211 jmc 1.27 CALL ZONAL_FILT_APPLY_TS( gT, gS, myThid )
212 heimbach 1.36 CALL TIMER_STOP ('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid)
213 jmc 1.27 ENDIF
214     #endif
215 heimbach 1.12
216 dimitri 1.45 C-- Step forward fields and calculate time tendency terms.
217     IF ( momStepping ) THEN
218 heimbach 1.36 CALL TIMER_START('DYNAMICS [FORWARD_STEP]',mythid)
219 heimbach 1.12 CALL DYNAMICS( myTime, myIter, myThid )
220 heimbach 1.36 CALL TIMER_STOP ('DYNAMICS [FORWARD_STEP]',mythid)
221 dimitri 1.45 ENDIF
222 heimbach 1.12
223 adcroft 1.1 #ifdef ALLOW_NONHYDROSTATIC
224     C-- Step forward W field in N-H algorithm
225 dimitri 1.45 IF ( momStepping .AND. nonHydrostatic ) THEN
226     CALL TIMER_START('CALC_GW [FORWARD_STEP]',myThid)
227     CALL CALC_GW(myThid)
228     CALL TIMER_STOP ('CALC_GW [FORWARD_STEP]',myThid)
229     ENDIF
230 adcroft 1.1 #endif
231 jmc 1.18
232     #ifdef NONLIN_FRSURF
233 jmc 1.21 C-- update hfacC,W,S and recip_hFac according to etaH(n+1) :
234 jmc 1.18 IF ( nonlinFreeSurf.GT.0) THEN
235 jmc 1.48 IF ( select_rStar.GT.0 ) THEN
236     CALL TIMER_START('UPDATE_R_STAR [FORWARD_STEP]',myThid)
237     CALL UPDATE_R_STAR( myTime, myIter, myThid )
238     CALL TIMER_STOP ('UPDATE_R_STAR [FORWARD_STEP]',myThid)
239     ELSE
240 dimitri 1.45 CALL TIMER_START('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
241 jmc 1.18 CALL UPDATE_SURF_DR( myTime, myIter, myThid )
242 heimbach 1.36 CALL TIMER_STOP ('UPDATE_SURF_DR [FORWARD_STEP]',myThid)
243 jmc 1.48 ENDIF
244 jmc 1.18 ENDIF
245     C- update also CG2D matrix (and preconditioner)
246 jmc 1.33 IF ( momStepping .AND. nonlinFreeSurf.GT.2 ) THEN
247 dimitri 1.45 CALL TIMER_START('UPDATE_CG2D [FORWARD_STEP]',myThid)
248 jmc 1.18 CALL UPDATE_CG2D( myTime, myIter, myThid )
249 jmc 1.47 CALL TIMER_STOP ('UPDATE_CG2D [FORWARD_STEP]',myThid)
250 adcroft 1.19 ENDIF
251 jmc 1.18 #endif
252 adcroft 1.1
253 jmc 1.27 C-- Apply Filters to u*,v* before SOLVE_FOR_PRESSURE
254     #ifdef ALLOW_SHAP_FILT
255     IF (useSHAP_FILT .AND. shap_filt_uvStar) THEN
256 heimbach 1.36 CALL TIMER_START('SHAP_FILT [FORWARD_STEP]',myThid)
257 jmc 1.27 CALL SHAP_FILT_APPLY_UV( gUnm1,gVnm1, myTime,myIter,myThid )
258     IF (implicDiv2Dflow.LT.1.) THEN
259     C-- Explicit+Implicit part of the Barotropic Flow Divergence
260     C => Filtering of uVel,vVel is necessary
261     CALL SHAP_FILT_APPLY_UV( uVel,vVel, myTime,myIter,myThid )
262     ENDIF
263 heimbach 1.36 CALL TIMER_STOP ('SHAP_FILT [FORWARD_STEP]',myThid)
264 jmc 1.27 ENDIF
265     #endif
266     #ifdef ALLOW_ZONAL_FILT
267     IF (useZONAL_FILT .AND. zonal_filt_uvStar) THEN
268 heimbach 1.36 CALL TIMER_START('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid)
269 jmc 1.27 CALL ZONAL_FILT_APPLY_UV( gUnm1, gVnm1, myThid )
270     IF (implicDiv2Dflow.LT.1.) THEN
271     C-- Explicit+Implicit part of the Barotropic Flow Divergence
272     C => Filtering of uVel,vVel is necessary
273     CALL ZONAL_FILT_APPLY_UV( uVel, vVel, myThid )
274     ENDIF
275 heimbach 1.36 CALL TIMER_STOP ('ZONAL_FILT_APPLY [FORWARD_STEP]',myThid)
276 jmc 1.27 ENDIF
277     #endif
278 heimbach 1.12
279 adcroft 1.1 C-- Solve elliptic equation(s).
280     C Two-dimensional only for conventional hydrostatic or
281     C three-dimensional for non-hydrostatic and/or IGW scheme.
282 adcroft 1.19 IF ( momStepping ) THEN
283 heimbach 1.36 CALL TIMER_START('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
284 jmc 1.31 CALL SOLVE_FOR_PRESSURE(myTime, myIter, myThid)
285 heimbach 1.36 CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [FORWARD_STEP]',myThid)
286 adcroft 1.19 ENDIF
287 adcroft 1.1
288 heimbach 1.37 #ifdef ALLOW_AUTODIFF_TAMC
289     cph This is needed because convective_adjustment calls
290     cph find_rho which may use pressure()
291     CADJ STORE pressure = comlev1, key = ikey_dynamics
292     #endif
293 adcroft 1.5 C-- Correct divergence in flow field and cycle time-stepping
294 jmc 1.7 C arrays (for all fields) ; update time-counter
295 heimbach 1.12 myIter = nIter0 + iLoop
296     myTime = startTime + deltaTClock * float(iLoop)
297 heimbach 1.36 CALL TIMER_START('THE_CORRECTION_STEP [FORWARD_STEP]',myThid)
298 heimbach 1.12 CALL THE_CORRECTION_STEP(myTime, myIter, myThid)
299 heimbach 1.36 CALL TIMER_STOP ('THE_CORRECTION_STEP [FORWARD_STEP]',myThid)
300 adcroft 1.5
301 adcroft 1.1 C-- Do "blocking" sends and receives for tendency "overlap" terms
302 heimbach 1.36 c CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
303 jmc 1.7 c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid )
304 heimbach 1.36 c CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
305 adcroft 1.5
306     C-- Do "blocking" sends and receives for field "overlap" terms
307 heimbach 1.36 CALL TIMER_START('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
308 adcroft 1.5 CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid )
309 heimbach 1.36 CALL TIMER_STOP ('BLOCKING_EXCHANGES [FORWARD_STEP]',myThid)
310 adcroft 1.20
311     #ifdef ALLOW_FLT
312     C-- Calculate float trajectories
313     IF (useFLT) THEN
314 heimbach 1.36 CALL TIMER_START('FLOATS [FORWARD_STEP]',myThid)
315 adcroft 1.20 CALL FLT_MAIN(myIter,myTime, myThid)
316 heimbach 1.36 CALL TIMER_STOP ('FLOATS [FORWARD_STEP]',myThid)
317 adcroft 1.20 ENDIF
318     #endif
319 heimbach 1.12
320     #ifndef EXCLUDE_MONITOR
321     C-- Check status of solution (statistics, cfl, etc...)
322 dimitri 1.45 CALL TIMER_START('MONITOR [FORWARD_STEP]',myThid)
323 heimbach 1.12 CALL MONITOR( myIter, myTime, myThid )
324 heimbach 1.36 CALL TIMER_STOP ('MONITOR [FORWARD_STEP]',myThid)
325 heimbach 1.12 #endif /* EXCLUDE_MONITOR */
326 adcroft 1.1
327 jmc 1.7 C-- Do IO if needed.
328 heimbach 1.36 CALL TIMER_START('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
329 heimbach 1.12 CALL DO_THE_MODEL_IO( myTime, myIter, myThid )
330 heimbach 1.36 CALL TIMER_STOP ('DO_THE_MODEL_IO [FORWARD_STEP]',myThid)
331 adcroft 1.1
332     C-- Save state for restarts
333 jmc 1.7 C Note: (jmc: is it still the case after ckp35 ?)
334 adcroft 1.1 C =====
335     C Because of the ordering of the timestepping code and
336     C tendency term code at end of loop model arrays hold
337     C U,V,T,S at "time-level" N but gu, gv, gs, gt, guNM1,...
338     C at "time-level" N+1/2 (guNM1 at "time-level" N+1/2 is
339 jmc 1.10 C gu at "time-level" N-1/2) and etaN at "time-level" N+1/2.
340 adcroft 1.1 C where N = I+timeLevBase-1
341     C Thus a checkpoint contains U.0000000000, GU.0000000001 and
342 jmc 1.10 C etaN.0000000001 in the indexing scheme used for the model
343 adcroft 1.1 C "state" files. This example is referred to as a checkpoint
344     C at time level 1
345 heimbach 1.36 CALL TIMER_START('WRITE_CHECKPOINT [FORWARD_STEP]',myThid)
346 adcroft 1.1 CALL WRITE_CHECKPOINT(
347 heimbach 1.12 & .FALSE., myTime, myIter, myThid )
348 heimbach 1.36 CALL TIMER_STOP ('WRITE_CHECKPOINT [FORWARD_STEP]',myThid)
349 adcroft 1.1
350     END

  ViewVC Help
Powered by ViewVC 1.1.22