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