1 |
C $Header: /u/gcmpack/models/MITgcmUV/model/src/forward_step.F,v 1.23 2001/09/27 20:12:10 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 |
#endif |
81 |
|
82 |
#ifdef EXACT_CONSERV |
83 |
IF (exactConserv) THEN |
84 |
C-- Update etaH(n+1) : |
85 |
DO bj=myByLo(myThid),myByHi(myThid) |
86 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
87 |
CALL CALC_EXACT_ETA( .FALSE., bi,bj, uVel,vVel, |
88 |
I startTime, nIter0, myThid ) |
89 |
ENDDO |
90 |
ENDDO |
91 |
IF (implicDiv2Dflow .NE. 1. _d 0 ) |
92 |
& _EXCH_XY_R8(etaH, myThid ) |
93 |
ENDIF |
94 |
#endif /* EXACT_CONSERV */ |
95 |
|
96 |
#ifdef NONLIN_FRSURF |
97 |
C-- compute the future surface level thickness |
98 |
C according to etaH(n+1) |
99 |
IF ( nonlinFreeSurf.GT.0) THEN |
100 |
CALL CALC_SURF_DR(etaH, myTime, myIter, myThid ) |
101 |
ENDIF |
102 |
#endif /* NONLIN_FRSURF */ |
103 |
|
104 |
C-- Load forcing/external data fields. |
105 |
#ifdef INCLUDE_EXTERNAL_FORCING_PACKAGE |
106 |
C NOTE, that although the exf package is part of the |
107 |
C distribution, it is not currently maintained, i.e. |
108 |
C exf is disabled by default in genmake. |
109 |
#ifdef ALLOW_BULKFORMULAE |
110 |
CADJ STORE theta = comlev1, key = ikey_dynamics |
111 |
#endif |
112 |
CALL EXF_GETFORCING( mytime, myiter, mythid ) |
113 |
#else |
114 |
CALL TIMER_START('EXTERNAL_FIELDS_LOAD[THE_MAIN_LOOP]',mythid) |
115 |
CALL EXTERNAL_FIELDS_LOAD( mytime, myiter, mythid ) |
116 |
CALL TIMER_STOP ('EXTERNAL_FIELDS_LOAD[THE_MAIN_LOOP]',mythid) |
117 |
#endif |
118 |
|
119 |
C-- Step forward fields and calculate time tendency terms. |
120 |
CALL TIMER_START('THERMODYNAMICS [THE_MAIN_LOOP]',mythid) |
121 |
CALL THERMODYNAMICS( myTime, myIter, myThid ) |
122 |
CALL TIMER_STOP ('THERMODYNAMICS [THE_MAIN_LOOP]',mythid) |
123 |
|
124 |
#ifdef ALLOW_SHAP_FILT |
125 |
IF (staggerTimeStep .AND. useSHAP_FILT) THEN |
126 |
c CALL TIMER_START('SHAP_FILT [THE_MAIN_LOOP]',myThid) |
127 |
CALL SHAP_FILT_APPLY_TS( gT, gS, myTime, myIter, myThid ) |
128 |
c CALL TIMER_STOP ('SHAP_FILT [THE_MAIN_LOOP]',myThid) |
129 |
ENDIF |
130 |
#endif |
131 |
|
132 |
C-- Step forward fields and calculate time tendency terms. |
133 |
IF ( momStepping ) THEN |
134 |
CALL TIMER_START('DYNAMICS [THE_MAIN_LOOP]',mythid) |
135 |
CALL DYNAMICS( myTime, myIter, myThid ) |
136 |
CALL TIMER_STOP ('DYNAMICS [THE_MAIN_LOOP]',mythid) |
137 |
ENDIF |
138 |
|
139 |
#ifndef ALLOW_AUTODIFF_TAMC |
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 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
149 |
|
150 |
#ifdef NONLIN_FRSURF |
151 |
C-- update hfacC,W,S and recip_hFac according to etaH(n+1) : |
152 |
IF ( momStepping ) THEN |
153 |
IF ( nonlinFreeSurf.GT.0) THEN |
154 |
CALL UPDATE_SURF_DR( myTime, myIter, myThid ) |
155 |
ENDIF |
156 |
C- update also CG2D matrix (and preconditioner) |
157 |
IF ( nonlinFreeSurf.GT.2) THEN |
158 |
CALL UPDATE_CG2D( myTime, myIter, myThid ) |
159 |
ENDIF |
160 |
ENDIF |
161 |
#endif |
162 |
|
163 |
C-- This block has been moved to the_correction_step(). We have |
164 |
C left this code here to indicate where the filters used to be |
165 |
C in the algorithm before JMC moved them to after the pressure |
166 |
C solver. |
167 |
C |
168 |
C#ifdef ALLOW_SHAP_FILT |
169 |
CC-- Step forward all tiles, filter and exchange. |
170 |
C CALL TIMER_START('SHAP_FILT [THE_MAIN_LOOP]',myThid) |
171 |
C CALL SHAP_FILT_APPLY( |
172 |
C I gUnm1, gVnm1, gTnm1, gSnm1, |
173 |
C I myTime, myIter, myThid ) |
174 |
C IF (implicDiv2Dflow.LT.1.) THEN |
175 |
CC-- Explicit+Implicit part of the Barotropic Flow Divergence |
176 |
CC => Filtering of uVel,vVel is necessary |
177 |
C CALL SHAP_FILT_UV( uVel, vVel, myTime, myThid ) |
178 |
C ENDIF |
179 |
C CALL TIMER_STOP ('SHAP_FILT [THE_MAIN_LOOP]',myThid) |
180 |
C#endif |
181 |
C |
182 |
C#ifdef ALLOW_ZONAL_FILT |
183 |
C IF (zonal_filt_lat.LT.90.) THEN |
184 |
C CALL TIMER_START('ZONAL_FILT_APPLY [THE_MAIN_LOOP]',myThid) |
185 |
C CALL ZONAL_FILT_APPLY( |
186 |
C U gUnm1, gVnm1, gTnm1, gSnm1, |
187 |
C I myThid ) |
188 |
C CALL TIMER_STOP ('ZONAL_FILT_APPLY [THE_MAIN_LOOP]',myThid) |
189 |
C ENDIF |
190 |
C#endif |
191 |
|
192 |
C-- Solve elliptic equation(s). |
193 |
C Two-dimensional only for conventional hydrostatic or |
194 |
C three-dimensional for non-hydrostatic and/or IGW scheme. |
195 |
IF ( momStepping ) THEN |
196 |
CALL TIMER_START('SOLVE_FOR_PRESSURE [THE_MAIN_LOOP]',myThid) |
197 |
CALL SOLVE_FOR_PRESSURE( myThid ) |
198 |
CALL TIMER_STOP ('SOLVE_FOR_PRESSURE [THE_MAIN_LOOP]',myThid) |
199 |
ENDIF |
200 |
|
201 |
C-- Correct divergence in flow field and cycle time-stepping |
202 |
C arrays (for all fields) ; update time-counter |
203 |
myIter = nIter0 + iLoop |
204 |
myTime = startTime + deltaTClock * float(iLoop) |
205 |
CALL TIMER_START('THE_CORRECTION_STEP [THE_MAIN_LOOP]',myThid) |
206 |
CALL THE_CORRECTION_STEP(myTime, myIter, myThid) |
207 |
CALL TIMER_STOP ('THE_CORRECTION_STEP [THE_MAIN_LOOP]',myThid) |
208 |
|
209 |
C-- Do "blocking" sends and receives for tendency "overlap" terms |
210 |
c CALL TIMER_START('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid) |
211 |
c CALL DO_GTERM_BLOCKING_EXCHANGES( myThid ) |
212 |
c CALL TIMER_STOP ('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid) |
213 |
|
214 |
C-- Do "blocking" sends and receives for field "overlap" terms |
215 |
CALL TIMER_START('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid) |
216 |
CALL DO_FIELDS_BLOCKING_EXCHANGES( myThid ) |
217 |
CALL TIMER_STOP ('BLOCKING_EXCHANGES [THE_MAIN_LOOP]',myThid) |
218 |
|
219 |
#ifdef ALLOW_FLT |
220 |
C-- Calculate float trajectories |
221 |
IF (useFLT) THEN |
222 |
CALL TIMER_START('FLOATS [THE_MAIN_LOOP]',myThid) |
223 |
CALL FLT_MAIN(myIter,myTime, myThid) |
224 |
CALL TIMER_STOP ('FLOATS [THE_MAIN_LOOP]',myThid) |
225 |
ENDIF |
226 |
#endif |
227 |
|
228 |
#ifndef EXCLUDE_MONITOR |
229 |
C-- Check status of solution (statistics, cfl, etc...) |
230 |
CALL MONITOR( myIter, myTime, myThid ) |
231 |
#endif /* EXCLUDE_MONITOR */ |
232 |
|
233 |
#ifndef ALLOW_AUTODIFF_TAMC |
234 |
C-- Do IO if needed. |
235 |
CALL TIMER_START('DO_THE_MODEL_IO [THE_MAIN_LOOP]',myThid) |
236 |
CALL DO_THE_MODEL_IO( myTime, myIter, myThid ) |
237 |
CALL TIMER_STOP ('DO_THE_MODEL_IO [THE_MAIN_LOOP]',myThid) |
238 |
|
239 |
C-- Save state for restarts |
240 |
C Note: (jmc: is it still the case after ckp35 ?) |
241 |
C ===== |
242 |
C Because of the ordering of the timestepping code and |
243 |
C tendency term code at end of loop model arrays hold |
244 |
C U,V,T,S at "time-level" N but gu, gv, gs, gt, guNM1,... |
245 |
C at "time-level" N+1/2 (guNM1 at "time-level" N+1/2 is |
246 |
C gu at "time-level" N-1/2) and etaN at "time-level" N+1/2. |
247 |
C where N = I+timeLevBase-1 |
248 |
C Thus a checkpoint contains U.0000000000, GU.0000000001 and |
249 |
C etaN.0000000001 in the indexing scheme used for the model |
250 |
C "state" files. This example is referred to as a checkpoint |
251 |
C at time level 1 |
252 |
CALL TIMER_START('WRITE_CHECKPOINT [THE_MAIN_LOOP]',myThid) |
253 |
CALL WRITE_CHECKPOINT( |
254 |
& .FALSE., myTime, myIter, myThid ) |
255 |
CALL TIMER_STOP ('WRITE_CHECKPOINT [THE_MAIN_LOOP]',myThid) |
256 |
|
257 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
258 |
|
259 |
END |