1 |
C $Header: /u/gcmpack/MITgcm/model/src/initialise_varia.F,v 1.79 2015/03/23 21:00:30 gforget Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "PACKAGES_CONFIG.h" |
5 |
#include "CPP_OPTIONS.h" |
6 |
#ifdef ALLOW_AUTODIFF |
7 |
# include "AUTODIFF_OPTIONS.h" |
8 |
#endif |
9 |
#ifdef ALLOW_CTRL |
10 |
# include "CTRL_OPTIONS.h" |
11 |
#endif |
12 |
|
13 |
CBOP |
14 |
C !ROUTINE: INITIALISE_VARIA |
15 |
C !INTERFACE: |
16 |
SUBROUTINE INITIALISE_VARIA( myThid ) |
17 |
C !DESCRIPTION: \bv |
18 |
C *==========================================================* |
19 |
C | SUBROUTINE INITIALISE_VARIA |
20 |
C | o Set the initial conditions for dynamics variables |
21 |
C | and time dependent arrays |
22 |
C *==========================================================* |
23 |
C | This routine reads/writes data from an input file and |
24 |
C | from various binary files. |
25 |
C | Each thread invokes an instance of this routine as does |
26 |
C | each process in a multi-process parallel environment like |
27 |
C | MPI. |
28 |
C *==========================================================* |
29 |
C \ev |
30 |
|
31 |
C !CALLING SEQUENCE: |
32 |
C INITIALISE_VARIA |
33 |
C | |
34 |
C #ifdef ALLOW_AUTODIFF |
35 |
C |-- INI_DEPTHS \ |
36 |
C |-- CTRL_DEPTH_INI \ |
37 |
C |-- UPDATE_MASKS_ETC } ALLOW_DEPTH_CONTROL case |
38 |
C |-- UPDATE_CG2D / |
39 |
C #endif |
40 |
C |-- INI_NLFS_VARS |
41 |
C |-- INI_DYNVARS |
42 |
C |-- INI_NH_VARS |
43 |
C |-- INI_FFIELDS |
44 |
C | |
45 |
C |-- INI_FIELDS |
46 |
C | |
47 |
C |-- INI_MIXING |
48 |
C | |
49 |
C |-- TAUEDDY_INIT_VARIA |
50 |
C | |
51 |
C |-- INI_FORCING |
52 |
C | |
53 |
C |-- AUTODIFF_INIT_VARIA |
54 |
C | |
55 |
C |-- PACKAGES_INIT_VARIABLES |
56 |
C | |
57 |
C |-- COST_INIT_VARIA |
58 |
C | |
59 |
C |-- CONVECTIVE_ADJUSTMENT_INI |
60 |
C | |
61 |
C |-- CALC_R_STAR |
62 |
C |-- UPDATE_R_STAR |
63 |
C |-- UPDATE_SIGMA |
64 |
C |-- CALC_SURF_DR |
65 |
C |-- UPDATE_SURF_DR |
66 |
C | |
67 |
C |-- UPDATE_CG2D |
68 |
C | |
69 |
C |-- INTEGR_CONTINUITY |
70 |
C | |
71 |
C |-- CALC_R_STAR |
72 |
C |-- CALC_SURF_DR |
73 |
C | |
74 |
C |-- STATE_SUMMARY |
75 |
C | |
76 |
C |-- MONITOR |
77 |
C | |
78 |
C |-- DO_STATEVARS_TAVE |
79 |
C | |
80 |
C |-- DO_THE_MODEL_IO |
81 |
|
82 |
C !USES: |
83 |
IMPLICIT NONE |
84 |
C == Global variables == |
85 |
#include "SIZE.h" |
86 |
#include "EEPARAMS.h" |
87 |
#include "PARAMS.h" |
88 |
#include "DYNVARS.h" |
89 |
#include "SURFACE.h" |
90 |
#ifdef ALLOW_AUTODIFF |
91 |
# include "GRID.h" |
92 |
# include "FFIELDS.h" |
93 |
# include "CTRL_FIELDS.h" |
94 |
#endif |
95 |
|
96 |
C !INPUT/OUTPUT PARAMETERS: |
97 |
C == Routine arguments == |
98 |
INTEGER myThid |
99 |
|
100 |
C !LOCAL VARIABLES: |
101 |
C == Local variables == |
102 |
INTEGER bi,bj |
103 |
CEOP |
104 |
|
105 |
#ifdef ALLOW_DEBUG |
106 |
IF (debugMode) CALL DEBUG_ENTER('INITIALISE_VARIA',myThid) |
107 |
#endif |
108 |
|
109 |
#ifdef ALLOW_AUTODIFF |
110 |
nIter0 = NINT( (startTime-baseTime)/deltaTClock ) |
111 |
#endif /* ALLOW_AUTODIFF */ |
112 |
|
113 |
#ifdef ALLOW_CTRL |
114 |
# ifdef ALLOW_DEPTH_CONTROL |
115 |
C-- Intialize the depth for TAF/TAMC |
116 |
CALL INI_DEPTHS( myThid ) |
117 |
C-- Get control parameter depth |
118 |
CALL CTRL_DEPTH_INI( myThid ) |
119 |
C-- Re-calculate hFacS/W and some other parameters from hFacC |
120 |
CALL UPDATE_MASKS_ETC( myThid ) |
121 |
C-- Update laplace operators for use in 2D conjugate gradient solver. |
122 |
CALL UPDATE_CG2D( startTime, nIter0, myThid ) |
123 |
# endif /* ALLOW_DEPTH_CONTROL */ |
124 |
#endif /* ALLOW_CTRL */ |
125 |
|
126 |
C-- Initialise Non-Lin FreeSurf variables: |
127 |
CALL INI_NLFS_VARS( myThid ) |
128 |
|
129 |
C-- Initialize DYNVARS arrays (state fields + G terms: Gu,Gv,...) to zero [always] |
130 |
#ifdef ALLOW_DEBUG |
131 |
IF (debugMode) CALL DEBUG_CALL('INI_DYNVARS',myThid) |
132 |
#endif |
133 |
CALL INI_DYNVARS( myThid ) |
134 |
|
135 |
C-- Initialize NH_VARS arrays to zero [always] |
136 |
#ifdef ALLOW_NONHYDROSTATIC |
137 |
CALL INI_NH_VARS( myThid ) |
138 |
#endif |
139 |
|
140 |
C-- Initialize FFIELDS arrays to zero [always] |
141 |
CALL INI_FFIELDS( myThid ) |
142 |
|
143 |
C-- Initialise model fields. |
144 |
C Starting values of U, V, W, temp., salt. and tendency terms |
145 |
C are set here. Fields are either set to default or read from |
146 |
C stored files. |
147 |
#ifdef ALLOW_DEBUG |
148 |
IF (debugMode) CALL DEBUG_CALL('INI_FIELDS',myThid) |
149 |
#endif |
150 |
CALL INI_FIELDS( myThid ) |
151 |
|
152 |
C-- Initialise 3-dim. diffusivities |
153 |
#ifdef ALLOW_DEBUG |
154 |
IF (debugMode) CALL DEBUG_CALL('INI_MIXING',myThid) |
155 |
#endif |
156 |
CALL INI_MIXING( myThid ) |
157 |
|
158 |
#ifdef ALLOW_EDDYPSI |
159 |
C-- Initialise eddy diffusivities |
160 |
CALL TAUEDDY_INIT_VARIA( myThid ) |
161 |
#endif |
162 |
|
163 |
C-- Initialise model forcing fields. |
164 |
#ifdef ALLOW_DEBUG |
165 |
IF (debugMode) CALL DEBUG_CALL('INI_FORCING',myThid) |
166 |
#endif |
167 |
CALL INI_FORCING( myThid ) |
168 |
|
169 |
#ifdef ALLOW_AUTODIFF |
170 |
C-- Initialise active fields to help TAMC |
171 |
if (useAUTODIFF) CALL AUTODIFF_INIT_VARIA( myThid ) |
172 |
#endif |
173 |
|
174 |
C-- Initialize variable data for packages |
175 |
#ifdef ALLOW_DEBUG |
176 |
IF (debugMode) CALL DEBUG_CALL('PACKAGES_INIT_VARIABLES',myThid) |
177 |
#endif |
178 |
#ifdef ALLOW_AUTODIFF_TAMC |
179 |
# ifdef NONLIN_FRSURF |
180 |
CADJ STORE recip_hFacC = tapelev_init, key = 1 |
181 |
# endif |
182 |
#endif |
183 |
CALL PACKAGES_INIT_VARIABLES( myThid ) |
184 |
|
185 |
#ifdef ALLOW_COST |
186 |
C-- Initialise the cost function (moved out of packages_init_variables to |
187 |
C here to prevent resetting cost-funct in adinitialise_varia recomput.) |
188 |
CALL COST_INIT_VARIA( myThid ) |
189 |
#endif /* ALLOW_COST */ |
190 |
|
191 |
c#ifndef ALLOW_AUTODIFF |
192 |
c IF ( usePickupBeforeC35 .AND. startTime .NE. baseTime ) THEN |
193 |
C-- IMPORTANT : Need to activate the following call to restart from a pickup |
194 |
C file written by MITgcmUV_checkpoint34 (Feb-08, 2001) or earlier. |
195 |
C- Disable this option on Jan-09, 2007. |
196 |
c CALL THE_CORRECTION_STEP(startTime, nIter0, myThid) |
197 |
c ENDIF |
198 |
c#endif |
199 |
|
200 |
#ifndef ALLOW_AUTODIFF_WHTAPEIO |
201 |
C-- Initial conditions are convectively adjusted (for historical reasons) |
202 |
IF ( startTime .EQ. baseTime .AND. cAdjFreq .NE. 0. ) THEN |
203 |
#ifdef ALLOW_DEBUG |
204 |
IF (debugMode) CALL DEBUG_CALL('CONVECTIVE_ADJUSTMENT_INI',myThid) |
205 |
#endif |
206 |
CADJ loop = parallel |
207 |
DO bj = myByLo(myThid), myByHi(myThid) |
208 |
CADJ loop = parallel |
209 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
210 |
CALL CONVECTIVE_ADJUSTMENT_INI( |
211 |
I bi, bj, startTime, nIter0, myThid ) |
212 |
ENDDO |
213 |
ENDDO |
214 |
ENDIF |
215 |
#endif /* ALLOW_AUTODIFF_WHTAPEIO */ |
216 |
|
217 |
#ifdef NONLIN_FRSURF |
218 |
C-- Compute the surface level thickness <-- function of etaH(n) |
219 |
C and modify hFac(C,W,S) accordingly : |
220 |
# ifndef DISABLE_RSTAR_CODE |
221 |
IF ( select_rStar.NE.0 ) |
222 |
& CALL CALC_R_STAR(etaH, startTime, -1 , myThid ) |
223 |
# endif /* DISABLE_RSTAR_CODE */ |
224 |
IF ( nonlinFreeSurf.GT.0 ) THEN |
225 |
IF ( select_rStar.GT.0 ) THEN |
226 |
# ifndef DISABLE_RSTAR_CODE |
227 |
CALL UPDATE_R_STAR( .TRUE., startTime, nIter0, myThid ) |
228 |
# endif /* DISABLE_RSTAR_CODE */ |
229 |
ELSEIF ( selectSigmaCoord.NE.0 ) THEN |
230 |
# ifndef DISABLE_SIGMA_CODE |
231 |
CALL UPDATE_SIGMA( etaH, startTime, nIter0, myThid ) |
232 |
# endif /* DISABLE_SIGMA_CODE */ |
233 |
ELSE |
234 |
CALL CALC_SURF_DR(etaH, startTime, -1 , myThid ) |
235 |
CALL UPDATE_SURF_DR( .TRUE., startTime, nIter0, myThid ) |
236 |
ENDIF |
237 |
ENDIF |
238 |
C- update also CG2D matrix (and preconditioner) |
239 |
IF ( nonlinFreeSurf.GT.2 ) THEN |
240 |
CALL UPDATE_CG2D( startTime, nIter0, myThid ) |
241 |
ENDIF |
242 |
#endif /* NONLIN_FRSURF */ |
243 |
|
244 |
#ifdef ALLOW_DEBUG |
245 |
IF (debugMode) CALL DEBUG_CALL('INTEGR_CONTINUITY',myThid) |
246 |
#endif |
247 |
C-- Integrate continuity vertically for vertical velocity |
248 |
CALL INTEGR_CONTINUITY( uVel, vVel, |
249 |
I startTime, nIter0, myThid ) |
250 |
|
251 |
#ifdef NONLIN_FRSURF |
252 |
IF ( select_rStar.NE.0 ) THEN |
253 |
#ifndef DISABLE_RSTAR_CODE |
254 |
C-- r* : compute the future level thickness according to etaH(n+1) |
255 |
CALL CALC_R_STAR(etaH, startTime, nIter0, myThid ) |
256 |
#endif |
257 |
ELSEIF ( nonlinFreeSurf.GT.0 .AND. selectSigmaCoord.EQ.0 ) THEN |
258 |
C-- compute the future surface level thickness according to etaH(n+1) |
259 |
CALL CALC_SURF_DR(etaH, startTime, nIter0, myThid ) |
260 |
ENDIF |
261 |
#endif /* NONLIN_FRSURF */ |
262 |
|
263 |
c IF ( nIter0.EQ.0 .AND. staggerTimeStep ) THEN |
264 |
C-- Filter initial T & S fields if staggerTimeStep |
265 |
C (only for backward compatibility ; to be removed later) |
266 |
#ifdef ALLOW_SHAP_FILT |
267 |
c IF ( useSHAP_FILT .AND. shap_filt_TrStagg ) THEN |
268 |
c CALL SHAP_FILT_APPLY_TS(theta,salt,startTime,nIter0,myThid) |
269 |
c ENDIF |
270 |
#endif |
271 |
#ifdef ALLOW_ZONAL_FILT |
272 |
c IF ( useZONAL_FILT .AND. zonal_filt_TrStagg ) THEN |
273 |
c CALL ZONAL_FILT_APPLY_TS( theta, salt, myThid ) |
274 |
c ENDIF |
275 |
#endif |
276 |
c ENDIF |
277 |
|
278 |
#ifdef ALLOW_GRIDALT |
279 |
IF (useGRIDALT) THEN |
280 |
CALL TIMER_START('GRIDALT_UPDATE [INITIALISE_VARIA]',myThid) |
281 |
CALL GRIDALT_UPDATE(myThid) |
282 |
CALL TIMER_STOP ('GRIDALT_UPDATE [INITIALISE_VARIA]',myThid) |
283 |
ENDIF |
284 |
#endif |
285 |
|
286 |
C-- Finally summarise the model state |
287 |
#ifdef ALLOW_DEBUG |
288 |
IF (debugMode) CALL DEBUG_CALL('STATE_SUMMARY',myThid) |
289 |
#endif |
290 |
CALL STATE_SUMMARY( myThid ) |
291 |
|
292 |
#ifdef ALLOW_MONITOR |
293 |
#ifdef ALLOW_DEBUG |
294 |
IF (debugMode) CALL DEBUG_CALL('MONITOR',myThid) |
295 |
#endif |
296 |
C-- Check status of initial state (statistics, cfl, etc...) |
297 |
CALL MONITOR( startTime, nIter0, myThid ) |
298 |
#endif /* ALLOW_MONITOR */ |
299 |
|
300 |
#ifdef ALLOW_TIMEAVE |
301 |
#ifdef ALLOW_DEBUG |
302 |
IF (debugMode) CALL DEBUG_CALL('DO_STATEVARS_TAVE',myThid) |
303 |
#endif |
304 |
C-- Initialise time-average arrays with initial state values |
305 |
CALL DO_STATEVARS_TAVE( startTime, nIter0, myThid ) |
306 |
#endif |
307 |
|
308 |
C-- Dump initial state to files |
309 |
#ifdef ALLOW_DEBUG |
310 |
IF (debugMode) CALL DEBUG_CALL('DO_THE_MODEL_IO',myThid) |
311 |
#endif |
312 |
CALL DO_THE_MODEL_IO( .FALSE., startTime, nIter0, myThid ) |
313 |
|
314 |
#ifdef ALLOW_DEBUG |
315 |
IF (debugMode) CALL DEBUG_LEAVE('INITIALISE_VARIA',myThid) |
316 |
#endif |
317 |
|
318 |
C-- Check barrier synchronization: |
319 |
CALL BAR_CHECK( 4, myThid ) |
320 |
|
321 |
RETURN |
322 |
END |