1 |
C $Header: /u/gcmpack/MITgcm/model/src/initialise_varia.F,v 1.39 2004/09/17 22:57:12 heimbach 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 |-- 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 |
C |-- INI_AUTODIFF |
43 |
C | |
44 |
C |-- PACKAGES_INIT_VARIABLES |
45 |
C | |
46 |
C |-- THE_CORRECTION_STEP (if restart from old pickup files) |
47 |
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 |
C |-- INTEGR_CONTINUITY |
57 |
C | |
58 |
C |-- TIMEAVE_STATVARS |
59 |
C | |
60 |
C |-- PTRACERS_STATVARS |
61 |
C | |
62 |
C |-- STATE_SUMMARY |
63 |
|
64 |
C !USES: |
65 |
IMPLICIT NONE |
66 |
C == Global variables == |
67 |
#include "SIZE.h" |
68 |
#include "EEPARAMS.h" |
69 |
#include "PARAMS.h" |
70 |
#include "DYNVARS.h" |
71 |
#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 |
|
79 |
C !INPUT/OUTPUT PARAMETERS: |
80 |
C == Routine arguments == |
81 |
INTEGER myThid |
82 |
|
83 |
C !LOCAL VARIABLES: |
84 |
C == Local variables == |
85 |
INTEGER bi,bj,K,iMin,iMax,jMin,jMax |
86 |
CEOP |
87 |
|
88 |
#ifdef ALLOW_AUTODIFF_TAMC |
89 |
|
90 |
nIter0 = INT( startTime/deltaTClock ) |
91 |
|
92 |
C-- Set Bo_surf => define the Linear Relation: Phi_surf(eta) |
93 |
CALL INI_LINEAR_PHISURF( myThid ) |
94 |
|
95 |
C-- Set coriolis operators |
96 |
CALL INI_CORI( myThid ) |
97 |
|
98 |
C-- Set laplace operators for use in 2D conjugate gradient solver. |
99 |
CALL INI_CG2D( myThid ) |
100 |
|
101 |
#ifdef ALLOW_NONHYDROSTATIC |
102 |
C-- Set laplace operators for use in 3D conjugate gradient solver. |
103 |
ceh3 should add an IF ( useNONHYDROSTATIC ) THEN |
104 |
CALL INI_CG3D( myThid ) |
105 |
#endif |
106 |
|
107 |
#endif /* ALLOW_AUTODIFF_TAMC */ |
108 |
_BARRIER |
109 |
|
110 |
C-- Initialise 3-dim. diffusivities |
111 |
CALL INI_MIXING( myThid ) |
112 |
_BARRIER |
113 |
|
114 |
C-- Initialize DYNVARS arrays (state fields + G terms: Gu,Gv,...) to zero [always] |
115 |
CALL INI_DYNVARS( myThid ) |
116 |
|
117 |
C-- Initialise model fields. |
118 |
C Starting values of U, V, W, temp., salt. and tendency terms |
119 |
C are set here. Fields are either set to default or read from |
120 |
C stored files. |
121 |
#ifndef ALLOW_OFFLINE |
122 |
CALL INI_FIELDS( myThid ) |
123 |
#endif |
124 |
_BARRIER |
125 |
|
126 |
#ifdef ALLOW_AUTODIFF_TAMC |
127 |
C-- Initialise active fields to help TAMC |
128 |
CALL INI_AUTODIFF( myThid ) |
129 |
_BARRIER |
130 |
#endif |
131 |
|
132 |
C-- Initialize variable data for packages |
133 |
CALL PACKAGES_INIT_VARIABLES( myThid ) |
134 |
|
135 |
#ifndef ALLOW_AUTODIFF_TAMC |
136 |
IF ( usePickupBeforeC35 ) THEN |
137 |
C-- IMPORTANT : Need to activate the following call to restart from |
138 |
C a pickup file written by MITgcmUV_checkpoint34 or earlier. |
139 |
IF ( startTime .NE. 0. ) THEN |
140 |
CALL THE_CORRECTION_STEP(startTime, nIter0, myThid) |
141 |
ENDIF |
142 |
ENDIF |
143 |
#endif |
144 |
|
145 |
C-- Initial conditions are convectively adjusted (for historical reasons) |
146 |
#ifndef ALLOW_OFFLINE |
147 |
IF ( startTime .EQ. 0. .AND. cAdjFreq .NE. 0. ) THEN |
148 |
CADJ loop = parallel |
149 |
DO bj = myByLo(myThid), myByHi(myThid) |
150 |
CADJ loop = parallel |
151 |
DO bi = myBxLo(myThid), myBxHi(myThid) |
152 |
iMin=1-Olx |
153 |
iMax=sNx+Olx |
154 |
jMin=1-Oly |
155 |
jMax=sNy+Oly |
156 |
CALL CONVECTIVE_ADJUSTMENT_INI( |
157 |
I bi, bj, iMin, iMax, jMin, jMax, |
158 |
I startTime, nIter0, myThid ) |
159 |
ENDDO |
160 |
ENDDO |
161 |
_BARRIER |
162 |
END IF |
163 |
#endif |
164 |
|
165 |
#ifdef NONLIN_FRSURF |
166 |
C-- Compute the surface level thickness <-- function of etaH(n) |
167 |
C and modify hFac(C,W,S) accordingly : |
168 |
IF ( select_rStar.NE.0 ) |
169 |
& CALL CALC_R_STAR(etaH, startTime, -1 , myThid ) |
170 |
IF (nonlinFreeSurf.GT.0) THEN |
171 |
IF ( select_rStar.GT.0 ) THEN |
172 |
CALL UPDATE_R_STAR( startTime, nIter0, myThid ) |
173 |
ELSE |
174 |
CALL CALC_SURF_DR(etaH, startTime, -1 , myThid ) |
175 |
CALL UPDATE_SURF_DR( startTime, nIter0, myThid ) |
176 |
ENDIF |
177 |
ENDIF |
178 |
C- update also CG2D matrix (and preconditioner) |
179 |
IF ( nonlinFreeSurf.GT.2) THEN |
180 |
CALL UPDATE_CG2D( startTime, nIter0, myThid ) |
181 |
ENDIF |
182 |
#endif /* NONLIN_FRSURF */ |
183 |
|
184 |
#ifndef ALLOW_OFFLINE |
185 |
DO bj=myByLo(myThid),myByHi(myThid) |
186 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
187 |
|
188 |
C-- Integrate continuity vertically for vertical velocity |
189 |
CALL INTEGR_CONTINUITY( bi, bj, uVel, vVel, |
190 |
I startTime, nIter0, myThid ) |
191 |
|
192 |
ENDDO |
193 |
ENDDO |
194 |
#endif |
195 |
|
196 |
#ifdef EXACT_CONSERV |
197 |
IF ( exactConserv ) THEN |
198 |
C-- Update etaH(n) : |
199 |
CALL UPDATE_ETAH( startTime, nIter0, myThid ) |
200 |
ENDIF |
201 |
#endif /* EXACT_CONSERV */ |
202 |
|
203 |
#ifdef NONLIN_FRSURF |
204 |
IF ( select_rStar.NE.0 ) THEN |
205 |
C-- r* : compute the future level thickness according to etaH(n+1) |
206 |
CALL CALC_R_STAR(etaH, startTime, nIter0, myThid ) |
207 |
ELSEIF ( nonlinFreeSurf.GT.0) THEN |
208 |
C-- compute the future surface level thickness according to etaH(n+1) |
209 |
CALL CALC_SURF_DR(etaH, startTime, nIter0, myThid ) |
210 |
ENDIF |
211 |
#endif /* NONLIN_FRSURF */ |
212 |
|
213 |
c IF ( nIter0.EQ.0 .AND. staggerTimeStep ) THEN |
214 |
C-- Filter initial T & S fields if staggerTimeStep |
215 |
C (only for backward compatibility ; to be removed later) |
216 |
#ifdef ALLOW_SHAP_FILT |
217 |
c IF ( useSHAP_FILT .AND. shap_filt_TrStagg ) THEN |
218 |
c CALL SHAP_FILT_APPLY_TS(theta,salt,startTime,nIter0,myThid) |
219 |
c ENDIF |
220 |
#endif |
221 |
#ifdef ALLOW_ZONAL_FILT |
222 |
c IF ( useZONAL_FILT .AND. zonal_filt_TrStagg ) THEN |
223 |
c CALL ZONAL_FILT_APPLY_TS( theta, salt, myThid ) |
224 |
c ENDIF |
225 |
#endif |
226 |
c ENDIF |
227 |
|
228 |
#ifdef ALLOW_TIMEAVE |
229 |
DO bj=myByLo(myThid),myByHi(myThid) |
230 |
DO bi=myBxLo(myThid),myBxHi(myThid) |
231 |
|
232 |
C-- initialise time-average arrays with initial state values |
233 |
IF (taveFreq.GT.0.) THEN |
234 |
CALL TIMEAVE_STATVARS(startTime, nIter0, bi, bj, myThid) |
235 |
ENDIF |
236 |
cswdptr -- add --- |
237 |
#ifdef ALLOW_PTRACERS |
238 |
che3 needs an IF ( usePTRACERS ) THEN |
239 |
CALL PTRACERS_STATVARS(startTime, nIter0, bi, bj, myThid) |
240 |
#endif |
241 |
cswdptr -- end add --- |
242 |
|
243 |
C-- end bi,bj loop. |
244 |
ENDDO |
245 |
ENDDO |
246 |
#endif /* ALLOW_TIMEAVE */ |
247 |
|
248 |
C AMM |
249 |
#ifdef ALLOW_GRIDALT |
250 |
if (useGRIDALT) then |
251 |
CALL TIMER_START('GRIDALT_UPDATE [INITIALISE_VARIA]',mythid) |
252 |
CALL GRIDALT_UPDATE(myThid) |
253 |
CALL TIMER_STOP ('GRIDALT_UPDATE [INITIALISE_VARIA]',mythid) |
254 |
endif |
255 |
#endif |
256 |
C AMM |
257 |
|
258 |
C-- Fill in overlap regions for wVel : |
259 |
#ifndef ALLOW_OFFLINE |
260 |
_EXCH_XYZ_R8(wVel,myThid) |
261 |
#endif |
262 |
|
263 |
C-- Finally summarise the model state |
264 |
CALL STATE_SUMMARY( myThid ) |
265 |
|
266 |
RETURN |
267 |
END |