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