1 |
dgoldberg |
1.1 |
C $Header: /u/gcmpack/MITgcm/model/src/ini_parms.F,v 1.271 2016/04/05 20:36:41 jmc Exp $ |
2 |
|
|
C $Name: $ |
3 |
|
|
|
4 |
|
|
#include "PACKAGES_CONFIG.h" |
5 |
|
|
#include "CPP_OPTIONS.h" |
6 |
|
|
#ifdef ALLOW_EXCH2 |
7 |
|
|
# include "W2_OPTIONS.h" |
8 |
|
|
#endif /* ALLOW_EXCH2 */ |
9 |
|
|
|
10 |
|
|
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
11 |
|
|
CBOP |
12 |
|
|
C !ROUTINE: INI_PARMS |
13 |
|
|
C !INTERFACE: |
14 |
|
|
SUBROUTINE INI_PARMS( myThid ) |
15 |
|
|
|
16 |
|
|
C !DESCRIPTION: |
17 |
|
|
C Routine to load model "parameters" from parameter file "data" |
18 |
|
|
|
19 |
|
|
C !USES: |
20 |
|
|
IMPLICIT NONE |
21 |
|
|
#include "SIZE.h" |
22 |
|
|
#include "EEPARAMS.h" |
23 |
|
|
#include "PARAMS.h" |
24 |
|
|
#ifdef ALLOW_EXCH2 |
25 |
|
|
# include "W2_EXCH2_SIZE.h" |
26 |
|
|
# include "W2_EXCH2_TOPOLOGY.h" |
27 |
|
|
#endif /* ALLOW_EXCH2 */ |
28 |
|
|
#include "EOS.h" |
29 |
|
|
C- need GRID.h to save, in rF(1), half retired param Ro_SeaLevel value: |
30 |
|
|
#include "GRID.h" |
31 |
|
|
#include "SET_GRID.h" |
32 |
|
|
|
33 |
|
|
C !INPUT/OUTPUT PARAMETERS: |
34 |
|
|
C myThid :: Number of this instance of INI_PARMS |
35 |
|
|
INTEGER myThid |
36 |
|
|
|
37 |
|
|
C !LOCAL VARIABLES: |
38 |
|
|
C dxSpacing, dySpacing :: Default spacing in X and Y. |
39 |
|
|
C :: Units are that of coordinate system |
40 |
|
|
C :: i.e. cartesian => metres |
41 |
|
|
C :: s. polar => degrees |
42 |
|
|
C SadournyCoriolis :: for backward compatibility |
43 |
|
|
C deltaTtracer :: Timestep for tracer equations ( s ) |
44 |
|
|
C forcing_In_AB :: flag to put all forcings (Temp,Salt,Tracers,Momentum) |
45 |
|
|
C :: contribution in (or out of) Adams-Bashforth time stepping. |
46 |
|
|
C goptCount :: Used to count the nuber of grid options (only one is allowed!) |
47 |
|
|
C msgBuf :: Informational/error message buffer |
48 |
|
|
C errIO :: IO error flag |
49 |
|
|
C errCount :: error counter (inconsitent params and other errors) |
50 |
|
|
C iUnit :: Work variable for IO unit number |
51 |
|
|
C k, i, j :: Loop counters |
52 |
|
|
C xxxDefault :: Default value for variable xxx |
53 |
|
|
_RL dxSpacing |
54 |
|
|
_RL dySpacing |
55 |
|
|
_RL deltaTtracer |
56 |
|
|
CHARACTER*(MAX_LEN_MBUF) msgBuf |
57 |
|
|
LOGICAL SadournyCoriolis |
58 |
|
|
LOGICAL forcing_In_AB |
59 |
|
|
INTEGER goptCount |
60 |
|
|
INTEGER gridNx, gridNy |
61 |
|
|
INTEGER k, i, j, iUnit |
62 |
|
|
INTEGER errIO, errCount |
63 |
|
|
C Default values for variables which have vertical coordinate system |
64 |
|
|
C dependency. |
65 |
|
|
_RL viscArDefault |
66 |
|
|
_RL diffKrTDefault |
67 |
|
|
_RL diffKrSDefault |
68 |
|
|
_RL hFacMinDrDefault |
69 |
|
|
_RL delRDefault(Nr) |
70 |
|
|
C zCoordInputData :: Variables used to select between different coordinate systems. |
71 |
|
|
C pCoordInputData :: The vertical coordinate system in the rest of the model is |
72 |
|
|
C rCoordInputData :: written in terms of r. In the model "data" file input data |
73 |
|
|
C coordsSet :: can be interms of z, p or r. |
74 |
|
|
C :: e.g. delZ or delP or delR for the vertical grid spacing. |
75 |
|
|
C :: The following rules apply: |
76 |
|
|
C :: All parameters must use the same vertical coordinate system. |
77 |
|
|
C :: e.g. delZ and viscAz is legal but |
78 |
|
|
C :: delZ and viscAr is an error. |
79 |
|
|
C :: Similarly specifying delZ and delP is an error. |
80 |
|
|
C :: zCoord..., pCoord..., rCoord... are used to flag when |
81 |
|
|
C :: z, p or r are used. |
82 |
|
|
C :: coordsSet counts how many vertical coordinate systems have |
83 |
|
|
C :: been used to specify variables. coordsSet > 1 is an error. |
84 |
|
|
C vertSetCount :: to count number of vertical array elements which are set. |
85 |
|
|
C specifiedDiffKrT :: internal flag, true when any diffK[r,z,p,Nr]T is specified |
86 |
|
|
C specifiedDiffKrS :: internal flag, true when any diffK[r,z,p,Nr]S is specified |
87 |
|
|
|
88 |
|
|
LOGICAL zCoordInputData |
89 |
|
|
LOGICAL pCoordInputData |
90 |
|
|
LOGICAL rCoordInputData |
91 |
|
|
INTEGER coordsSet |
92 |
|
|
INTEGER vertSetCount |
93 |
|
|
LOGICAL specifiedDiffKrT, specifiedDiffKrS |
94 |
|
|
|
95 |
|
|
C Variables which have vertical coordinate system dependency. |
96 |
|
|
C delZ :: Vertical grid spacing ( m ). |
97 |
|
|
C delP :: Vertical grid spacing ( Pa ). |
98 |
|
|
C viscAz :: Eddy viscosity coeff. for mixing of momentum vertically ( m^2/s ) |
99 |
|
|
C viscAp :: Eddy viscosity coeff. for mixing of momentum vertically ( Pa^2/s ) |
100 |
|
|
C diffKzT :: Laplacian diffusion coeff. for mixing of heat vertically ( m^2/s ) |
101 |
|
|
C diffKpT :: Laplacian diffusion coeff. for mixing of heat vertically ( Pa^2/s ) |
102 |
|
|
C diffKzS :: Laplacian diffusion coeff. for mixing of salt vertically ( m^2/s ) |
103 |
|
|
C diffKpS :: Laplacian diffusion coeff. for mixing of salt vertically ( Pa^2/s ) |
104 |
|
|
_RL delZ(Nr) |
105 |
|
|
_RL delP(Nr) |
106 |
|
|
_RL viscAz |
107 |
|
|
_RL viscAp |
108 |
|
|
_RL viscAr |
109 |
|
|
_RL diffKzT |
110 |
|
|
_RL diffKpT |
111 |
|
|
_RL diffKrT |
112 |
|
|
_RL diffKzS |
113 |
|
|
_RL diffKpS |
114 |
|
|
_RL diffKrS |
115 |
|
|
|
116 |
|
|
C Retired main data file parameters. Kept here to trap use of old data files. |
117 |
|
|
C nRetired :: Counter used to trap gracefully namelists containing "retired" |
118 |
|
|
C :: parameters. These are parameters that are either no-longer used |
119 |
|
|
C or that have moved to a different input file and/or namelist. |
120 |
|
|
C Namelist PARM01: |
121 |
|
|
C useConstantF :: Coriolis coeff set to f0 (replaced by selectCoriMap=0) |
122 |
|
|
C useBetaPlaneF :: Coriolis coeff = f0 + beta.y (replaced by selectCoriMap=1) |
123 |
|
|
C useSphereF :: Coriolis = 2.omega.sin(phi) (replaced by selectCoriMap=2) |
124 |
|
|
C tracerAdvScheme :: tracer advection scheme (old passive tracer code) |
125 |
|
|
C trac_EvPrRn :: tracer conc. in Rain & Evap (old passive tracer code) |
126 |
|
|
C saltDiffusion :: diffusion of salinity on/off (flag not used) |
127 |
|
|
C tempDiffusion :: diffusion of temperature on/off (flag not used) |
128 |
|
|
C zonal_filt_lat :: Moved to package "zonal_filt" |
129 |
|
|
C gravitySign :: direction of gravity relative to R direction |
130 |
|
|
C :: (removed from namelist and set according to z/p coordinate) |
131 |
|
|
C viscAstrain :: replaced by standard viscosity coeff & useStrainTensionVisc |
132 |
|
|
C viscAtension :: replaced by standard viscosity coeff & useStrainTensionVisc |
133 |
|
|
C useAnisotropicViscAgridMax :: Changed to be default behavior. Can |
134 |
|
|
C use old method by setting useAreaViscLength=.true. |
135 |
|
|
C usePickupBeforeC35 :: to restart from old-pickup files (generated with code |
136 |
|
|
C from before checkpoint-35, Feb 08, 2001): disabled (Jan 2007) |
137 |
|
|
C debugMode :: to print debug msg. now read from parameter file eedata |
138 |
|
|
C allowInteriorFreezing :: Allow water at depth to freeze and rise to the surface |
139 |
|
|
C (replaced by pkg/frazil) |
140 |
|
|
C useOldFreezing :: use the old version (before checkpoint52a_pre, 2003-11-12) |
141 |
|
|
C Namelist PARM03: |
142 |
|
|
C tauThetaClimRelax3Dim :: replaced by pkg/rbcs (3.D Relaxation B.Cs) |
143 |
|
|
C tauSaltClimRelax3Dim :: replaced by pkg/rbcs (3.D Relaxation B.Cs) |
144 |
|
|
C calendarDumps :: moved to package "cal" (calendar) |
145 |
|
|
C Namelist PARM04: |
146 |
|
|
C Ro_SeaLevel :: origin of the vertical R-coords axis ; |
147 |
|
|
C :: replaced by top_Pres or seaLev_Z setting |
148 |
|
|
C groundAtK1 :: put the surface(k=1) at the ground (replaced by usingPCoords) |
149 |
|
|
C rkFac :: removed from namelist ; replaced by -rkSign |
150 |
|
|
C thetaMin :: unfortunate variable name ; replaced by xgOrigin |
151 |
|
|
C phiMin :: unfortunate variable name ; replaced by ygOrigin |
152 |
|
|
C Namelist PARM05: |
153 |
|
|
C shelfIceFile :: File containing the topography of the shelfice draught |
154 |
|
|
C (replaced by SHELFICEtopoFile in SHELFICE.h) |
155 |
|
|
C dQdTfile :: File containing thermal relaxation coefficient |
156 |
|
|
|
157 |
|
|
INTEGER nRetired |
158 |
|
|
LOGICAL useConstantF, useBetaPlaneF, useSphereF |
159 |
|
|
LOGICAL tempDiffusion, saltDiffusion |
160 |
|
|
INTEGER tracerAdvScheme |
161 |
|
|
_RL trac_EvPrRn |
162 |
|
|
_RL zonal_filt_lat |
163 |
|
|
c _RL gravitySign |
164 |
|
|
_RL viscAstrain, viscAtension |
165 |
|
|
LOGICAL useAnisotropicViscAgridMax |
166 |
|
|
LOGICAL usePickupBeforeC35 |
167 |
|
|
LOGICAL saveDebugMode |
168 |
|
|
LOGICAL allowInteriorFreezing, useOldFreezing |
169 |
|
|
C- |
170 |
|
|
_RL tauThetaClimRelax3Dim, tauSaltClimRelax3Dim |
171 |
|
|
LOGICAL calendarDumps |
172 |
|
|
C- |
173 |
|
|
LOGICAL groundAtK1 |
174 |
|
|
_RL Ro_SeaLevel |
175 |
|
|
_RL rkFac |
176 |
|
|
_RL thetaMin, phiMin |
177 |
|
|
CHARACTER*(MAX_LEN_FNAM) shelfIceFile |
178 |
|
|
CHARACTER*(MAX_LEN_FNAM) dQdTfile |
179 |
|
|
|
180 |
|
|
C-- Continuous equation parameters |
181 |
|
|
NAMELIST /PARM01/ |
182 |
|
|
& gravitySign, nh_Am2, |
183 |
|
|
& gravity, gBaro, gravityFile, rhonil, tAlpha, sBeta, |
184 |
|
|
& selectCoriMap, f0, beta, fPrime, omega, rotationPeriod, |
185 |
|
|
& viscAh, viscAhW, viscAhMax, |
186 |
|
|
& viscAhGrid, viscAhGridMax, viscAhGridMin, |
187 |
|
|
& viscC2leith, viscC4leith, smag3D_coeff, useSmag3D, |
188 |
|
|
& useFullLeith, useAnisotropicViscAgridMax, useStrainTensionVisc, |
189 |
|
|
& useAreaViscLength, |
190 |
|
|
& viscC2leithD, viscC4leithD, viscC2smag, viscC4smag, |
191 |
|
|
& viscAhD, viscAhZ, viscA4D, viscA4Z, |
192 |
|
|
& viscA4, viscA4W, |
193 |
|
|
& viscA4Max, viscA4Grid, viscA4GridMax, viscA4GridMin, |
194 |
|
|
& viscA4ReMax, viscAhReMax, |
195 |
|
|
& cosPower, viscAstrain, viscAtension, |
196 |
|
|
& diffKhT, diffK4T, diffKhS, diffK4S, |
197 |
|
|
& tRef, sRef, tRefFile, sRefFile, rhoRefFile, |
198 |
|
|
& eosType, selectP_inEOS_Zc, integr_GeoPot, selectFindRoSurf, |
199 |
|
|
& HeatCapacity_Cp, celsius2K, atm_Cp, atm_Rd, atm_Rq, atm_Po, |
200 |
|
|
& no_slip_sides, sideDragFactor, |
201 |
|
|
& no_slip_bottom, bottomVisc_pCell, |
202 |
|
|
& bottomDragLinear, bottomDragQuadratic, selectBotDragQuadr, |
203 |
|
|
& momViscosity, momAdvection, momForcing, useCoriolis, |
204 |
|
|
& useConstantF, useBetaPlaneF, useSphereF, use3dCoriolis, |
205 |
|
|
& momPressureForcing, metricTerms, vectorInvariantMomentum, |
206 |
|
|
& tempDiffusion, tempAdvection, tempForcing, addFrictionHeating, |
207 |
|
|
& saltDiffusion, saltAdvection, saltForcing, |
208 |
|
|
& implicSurfPress, implicDiv2Dflow, implicitNHPress, |
209 |
|
|
& implicitFreeSurface, rigidLid, freeSurfFac, |
210 |
|
|
& hFacMin, hFacMinDz, hFacMinDp, hFacMinDr, |
211 |
|
|
& exactConserv, linFSConserveTr, uniformLin_PhiSurf, |
212 |
|
|
& nonlinFreeSurf, hFacInf, hFacSup, select_rStar, |
213 |
|
|
& nonHydrostatic, selectNHfreeSurf, quasiHydrostatic, |
214 |
|
|
& implicitIntGravWave, staggerTimeStep, doResetHFactors, |
215 |
|
|
& tempStepping, saltStepping, momStepping, |
216 |
|
|
& implicitDiffusion, implicitViscosity, implBottomFriction, |
217 |
|
|
& tempImplVertAdv, saltImplVertAdv, momImplVertAdv, |
218 |
|
|
& viscAz, diffKzT, diffKzS, viscAp, diffKpT, diffKpS, |
219 |
|
|
& viscAr, diffKrT, diffKrS, viscArNr, diffKrNrT, diffKrNrS, |
220 |
|
|
& diffKr4T, diffKr4S, BL79LatVary, |
221 |
|
|
& diffKrBL79surf, diffKrBL79deep, diffKrBL79scl, diffKrBL79Ho, |
222 |
|
|
& diffKrBLEQsurf, diffKrBLEQdeep, diffKrBLEQscl, diffKrBLEQHo, |
223 |
|
|
& rhoConst, thetaConst, rhoConstFresh, buoyancyRelation, |
224 |
|
|
& writeBinaryPrec, readBinaryPrec, writeStatePrec, |
225 |
|
|
& globalFiles, useSingleCpuIO, useSingleCpuInput, |
226 |
|
|
& allowFreezing, allowInteriorFreezing, useOldFreezing, ivdc_kappa, |
227 |
|
|
& hMixCriteria, dRhoSmall, hMixSmooth, |
228 |
|
|
& usePickupBeforeC35, usePickupBeforeC54, debugMode, debugLevel, |
229 |
|
|
& tempAdvScheme, tempVertAdvScheme, |
230 |
|
|
& saltAdvScheme, saltVertAdvScheme, tracerAdvScheme, |
231 |
|
|
& multiDimAdvection, useEnergyConservingCoriolis, |
232 |
|
|
& useCDscheme, useJamartWetPoints, useJamartMomAdv, useNHMTerms, |
233 |
|
|
& selectVortScheme, upwindVorticity, highOrderVorticity, |
234 |
|
|
& SadournyCoriolis, useAbsVorticity, upwindShear, selectKEscheme, |
235 |
|
|
& selectAddFluid, useRealFreshWaterFlux, convertFW2Salt, |
236 |
|
|
& temp_EvPrRn, salt_EvPrRn, trac_EvPrRn, |
237 |
|
|
& temp_addMass, salt_addMass, zonal_filt_lat, |
238 |
|
|
& smoothAbsFuncRange, |
239 |
|
|
& balanceEmPmR, balanceQnet, balancePrintMean, |
240 |
|
|
& balanceThetaClimRelax, balanceSaltClimRelax |
241 |
|
|
|
242 |
|
|
C-- Elliptic solver parameters |
243 |
|
|
NAMELIST /PARM02/ |
244 |
|
|
& cg2dMaxIters, cg2dChkResFreq, cg2dUseMinResSol, |
245 |
|
|
& cg2dTargetResidual, cg2dTargetResWunit, |
246 |
|
|
& cg2dpcOffDFac, cg2dPreCondFreq, |
247 |
|
|
& cg3dMaxIters, cg3dChkResFreq, cg3dTargetResidual, |
248 |
|
|
& useSRCGSolver, printResidualFreq |
249 |
|
|
#ifdef ALLOW_PETSC |
250 |
|
|
#ifdef ALLOW_NONHYDROSTATIC |
251 |
|
|
& , CG3D_PETSC_SOLVER_TYPE, |
252 |
|
|
& CG3D_PETSC_PRECOND_TYPE, |
253 |
|
|
& use_cg3d_petsc, |
254 |
|
|
& cg3d_petsc_reuse_mat, |
255 |
|
|
& cg3d_petsc_cpuInVert |
256 |
|
|
#endif |
257 |
|
|
#endif |
258 |
|
|
|
259 |
|
|
C-- Time stepping parammeters |
260 |
|
|
NAMELIST /PARM03/ |
261 |
|
|
& nIter0, nTimeSteps, nTimeSteps_l2, nEndIter, |
262 |
|
|
& baseTime, startTime, endTime, |
263 |
|
|
& deltaT, deltaTClock, deltaTMom, |
264 |
|
|
& deltaTtracer, dTtracerLev, deltaTFreeSurf, |
265 |
|
|
& forcing_In_AB, momForcingOutAB, tracForcingOutAB, |
266 |
|
|
& momDissip_In_AB, doAB_onGtGs, |
267 |
|
|
& abEps, alph_AB, beta_AB, startFromPickupAB2, applyExchUV_early, |
268 |
|
|
& tauCD, rCD, epsAB_CD, cAdjFreq, |
269 |
|
|
& chkPtFreq, pChkPtFreq, pickupSuff, pickupStrictlyMatch, |
270 |
|
|
& writePickupAtEnd, |
271 |
|
|
& dumpFreq, dumpInitAndLast, adjDumpFreq, taveFreq, tave_lastIter, |
272 |
|
|
& diagFreq, monitorFreq, adjMonitorFreq, monitorSelect, |
273 |
|
|
& outputTypesInclusive, |
274 |
|
|
& tauThetaClimRelax, tauSaltClimRelax, latBandClimRelax, |
275 |
|
|
& tauThetaClimRelax3Dim, tauSaltClimRelax3Dim, |
276 |
|
|
& periodicExternalForcing, externForcingPeriod, externForcingCycle, |
277 |
|
|
& calendarDumps |
278 |
|
|
|
279 |
|
|
C-- Gridding parameters |
280 |
|
|
NAMELIST /PARM04/ |
281 |
|
|
& usingCartesianGrid, usingCylindricalGrid, |
282 |
|
|
& usingSphericalPolarGrid, usingCurvilinearGrid, |
283 |
|
|
& xgOrigin, ygOrigin, dxSpacing, dySpacing, |
284 |
|
|
& delX, delY, delXFile, delYFile, horizGridFile, |
285 |
|
|
& phiEuler, thetaEuler, psiEuler, |
286 |
|
|
& rSphere, radius_fromHorizGrid, deepAtmosphere, seaLev_Z, |
287 |
|
|
& top_Pres, delZ, delP, delR, delRc, delRFile, delRcFile, |
288 |
|
|
& selectSigmaCoord, rSigmaBnd, hybSigmFile, |
289 |
|
|
& Ro_SeaLevel, rkFac, groundAtK1, thetaMin, phiMin |
290 |
|
|
|
291 |
|
|
C-- Input files |
292 |
|
|
NAMELIST /PARM05/ |
293 |
|
|
& bathyFile, topoFile, addWwallFile, addSwallFile, shelfIceFile, |
294 |
|
|
& diffKrFile, viscAhDfile, viscAhZfile, viscA4Dfile, viscA4Zfile, |
295 |
|
|
& hydrogThetaFile, hydrogSaltFile, |
296 |
|
|
& maskIniTemp, maskIniSalt, checkIniTemp, checkIniSalt, |
297 |
|
|
& zonalWindFile, meridWindFile, |
298 |
|
|
& thetaClimFile, saltClimFile, |
299 |
|
|
& surfQfile, surfQnetFile, surfQswFile, EmPmRfile, saltFluxFile, |
300 |
|
|
& lambdaThetaFile, lambdaSaltFile, |
301 |
|
|
& uVelInitFile, vVelInitFile, pSurfInitFile, |
302 |
|
|
& dQdTFile, ploadFile, addMassFile, tCylIn, tCylOut, |
303 |
|
|
& eddyPsiXFile, eddyPsiYFile, geothermalFile, |
304 |
|
|
& mdsioLocalDir, adTapeDir, |
305 |
|
|
& the_run_name |
306 |
|
|
CEOP |
307 |
|
|
|
308 |
|
|
#ifdef ALLOW_EXCH2 |
309 |
|
|
gridNx = exch2_mydNx(1) |
310 |
|
|
gridNy = exch2_mydNy(1) |
311 |
|
|
#else /* ALLOW_EXCH2 */ |
312 |
|
|
gridNx = Nx |
313 |
|
|
gridNy = Ny |
314 |
|
|
#endif /* ALLOW_EXCH2 */ |
315 |
|
|
|
316 |
|
|
_BEGIN_MASTER(myThid) |
317 |
|
|
|
318 |
|
|
C Defaults values for input parameters |
319 |
|
|
CALL SET_DEFAULTS( |
320 |
|
|
O viscArDefault, diffKrTDefault, diffKrSDefault, |
321 |
|
|
O hFacMinDrDefault, delRDefault, |
322 |
|
|
I myThid ) |
323 |
|
|
SadournyCoriolis = .FALSE. |
324 |
|
|
|
325 |
|
|
C-- Initialise "which vertical coordinate system used" flags. |
326 |
|
|
zCoordInputData = .FALSE. |
327 |
|
|
pCoordInputData = .FALSE. |
328 |
|
|
rCoordInputData = .FALSE. |
329 |
|
|
coordsSet = 0 |
330 |
|
|
|
331 |
|
|
C-- Initialise retired parameters to unlikely value |
332 |
|
|
nRetired = 0 |
333 |
|
|
useConstantF = .FALSE. |
334 |
|
|
useBetaPlaneF = .FALSE. |
335 |
|
|
useSphereF = .TRUE. |
336 |
|
|
tempDiffusion = .TRUE. |
337 |
|
|
saltDiffusion = .TRUE. |
338 |
|
|
tracerAdvScheme = UNSET_I |
339 |
|
|
trac_EvPrRn = UNSET_RL |
340 |
|
|
zonal_filt_lat = UNSET_RL |
341 |
|
|
gravitySign = UNSET_RL |
342 |
|
|
viscAstrain = UNSET_RL |
343 |
|
|
viscAtension = UNSET_RL |
344 |
|
|
useAnisotropicViscAgridMax=.TRUE. |
345 |
|
|
usePickupBeforeC35 = .FALSE. |
346 |
|
|
saveDebugMode = debugMode |
347 |
|
|
allowInteriorFreezing = .FALSE. |
348 |
|
|
useOldFreezing = .FALSE. |
349 |
|
|
tauThetaClimRelax3Dim = UNSET_RL |
350 |
|
|
tauSaltClimRelax3Dim = UNSET_RL |
351 |
|
|
calendarDumps = .FALSE. |
352 |
|
|
Ro_SeaLevel = UNSET_RL |
353 |
|
|
rkFac = UNSET_RL |
354 |
|
|
groundAtK1 = .FALSE. |
355 |
|
|
thetaMin = UNSET_RL |
356 |
|
|
phiMin = UNSET_RL |
357 |
|
|
shelfIceFile = ' ' |
358 |
|
|
dQdTFile = ' ' |
359 |
|
|
|
360 |
|
|
C-- Open the parameter file |
361 |
|
|
WRITE(msgBuf,'(A)') |
362 |
|
|
& ' INI_PARMS: opening model parameter file "data"' |
363 |
|
|
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
364 |
|
|
& SQUEEZE_RIGHT, myThid ) |
365 |
|
|
|
366 |
|
|
CALL OPEN_COPY_DATA_FILE( 'data', 'INI_PARMS', |
367 |
|
|
O iUnit, myThid ) |
368 |
|
|
|
369 |
|
|
C-- Read settings from iUnit (= a copy of model parameter file "data"). |
370 |
|
|
errIO = 0 |
371 |
|
|
errCount = 0 |
372 |
|
|
|
373 |
|
|
C-- Set default "physical" parameters |
374 |
|
|
viscAhW = UNSET_RL |
375 |
|
|
viscA4W = UNSET_RL |
376 |
|
|
viscAhD = UNSET_RL |
377 |
|
|
viscAhZ = UNSET_RL |
378 |
|
|
viscA4D = UNSET_RL |
379 |
|
|
viscA4Z = UNSET_RL |
380 |
|
|
viscAz = UNSET_RL |
381 |
|
|
viscAp = UNSET_RL |
382 |
|
|
viscAr = UNSET_RL |
383 |
|
|
diffKzT = UNSET_RL |
384 |
|
|
diffKpT = UNSET_RL |
385 |
|
|
diffKrT = UNSET_RL |
386 |
|
|
diffKzS = UNSET_RL |
387 |
|
|
diffKpS = UNSET_RL |
388 |
|
|
diffKrS = UNSET_RL |
389 |
|
|
hFacMinDr = UNSET_RL |
390 |
|
|
hFacMinDz = UNSET_RL |
391 |
|
|
hFacMinDp = UNSET_RL |
392 |
|
|
tAlpha = UNSET_RL |
393 |
|
|
sBeta = UNSET_RL |
394 |
|
|
implicitNHPress = UNSET_RL |
395 |
|
|
tempVertAdvScheme = 0 |
396 |
|
|
saltVertAdvScheme = 0 |
397 |
|
|
C-- z,p,r coord input switching. |
398 |
|
|
WRITE(msgBuf,'(A)') ' INI_PARMS ; starts to read PARM01' |
399 |
|
|
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
400 |
|
|
& SQUEEZE_RIGHT, myThid ) |
401 |
|
|
READ(UNIT=iUnit,NML=PARM01) !,IOSTAT=errIO) |
402 |
|
|
IF ( errIO .LT. 0 ) THEN |
403 |
|
|
WRITE(msgBuf,'(A)') |
404 |
|
|
& 'S/R INI_PARMS: Error reading model parameter file "data"' |
405 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
406 |
|
|
WRITE(msgBuf,'(A)') 'S/R INI_PARMS: Problem in namelist PARM01' |
407 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
408 |
|
|
STOP 'ABNORMAL END: S/R INI_PARMS' |
409 |
|
|
ELSE |
410 |
|
|
WRITE(msgBuf,'(A)') ' INI_PARMS ; read PARM01 : OK' |
411 |
|
|
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
412 |
|
|
& SQUEEZE_RIGHT, myThid ) |
413 |
|
|
ENDIF |
414 |
|
|
|
415 |
|
|
C- set the type of vertical coordinate and type of fluid |
416 |
|
|
C according to buoyancyRelation |
417 |
|
|
usingPCoords = .FALSE. |
418 |
|
|
usingZCoords = .FALSE. |
419 |
|
|
fluidIsAir = .FALSE. |
420 |
|
|
fluidIsWater = .FALSE. |
421 |
|
|
IF ( buoyancyRelation.EQ.'ATMOSPHERIC' ) THEN |
422 |
|
|
usingPCoords = .TRUE. |
423 |
|
|
fluidIsAir = .TRUE. |
424 |
|
|
ELSEIF ( buoyancyRelation.EQ.'OCEANICP') THEN |
425 |
|
|
usingPCoords = .TRUE. |
426 |
|
|
fluidIsWater = .TRUE. |
427 |
|
|
ELSEIF ( buoyancyRelation.EQ.'OCEANIC' ) THEN |
428 |
|
|
usingZCoords = .TRUE. |
429 |
|
|
fluidIsWater = .TRUE. |
430 |
|
|
ELSE |
431 |
|
|
WRITE(msgBuf,'(2A)') 'S/R INI_PARMS:', |
432 |
|
|
& ' Bad value of buoyancyRelation ' |
433 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
434 |
|
|
errCount = errCount + 1 |
435 |
|
|
ENDIF |
436 |
|
|
|
437 |
|
|
IF ( .NOT.rigidLid .AND. |
438 |
|
|
& .NOT.implicitFreeSurface ) THEN |
439 |
|
|
C- No barotropic solver selected => use implicitFreeSurface as default |
440 |
|
|
WRITE(msgBuf,'(A)') |
441 |
|
|
& 'S/R INI_PARMS: No request for barotropic solver' |
442 |
|
|
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
443 |
|
|
& SQUEEZE_RIGHT, myThid ) |
444 |
|
|
WRITE(msgBuf,'(A)') |
445 |
|
|
& 'S/R INI_PARMS: => Use implicitFreeSurface as default' |
446 |
|
|
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
447 |
|
|
& SQUEEZE_RIGHT, myThid ) |
448 |
|
|
implicitFreeSurface = .TRUE. |
449 |
|
|
ENDIF |
450 |
|
|
IF ( implicitFreeSurface ) freeSurfFac = 1. _d 0 |
451 |
|
|
IF ( rigidLid ) freeSurfFac = 0. _d 0 |
452 |
|
|
IF ( gBaro .EQ. UNSET_RL ) gBaro=gravity |
453 |
|
|
IF ( rhoConst .EQ. UNSET_RL ) rhoConst=rhoNil |
454 |
|
|
IF ( rhoConstFresh .EQ. UNSET_RL ) rhoConstFresh=rhoConst |
455 |
|
|
IF ( implicitNHPress.EQ.UNSET_RL ) |
456 |
|
|
& implicitNHPress = implicSurfPress |
457 |
|
|
IF ( omega .EQ. UNSET_RL ) THEN |
458 |
|
|
omega = 0. _d 0 |
459 |
|
|
IF ( rotationPeriod .NE. 0. _d 0 ) |
460 |
|
|
& omega = 2. _d 0 * PI / rotationPeriod |
461 |
|
|
ELSEIF ( omega .EQ. 0. _d 0 ) THEN |
462 |
|
|
rotationPeriod = 0. _d 0 |
463 |
|
|
ELSE |
464 |
|
|
rotationPeriod = 2. _d 0 * PI / omega |
465 |
|
|
ENDIF |
466 |
|
|
IF ( atm_Rd .EQ. UNSET_RL ) THEN |
467 |
|
|
atm_Rd = atm_Cp * atm_kappa |
468 |
|
|
ELSE |
469 |
|
|
atm_kappa = atm_Rd / atm_Cp |
470 |
|
|
ENDIF |
471 |
|
|
C-- Non-hydrostatic/quasi-hydrostatic |
472 |
|
|
IF (nonHydrostatic.AND.quasiHydrostatic) THEN |
473 |
|
|
WRITE(msgBuf,'(A)') |
474 |
|
|
& 'Illegal: both nonHydrostatic = quasiHydrostatic = TRUE' |
475 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
476 |
|
|
errCount = errCount + 1 |
477 |
|
|
ENDIF |
478 |
|
|
C-- Advection and Forcing for Temp and salt |
479 |
|
|
IF (tempVertAdvScheme.EQ.0) tempVertAdvScheme = tempAdvScheme |
480 |
|
|
IF (saltVertAdvScheme.EQ.0) saltVertAdvScheme = saltAdvScheme |
481 |
|
|
C-- horizontal viscosity (acting on Divergence or Vorticity) |
482 |
|
|
IF ( viscAhD .EQ. UNSET_RL ) viscAhD = viscAh |
483 |
|
|
IF ( viscAhZ .EQ. UNSET_RL ) viscAhZ = viscAh |
484 |
|
|
IF ( viscA4D .EQ. UNSET_RL ) viscA4D = viscA4 |
485 |
|
|
IF ( viscA4Z .EQ. UNSET_RL ) viscA4Z = viscA4 |
486 |
|
|
C-- horizontal viscosity for vertical moments |
487 |
|
|
IF ( viscAhW .EQ. UNSET_RL ) viscAhW = viscAhD |
488 |
|
|
IF ( viscA4W .EQ. UNSET_RL ) viscA4W = viscA4D |
489 |
|
|
C-- z,p,r coord input switching. |
490 |
|
|
IF ( viscAz .NE. UNSET_RL ) zCoordInputData = .TRUE. |
491 |
|
|
IF ( viscAp .NE. UNSET_RL ) pCoordInputData = .TRUE. |
492 |
|
|
IF ( viscAr .NE. UNSET_RL ) rCoordInputData = .TRUE. |
493 |
|
|
IF ( viscAr .EQ. UNSET_RL ) viscAr = viscAz |
494 |
|
|
IF ( viscAr .EQ. UNSET_RL ) viscAr = viscAp |
495 |
|
|
vertSetCount = 0 |
496 |
|
|
DO k=1,Nr |
497 |
|
|
IF ( viscArNr(k).NE.UNSET_RL ) vertSetCount = vertSetCount + 1 |
498 |
|
|
ENDDO |
499 |
|
|
IF ( vertSetCount.GT.0 .AND. vertSetCount.LT.Nr ) THEN |
500 |
|
|
WRITE(msgBuf,'(A,2(I5,A))') 'S/R INI_PARMS: Partial setting (', |
501 |
|
|
& vertSetCount, ' /', Nr, ') of viscArNr is not allowed' |
502 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
503 |
|
|
errCount = errCount + 1 |
504 |
|
|
ENDIF |
505 |
|
|
IF ( viscAr .EQ. UNSET_RL ) THEN |
506 |
|
|
viscAr = viscArDefault |
507 |
|
|
ELSEIF ( vertSetCount.GT.0 ) THEN |
508 |
|
|
WRITE(msgBuf,'(2A)') 'S/R INI_PARMS: Cannot set both ', |
509 |
|
|
& 'viscArNr and viscAr (or Ap,Az) in param file data' |
510 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
511 |
|
|
errCount = errCount + 1 |
512 |
|
|
ENDIF |
513 |
|
|
IF ( vertSetCount.EQ.0 ) THEN |
514 |
|
|
DO k=1,Nr |
515 |
|
|
viscArNr(k) = viscAr |
516 |
|
|
ENDDO |
517 |
|
|
ENDIF |
518 |
|
|
#ifdef ALLOW_MOM_COMMON |
519 |
|
|
C- set default scheme for quadratic bottom-drag |
520 |
|
|
IF ( selectBotDragQuadr.EQ.-1 .AND. bottomDragQuadratic.NE.0. ) |
521 |
|
|
& selectBotDragQuadr = 0 |
522 |
|
|
#endif /* ALLOW_MOM_COMMON */ |
523 |
|
|
|
524 |
|
|
IF ( diffKzT .NE. UNSET_RL ) zCoordInputData = .TRUE. |
525 |
|
|
IF ( diffKpT .NE. UNSET_RL ) pCoordInputData = .TRUE. |
526 |
|
|
IF ( diffKrT .NE. UNSET_RL ) rCoordInputData = .TRUE. |
527 |
|
|
IF ( diffKrT .EQ. UNSET_RL ) diffKrT = diffKzT |
528 |
|
|
IF ( diffKrT .EQ. UNSET_RL ) diffKrT = diffKpT |
529 |
|
|
vertSetCount = 0 |
530 |
|
|
DO k=1,Nr |
531 |
|
|
IF ( diffKrNrT(k).NE.UNSET_RL ) vertSetCount = vertSetCount + 1 |
532 |
|
|
ENDDO |
533 |
|
|
IF ( vertSetCount.GT.0 .AND. vertSetCount.LT.Nr ) THEN |
534 |
|
|
WRITE(msgBuf,'(A,2(I5,A))') 'S/R INI_PARMS: Partial setting (', |
535 |
|
|
& vertSetCount, ' /', Nr, ') of diffKrNrT is not allowed' |
536 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
537 |
|
|
errCount = errCount + 1 |
538 |
|
|
ENDIF |
539 |
|
|
specifiedDiffKrT = vertSetCount.EQ.Nr |
540 |
|
|
IF ( diffKrT .EQ. UNSET_RL ) THEN |
541 |
|
|
diffKrT = diffKrTDefault |
542 |
|
|
ELSEIF ( vertSetCount.GT.0 ) THEN |
543 |
|
|
WRITE(msgBuf,'(2A)') 'S/R INI_PARMS: Cannot set both ', |
544 |
|
|
& 'diffKrNrT and diffKrT (or Kp,Kz) in param file data' |
545 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
546 |
|
|
errCount = errCount + 1 |
547 |
|
|
ELSE |
548 |
|
|
specifiedDiffKrT = .TRUE. |
549 |
|
|
ENDIF |
550 |
|
|
IF ( vertSetCount.EQ.0 ) THEN |
551 |
|
|
DO k=1,Nr |
552 |
|
|
diffKrNrT(k) = diffKrT |
553 |
|
|
ENDDO |
554 |
|
|
ENDIF |
555 |
|
|
|
556 |
|
|
IF ( diffKzS .NE. UNSET_RL ) zCoordInputData = .TRUE. |
557 |
|
|
IF ( diffKpS .NE. UNSET_RL ) pCoordInputData = .TRUE. |
558 |
|
|
IF ( diffKrS .NE. UNSET_RL ) rCoordInputData = .TRUE. |
559 |
|
|
IF ( diffKrS .EQ. UNSET_RL ) diffKrS = diffKzS |
560 |
|
|
IF ( diffKrS .EQ. UNSET_RL ) diffKrS = diffKpS |
561 |
|
|
vertSetCount = 0 |
562 |
|
|
DO k=1,Nr |
563 |
|
|
IF ( diffKrNrS(k).NE.UNSET_RL ) vertSetCount = vertSetCount + 1 |
564 |
|
|
ENDDO |
565 |
|
|
IF ( vertSetCount.GT.0 .AND. vertSetCount.LT.Nr ) THEN |
566 |
|
|
WRITE(msgBuf,'(A,2(I5,A))') 'S/R INI_PARMS: Partial setting (', |
567 |
|
|
& vertSetCount, ' /', Nr, ') of diffKrNrS is not allowed' |
568 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
569 |
|
|
errCount = errCount + 1 |
570 |
|
|
ENDIF |
571 |
|
|
IF ( vertSetCount.EQ.Nr ) THEN |
572 |
|
|
specifiedDiffKrS = .TRUE. |
573 |
|
|
IF ( diffKrS.NE.UNSET_RL ) THEN |
574 |
|
|
WRITE(msgBuf,'(2A)') 'S/R INI_PARMS: Cannot set both ', |
575 |
|
|
& 'diffKrNrS and diffKrS (or Kp,Kz) in param file data' |
576 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
577 |
|
|
errCount = errCount + 1 |
578 |
|
|
ENDIF |
579 |
|
|
ELSEIF ( diffKrS.NE.UNSET_RL ) THEN |
580 |
|
|
specifiedDiffKrS = .TRUE. |
581 |
|
|
DO k=1,Nr |
582 |
|
|
diffKrNrS(k) = diffKrS |
583 |
|
|
ENDDO |
584 |
|
|
ELSE |
585 |
|
|
specifiedDiffKrS = .FALSE. |
586 |
|
|
diffKrS = diffKrSDefault |
587 |
|
|
C- use temp diffusivity as default salt diffusivity: |
588 |
|
|
DO k=1,Nr |
589 |
|
|
diffKrNrS(k) = diffKrNrT(k) |
590 |
|
|
ENDDO |
591 |
|
|
ENDIF |
592 |
|
|
|
593 |
|
|
IF (diffKrBLEQsurf .EQ. UNSET_RL) diffKrBLEQsurf = diffKrBL79surf |
594 |
|
|
IF (diffKrBLEQdeep .EQ. UNSET_RL) diffKrBLEQdeep = diffKrBL79deep |
595 |
|
|
IF (diffKrBLEQscl .EQ. UNSET_RL) diffKrBLEQscl = diffKrBL79scl |
596 |
|
|
IF (diffKrBLEQHo .EQ. UNSET_RL) diffKrBLEQHo = diffKrBL79Ho |
597 |
|
|
|
598 |
|
|
IF ( hFacMinDz .NE. UNSET_RL ) zCoordInputData = .TRUE. |
599 |
|
|
IF ( hFacMinDp .NE. UNSET_RL ) pCoordInputData = .TRUE. |
600 |
|
|
IF ( hFacMinDr .NE. UNSET_RL ) rCoordInputData = .TRUE. |
601 |
|
|
IF ( hFacMinDr .EQ. UNSET_RL ) hFacMinDr = hFacMinDz |
602 |
|
|
IF ( hFacMinDr .EQ. UNSET_RL ) hFacMinDr = hFacMinDp |
603 |
|
|
IF ( hFacMinDr .EQ. UNSET_RL ) hFacMinDr = hFacMinDrDefault |
604 |
|
|
|
605 |
|
|
IF (convertFW2Salt.EQ.UNSET_RL) THEN |
606 |
|
|
convertFW2Salt = 35. |
607 |
|
|
IF (useRealFreshWaterFlux) convertFW2Salt=-1 |
608 |
|
|
IF ( selectAddFluid.GE.1 ) convertFW2Salt=-1 |
609 |
|
|
ENDIF |
610 |
|
|
|
611 |
|
|
IF ( SadournyCoriolis ) THEN |
612 |
|
|
C-- for backward compatibility : |
613 |
|
|
IF ( selectVortScheme.EQ.UNSET_I ) selectVortScheme = 2 |
614 |
|
|
IF ( selectVortScheme.NE.2 ) THEN |
615 |
|
|
WRITE(msgBuf,'(A,I5,A)') |
616 |
|
|
& 'S/R INI_PARMS: selectVortScheme=', selectVortScheme, |
617 |
|
|
& ' conflicts with "SadournyCoriolis"' |
618 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
619 |
|
|
errCount = errCount + 1 |
620 |
|
|
ENDIF |
621 |
|
|
ENDIF |
622 |
|
|
|
623 |
|
|
IF ( ivdc_kappa.NE.zeroRL .AND. .NOT.implicitDiffusion ) THEN |
624 |
|
|
WRITE(msgBuf,'(A,A)') |
625 |
|
|
& 'S/R INI_PARMS: To use ivdc_kappa you must enable implicit', |
626 |
|
|
& ' vertical diffusion.' |
627 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
628 |
|
|
errCount = errCount + 1 |
629 |
|
|
ENDIF |
630 |
|
|
|
631 |
|
|
coordsSet = 0 |
632 |
|
|
IF ( zCoordInputData ) coordsSet = coordsSet + 1 |
633 |
|
|
IF ( pCoordInputData ) coordsSet = coordsSet + 1 |
634 |
|
|
IF ( rCoordInputData ) coordsSet = coordsSet + 1 |
635 |
|
|
IF ( coordsSet .GT. 1 ) THEN |
636 |
|
|
WRITE(msgBuf,'(A)') |
637 |
|
|
& 'S/R INI_PARMS: Cannot mix z, p and r in the input data.' |
638 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
639 |
|
|
errCount = errCount + 1 |
640 |
|
|
ENDIF |
641 |
|
|
IF ( rhoConst .LE. 0. ) THEN |
642 |
|
|
WRITE(msgBuf,'(A)') |
643 |
|
|
& 'S/R INI_PARMS: rhoConst must be greater than 0.' |
644 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
645 |
|
|
errCount = errCount + 1 |
646 |
|
|
recip_rhoConst = 1. _d 0 |
647 |
|
|
ELSE |
648 |
|
|
recip_rhoConst = 1. _d 0 / rhoConst |
649 |
|
|
ENDIF |
650 |
|
|
IF ( eosType.EQ.'LINEAR' .AND. rhoNil.LE.0. ) THEN |
651 |
|
|
WRITE(msgBuf,'(A)') |
652 |
|
|
& 'S/R INI_PARMS: rhoNil must be greater than 0.' |
653 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
654 |
|
|
errCount = errCount + 1 |
655 |
|
|
ENDIF |
656 |
|
|
IF ( HeatCapacity_Cp .LE. 0. ) THEN |
657 |
|
|
WRITE(msgBuf,'(A)') |
658 |
|
|
& 'S/R INI_PARMS: HeatCapacity_Cp must be greater than 0.' |
659 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
660 |
|
|
errCount = errCount + 1 |
661 |
|
|
ENDIF |
662 |
|
|
IF ( gravity .LE. 0. ) THEN |
663 |
|
|
WRITE(msgBuf,'(A)') |
664 |
|
|
& 'S/R INI_PARMS: gravity must be greater than 0.' |
665 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
666 |
|
|
errCount = errCount + 1 |
667 |
|
|
recip_gravity = 1. _d 0 |
668 |
|
|
ELSE |
669 |
|
|
recip_gravity = 1. _d 0 / gravity |
670 |
|
|
ENDIF |
671 |
|
|
|
672 |
|
|
C- set default printResidualFreq according to debugLevel |
673 |
|
|
printResidualFreq = 0 |
674 |
|
|
IF ( debugLevel.GE.debLevE ) printResidualFreq = 1 |
675 |
|
|
|
676 |
|
|
C- set useSingleCpuInput=.TRUE. if useSingleCpuIO==.TRUE. |
677 |
|
|
IF ( useSingleCpuIO ) useSingleCpuInput=.TRUE. |
678 |
|
|
|
679 |
|
|
C Check for retired parameters still being used |
680 |
|
|
nRetired = 0 |
681 |
|
|
IF ( useConstantF ) THEN |
682 |
|
|
nRetired = nRetired+1 |
683 |
|
|
WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "useConstantF" ', |
684 |
|
|
& 'is no longer allowed in file "data"' |
685 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
686 |
|
|
WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: set "selectCoriMap"', |
687 |
|
|
& ' [0,1,2,3] to impose a setting over grid default' |
688 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
689 |
|
|
ENDIF |
690 |
|
|
IF ( useBetaPlaneF ) THEN |
691 |
|
|
nRetired = nRetired+1 |
692 |
|
|
WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "useBetaPlaneF" ', |
693 |
|
|
& 'is no longer allowed in file "data"' |
694 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
695 |
|
|
WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: set "selectCoriMap"', |
696 |
|
|
& ' [0,1,2,3] to impose a setting over grid default' |
697 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
698 |
|
|
ENDIF |
699 |
|
|
IF ( .NOT. useSphereF ) THEN |
700 |
|
|
nRetired = nRetired+1 |
701 |
|
|
WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "useSphereF" ', |
702 |
|
|
& 'is no longer allowed in file "data"' |
703 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
704 |
|
|
WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: set "selectCoriMap"', |
705 |
|
|
& ' [0,1,2,3] to impose a setting over grid default' |
706 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
707 |
|
|
ENDIF |
708 |
|
|
IF ( zonal_filt_lat .NE. UNSET_RL ) THEN |
709 |
|
|
nRetired = nRetired+1 |
710 |
|
|
WRITE(msgBuf,'(A,A)') |
711 |
|
|
& 'S/R INI_PARMS: Paramater "zonal_filt_lat" is', |
712 |
|
|
& ' no longer allowed in file "data".' |
713 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
714 |
|
|
WRITE(msgBuf,'(A,A)') |
715 |
|
|
& 'S/R INI_PARMS: Paramater "zonal_filt_lat" is', |
716 |
|
|
& ' now read from file "data.zonfilt".' |
717 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
718 |
|
|
ENDIF |
719 |
|
|
IF ( gravitySign .NE. UNSET_RL ) THEN |
720 |
|
|
nRetired = nRetired+1 |
721 |
|
|
WRITE(msgBuf,'(A,A)') |
722 |
|
|
& 'S/R INI_PARMS: "gravitySign" is set according to vertical ', |
723 |
|
|
& ' coordinate and is no longer allowed in file "data".' |
724 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
725 |
|
|
ENDIF |
726 |
|
|
IF ( tracerAdvScheme .NE. UNSET_I ) THEN |
727 |
|
|
nRetired = nRetired+1 |
728 |
|
|
WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "tracerAdvScheme" ', |
729 |
|
|
& '(old passive tracer code) is no longer allowed in file "data"' |
730 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
731 |
|
|
ENDIF |
732 |
|
|
IF ( trac_EvPrRn .NE. UNSET_RL ) THEN |
733 |
|
|
nRetired = nRetired+1 |
734 |
|
|
WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "trac_EvPrRn" ', |
735 |
|
|
& '(old passive tracer code) is no longer allowed in file "data"' |
736 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
737 |
|
|
ENDIF |
738 |
|
|
IF ( .NOT. tempDiffusion ) THEN |
739 |
|
|
nRetired = nRetired+1 |
740 |
|
|
WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "tempDiffusion" ', |
741 |
|
|
& 'is no longer allowed in file "data"' |
742 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
743 |
|
|
WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: to turn off diffusion', |
744 |
|
|
& ' => set diffusivity to zero' |
745 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
746 |
|
|
ENDIF |
747 |
|
|
IF ( .NOT. saltDiffusion ) THEN |
748 |
|
|
nRetired = nRetired+1 |
749 |
|
|
WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "saltDiffusion" ', |
750 |
|
|
& 'is no longer allowed in file "data"' |
751 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
752 |
|
|
WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: to turn off diffusion', |
753 |
|
|
& ' => set diffusivity to zero' |
754 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
755 |
|
|
ENDIF |
756 |
|
|
IF ( viscAstrain .NE. UNSET_RL ) THEN |
757 |
|
|
nRetired = nRetired+1 |
758 |
|
|
WRITE(msgBuf,'(A,A)') |
759 |
|
|
& 'S/R INI_PARMS: "viscAstrain" ', |
760 |
|
|
& 'is no longer allowed in file "data"' |
761 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
762 |
|
|
WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: to use Strain & Tension', |
763 |
|
|
& ' formulation => set useStrainTensionVisc to TRUE' |
764 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
765 |
|
|
ENDIF |
766 |
|
|
IF ( viscAtension .NE. UNSET_RL ) THEN |
767 |
|
|
nRetired = nRetired+1 |
768 |
|
|
WRITE(msgBuf,'(A,A)') |
769 |
|
|
& 'S/R INI_PARMS: "viscAtension" ', |
770 |
|
|
& 'is no longer allowed in file "data"' |
771 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
772 |
|
|
WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: to use Strain & Tension', |
773 |
|
|
& ' formulation => set useStrainTensionVisc to TRUE' |
774 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
775 |
|
|
ENDIF |
776 |
|
|
IF ( .NOT.useAnisotropicViscAgridMax ) THEN |
777 |
|
|
nRetired = nRetired+1 |
778 |
|
|
WRITE(msgBuf,'(A,A)') |
779 |
|
|
& 'S/R INI_PARMS: "useAnisotropicViscAgridMax" ', |
780 |
|
|
& 'is not allowed in "data" substitute useAreaViscLength=true' |
781 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
782 |
|
|
ENDIF |
783 |
|
|
IF ( usePickupBeforeC35 ) THEN |
784 |
|
|
nRetired = nRetired+1 |
785 |
|
|
WRITE(msgBuf,'(A,A)') |
786 |
|
|
& 'S/R INI_PARMS: "usePickupBeforeC35" ', |
787 |
|
|
& 'is no longer supported & not longer allowed in file "data"' |
788 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
789 |
|
|
ENDIF |
790 |
|
|
IF ( debugMode.NEQV.saveDebugMode ) THEN |
791 |
|
|
nRetired = nRetired+1 |
792 |
|
|
WRITE(msgBuf,'(A,A)') |
793 |
|
|
& 'S/R INI_PARMS: "debugMode" has been moved to "eedata"', |
794 |
|
|
& ' and is no longer allowed in file "data"' |
795 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
796 |
|
|
ENDIF |
797 |
|
|
IF ( allowInteriorFreezing ) THEN |
798 |
|
|
nRetired = nRetired+1 |
799 |
|
|
WRITE(msgBuf,'(A,A)') |
800 |
|
|
& 'S/R INI_PARMS: "allowInteriorFreezing" has been replaced', |
801 |
|
|
& ' by pkg/frazil and is no longer allowed in file "data"' |
802 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
803 |
|
|
ENDIF |
804 |
|
|
IF ( useOldFreezing ) THEN |
805 |
|
|
nRetired = nRetired+1 |
806 |
|
|
WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "useOldFreezing" ', |
807 |
|
|
& 'is no longer supported & not longer allowed in file "data"' |
808 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
809 |
|
|
ENDIF |
810 |
|
|
|
811 |
|
|
C-- Elliptic solver parameters |
812 |
|
|
WRITE(msgBuf,'(A)') ' INI_PARMS ; starts to read PARM02' |
813 |
|
|
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
814 |
|
|
& SQUEEZE_RIGHT, myThid ) |
815 |
|
|
READ(UNIT=iUnit,NML=PARM02) !,IOSTAT=errIO) |
816 |
|
|
IF ( errIO .LT. 0 ) THEN |
817 |
|
|
WRITE(msgBuf,'(A)') |
818 |
|
|
& 'S/R INI_PARMS: Error reading model parameter file "data"' |
819 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
820 |
|
|
WRITE(msgBuf,'(A)') 'S/R INI_PARMS: Problem in namelist PARM02' |
821 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
822 |
|
|
STOP 'ABNORMAL END: S/R INI_PARMS' |
823 |
|
|
ELSE |
824 |
|
|
WRITE(msgBuf,'(A)') ' INI_PARMS ; read PARM02 : OK' |
825 |
|
|
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
826 |
|
|
& SQUEEZE_RIGHT, myThid ) |
827 |
|
|
ENDIF |
828 |
|
|
|
829 |
|
|
C-- Time stepping parameters |
830 |
|
|
rCD = -1. _d 0 |
831 |
|
|
epsAB_CD = UNSET_RL |
832 |
|
|
latBandClimRelax = UNSET_RL |
833 |
|
|
deltaTtracer = 0. _d 0 |
834 |
|
|
forcing_In_AB = .TRUE. |
835 |
|
|
WRITE(msgBuf,'(A)') ' INI_PARMS ; starts to read PARM03' |
836 |
|
|
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
837 |
|
|
& SQUEEZE_RIGHT, myThid ) |
838 |
|
|
READ(UNIT=iUnit,NML=PARM03) !,IOSTAT=errIO) |
839 |
|
|
IF ( errIO .LT. 0 ) THEN |
840 |
|
|
WRITE(msgBuf,'(A)') |
841 |
|
|
& 'S/R INI_PARMS: Error reading model parameter file "data"' |
842 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
843 |
|
|
WRITE(msgBuf,'(A)') 'S/R INI_PARMS: Problem in namelist PARM03' |
844 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
845 |
|
|
STOP 'ABNORMAL END: S/R INI_PARMS' |
846 |
|
|
ELSE |
847 |
|
|
WRITE(msgBuf,'(A)') ' INI_PARMS ; read PARM03 : OK' |
848 |
|
|
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
849 |
|
|
& SQUEEZE_RIGHT, myThid ) |
850 |
|
|
ENDIF |
851 |
|
|
C Check for retired parameters still being used |
852 |
|
|
IF ( tauThetaClimRelax3Dim .NE. UNSET_RL ) THEN |
853 |
|
|
nRetired = nRetired+1 |
854 |
|
|
WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "tauThetaClimRelax3Dim" ', |
855 |
|
|
& 'is no longer allowed in file "data"' |
856 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
857 |
|
|
WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: 3-dim. relaxation code', |
858 |
|
|
& ' has moved to separate pkg/rbcs.' |
859 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
860 |
|
|
ENDIF |
861 |
|
|
IF ( tauSaltClimRelax3Dim .NE. UNSET_RL ) THEN |
862 |
|
|
nRetired = nRetired+1 |
863 |
|
|
WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "tauSaltClimRelax3Dim" ', |
864 |
|
|
& 'is no longer allowed in file "data"' |
865 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
866 |
|
|
WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: 3-dim. relaxation code', |
867 |
|
|
& ' has moved to separate pkg/rbcs.' |
868 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
869 |
|
|
ENDIF |
870 |
|
|
IF ( calendarDumps ) THEN |
871 |
|
|
nRetired = nRetired+1 |
872 |
|
|
WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "calendarDumps" ', |
873 |
|
|
& 'is no longer allowed in file "data"' |
874 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
875 |
|
|
WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: calendarDumps', |
876 |
|
|
& ' has moved to "data.cal"' |
877 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
878 |
|
|
ENDIF |
879 |
|
|
|
880 |
|
|
C Process "timestepping" params |
881 |
|
|
C o Time step size |
882 |
|
|
IF ( deltaTtracer .NE. dTtracerLev(1) .AND. |
883 |
|
|
& deltaTtracer .NE. 0. .AND. dTtracerLev(1) .NE. 0. ) THEN |
884 |
|
|
WRITE(msgBuf,'(A)') |
885 |
|
|
& 'S/R INI_PARMS: deltaTtracer & dTtracerLev(1) not equal' |
886 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
887 |
|
|
errCount = errCount + 1 |
888 |
|
|
ELSEIF ( dTtracerLev(1) .NE. 0. ) THEN |
889 |
|
|
deltaTtracer = dTtracerLev(1) |
890 |
|
|
ENDIF |
891 |
|
|
IF ( deltaT .EQ. 0. ) deltaT = deltaTClock |
892 |
|
|
IF ( deltaT .EQ. 0. ) deltaT = deltaTtracer |
893 |
|
|
IF ( deltaT .EQ. 0. ) deltaT = deltaTMom |
894 |
|
|
IF ( deltaT .EQ. 0. ) deltaT = deltaTFreeSurf |
895 |
|
|
IF ( deltaT .EQ. 0. ) THEN |
896 |
|
|
WRITE(msgBuf,'(2A)') 'S/R INI_PARMS: ', |
897 |
|
|
& 'need to specify in file "data", namelist "PARM03"' |
898 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
899 |
|
|
WRITE(msgBuf,'(2A)') 'S/R INI_PARMS: ', |
900 |
|
|
& ' a model timestep (in s) deltaT or deltaTClock= ?' |
901 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
902 |
|
|
errCount = errCount + 1 |
903 |
|
|
deltaT = 1. |
904 |
|
|
ENDIF |
905 |
|
|
IF ( deltaTMom .EQ. 0. ) deltaTMom = deltaT |
906 |
|
|
IF ( deltaTtracer .EQ. 0. ) deltaTtracer = deltaT |
907 |
|
|
IF ( deltaTClock .EQ. 0. ) deltaTClock = deltaT |
908 |
|
|
DO k=1,Nr |
909 |
|
|
IF (dTtracerLev(k).EQ.0.) dTtracerLev(k)= deltaTtracer |
910 |
|
|
ENDDO |
911 |
|
|
C Note that this line should set deltaFreesurf=deltaTtracer |
912 |
|
|
C but this would change a lot of existing set-ups so we are |
913 |
|
|
C obliged to set the default inappropriately. |
914 |
|
|
C Be advised that when using asynchronous time stepping |
915 |
|
|
C it is better to set deltaTreesurf=deltaTtracer |
916 |
|
|
IF ( deltaTFreeSurf .EQ. 0. ) deltaTFreeSurf = deltaTMom |
917 |
|
|
IF ( periodicExternalForcing ) THEN |
918 |
|
|
IF ( externForcingCycle*externForcingPeriod .EQ. 0. ) THEN |
919 |
|
|
WRITE(msgBuf,'(A)') |
920 |
|
|
& 'S/R INI_PARMS: externForcingCycle,externForcingPeriod =0' |
921 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
922 |
|
|
errCount = errCount + 1 |
923 |
|
|
ELSEIF ( INT(externForcingCycle/externForcingPeriod) .NE. |
924 |
|
|
& externForcingCycle/externForcingPeriod ) THEN |
925 |
|
|
WRITE(msgBuf,'(A)') |
926 |
|
|
& 'S/R INI_PARMS: externForcingCycle <> N*externForcingPeriod' |
927 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
928 |
|
|
errCount = errCount + 1 |
929 |
|
|
ELSEIF ( externForcingCycle.LT.externForcingPeriod ) THEN |
930 |
|
|
WRITE(msgBuf,'(A)') |
931 |
|
|
& 'S/R INI_PARMS: externForcingCycle < externForcingPeriod' |
932 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
933 |
|
|
errCount = errCount + 1 |
934 |
|
|
ENDIF |
935 |
|
|
IF ( externForcingPeriod.LT.deltaTClock ) THEN |
936 |
|
|
WRITE(msgBuf,'(A)') |
937 |
|
|
& 'S/R INI_PARMS: externForcingPeriod < deltaTClock' |
938 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
939 |
|
|
errCount = errCount + 1 |
940 |
|
|
ENDIF |
941 |
|
|
ENDIF |
942 |
|
|
C o Adams-Bashforth time stepping: |
943 |
|
|
IF ( momForcingOutAB .EQ. UNSET_I ) THEN |
944 |
|
|
momForcingOutAB = 1 |
945 |
|
|
IF ( forcing_In_AB ) momForcingOutAB = 0 |
946 |
|
|
ENDIF |
947 |
|
|
IF ( tracForcingOutAB .EQ. UNSET_I ) THEN |
948 |
|
|
tracForcingOutAB = 1 |
949 |
|
|
IF ( forcing_In_AB ) tracForcingOutAB = 0 |
950 |
|
|
ENDIF |
951 |
|
|
C o Convection frequency |
952 |
|
|
IF ( cAdjFreq .LT. 0. ) THEN |
953 |
|
|
cAdjFreq = deltaTClock |
954 |
|
|
ENDIF |
955 |
|
|
IF ( ivdc_kappa .NE. 0. .AND. cAdjFreq .NE. 0. ) THEN |
956 |
|
|
WRITE(msgBuf,'(A,A)') |
957 |
|
|
& 'S/R INI_PARMS: You have enabled both ivdc_kappa and', |
958 |
|
|
& ' convective adjustment.' |
959 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
960 |
|
|
errCount = errCount + 1 |
961 |
|
|
ENDIF |
962 |
|
|
IF (useCDscheme) THEN |
963 |
|
|
C o CD coupling (CD scheme): |
964 |
|
|
IF ( tauCD .EQ. 0. _d 0 ) tauCD = deltaTMom |
965 |
|
|
IF ( rCD .LT. 0. ) rCD = 1. _d 0 - deltaTMom/tauCD |
966 |
|
|
IF ( epsAB_CD .EQ. UNSET_RL ) epsAB_CD = abEps |
967 |
|
|
ENDIF |
968 |
|
|
|
969 |
|
|
IF ( startTime.EQ.UNSET_RL .AND. nIter0.EQ.-1 ) THEN |
970 |
|
|
C o set default value for start time & nIter0 |
971 |
|
|
startTime = baseTime |
972 |
|
|
nIter0 = 0 |
973 |
|
|
ELSEIF ( startTime.EQ.UNSET_RL ) THEN |
974 |
|
|
C o set start time from nIter0 |
975 |
|
|
startTime = baseTime + deltaTClock*DFLOAT(nIter0) |
976 |
|
|
ELSEIF ( nIter0.EQ.-1 ) THEN |
977 |
|
|
C o set nIter0 from start time |
978 |
|
|
nIter0 = NINT( (startTime-baseTime)/deltaTClock ) |
979 |
|
|
ELSEIF ( baseTime.EQ.0. ) THEN |
980 |
|
|
C o set base time from the 2 others |
981 |
|
|
baseTime = startTime - deltaTClock*DFLOAT(nIter0) |
982 |
|
|
ENDIF |
983 |
|
|
|
984 |
|
|
nTimeSteps_l2 = 4 |
985 |
|
|
C o nTimeSteps 1 |
986 |
|
|
IF ( nTimeSteps .EQ. 0 .AND. nEndIter .NE. 0 ) |
987 |
|
|
& nTimeSteps = nEndIter-nIter0 |
988 |
|
|
C o nTimeSteps 2 |
989 |
|
|
IF ( nTimeSteps .EQ. 0 .AND. endTime .NE. 0. ) |
990 |
|
|
& nTimeSteps = NINT((endTime-startTime)/deltaTClock) |
991 |
|
|
C o nEndIter 1 |
992 |
|
|
IF ( nEndIter .EQ. 0 .AND. nTimeSteps .NE. 0 ) |
993 |
|
|
& nEndIter = nIter0+nTimeSteps |
994 |
|
|
C o nEndIter 2 |
995 |
|
|
IF ( nEndIter .EQ. 0 .AND. endTime .NE. 0. ) |
996 |
|
|
& nEndIter = NINT((endTime-baseTime)/deltaTClock) |
997 |
|
|
C o End Time 1 |
998 |
|
|
IF ( endTime .EQ. 0. .AND. nTimeSteps .NE. 0 ) |
999 |
|
|
& endTime = startTime + deltaTClock*DFLOAT(nTimeSteps) |
1000 |
|
|
C o End Time 2 |
1001 |
|
|
IF ( endTime .EQ. 0. .AND. nEndIter .NE. 0 ) |
1002 |
|
|
& endTime = baseTime + deltaTClock*DFLOAT(nEndIter) |
1003 |
|
|
|
1004 |
|
|
C o Consistent? |
1005 |
|
|
IF ( startTime .NE. baseTime+deltaTClock*DFLOAT(nIter0) ) THEN |
1006 |
|
|
WRITE(msgBuf,'(A)') |
1007 |
|
|
& 'S/R INI_PARMS: startTime, baseTime and nIter0 are inconsistent' |
1008 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1009 |
|
|
WRITE(msgBuf,'(A)') |
1010 |
|
|
& 'S/R INI_PARMS: Perhaps more than two were set at once' |
1011 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1012 |
|
|
errCount = errCount + 1 |
1013 |
|
|
ENDIF |
1014 |
|
|
IF ( nEndIter .NE. nIter0+nTimeSteps ) THEN |
1015 |
|
|
WRITE(msgBuf,'(A)') |
1016 |
|
|
& 'S/R INI_PARMS: nIter0, nTimeSteps and nEndIter are inconsistent' |
1017 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1018 |
|
|
WRITE(msgBuf,'(A)') |
1019 |
|
|
& 'S/R INI_PARMS: Perhaps more than two were set at once' |
1020 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1021 |
|
|
errCount = errCount + 1 |
1022 |
|
|
ENDIF |
1023 |
|
|
IF ( nTimeSteps .NE. NINT((endTime-startTime)/deltaTClock) |
1024 |
|
|
& ) THEN |
1025 |
|
|
WRITE(msgBuf,'(A)') |
1026 |
|
|
& 'S/R INI_PARMS: both endTime and nTimeSteps have been set' |
1027 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1028 |
|
|
WRITE(msgBuf,'(A)') |
1029 |
|
|
& 'S/R INI_PARMS: but are inconsistent' |
1030 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1031 |
|
|
errCount = errCount + 1 |
1032 |
|
|
ENDIF |
1033 |
|
|
|
1034 |
|
|
C o Monitor (should also add CPP flag for monitor?) |
1035 |
|
|
IF (monitorFreq.LT.0.) THEN |
1036 |
|
|
monitorFreq=0. |
1037 |
|
|
IF (dumpFreq.NE.0.) monitorFreq=dumpFreq |
1038 |
|
|
IF (diagFreq.NE.0..AND.diagFreq.LT.monitorFreq) |
1039 |
|
|
& monitorFreq=diagFreq |
1040 |
|
|
IF (taveFreq.NE.0..AND.taveFreq.LT.monitorFreq) |
1041 |
|
|
& monitorFreq=taveFreq |
1042 |
|
|
IF (chkPtFreq.NE.0..AND.chkPtFreq.LT.monitorFreq) |
1043 |
|
|
& monitorFreq=chkPtFreq |
1044 |
|
|
IF (pChkPtFreq.NE.0..AND.pChkPtFreq.LT.monitorFreq) |
1045 |
|
|
& monitorFreq=pChkPtFreq |
1046 |
|
|
IF (monitorFreq.EQ.0.) monitorFreq=deltaTClock |
1047 |
|
|
ENDIF |
1048 |
|
|
IF ( monitorSelect.EQ.UNSET_I ) THEN |
1049 |
|
|
monitorSelect = 2 |
1050 |
|
|
IF ( fluidIsWater ) monitorSelect = 3 |
1051 |
|
|
ENDIF |
1052 |
|
|
|
1053 |
|
|
C-- Grid parameters |
1054 |
|
|
C In cartesian coords distances are in metres |
1055 |
|
|
DO k =1,Nr |
1056 |
|
|
delZ(k) = UNSET_RL |
1057 |
|
|
delP(k) = UNSET_RL |
1058 |
|
|
delR(k) = UNSET_RL |
1059 |
|
|
ENDDO |
1060 |
|
|
C In spherical polar distances are in degrees |
1061 |
|
|
dxSpacing = UNSET_RL |
1062 |
|
|
dySpacing = UNSET_RL |
1063 |
|
|
WRITE(msgBuf,'(A)') ' INI_PARMS ; starts to read PARM04' |
1064 |
|
|
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
1065 |
|
|
& SQUEEZE_RIGHT, myThid ) |
1066 |
|
|
READ(UNIT=iUnit,NML=PARM04,IOSTAT=errIO) |
1067 |
|
|
IF ( errIO .LT. 0 ) THEN |
1068 |
|
|
WRITE(msgBuf,'(A)') |
1069 |
|
|
& 'S/R INI_PARMS: Error reading model parameter file "data"' |
1070 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1071 |
|
|
WRITE(msgBuf,'(A)') 'S/R INI_PARMS: Problem in namelist PARM04' |
1072 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1073 |
|
|
STOP 'ABNORMAL END: S/R INI_PARMS' |
1074 |
|
|
ELSE |
1075 |
|
|
WRITE(msgBuf,'(A)') ' INI_PARMS ; read PARM04 : OK' |
1076 |
|
|
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
1077 |
|
|
& SQUEEZE_RIGHT, myThid ) |
1078 |
|
|
ENDIF |
1079 |
|
|
|
1080 |
|
|
C Check for retired parameters still being used |
1081 |
|
|
IF ( Ro_SeaLevel .NE. UNSET_RL ) THEN |
1082 |
|
|
c nRetired = nRetired+1 |
1083 |
|
|
IF ( usingPCoords ) THEN |
1084 |
|
|
WRITE(msgBuf,'(2A)') '** WARNING ** INI_PARMS: ', |
1085 |
|
|
& '"Ro_SeaLevel" (P @ bottom) depreciated (backward compat' |
1086 |
|
|
ELSEIF ( usingZCoords ) THEN |
1087 |
|
|
WRITE(msgBuf,'(2A)') '** WARNING ** INI_PARMS: ', |
1088 |
|
|
& '"Ro_SeaLevel" (Z @ top) depreciated (backward compat' |
1089 |
|
|
ENDIF |
1090 |
|
|
CALL PRINT_MESSAGE( msgBuf, errorMessageUnit, |
1091 |
|
|
& SQUEEZE_RIGHT, myThid ) |
1092 |
|
|
IF ( usingPCoords ) THEN |
1093 |
|
|
WRITE(msgBuf,'(2A)') '** WARNING ** INI_PARMS: ', |
1094 |
|
|
& ' only). To set vert. axis, use instead "top_Pres".' |
1095 |
|
|
ELSEIF ( usingZCoords ) THEN |
1096 |
|
|
WRITE(msgBuf,'(2A)') '** WARNING ** INI_PARMS: ', |
1097 |
|
|
& ' only). To set vert. axis, use instead "seaLev_Z".' |
1098 |
|
|
ENDIF |
1099 |
|
|
CALL PRINT_MESSAGE( msgBuf, errorMessageUnit, |
1100 |
|
|
& SQUEEZE_RIGHT, myThid ) |
1101 |
|
|
ENDIF |
1102 |
|
|
IF ( rkFac .NE. UNSET_RL ) THEN |
1103 |
|
|
nRetired = nRetired+1 |
1104 |
|
|
WRITE(msgBuf,'(A,A)') |
1105 |
|
|
& 'S/R INI_PARMS: "rkFac" has been replaced by -rkSign ', |
1106 |
|
|
& ' and is no longer allowed in file "data".' |
1107 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1108 |
|
|
ENDIF |
1109 |
|
|
IF ( groundAtK1 ) THEN |
1110 |
|
|
c nRetired = nRetired+1 |
1111 |
|
|
WRITE(msgBuf,'(A,A)') |
1112 |
|
|
& 'S/R INI_PARMS: "groundAtK1" is set according to vertical ', |
1113 |
|
|
& ' coordinate and is no longer allowed in file "data".' |
1114 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1115 |
|
|
ENDIF |
1116 |
|
|
IF ( thetaMin .NE. UNSET_RL ) THEN |
1117 |
|
|
nRetired = nRetired+1 |
1118 |
|
|
WRITE(msgBuf,'(A,A)') |
1119 |
|
|
& 'S/R INI_PARMS: "thetaMin" no longer allowed,', |
1120 |
|
|
& ' has been replaced by "xgOrigin"' |
1121 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1122 |
|
|
ENDIF |
1123 |
|
|
IF ( phiMin .NE. UNSET_RL ) THEN |
1124 |
|
|
nRetired = nRetired+1 |
1125 |
|
|
WRITE(msgBuf,'(A,A)') |
1126 |
|
|
& 'S/R INI_PARMS: "phiMin" no longer allowed,', |
1127 |
|
|
& ' has been replaced by "ygOrigin"' |
1128 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1129 |
|
|
ENDIF |
1130 |
|
|
|
1131 |
|
|
C X coordinate : Check for multiple definitions |
1132 |
|
|
goptCount = 0 |
1133 |
|
|
IF ( delX(1) .NE. UNSET_RL ) goptCount = goptCount + 1 |
1134 |
|
|
IF ( dxSpacing .NE. UNSET_RL ) goptCount = goptCount + 1 |
1135 |
|
|
IF ( delXFile .NE. ' ' ) goptCount = goptCount + 1 |
1136 |
|
|
IF ( goptCount.GT.1 ) THEN |
1137 |
|
|
WRITE(msgBuf,'(A,A)') 'Too many specifications for delX:', |
1138 |
|
|
& 'Specify only one of delX, dxSpacing or delXfile' |
1139 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1140 |
|
|
errCount = errCount + 1 |
1141 |
|
|
ENDIF |
1142 |
|
|
IF ( dxSpacing .NE. UNSET_RL ) THEN |
1143 |
|
|
DO i=1,gridNx |
1144 |
|
|
delX(i) = dxSpacing |
1145 |
|
|
ENDDO |
1146 |
|
|
ENDIF |
1147 |
|
|
C Y coordinate : Check for multiple definitions |
1148 |
|
|
goptCount = 0 |
1149 |
|
|
IF ( delY(1) .NE. UNSET_RL ) goptCount = goptCount + 1 |
1150 |
|
|
IF ( dySpacing .NE. UNSET_RL ) goptCount = goptCount + 1 |
1151 |
|
|
IF ( delYFile .NE. ' ' ) goptCount = goptCount + 1 |
1152 |
|
|
IF ( goptCount.GT.1 ) THEN |
1153 |
|
|
WRITE(msgBuf,'(A,A)') 'Too many specifications for delY:', |
1154 |
|
|
& 'Specify only one of delY, dySpacing or delYfile' |
1155 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1156 |
|
|
errCount = errCount + 1 |
1157 |
|
|
ENDIF |
1158 |
|
|
IF ( dySpacing .NE. UNSET_RL ) THEN |
1159 |
|
|
DO j=1,gridNy |
1160 |
|
|
delY(j) = dySpacing |
1161 |
|
|
ENDDO |
1162 |
|
|
ENDIF |
1163 |
|
|
|
1164 |
|
|
C-- Check for conflicting grid definitions. |
1165 |
|
|
goptCount = 0 |
1166 |
|
|
IF ( usingCartesianGrid ) goptCount = goptCount+1 |
1167 |
|
|
IF ( usingSphericalPolarGrid ) goptCount = goptCount+1 |
1168 |
|
|
IF ( usingCurvilinearGrid ) goptCount = goptCount+1 |
1169 |
|
|
IF ( usingCylindricalGrid ) goptCount = goptCount+1 |
1170 |
|
|
IF ( goptCount .GT. 1 ) THEN |
1171 |
|
|
WRITE(msgBuf,'(A)') |
1172 |
|
|
& 'S/R INI_PARMS: More than one coordinate system requested' |
1173 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1174 |
|
|
errCount = errCount + 1 |
1175 |
|
|
ENDIF |
1176 |
|
|
IF ( goptCount .LT. 1 ) THEN |
1177 |
|
|
C- No horizontal grid is specified => use Cartesian grid as default: |
1178 |
|
|
WRITE(msgBuf,'(A)') |
1179 |
|
|
& 'S/R INI_PARMS: No horizontal grid requested' |
1180 |
|
|
CALL PRINT_MESSAGE( msgBuf, errorMessageUnit, |
1181 |
|
|
& SQUEEZE_RIGHT, myThid ) |
1182 |
|
|
WRITE(msgBuf,'(A)') |
1183 |
|
|
& 'S/R INI_PARMS: => Use Cartesian Grid as default' |
1184 |
|
|
CALL PRINT_MESSAGE( msgBuf, errorMessageUnit, |
1185 |
|
|
& SQUEEZE_RIGHT, myThid ) |
1186 |
|
|
usingCartesianGrid = .TRUE. |
1187 |
|
|
ENDIF |
1188 |
|
|
C- Radius of the Planet: |
1189 |
|
|
IF ( rSphere.EQ.UNSET_RL ) THEN |
1190 |
|
|
IF ( usingCurvilinearGrid .AND. |
1191 |
|
|
& radius_fromHorizGrid.NE.UNSET_RL ) THEN |
1192 |
|
|
rSphere = radius_fromHorizGrid |
1193 |
|
|
ELSE |
1194 |
|
|
rSphere = 6370. _d 3 |
1195 |
|
|
ENDIF |
1196 |
|
|
ENDIF |
1197 |
|
|
IF ( radius_fromHorizGrid.EQ.UNSET_RL ) THEN |
1198 |
|
|
radius_fromHorizGrid = rSphere |
1199 |
|
|
ENDIF |
1200 |
|
|
IF ( rSphere .NE. 0. ) THEN |
1201 |
|
|
recip_rSphere = 1. _d 0/rSphere |
1202 |
|
|
ELSE |
1203 |
|
|
recip_rSphere = 0. |
1204 |
|
|
ENDIF |
1205 |
|
|
C-- Default vertical axis origin |
1206 |
|
|
IF ( Ro_SeaLevel .NE. UNSET_RL ) THEN |
1207 |
|
|
IF ( usingPCoords .AND. top_Pres.NE.UNSET_RL ) THEN |
1208 |
|
|
WRITE(msgBuf,'(2A)') 'S/R INI_PARMS: ', |
1209 |
|
|
& 'Cannot set both "Ro_SeaLevel" and "top_Pres"' |
1210 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1211 |
|
|
errCount = errCount + 1 |
1212 |
|
|
ENDIF |
1213 |
|
|
IF ( usingZCoords .AND. seaLev_Z.NE.UNSET_RL ) THEN |
1214 |
|
|
WRITE(msgBuf,'(2A)') 'S/R INI_PARMS: ', |
1215 |
|
|
& 'Cannot set both "Ro_SeaLevel" and "seaLev_Z"' |
1216 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1217 |
|
|
errCount = errCount + 1 |
1218 |
|
|
ENDIF |
1219 |
|
|
rF(1) = Ro_SeaLevel |
1220 |
|
|
ELSE |
1221 |
|
|
rF(1) = UNSET_RS |
1222 |
|
|
ENDIF |
1223 |
|
|
IF ( top_Pres.EQ.UNSET_RL ) top_Pres = 0. |
1224 |
|
|
IF ( seaLev_Z.EQ.UNSET_RL ) seaLev_Z = 0. |
1225 |
|
|
C-- Default origin for X & Y axis (longitude & latitude origin): |
1226 |
|
|
IF ( xgOrigin .EQ. UNSET_RL ) xgOrigin = 0. |
1227 |
|
|
IF ( ygOrigin .EQ. UNSET_RL ) THEN |
1228 |
|
|
IF ( usingSphericalPolarGrid ) THEN |
1229 |
|
|
ygOrigin = 0. |
1230 |
|
|
ELSEIF ( usingCartesianGrid ) THEN |
1231 |
|
|
ygOrigin = 0. |
1232 |
|
|
ELSEIF ( usingCylindricalGrid ) THEN |
1233 |
|
|
ygOrigin = 0. |
1234 |
|
|
ELSEIF ( usingCurvilinearGrid ) THEN |
1235 |
|
|
ygOrigin = 0. |
1236 |
|
|
ELSE |
1237 |
|
|
WRITE(msgBuf,'(A)') |
1238 |
|
|
& 'S/R INI_PARMS: found no coordinate system to set ygOrigin' |
1239 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1240 |
|
|
errCount = errCount + 1 |
1241 |
|
|
ENDIF |
1242 |
|
|
ENDIF |
1243 |
|
|
C-- Make metric term & Coriolis settings consistent with underlying grid. |
1244 |
|
|
IF ( usingCartesianGrid ) THEN |
1245 |
|
|
metricTerms = .FALSE. |
1246 |
|
|
useNHMTerms = .FALSE. |
1247 |
|
|
ENDIF |
1248 |
|
|
IF ( usingCylindricalGrid ) THEN |
1249 |
|
|
useNHMTerms = .FALSE. |
1250 |
|
|
WRITE(msgBuf,'(A)') 'S/R INI_PARMS ; Cylinder OK' |
1251 |
|
|
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
1252 |
|
|
& SQUEEZE_RIGHT, myThid ) |
1253 |
|
|
ENDIF |
1254 |
|
|
IF ( usingCurvilinearGrid ) THEN |
1255 |
|
|
metricTerms = .FALSE. |
1256 |
|
|
ENDIF |
1257 |
|
|
IF ( selectCoriMap.EQ.-1 ) THEN |
1258 |
|
|
IF ( usingCartesianGrid.OR.usingCylindricalGrid ) THEN |
1259 |
|
|
C default is to use Beta-Plane Coriolis |
1260 |
|
|
selectCoriMap = 1 |
1261 |
|
|
ELSE |
1262 |
|
|
C default for other grids is to use Spherical Coriolis |
1263 |
|
|
selectCoriMap = 2 |
1264 |
|
|
ENDIF |
1265 |
|
|
ENDIF |
1266 |
|
|
IF ( .NOT.(nonHydrostatic.OR.quasiHydrostatic) ) |
1267 |
|
|
& use3dCoriolis = .FALSE. |
1268 |
|
|
IF ( (selectCoriMap.EQ.0 .OR.selectCoriMap.EQ.1) |
1269 |
|
|
& .AND. fPrime.EQ.0. ) use3dCoriolis = .FALSE. |
1270 |
|
|
|
1271 |
|
|
C-- Grid rotation |
1272 |
|
|
IF ( phiEuler .NE. 0. _d 0 .OR. thetaEuler .NE. 0. _d 0 |
1273 |
|
|
& .OR. psiEuler .NE. 0. _d 0 ) rotateGrid = .TRUE. |
1274 |
|
|
|
1275 |
|
|
C-- Set default for latitude-band where relaxation to climatology applies |
1276 |
|
|
C note: done later (once domain size is known) if using CartesianGrid |
1277 |
|
|
IF ( latBandClimRelax .EQ. UNSET_RL) THEN |
1278 |
|
|
IF ( usingSphericalPolarGrid ) latBandClimRelax= 180. _d 0 |
1279 |
|
|
IF ( usingCurvilinearGrid ) latBandClimRelax= 180. _d 0 |
1280 |
|
|
ENDIF |
1281 |
|
|
|
1282 |
|
|
C-- set cell Center depth and put Interface at the middle between 2 C |
1283 |
|
|
setCenterDr = .FALSE. |
1284 |
|
|
DO k=1,Nr+1 |
1285 |
|
|
IF ( delRc(k).EQ.UNSET_RL ) THEN |
1286 |
|
|
IF ( setCenterDr ) THEN |
1287 |
|
|
WRITE(msgBuf,'(A,I4)') |
1288 |
|
|
& 'S/R INI_PARMS: No value for delRc at k =', k |
1289 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1290 |
|
|
errCount = errCount + 1 |
1291 |
|
|
ENDIF |
1292 |
|
|
ELSE |
1293 |
|
|
IF ( k.EQ.1 ) setCenterDr = .TRUE. |
1294 |
|
|
IF ( .NOT.setCenterDr ) THEN |
1295 |
|
|
WRITE(msgBuf,'(A,I4)') |
1296 |
|
|
& 'S/R INI_PARMS: No value for delRc at k <', k |
1297 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1298 |
|
|
errCount = errCount + 1 |
1299 |
|
|
ENDIF |
1300 |
|
|
ENDIF |
1301 |
|
|
ENDDO |
1302 |
|
|
IF ( setCenterDr ) rCoordInputData = .TRUE. |
1303 |
|
|
C-- p, z, r coord parameters |
1304 |
|
|
setInterFDr = .FALSE. |
1305 |
|
|
DO k = 1, Nr |
1306 |
|
|
IF ( delZ(k) .NE. UNSET_RL ) zCoordInputData = .TRUE. |
1307 |
|
|
IF ( delP(k) .NE. UNSET_RL ) pCoordInputData = .TRUE. |
1308 |
|
|
IF ( delR(k) .NE. UNSET_RL ) rCoordInputData = .TRUE. |
1309 |
|
|
IF ( delR(k) .EQ. UNSET_RL ) delR(k) = delZ(k) |
1310 |
|
|
IF ( delR(k) .EQ. UNSET_RL ) delR(k) = delP(k) |
1311 |
|
|
IF ( delR(k) .EQ. UNSET_RL ) THEN |
1312 |
|
|
IF ( setInterFDr ) THEN |
1313 |
|
|
WRITE(msgBuf,'(A,I4)') |
1314 |
|
|
& 'S/R INI_PARMS: No value for delZ/delP/delR at k =', k |
1315 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1316 |
|
|
errCount = errCount + 1 |
1317 |
|
|
ENDIF |
1318 |
|
|
ELSE |
1319 |
|
|
IF ( k.EQ.1 ) setInterFDr = .TRUE. |
1320 |
|
|
IF ( .NOT.setInterFDr ) THEN |
1321 |
|
|
WRITE(msgBuf,'(A,I4)') |
1322 |
|
|
& 'S/R INI_PARMS: No value for delZ/delP/delR at k <', k |
1323 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1324 |
|
|
errCount = errCount + 1 |
1325 |
|
|
ENDIF |
1326 |
|
|
ENDIF |
1327 |
|
|
ENDDO |
1328 |
|
|
C Check for multiple coordinate systems |
1329 |
|
|
coordsSet = 0 |
1330 |
|
|
IF ( zCoordInputData ) coordsSet = coordsSet + 1 |
1331 |
|
|
IF ( pCoordInputData ) coordsSet = coordsSet + 1 |
1332 |
|
|
IF ( rCoordInputData ) coordsSet = coordsSet + 1 |
1333 |
|
|
IF ( coordsSet .GT. 1 ) THEN |
1334 |
|
|
WRITE(msgBuf,'(A)') |
1335 |
|
|
& 'S/R INI_PARMS: Cannot mix z, p and r in the input data.' |
1336 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1337 |
|
|
errCount = errCount + 1 |
1338 |
|
|
ENDIF |
1339 |
|
|
C- Check for double definition (file & namelist) |
1340 |
|
|
IF ( delRcFile.NE.' ' ) THEN |
1341 |
|
|
IF ( setCenterDr ) THEN |
1342 |
|
|
WRITE(msgBuf,'(A)') |
1343 |
|
|
& 'S/R INI_PARMS: Cannot set both delRc and delRcFile' |
1344 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1345 |
|
|
errCount = errCount + 1 |
1346 |
|
|
ENDIF |
1347 |
|
|
setCenterDr = .TRUE. |
1348 |
|
|
ENDIF |
1349 |
|
|
IF ( delRFile.NE.' ' ) THEN |
1350 |
|
|
IF ( setInterFDr ) THEN |
1351 |
|
|
WRITE(msgBuf,'(A)') |
1352 |
|
|
& 'S/R INI_PARMS: Cannot set both delR and delRFile' |
1353 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1354 |
|
|
errCount = errCount + 1 |
1355 |
|
|
ENDIF |
1356 |
|
|
setInterFDr = .TRUE. |
1357 |
|
|
ENDIF |
1358 |
|
|
|
1359 |
|
|
C-- Input files |
1360 |
|
|
WRITE(msgBuf,'(A)') ' INI_PARMS ; starts to read PARM05' |
1361 |
|
|
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
1362 |
|
|
& SQUEEZE_RIGHT, myThid ) |
1363 |
|
|
READ(UNIT=iUnit,NML=PARM05) !,IOSTAT=errIO) |
1364 |
|
|
IF ( errIO .LT. 0 ) THEN |
1365 |
|
|
WRITE(msgBuf,'(A)') |
1366 |
|
|
& 'S/R INI_PARMS: Error reading model parameter file "data"' |
1367 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1368 |
|
|
WRITE(msgBuf,'(A)') 'S/R INI_PARMS: Problem in namelist PARM05' |
1369 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1370 |
|
|
STOP 'ABNORMAL END: S/R INI_PARMS' |
1371 |
|
|
ELSE |
1372 |
|
|
WRITE(msgBuf,'(A)') ' INI_PARMS ; read PARM05 : OK' |
1373 |
|
|
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
1374 |
|
|
& SQUEEZE_RIGHT, myThid ) |
1375 |
|
|
ENDIF |
1376 |
|
|
C Check for retired parameters still being used |
1377 |
|
|
IF ( shelfIceFile .NE. ' ' ) THEN |
1378 |
|
|
nRetired = nRetired+1 |
1379 |
|
|
WRITE(msgBuf,'(A,A)') |
1380 |
|
|
& 'S/R INI_PARMS: "shelfIceFile" is not allowed in "data", ', |
1381 |
|
|
& 'substitute "SHELFICEtopoFile" in data.shelfice' |
1382 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1383 |
|
|
ENDIF |
1384 |
|
|
IF ( dQdTFile .NE. ' ' ) THEN |
1385 |
|
|
nRetired = nRetired+1 |
1386 |
|
|
WRITE(msgBuf,'(A,A)') |
1387 |
|
|
& 'S/R INI_PARMS: "dQdTFile" has been retired from file "data"' |
1388 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1389 |
|
|
ENDIF |
1390 |
|
|
C Check if two conflicting I/O directories have been specified |
1391 |
|
|
IF (mdsioLocalDir .NE. ' ' .AND. adTapeDir .NE. ' ') THEN |
1392 |
|
|
WRITE(msgBuf,'(A)') |
1393 |
|
|
& 'S/R INI_PARMS: mdsioLocalDir and adTapeDir cannot be' |
1394 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1395 |
|
|
WRITE(msgBuf,'(A)') |
1396 |
|
|
& 'S/R INI_PARMS: specified at the same time' |
1397 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1398 |
|
|
errCount = errCount + 1 |
1399 |
|
|
ENDIF |
1400 |
|
|
|
1401 |
|
|
C Check for relaxation term: |
1402 |
|
|
IF ( (tauThetaClimRelax.GT.0.).AND. |
1403 |
|
|
& (thetaClimFile.EQ.' ') ) THEN |
1404 |
|
|
WRITE(msgBuf,'(A)') |
1405 |
|
|
& 'S/R INI_PARMS: tauThetaClimRelax > 0 but' |
1406 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1407 |
|
|
WRITE(msgBuf,'(A)') |
1408 |
|
|
& 'S/R INI_PARMS: thetaClimFile is undefined' |
1409 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1410 |
|
|
errCount = errCount + 1 |
1411 |
|
|
ENDIF |
1412 |
|
|
IF ( (tauSaltClimRelax.GT.0.).AND. |
1413 |
|
|
& (saltClimFile.EQ.' ') ) THEN |
1414 |
|
|
WRITE(msgBuf,'(A)') |
1415 |
|
|
& 'S/R INI_PARMS: tauSaltClimRelax > 0 but' |
1416 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1417 |
|
|
WRITE(msgBuf,'(A)') |
1418 |
|
|
& 'S/R INI_PARMS: saltClimFile is undefined' |
1419 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1420 |
|
|
errCount = errCount + 1 |
1421 |
|
|
ENDIF |
1422 |
|
|
C Check vertical diffusivity setting: |
1423 |
|
|
#ifdef ALLOW_3D_DIFFKR |
1424 |
|
|
IF ( specifiedDiffKrT ) THEN |
1425 |
|
|
WRITE(msgBuf,'(2A)') '** WARNING ** INI_PARMS: Ignores diffKr', |
1426 |
|
|
& 'T (or Kp,Kz) setting in file "data" with ALLOW_3D_DIFFKR' |
1427 |
|
|
CALL PRINT_MESSAGE( msgBuf, errorMessageUnit, |
1428 |
|
|
& SQUEEZE_RIGHT, myThid ) |
1429 |
|
|
ENDIF |
1430 |
|
|
IF ( specifiedDiffKrS .AND. diffKrFile.NE.' ' ) THEN |
1431 |
|
|
WRITE(msgBuf,'(2A)') '** WARNING ** INI_PARMS: Ignores diffKr', |
1432 |
|
|
& 'S (or Kp,Kz) setting in file "data" and uses diffKrFile' |
1433 |
|
|
CALL PRINT_MESSAGE( msgBuf, errorMessageUnit, |
1434 |
|
|
& SQUEEZE_RIGHT, myThid ) |
1435 |
|
|
ENDIF |
1436 |
|
|
#endif |
1437 |
|
|
|
1438 |
|
|
C-- Set Units conversion factor required to incorporate |
1439 |
|
|
C surface forcing into z-p isomorphic equations: |
1440 |
|
|
C mass2rUnit: from mass per unit area [kg/m2] to r-coordinate (z:=1/rho;p:=g) |
1441 |
|
|
C rUnit2mass: from r-coordinate to mass per unit area [kg/m2] (z:=rho;p:=1/g) |
1442 |
|
|
IF ( usingPCoords ) THEN |
1443 |
|
|
mass2rUnit = gravity |
1444 |
|
|
rUnit2mass = recip_gravity |
1445 |
|
|
ELSE |
1446 |
|
|
mass2rUnit = recip_rhoConst |
1447 |
|
|
rUnit2mass = rhoConst |
1448 |
|
|
ENDIF |
1449 |
|
|
|
1450 |
|
|
C-- For backward compatibility, set temp_addMass and salt_addMass |
1451 |
|
|
C to temp_EvPrRn and salt_EvPrRn if not set in parameter file "data" |
1452 |
|
|
IF (temp_addMass .EQ. UNSET_RL) temp_addMass = temp_EvPrRn |
1453 |
|
|
IF (salt_addMass .EQ. UNSET_RL) salt_addMass = salt_EvPrRn |
1454 |
|
|
|
1455 |
|
|
C-- Make a choice (if unset) regarding using CG2D solver min.-residual sol. |
1456 |
|
|
C for simple set-up (cartesian grid + flat bottom), default is to |
1457 |
|
|
C use the solver minimum residual solution (cg2dUseMinResSol=1). |
1458 |
|
|
IF ( cg2dUseMinResSol.EQ.UNSET_I ) THEN |
1459 |
|
|
cg2dUseMinResSol = 0 |
1460 |
|
|
IF ( topoFile.EQ.' ' .AND. bathyFile.EQ.' ' |
1461 |
|
|
& .AND. usingCartesianGrid ) cg2dUseMinResSol = 1 |
1462 |
|
|
ENDIF |
1463 |
|
|
|
1464 |
|
|
C-- Close the open data file |
1465 |
|
|
CLOSE(iUnit) |
1466 |
|
|
|
1467 |
|
|
WRITE(msgBuf,'(A)') ' INI_PARMS: finished reading file "data"' |
1468 |
|
|
CALL PRINT_MESSAGE( msgBuf, standardMessageUnit, |
1469 |
|
|
& SQUEEZE_RIGHT , myThid ) |
1470 |
|
|
|
1471 |
|
|
C-- Check whether any retired parameters were found. |
1472 |
|
|
IF ( nRetired .GT. 0 ) THEN |
1473 |
|
|
WRITE(msgBuf,'(A)') |
1474 |
|
|
& 'S/R INI_PARMS: Error reading parameter file "data":' |
1475 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1476 |
|
|
WRITE(msgBuf,'(I4,A)') nRetired, |
1477 |
|
|
& ' out of date parameters were found in the namelist' |
1478 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1479 |
|
|
errCount = errCount + 1 |
1480 |
|
|
ENDIF |
1481 |
|
|
C-- Stop if any error was found (including retired params): |
1482 |
|
|
IF ( errCount .GE. 1 ) THEN |
1483 |
|
|
WRITE(msgBuf,'(A,I3,A)') |
1484 |
|
|
& 'S/R INI_PARMS: detected', errCount,' fatal error(s)' |
1485 |
|
|
CALL PRINT_ERROR( msgBuf, myThid ) |
1486 |
|
|
CALL ALL_PROC_DIE( 0 ) |
1487 |
|
|
STOP 'ABNORMAL END: S/R INI_PARMS' |
1488 |
|
|
ENDIF |
1489 |
|
|
|
1490 |
|
|
_END_MASTER(myThid) |
1491 |
|
|
|
1492 |
|
|
C-- Everyone else must wait for the parameters to be loaded |
1493 |
|
|
_BARRIER |
1494 |
|
|
|
1495 |
|
|
RETURN |
1496 |
|
|
END |