1 |
jmc |
1.45 |
C $Header: /u/gcmpack/MITgcm/model/src/initialise_varia.F,v 1.44 2005/06/19 21:36:33 jmc Exp $ |
2 |
jmc |
1.22 |
C $Name: $ |
3 |
adcroft |
1.1 |
|
4 |
edhill |
1.32 |
#include "PACKAGES_CONFIG.h" |
5 |
adcroft |
1.1 |
#include "CPP_OPTIONS.h" |
6 |
|
|
|
7 |
cnh |
1.20 |
CBOP |
8 |
|
|
C !ROUTINE: INITIALISE_VARIA |
9 |
|
|
C !INTERFACE: |
10 |
adcroft |
1.1 |
SUBROUTINE INITIALISE_VARIA(myThid) |
11 |
cnh |
1.20 |
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 |-- INI_LINEAR_PHISURF |
29 |
|
|
C | |
30 |
|
|
C |-- INI_CORI |
31 |
|
|
C | |
32 |
|
|
C |-- INI_CG2D |
33 |
|
|
C | |
34 |
|
|
C |-- INI_CG3D |
35 |
|
|
C | |
36 |
|
|
C |-- INI_MIXING |
37 |
|
|
C | |
38 |
|
|
C |-- INI_DYNVARS |
39 |
|
|
C | |
40 |
|
|
C |-- INI_FIELDS |
41 |
|
|
C | |
42 |
dimitri |
1.28 |
C |-- INI_AUTODIFF |
43 |
|
|
C | |
44 |
adcroft |
1.23 |
C |-- PACKAGES_INIT_VARIABLES |
45 |
|
|
C | |
46 |
jmc |
1.30 |
C |-- THE_CORRECTION_STEP (if restart from old pickup files) |
47 |
cnh |
1.20 |
C | |
48 |
|
|
C |-- CALL CONVECTIVE_ADJUSTMENT_INI |
49 |
|
|
C | |
50 |
|
|
C |-- CALC_SURF_DR |
51 |
|
|
C | |
52 |
|
|
C |-- UPDATE_SURF_DR |
53 |
|
|
C | |
54 |
|
|
C |-- UPDATE_CG2D |
55 |
|
|
C | |
56 |
jmc |
1.26 |
C |-- INTEGR_CONTINUITY |
57 |
cnh |
1.20 |
C | |
58 |
|
|
C |-- TIMEAVE_STATVARS |
59 |
dimitri |
1.33 |
C | |
60 |
|
|
C |-- PTRACERS_STATVARS |
61 |
jmc |
1.26 |
C | |
62 |
|
|
C |-- STATE_SUMMARY |
63 |
cnh |
1.20 |
|
64 |
|
|
C !USES: |
65 |
adcroft |
1.1 |
IMPLICIT NONE |
66 |
|
|
C == Global variables == |
67 |
|
|
#include "SIZE.h" |
68 |
|
|
#include "EEPARAMS.h" |
69 |
|
|
#include "PARAMS.h" |
70 |
adcroft |
1.5 |
#include "DYNVARS.h" |
71 |
jmc |
1.36 |
#include "SURFACE.h" |
72 |
|
|
#ifdef ALLOW_SHAP_FILT |
73 |
|
|
# include "SHAP_FILT.h" |
74 |
|
|
#endif |
75 |
|
|
#ifdef ALLOW_ZONAL_FILT |
76 |
|
|
# include "ZONAL_FILT.h" |
77 |
|
|
#endif |
78 |
adcroft |
1.1 |
|
79 |
cnh |
1.20 |
C !INPUT/OUTPUT PARAMETERS: |
80 |
adcroft |
1.1 |
C == Routine arguments == |
81 |
|
|
INTEGER myThid |
82 |
|
|
|
83 |
cnh |
1.20 |
C !LOCAL VARIABLES: |
84 |
adcroft |
1.1 |
C == Local variables == |
85 |
jmc |
1.45 |
INTEGER bi,bj,iMin,iMax,jMin,jMax |
86 |
cnh |
1.20 |
CEOP |
87 |
heimbach |
1.14 |
|
88 |
adcroft |
1.41 |
#ifdef ALLOW_DEBUG |
89 |
|
|
IF (debugMode) CALL DEBUG_ENTER('INITIALISE_VARIA',myThid) |
90 |
|
|
#endif |
91 |
|
|
|
92 |
heimbach |
1.40 |
#ifdef ALLOW_AUTODIFF_TAMC |
93 |
heimbach |
1.15 |
|
94 |
jmc |
1.43 |
nIter0 = NINT( (startTime-baseTime)/deltaTClock ) |
95 |
heimbach |
1.15 |
|
96 |
|
|
C-- Set Bo_surf => define the Linear Relation: Phi_surf(eta) |
97 |
adcroft |
1.41 |
#ifdef ALLOW_DEBUG |
98 |
|
|
IF (debugMode) CALL DEBUG_CALL('INI_LINEAR_PHISURF',myThid) |
99 |
|
|
#endif |
100 |
heimbach |
1.15 |
CALL INI_LINEAR_PHISURF( myThid ) |
101 |
|
|
|
102 |
|
|
C-- Set coriolis operators |
103 |
adcroft |
1.41 |
#ifdef ALLOW_DEBUG |
104 |
|
|
IF (debugMode) CALL DEBUG_CALL('INI_CORI',myThid) |
105 |
|
|
#endif |
106 |
heimbach |
1.15 |
CALL INI_CORI( myThid ) |
107 |
|
|
|
108 |
|
|
C-- Set laplace operators for use in 2D conjugate gradient solver. |
109 |
adcroft |
1.41 |
#ifdef ALLOW_DEBUG |
110 |
|
|
IF (debugMode) CALL DEBUG_CALL('INI_CG2D',myThid) |
111 |
|
|
#endif |
112 |
heimbach |
1.15 |
CALL INI_CG2D( myThid ) |
113 |
|
|
|
114 |
|
|
#ifdef ALLOW_NONHYDROSTATIC |
115 |
|
|
C-- Set laplace operators for use in 3D conjugate gradient solver. |
116 |
edhill |
1.32 |
ceh3 should add an IF ( useNONHYDROSTATIC ) THEN |
117 |
adcroft |
1.41 |
#ifdef ALLOW_DEBUG |
118 |
|
|
IF (debugMode) CALL DEBUG_CALL('INI_CG3D',myThid) |
119 |
|
|
#endif |
120 |
heimbach |
1.15 |
CALL INI_CG3D( myThid ) |
121 |
heimbach |
1.14 |
#endif |
122 |
adcroft |
1.1 |
|
123 |
heimbach |
1.40 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
124 |
heimbach |
1.16 |
_BARRIER |
125 |
heimbach |
1.15 |
|
126 |
heimbach |
1.16 |
C-- Initialise 3-dim. diffusivities |
127 |
adcroft |
1.41 |
#ifdef ALLOW_DEBUG |
128 |
|
|
IF (debugMode) CALL DEBUG_CALL('INI_MIXING',myThid) |
129 |
|
|
#endif |
130 |
heimbach |
1.16 |
CALL INI_MIXING( myThid ) |
131 |
adcroft |
1.5 |
_BARRIER |
132 |
heimbach |
1.16 |
|
133 |
heimbach |
1.42 |
#ifdef ALLOW_TAU_EDDY |
134 |
|
|
C-- Initialise eddy diffusivities |
135 |
|
|
CALL TAUEDDY_INIT_VARIA( myThid ) |
136 |
|
|
#endif |
137 |
|
|
|
138 |
adcroft |
1.12 |
C-- Initialize DYNVARS arrays (state fields + G terms: Gu,Gv,...) to zero [always] |
139 |
adcroft |
1.41 |
#ifdef ALLOW_DEBUG |
140 |
|
|
IF (debugMode) CALL DEBUG_CALL('INI_DYNVARS',myThid) |
141 |
|
|
#endif |
142 |
adcroft |
1.12 |
CALL INI_DYNVARS( myThid ) |
143 |
|
|
|
144 |
adcroft |
1.1 |
C-- Initialise model fields. |
145 |
|
|
C Starting values of U, V, W, temp., salt. and tendency terms |
146 |
|
|
C are set here. Fields are either set to default or read from |
147 |
|
|
C stored files. |
148 |
stephd |
1.38 |
#ifndef ALLOW_OFFLINE |
149 |
adcroft |
1.41 |
#ifdef ALLOW_DEBUG |
150 |
|
|
IF (debugMode) CALL DEBUG_CALL('INI_FIELDS',myThid) |
151 |
|
|
#endif |
152 |
adcroft |
1.1 |
CALL INI_FIELDS( myThid ) |
153 |
stephd |
1.38 |
#endif |
154 |
adcroft |
1.1 |
_BARRIER |
155 |
jmc |
1.8 |
|
156 |
heimbach |
1.21 |
#ifdef ALLOW_AUTODIFF_TAMC |
157 |
|
|
C-- Initialise active fields to help TAMC |
158 |
|
|
CALL INI_AUTODIFF( myThid ) |
159 |
heimbach |
1.13 |
_BARRIER |
160 |
heimbach |
1.15 |
#endif |
161 |
heimbach |
1.13 |
|
162 |
heimbach |
1.25 |
C-- Initialize variable data for packages |
163 |
adcroft |
1.41 |
#ifdef ALLOW_DEBUG |
164 |
|
|
IF (debugMode) CALL DEBUG_CALL('PACKAGES_INIT_VARIABLES',myThid) |
165 |
|
|
#endif |
166 |
heimbach |
1.25 |
CALL PACKAGES_INIT_VARIABLES( myThid ) |
167 |
|
|
|
168 |
heimbach |
1.18 |
#ifndef ALLOW_AUTODIFF_TAMC |
169 |
heimbach |
1.11 |
IF ( usePickupBeforeC35 ) THEN |
170 |
jmc |
1.8 |
C-- IMPORTANT : Need to activate the following call to restart from |
171 |
|
|
C a pickup file written by MITgcmUV_checkpoint34 or earlier. |
172 |
jmc |
1.43 |
IF ( startTime .NE. baseTime ) THEN |
173 |
adcroft |
1.41 |
#ifdef ALLOW_DEBUG |
174 |
|
|
IF (debugMode) CALL DEBUG_CALL('THE_CORRECTION_STEP',myThid) |
175 |
|
|
#endif |
176 |
heimbach |
1.11 |
CALL THE_CORRECTION_STEP(startTime, nIter0, myThid) |
177 |
|
|
ENDIF |
178 |
|
|
ENDIF |
179 |
heimbach |
1.18 |
#endif |
180 |
adcroft |
1.1 |
|
181 |
adcroft |
1.5 |
C-- Initial conditions are convectively adjusted (for historical reasons) |
182 |
stephd |
1.38 |
#ifndef ALLOW_OFFLINE |
183 |
jmc |
1.43 |
IF ( startTime .EQ. baseTime .AND. cAdjFreq .NE. 0. ) THEN |
184 |
adcroft |
1.41 |
#ifdef ALLOW_DEBUG |
185 |
|
|
IF (debugMode) CALL DEBUG_CALL('CONVECTIVE_ADJUSTMENT_INI',myThid) |
186 |
|
|
#endif |
187 |
heimbach |
1.13 |
CADJ loop = parallel |
188 |
adcroft |
1.1 |
DO bj = myByLo(myThid), myByHi(myThid) |
189 |
heimbach |
1.13 |
CADJ loop = parallel |
190 |
adcroft |
1.1 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
191 |
adcroft |
1.5 |
iMin=1-Olx |
192 |
|
|
iMax=sNx+Olx |
193 |
|
|
jMin=1-Oly |
194 |
|
|
jMax=sNy+Oly |
195 |
heimbach |
1.11 |
CALL CONVECTIVE_ADJUSTMENT_INI( |
196 |
adcroft |
1.5 |
I bi, bj, iMin, iMax, jMin, jMax, |
197 |
jmc |
1.7 |
I startTime, nIter0, myThid ) |
198 |
adcroft |
1.1 |
ENDDO |
199 |
|
|
ENDDO |
200 |
adcroft |
1.5 |
_BARRIER |
201 |
|
|
END IF |
202 |
stephd |
1.38 |
#endif |
203 |
jmc |
1.17 |
|
204 |
|
|
#ifdef NONLIN_FRSURF |
205 |
jmc |
1.19 |
C-- Compute the surface level thickness <-- function of etaH(n) |
206 |
jmc |
1.17 |
C and modify hFac(C,W,S) accordingly : |
207 |
jmc |
1.36 |
IF ( select_rStar.NE.0 ) |
208 |
|
|
& CALL CALC_R_STAR(etaH, startTime, -1 , myThid ) |
209 |
jmc |
1.17 |
IF (nonlinFreeSurf.GT.0) THEN |
210 |
jmc |
1.29 |
IF ( select_rStar.GT.0 ) THEN |
211 |
|
|
CALL UPDATE_R_STAR( startTime, nIter0, myThid ) |
212 |
|
|
ELSE |
213 |
jmc |
1.22 |
CALL CALC_SURF_DR(etaH, startTime, -1 , myThid ) |
214 |
jmc |
1.17 |
CALL UPDATE_SURF_DR( startTime, nIter0, myThid ) |
215 |
jmc |
1.29 |
ENDIF |
216 |
jmc |
1.17 |
ENDIF |
217 |
|
|
C- update also CG2D matrix (and preconditioner) |
218 |
|
|
IF ( nonlinFreeSurf.GT.2) THEN |
219 |
|
|
CALL UPDATE_CG2D( startTime, nIter0, myThid ) |
220 |
|
|
ENDIF |
221 |
jmc |
1.29 |
#endif /* NONLIN_FRSURF */ |
222 |
adcroft |
1.1 |
|
223 |
stephd |
1.38 |
#ifndef ALLOW_OFFLINE |
224 |
adcroft |
1.41 |
#ifdef ALLOW_DEBUG |
225 |
|
|
IF (debugMode) CALL DEBUG_CALL('INTEGR_CONTINUITY',myThid) |
226 |
|
|
#endif |
227 |
jmc |
1.26 |
DO bj=myByLo(myThid),myByHi(myThid) |
228 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
229 |
adcroft |
1.1 |
|
230 |
jmc |
1.26 |
C-- Integrate continuity vertically for vertical velocity |
231 |
|
|
CALL INTEGR_CONTINUITY( bi, bj, uVel, vVel, |
232 |
|
|
I startTime, nIter0, myThid ) |
233 |
|
|
|
234 |
jmc |
1.36 |
ENDDO |
235 |
|
|
ENDDO |
236 |
stephd |
1.38 |
#endif |
237 |
jmc |
1.36 |
|
238 |
|
|
#ifdef EXACT_CONSERV |
239 |
|
|
IF ( exactConserv ) THEN |
240 |
jmc |
1.44 |
#ifdef ALLOW_DEBUG |
241 |
|
|
IF (debugMode) CALL DEBUG_CALL('UPDATE_ETAH',myThid) |
242 |
|
|
#endif |
243 |
jmc |
1.36 |
C-- Update etaH(n) : |
244 |
|
|
CALL UPDATE_ETAH( startTime, nIter0, myThid ) |
245 |
|
|
ENDIF |
246 |
|
|
#endif /* EXACT_CONSERV */ |
247 |
|
|
|
248 |
|
|
#ifdef NONLIN_FRSURF |
249 |
|
|
IF ( select_rStar.NE.0 ) THEN |
250 |
|
|
C-- r* : compute the future level thickness according to etaH(n+1) |
251 |
|
|
CALL CALC_R_STAR(etaH, startTime, nIter0, myThid ) |
252 |
|
|
ELSEIF ( nonlinFreeSurf.GT.0) THEN |
253 |
|
|
C-- compute the future surface level thickness according to etaH(n+1) |
254 |
|
|
CALL CALC_SURF_DR(etaH, startTime, nIter0, myThid ) |
255 |
|
|
ENDIF |
256 |
|
|
#endif /* NONLIN_FRSURF */ |
257 |
|
|
|
258 |
jmc |
1.37 |
c IF ( nIter0.EQ.0 .AND. staggerTimeStep ) THEN |
259 |
jmc |
1.36 |
C-- Filter initial T & S fields if staggerTimeStep |
260 |
|
|
C (only for backward compatibility ; to be removed later) |
261 |
|
|
#ifdef ALLOW_SHAP_FILT |
262 |
jmc |
1.37 |
c IF ( useSHAP_FILT .AND. shap_filt_TrStagg ) THEN |
263 |
|
|
c CALL SHAP_FILT_APPLY_TS(theta,salt,startTime,nIter0,myThid) |
264 |
|
|
c ENDIF |
265 |
jmc |
1.36 |
#endif |
266 |
|
|
#ifdef ALLOW_ZONAL_FILT |
267 |
jmc |
1.37 |
c IF ( useZONAL_FILT .AND. zonal_filt_TrStagg ) THEN |
268 |
|
|
c CALL ZONAL_FILT_APPLY_TS( theta, salt, myThid ) |
269 |
|
|
c ENDIF |
270 |
jmc |
1.36 |
#endif |
271 |
jmc |
1.37 |
c ENDIF |
272 |
jmc |
1.36 |
|
273 |
jmc |
1.9 |
#ifdef ALLOW_TIMEAVE |
274 |
adcroft |
1.41 |
#ifdef ALLOW_DEBUG |
275 |
|
|
IF (debugMode) CALL DEBUG_CALL('TIMEAVE_STATVARS',myThid) |
276 |
|
|
#endif |
277 |
jmc |
1.36 |
DO bj=myByLo(myThid),myByHi(myThid) |
278 |
|
|
DO bi=myBxLo(myThid),myBxHi(myThid) |
279 |
|
|
|
280 |
jmc |
1.7 |
C-- initialise time-average arrays with initial state values |
281 |
jmc |
1.26 |
IF (taveFreq.GT.0.) THEN |
282 |
jmc |
1.9 |
CALL TIMEAVE_STATVARS(startTime, nIter0, bi, bj, myThid) |
283 |
jmc |
1.26 |
ENDIF |
284 |
stephd |
1.31 |
cswdptr -- add --- |
285 |
|
|
#ifdef ALLOW_PTRACERS |
286 |
edhill |
1.32 |
che3 needs an IF ( usePTRACERS ) THEN |
287 |
adcroft |
1.41 |
#ifdef ALLOW_DEBUG |
288 |
|
|
IF (debugMode) CALL DEBUG_CALL('PTRACERS_STATVARS',myThid) |
289 |
|
|
#endif |
290 |
stephd |
1.31 |
CALL PTRACERS_STATVARS(startTime, nIter0, bi, bj, myThid) |
291 |
|
|
#endif |
292 |
|
|
cswdptr -- end add --- |
293 |
jmc |
1.26 |
|
294 |
|
|
C-- end bi,bj loop. |
295 |
|
|
ENDDO |
296 |
|
|
ENDDO |
297 |
jmc |
1.36 |
#endif /* ALLOW_TIMEAVE */ |
298 |
molod |
1.34 |
|
299 |
|
|
C AMM |
300 |
|
|
#ifdef ALLOW_GRIDALT |
301 |
molod |
1.35 |
if (useGRIDALT) then |
302 |
|
|
CALL TIMER_START('GRIDALT_UPDATE [INITIALISE_VARIA]',mythid) |
303 |
|
|
CALL GRIDALT_UPDATE(myThid) |
304 |
|
|
CALL TIMER_STOP ('GRIDALT_UPDATE [INITIALISE_VARIA]',mythid) |
305 |
|
|
endif |
306 |
molod |
1.34 |
#endif |
307 |
|
|
C AMM |
308 |
jmc |
1.26 |
|
309 |
|
|
C-- Fill in overlap regions for wVel : |
310 |
stephd |
1.38 |
#ifndef ALLOW_OFFLINE |
311 |
jmc |
1.26 |
_EXCH_XYZ_R8(wVel,myThid) |
312 |
stephd |
1.38 |
#endif |
313 |
jmc |
1.26 |
|
314 |
|
|
C-- Finally summarise the model state |
315 |
adcroft |
1.41 |
#ifdef ALLOW_DEBUG |
316 |
|
|
IF (debugMode) CALL DEBUG_CALL('STATE_SUMMARY',myThid) |
317 |
|
|
#endif |
318 |
jmc |
1.26 |
CALL STATE_SUMMARY( myThid ) |
319 |
jmc |
1.7 |
|
320 |
adcroft |
1.41 |
#ifdef ALLOW_DEBUG |
321 |
|
|
IF (debugMode) CALL DEBUG_LEAVE('INITIALISE_VARIA',myThid) |
322 |
|
|
#endif |
323 |
jmc |
1.7 |
RETURN |
324 |
adcroft |
1.1 |
END |