/[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.14 - (show annotations) (download)
Wed Aug 1 22:12:12 2001 UTC (22 years, 10 months ago) by heimbach
Branch: MAIN
Changes since 1.13: +12 -12 lines
Moved dummy_in_stepping routine to different place.

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

  ViewVC Help
Powered by ViewVC 1.1.22