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