/[MITgcm]/MITgcm_contrib/ksnow/press_release/code/ini_parms.F
ViewVC logotype

Contents of /MITgcm_contrib/ksnow/press_release/code/ini_parms.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.3 - (show annotations) (download)
Sat Mar 4 11:57:16 2017 UTC (8 years, 4 months ago) by dgoldberg
Branch: MAIN
CVS Tags: HEAD
Changes since 1.2: +2 -2 lines
further changes

1 C $Header: /u/gcmpack/MITgcm_contrib/ksnow/press_release/code/ini_parms.F,v 1.2 2017/01/30 16:32:18 ksnow 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, selectImplicitDrag,
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_PRESSURE_RELEASE_CODE
250 & ,cg2dMinColumnEps, pReleaseVisc, pReleaseDamp
251 #endif
252
253 C-- Time stepping parammeters
254 NAMELIST /PARM03/
255 & nIter0, nTimeSteps, nTimeSteps_l2, nEndIter,
256 & baseTime, startTime, endTime,
257 & deltaT, deltaTClock, deltaTMom,
258 & deltaTtracer, dTtracerLev, deltaTFreeSurf,
259 & forcing_In_AB, momForcingOutAB, tracForcingOutAB,
260 & momDissip_In_AB, doAB_onGtGs,
261 & abEps, alph_AB, beta_AB, startFromPickupAB2, applyExchUV_early,
262 & tauCD, rCD, epsAB_CD, cAdjFreq,
263 & chkPtFreq, pChkPtFreq, pickupSuff, pickupStrictlyMatch,
264 & writePickupAtEnd,
265 & dumpFreq, dumpInitAndLast, adjDumpFreq, taveFreq, tave_lastIter,
266 & diagFreq, monitorFreq, adjMonitorFreq, monitorSelect,
267 & outputTypesInclusive,
268 & tauThetaClimRelax, tauSaltClimRelax, latBandClimRelax,
269 & tauThetaClimRelax3Dim, tauSaltClimRelax3Dim,
270 & periodicExternalForcing, externForcingPeriod, externForcingCycle,
271 & calendarDumps
272
273 C-- Gridding parameters
274 NAMELIST /PARM04/
275 & usingCartesianGrid, usingCylindricalGrid,
276 & usingSphericalPolarGrid, usingCurvilinearGrid,
277 & xgOrigin, ygOrigin, dxSpacing, dySpacing,
278 & delX, delY, delXFile, delYFile, horizGridFile,
279 & phiEuler, thetaEuler, psiEuler,
280 & rSphere, radius_fromHorizGrid, deepAtmosphere, seaLev_Z,
281 & top_Pres, delZ, delP, delR, delRc, delRFile, delRcFile,
282 & selectSigmaCoord, rSigmaBnd, hybSigmFile,
283 & Ro_SeaLevel, rkFac, groundAtK1, thetaMin, phiMin
284
285 C-- Input files
286 NAMELIST /PARM05/
287 & bathyFile, topoFile, addWwallFile, addSwallFile, shelfIceFile,
288 & diffKrFile, viscAhDfile, viscAhZfile, viscA4Dfile, viscA4Zfile,
289 & hydrogThetaFile, hydrogSaltFile,
290 & maskIniTemp, maskIniSalt, checkIniTemp, checkIniSalt,
291 & zonalWindFile, meridWindFile,
292 & thetaClimFile, saltClimFile,
293 & surfQfile, surfQnetFile, surfQswFile, EmPmRfile, saltFluxFile,
294 & lambdaThetaFile, lambdaSaltFile,
295 & uVelInitFile, vVelInitFile, pSurfInitFile,
296 & dQdTFile, ploadFile, addMassFile, tCylIn, tCylOut,
297 & eddyPsiXFile, eddyPsiYFile, geothermalFile,
298 & mdsioLocalDir, adTapeDir,
299 & the_run_name
300 CEOP
301
302 #ifdef ALLOW_EXCH2
303 gridNx = exch2_mydNx(1)
304 gridNy = exch2_mydNy(1)
305 #else /* ALLOW_EXCH2 */
306 gridNx = Nx
307 gridNy = Ny
308 #endif /* ALLOW_EXCH2 */
309
310 _BEGIN_MASTER(myThid)
311
312 C Defaults values for input parameters
313 CALL SET_DEFAULTS(
314 O viscArDefault, diffKrTDefault, diffKrSDefault,
315 O hFacMinDrDefault, delRDefault,
316 I myThid )
317 SadournyCoriolis = .FALSE.
318
319 C-- Initialise "which vertical coordinate system used" flags.
320 zCoordInputData = .FALSE.
321 pCoordInputData = .FALSE.
322 rCoordInputData = .FALSE.
323 coordsSet = 0
324
325 C-- Initialise retired parameters to unlikely value
326 nRetired = 0
327 useConstantF = .FALSE.
328 useBetaPlaneF = .FALSE.
329 useSphereF = .TRUE.
330 tempDiffusion = .TRUE.
331 saltDiffusion = .TRUE.
332 tracerAdvScheme = UNSET_I
333 trac_EvPrRn = UNSET_RL
334 zonal_filt_lat = UNSET_RL
335 gravitySign = UNSET_RL
336 viscAstrain = UNSET_RL
337 viscAtension = UNSET_RL
338 useAnisotropicViscAgridMax=.TRUE.
339 usePickupBeforeC35 = .FALSE.
340 saveDebugMode = debugMode
341 allowInteriorFreezing = .FALSE.
342 useOldFreezing = .FALSE.
343 tauThetaClimRelax3Dim = UNSET_RL
344 tauSaltClimRelax3Dim = UNSET_RL
345 calendarDumps = .FALSE.
346 Ro_SeaLevel = UNSET_RL
347 rkFac = UNSET_RL
348 groundAtK1 = .FALSE.
349 thetaMin = UNSET_RL
350 phiMin = UNSET_RL
351 shelfIceFile = ' '
352 dQdTFile = ' '
353
354 C-- Open the parameter file
355 WRITE(msgBuf,'(A)')
356 & ' INI_PARMS: opening model parameter file "data"'
357 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
358 & SQUEEZE_RIGHT, myThid )
359
360 CALL OPEN_COPY_DATA_FILE( 'data', 'INI_PARMS',
361 O iUnit, myThid )
362
363 C-- Read settings from iUnit (= a copy of model parameter file "data").
364 errIO = 0
365 errCount = 0
366
367 C-- Set default "physical" parameters
368 viscAhW = UNSET_RL
369 viscA4W = UNSET_RL
370 viscAhD = UNSET_RL
371 viscAhZ = UNSET_RL
372 viscA4D = UNSET_RL
373 viscA4Z = UNSET_RL
374 viscAz = UNSET_RL
375 viscAp = UNSET_RL
376 viscAr = UNSET_RL
377 diffKzT = UNSET_RL
378 diffKpT = UNSET_RL
379 diffKrT = UNSET_RL
380 diffKzS = UNSET_RL
381 diffKpS = UNSET_RL
382 diffKrS = UNSET_RL
383 hFacMinDr = UNSET_RL
384 hFacMinDz = UNSET_RL
385 hFacMinDp = UNSET_RL
386 tAlpha = UNSET_RL
387 sBeta = UNSET_RL
388 implicitNHPress = UNSET_RL
389 tempVertAdvScheme = 0
390 saltVertAdvScheme = 0
391 C-- z,p,r coord input switching.
392 WRITE(msgBuf,'(A)') ' INI_PARMS ; starts to read PARM01'
393 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
394 & SQUEEZE_RIGHT, myThid )
395 READ(UNIT=iUnit,NML=PARM01) !,IOSTAT=errIO)
396 IF ( errIO .LT. 0 ) THEN
397 WRITE(msgBuf,'(A)')
398 & 'S/R INI_PARMS: Error reading model parameter file "data"'
399 CALL PRINT_ERROR( msgBuf, myThid )
400 WRITE(msgBuf,'(A)') 'S/R INI_PARMS: Problem in namelist PARM01'
401 CALL PRINT_ERROR( msgBuf, myThid )
402 STOP 'ABNORMAL END: S/R INI_PARMS'
403 ELSE
404 WRITE(msgBuf,'(A)') ' INI_PARMS ; read PARM01 : OK'
405 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
406 & SQUEEZE_RIGHT, myThid )
407 ENDIF
408
409 C- set the type of vertical coordinate and type of fluid
410 C according to buoyancyRelation
411 usingPCoords = .FALSE.
412 usingZCoords = .FALSE.
413 fluidIsAir = .FALSE.
414 fluidIsWater = .FALSE.
415 IF ( buoyancyRelation.EQ.'ATMOSPHERIC' ) THEN
416 usingPCoords = .TRUE.
417 fluidIsAir = .TRUE.
418 ELSEIF ( buoyancyRelation.EQ.'OCEANICP') THEN
419 usingPCoords = .TRUE.
420 fluidIsWater = .TRUE.
421 ELSEIF ( buoyancyRelation.EQ.'OCEANIC' ) THEN
422 usingZCoords = .TRUE.
423 fluidIsWater = .TRUE.
424 ELSE
425 WRITE(msgBuf,'(2A)') 'S/R INI_PARMS:',
426 & ' Bad value of buoyancyRelation '
427 CALL PRINT_ERROR( msgBuf, myThid )
428 errCount = errCount + 1
429 ENDIF
430
431 IF ( .NOT.rigidLid .AND.
432 & .NOT.implicitFreeSurface ) THEN
433 C- No barotropic solver selected => use implicitFreeSurface as default
434 WRITE(msgBuf,'(A)')
435 & 'S/R INI_PARMS: No request for barotropic solver'
436 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
437 & SQUEEZE_RIGHT, myThid )
438 WRITE(msgBuf,'(A)')
439 & 'S/R INI_PARMS: => Use implicitFreeSurface as default'
440 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
441 & SQUEEZE_RIGHT, myThid )
442 implicitFreeSurface = .TRUE.
443 ENDIF
444 IF ( implicitFreeSurface ) freeSurfFac = 1. _d 0
445 IF ( rigidLid ) freeSurfFac = 0. _d 0
446 IF ( gBaro .EQ. UNSET_RL ) gBaro=gravity
447 IF ( rhoConst .EQ. UNSET_RL ) rhoConst=rhoNil
448 IF ( rhoConstFresh .EQ. UNSET_RL ) rhoConstFresh=rhoConst
449 IF ( implicitNHPress.EQ.UNSET_RL )
450 & implicitNHPress = implicSurfPress
451 IF ( omega .EQ. UNSET_RL ) THEN
452 omega = 0. _d 0
453 IF ( rotationPeriod .NE. 0. _d 0 )
454 & omega = 2. _d 0 * PI / rotationPeriod
455 ELSEIF ( omega .EQ. 0. _d 0 ) THEN
456 rotationPeriod = 0. _d 0
457 ELSE
458 rotationPeriod = 2. _d 0 * PI / omega
459 ENDIF
460 IF ( atm_Rd .EQ. UNSET_RL ) THEN
461 atm_Rd = atm_Cp * atm_kappa
462 ELSE
463 atm_kappa = atm_Rd / atm_Cp
464 ENDIF
465 C-- Non-hydrostatic/quasi-hydrostatic
466 IF (nonHydrostatic.AND.quasiHydrostatic) THEN
467 WRITE(msgBuf,'(A)')
468 & 'Illegal: both nonHydrostatic = quasiHydrostatic = TRUE'
469 CALL PRINT_ERROR( msgBuf, myThid )
470 errCount = errCount + 1
471 ENDIF
472 C-- Advection and Forcing for Temp and salt
473 IF (tempVertAdvScheme.EQ.0) tempVertAdvScheme = tempAdvScheme
474 IF (saltVertAdvScheme.EQ.0) saltVertAdvScheme = saltAdvScheme
475 C-- horizontal viscosity (acting on Divergence or Vorticity)
476 IF ( viscAhD .EQ. UNSET_RL ) viscAhD = viscAh
477 IF ( viscAhZ .EQ. UNSET_RL ) viscAhZ = viscAh
478 IF ( viscA4D .EQ. UNSET_RL ) viscA4D = viscA4
479 IF ( viscA4Z .EQ. UNSET_RL ) viscA4Z = viscA4
480 C-- horizontal viscosity for vertical moments
481 IF ( viscAhW .EQ. UNSET_RL ) viscAhW = viscAhD
482 IF ( viscA4W .EQ. UNSET_RL ) viscA4W = viscA4D
483 C-- z,p,r coord input switching.
484 IF ( viscAz .NE. UNSET_RL ) zCoordInputData = .TRUE.
485 IF ( viscAp .NE. UNSET_RL ) pCoordInputData = .TRUE.
486 IF ( viscAr .NE. UNSET_RL ) rCoordInputData = .TRUE.
487 IF ( viscAr .EQ. UNSET_RL ) viscAr = viscAz
488 IF ( viscAr .EQ. UNSET_RL ) viscAr = viscAp
489 vertSetCount = 0
490 DO k=1,Nr
491 IF ( viscArNr(k).NE.UNSET_RL ) vertSetCount = vertSetCount + 1
492 ENDDO
493 IF ( vertSetCount.GT.0 .AND. vertSetCount.LT.Nr ) THEN
494 WRITE(msgBuf,'(A,2(I5,A))') 'S/R INI_PARMS: Partial setting (',
495 & vertSetCount, ' /', Nr, ') of viscArNr is not allowed'
496 CALL PRINT_ERROR( msgBuf, myThid )
497 errCount = errCount + 1
498 ENDIF
499 IF ( viscAr .EQ. UNSET_RL ) THEN
500 viscAr = viscArDefault
501 ELSEIF ( vertSetCount.GT.0 ) THEN
502 WRITE(msgBuf,'(2A)') 'S/R INI_PARMS: Cannot set both ',
503 & 'viscArNr and viscAr (or Ap,Az) in param file data'
504 CALL PRINT_ERROR( msgBuf, myThid )
505 errCount = errCount + 1
506 ENDIF
507 IF ( vertSetCount.EQ.0 ) THEN
508 DO k=1,Nr
509 viscArNr(k) = viscAr
510 ENDDO
511 ENDIF
512 #ifdef ALLOW_MOM_COMMON
513 C- set default scheme for quadratic bottom-drag
514 IF ( selectBotDragQuadr.EQ.-1 .AND. bottomDragQuadratic.NE.0. )
515 & selectBotDragQuadr = 0
516 #endif /* ALLOW_MOM_COMMON */
517
518 IF ( diffKzT .NE. UNSET_RL ) zCoordInputData = .TRUE.
519 IF ( diffKpT .NE. UNSET_RL ) pCoordInputData = .TRUE.
520 IF ( diffKrT .NE. UNSET_RL ) rCoordInputData = .TRUE.
521 IF ( diffKrT .EQ. UNSET_RL ) diffKrT = diffKzT
522 IF ( diffKrT .EQ. UNSET_RL ) diffKrT = diffKpT
523 vertSetCount = 0
524 DO k=1,Nr
525 IF ( diffKrNrT(k).NE.UNSET_RL ) vertSetCount = vertSetCount + 1
526 ENDDO
527 IF ( vertSetCount.GT.0 .AND. vertSetCount.LT.Nr ) THEN
528 WRITE(msgBuf,'(A,2(I5,A))') 'S/R INI_PARMS: Partial setting (',
529 & vertSetCount, ' /', Nr, ') of diffKrNrT is not allowed'
530 CALL PRINT_ERROR( msgBuf, myThid )
531 errCount = errCount + 1
532 ENDIF
533 specifiedDiffKrT = vertSetCount.EQ.Nr
534 IF ( diffKrT .EQ. UNSET_RL ) THEN
535 diffKrT = diffKrTDefault
536 ELSEIF ( vertSetCount.GT.0 ) THEN
537 WRITE(msgBuf,'(2A)') 'S/R INI_PARMS: Cannot set both ',
538 & 'diffKrNrT and diffKrT (or Kp,Kz) in param file data'
539 CALL PRINT_ERROR( msgBuf, myThid )
540 errCount = errCount + 1
541 ELSE
542 specifiedDiffKrT = .TRUE.
543 ENDIF
544 IF ( vertSetCount.EQ.0 ) THEN
545 DO k=1,Nr
546 diffKrNrT(k) = diffKrT
547 ENDDO
548 ENDIF
549
550 IF ( diffKzS .NE. UNSET_RL ) zCoordInputData = .TRUE.
551 IF ( diffKpS .NE. UNSET_RL ) pCoordInputData = .TRUE.
552 IF ( diffKrS .NE. UNSET_RL ) rCoordInputData = .TRUE.
553 IF ( diffKrS .EQ. UNSET_RL ) diffKrS = diffKzS
554 IF ( diffKrS .EQ. UNSET_RL ) diffKrS = diffKpS
555 vertSetCount = 0
556 DO k=1,Nr
557 IF ( diffKrNrS(k).NE.UNSET_RL ) vertSetCount = vertSetCount + 1
558 ENDDO
559 IF ( vertSetCount.GT.0 .AND. vertSetCount.LT.Nr ) THEN
560 WRITE(msgBuf,'(A,2(I5,A))') 'S/R INI_PARMS: Partial setting (',
561 & vertSetCount, ' /', Nr, ') of diffKrNrS is not allowed'
562 CALL PRINT_ERROR( msgBuf, myThid )
563 errCount = errCount + 1
564 ENDIF
565 IF ( vertSetCount.EQ.Nr ) THEN
566 specifiedDiffKrS = .TRUE.
567 IF ( diffKrS.NE.UNSET_RL ) THEN
568 WRITE(msgBuf,'(2A)') 'S/R INI_PARMS: Cannot set both ',
569 & 'diffKrNrS and diffKrS (or Kp,Kz) in param file data'
570 CALL PRINT_ERROR( msgBuf, myThid )
571 errCount = errCount + 1
572 ENDIF
573 ELSEIF ( diffKrS.NE.UNSET_RL ) THEN
574 specifiedDiffKrS = .TRUE.
575 DO k=1,Nr
576 diffKrNrS(k) = diffKrS
577 ENDDO
578 ELSE
579 specifiedDiffKrS = .FALSE.
580 diffKrS = diffKrSDefault
581 C- use temp diffusivity as default salt diffusivity:
582 DO k=1,Nr
583 diffKrNrS(k) = diffKrNrT(k)
584 ENDDO
585 ENDIF
586
587 IF (diffKrBLEQsurf .EQ. UNSET_RL) diffKrBLEQsurf = diffKrBL79surf
588 IF (diffKrBLEQdeep .EQ. UNSET_RL) diffKrBLEQdeep = diffKrBL79deep
589 IF (diffKrBLEQscl .EQ. UNSET_RL) diffKrBLEQscl = diffKrBL79scl
590 IF (diffKrBLEQHo .EQ. UNSET_RL) diffKrBLEQHo = diffKrBL79Ho
591
592 IF ( hFacMinDz .NE. UNSET_RL ) zCoordInputData = .TRUE.
593 IF ( hFacMinDp .NE. UNSET_RL ) pCoordInputData = .TRUE.
594 IF ( hFacMinDr .NE. UNSET_RL ) rCoordInputData = .TRUE.
595 IF ( hFacMinDr .EQ. UNSET_RL ) hFacMinDr = hFacMinDz
596 IF ( hFacMinDr .EQ. UNSET_RL ) hFacMinDr = hFacMinDp
597 IF ( hFacMinDr .EQ. UNSET_RL ) hFacMinDr = hFacMinDrDefault
598
599 IF (convertFW2Salt.EQ.UNSET_RL) THEN
600 convertFW2Salt = 35.
601 IF (useRealFreshWaterFlux) convertFW2Salt=-1
602 IF ( selectAddFluid.GE.1 ) convertFW2Salt=-1
603 ENDIF
604
605 IF ( SadournyCoriolis ) THEN
606 C-- for backward compatibility :
607 IF ( selectVortScheme.EQ.UNSET_I ) selectVortScheme = 2
608 IF ( selectVortScheme.NE.2 ) THEN
609 WRITE(msgBuf,'(A,I5,A)')
610 & 'S/R INI_PARMS: selectVortScheme=', selectVortScheme,
611 & ' conflicts with "SadournyCoriolis"'
612 CALL PRINT_ERROR( msgBuf, myThid )
613 errCount = errCount + 1
614 ENDIF
615 ENDIF
616
617 IF ( ivdc_kappa.NE.zeroRL .AND. .NOT.implicitDiffusion ) THEN
618 WRITE(msgBuf,'(A,A)')
619 & 'S/R INI_PARMS: To use ivdc_kappa you must enable implicit',
620 & ' vertical diffusion.'
621 CALL PRINT_ERROR( msgBuf, myThid )
622 errCount = errCount + 1
623 ENDIF
624
625 coordsSet = 0
626 IF ( zCoordInputData ) coordsSet = coordsSet + 1
627 IF ( pCoordInputData ) coordsSet = coordsSet + 1
628 IF ( rCoordInputData ) coordsSet = coordsSet + 1
629 IF ( coordsSet .GT. 1 ) THEN
630 WRITE(msgBuf,'(A)')
631 & 'S/R INI_PARMS: Cannot mix z, p and r in the input data.'
632 CALL PRINT_ERROR( msgBuf, myThid )
633 errCount = errCount + 1
634 ENDIF
635 IF ( rhoConst .LE. 0. ) THEN
636 WRITE(msgBuf,'(A)')
637 & 'S/R INI_PARMS: rhoConst must be greater than 0.'
638 CALL PRINT_ERROR( msgBuf, myThid )
639 errCount = errCount + 1
640 recip_rhoConst = 1. _d 0
641 ELSE
642 recip_rhoConst = 1. _d 0 / rhoConst
643 ENDIF
644 IF ( eosType.EQ.'LINEAR' .AND. rhoNil.LE.0. ) THEN
645 WRITE(msgBuf,'(A)')
646 & 'S/R INI_PARMS: rhoNil must be greater than 0.'
647 CALL PRINT_ERROR( msgBuf, myThid )
648 errCount = errCount + 1
649 ENDIF
650 IF ( HeatCapacity_Cp .LE. 0. ) THEN
651 WRITE(msgBuf,'(A)')
652 & 'S/R INI_PARMS: HeatCapacity_Cp must be greater than 0.'
653 CALL PRINT_ERROR( msgBuf, myThid )
654 errCount = errCount + 1
655 ENDIF
656 IF ( gravity .LE. 0. ) THEN
657 WRITE(msgBuf,'(A)')
658 & 'S/R INI_PARMS: gravity must be greater than 0.'
659 CALL PRINT_ERROR( msgBuf, myThid )
660 errCount = errCount + 1
661 recip_gravity = 1. _d 0
662 ELSE
663 recip_gravity = 1. _d 0 / gravity
664 ENDIF
665
666 C- set default printResidualFreq according to debugLevel
667 printResidualFreq = 0
668 IF ( debugLevel.GE.debLevE ) printResidualFreq = 1
669
670 C- set useSingleCpuInput=.TRUE. if useSingleCpuIO==.TRUE.
671 IF ( useSingleCpuIO ) useSingleCpuInput=.TRUE.
672
673 C Check for retired parameters still being used
674 nRetired = 0
675 IF ( useConstantF ) THEN
676 nRetired = nRetired+1
677 WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "useConstantF" ',
678 & 'is no longer allowed in file "data"'
679 CALL PRINT_ERROR( msgBuf, myThid )
680 WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: set "selectCoriMap"',
681 & ' [0,1,2,3] to impose a setting over grid default'
682 CALL PRINT_ERROR( msgBuf, myThid )
683 ENDIF
684 IF ( useBetaPlaneF ) THEN
685 nRetired = nRetired+1
686 WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "useBetaPlaneF" ',
687 & 'is no longer allowed in file "data"'
688 CALL PRINT_ERROR( msgBuf, myThid )
689 WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: set "selectCoriMap"',
690 & ' [0,1,2,3] to impose a setting over grid default'
691 CALL PRINT_ERROR( msgBuf, myThid )
692 ENDIF
693 IF ( .NOT. useSphereF ) THEN
694 nRetired = nRetired+1
695 WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "useSphereF" ',
696 & 'is no longer allowed in file "data"'
697 CALL PRINT_ERROR( msgBuf, myThid )
698 WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: set "selectCoriMap"',
699 & ' [0,1,2,3] to impose a setting over grid default'
700 CALL PRINT_ERROR( msgBuf, myThid )
701 ENDIF
702 IF ( zonal_filt_lat .NE. UNSET_RL ) THEN
703 nRetired = nRetired+1
704 WRITE(msgBuf,'(A,A)')
705 & 'S/R INI_PARMS: Paramater "zonal_filt_lat" is',
706 & ' no longer allowed in file "data".'
707 CALL PRINT_ERROR( msgBuf, myThid )
708 WRITE(msgBuf,'(A,A)')
709 & 'S/R INI_PARMS: Paramater "zonal_filt_lat" is',
710 & ' now read from file "data.zonfilt".'
711 CALL PRINT_ERROR( msgBuf, myThid )
712 ENDIF
713 IF ( gravitySign .NE. UNSET_RL ) THEN
714 nRetired = nRetired+1
715 WRITE(msgBuf,'(A,A)')
716 & 'S/R INI_PARMS: "gravitySign" is set according to vertical ',
717 & ' coordinate and is no longer allowed in file "data".'
718 CALL PRINT_ERROR( msgBuf, myThid )
719 ENDIF
720 IF ( tracerAdvScheme .NE. UNSET_I ) THEN
721 nRetired = nRetired+1
722 WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "tracerAdvScheme" ',
723 & '(old passive tracer code) is no longer allowed in file "data"'
724 CALL PRINT_ERROR( msgBuf, myThid )
725 ENDIF
726 IF ( trac_EvPrRn .NE. UNSET_RL ) THEN
727 nRetired = nRetired+1
728 WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "trac_EvPrRn" ',
729 & '(old passive tracer code) is no longer allowed in file "data"'
730 CALL PRINT_ERROR( msgBuf, myThid )
731 ENDIF
732 IF ( .NOT. tempDiffusion ) THEN
733 nRetired = nRetired+1
734 WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "tempDiffusion" ',
735 & 'is no longer allowed in file "data"'
736 CALL PRINT_ERROR( msgBuf, myThid )
737 WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: to turn off diffusion',
738 & ' => set diffusivity to zero'
739 CALL PRINT_ERROR( msgBuf, myThid )
740 ENDIF
741 IF ( .NOT. saltDiffusion ) THEN
742 nRetired = nRetired+1
743 WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "saltDiffusion" ',
744 & 'is no longer allowed in file "data"'
745 CALL PRINT_ERROR( msgBuf, myThid )
746 WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: to turn off diffusion',
747 & ' => set diffusivity to zero'
748 CALL PRINT_ERROR( msgBuf, myThid )
749 ENDIF
750 IF ( viscAstrain .NE. UNSET_RL ) THEN
751 nRetired = nRetired+1
752 WRITE(msgBuf,'(A,A)')
753 & 'S/R INI_PARMS: "viscAstrain" ',
754 & 'is no longer allowed in file "data"'
755 CALL PRINT_ERROR( msgBuf, myThid )
756 WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: to use Strain & Tension',
757 & ' formulation => set useStrainTensionVisc to TRUE'
758 CALL PRINT_ERROR( msgBuf, myThid )
759 ENDIF
760 IF ( viscAtension .NE. UNSET_RL ) THEN
761 nRetired = nRetired+1
762 WRITE(msgBuf,'(A,A)')
763 & 'S/R INI_PARMS: "viscAtension" ',
764 & 'is no longer allowed in file "data"'
765 CALL PRINT_ERROR( msgBuf, myThid )
766 WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: to use Strain & Tension',
767 & ' formulation => set useStrainTensionVisc to TRUE'
768 CALL PRINT_ERROR( msgBuf, myThid )
769 ENDIF
770 IF ( .NOT.useAnisotropicViscAgridMax ) THEN
771 nRetired = nRetired+1
772 WRITE(msgBuf,'(A,A)')
773 & 'S/R INI_PARMS: "useAnisotropicViscAgridMax" ',
774 & 'is not allowed in "data" substitute useAreaViscLength=true'
775 CALL PRINT_ERROR( msgBuf, myThid )
776 ENDIF
777 IF ( usePickupBeforeC35 ) THEN
778 nRetired = nRetired+1
779 WRITE(msgBuf,'(A,A)')
780 & 'S/R INI_PARMS: "usePickupBeforeC35" ',
781 & 'is no longer supported & not longer allowed in file "data"'
782 CALL PRINT_ERROR( msgBuf, myThid )
783 ENDIF
784 IF ( debugMode.NEQV.saveDebugMode ) THEN
785 nRetired = nRetired+1
786 WRITE(msgBuf,'(A,A)')
787 & 'S/R INI_PARMS: "debugMode" has been moved to "eedata"',
788 & ' and is no longer allowed in file "data"'
789 CALL PRINT_ERROR( msgBuf, myThid )
790 ENDIF
791 IF ( allowInteriorFreezing ) THEN
792 nRetired = nRetired+1
793 WRITE(msgBuf,'(A,A)')
794 & 'S/R INI_PARMS: "allowInteriorFreezing" has been replaced',
795 & ' by pkg/frazil and is no longer allowed in file "data"'
796 CALL PRINT_ERROR( msgBuf, myThid )
797 ENDIF
798 IF ( useOldFreezing ) THEN
799 nRetired = nRetired+1
800 WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "useOldFreezing" ',
801 & 'is no longer supported & not longer allowed in file "data"'
802 CALL PRINT_ERROR( msgBuf, myThid )
803 ENDIF
804
805 C-- Elliptic solver parameters
806 WRITE(msgBuf,'(A)') ' INI_PARMS ; starts to read PARM02'
807 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
808 & SQUEEZE_RIGHT, myThid )
809 READ(UNIT=iUnit,NML=PARM02) !,IOSTAT=errIO)
810 IF ( errIO .LT. 0 ) THEN
811 WRITE(msgBuf,'(A)')
812 & 'S/R INI_PARMS: Error reading model parameter file "data"'
813 CALL PRINT_ERROR( msgBuf, myThid )
814 WRITE(msgBuf,'(A)') 'S/R INI_PARMS: Problem in namelist PARM02'
815 CALL PRINT_ERROR( msgBuf, myThid )
816 STOP 'ABNORMAL END: S/R INI_PARMS'
817 ELSE
818 WRITE(msgBuf,'(A)') ' INI_PARMS ; read PARM02 : OK'
819 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
820 & SQUEEZE_RIGHT, myThid )
821 ENDIF
822
823 C-- Time stepping parameters
824 rCD = -1. _d 0
825 epsAB_CD = UNSET_RL
826 latBandClimRelax = UNSET_RL
827 deltaTtracer = 0. _d 0
828 forcing_In_AB = .TRUE.
829 WRITE(msgBuf,'(A)') ' INI_PARMS ; starts to read PARM03'
830 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
831 & SQUEEZE_RIGHT, myThid )
832 READ(UNIT=iUnit,NML=PARM03) !,IOSTAT=errIO)
833 IF ( errIO .LT. 0 ) THEN
834 WRITE(msgBuf,'(A)')
835 & 'S/R INI_PARMS: Error reading model parameter file "data"'
836 CALL PRINT_ERROR( msgBuf, myThid )
837 WRITE(msgBuf,'(A)') 'S/R INI_PARMS: Problem in namelist PARM03'
838 CALL PRINT_ERROR( msgBuf, myThid )
839 STOP 'ABNORMAL END: S/R INI_PARMS'
840 ELSE
841 WRITE(msgBuf,'(A)') ' INI_PARMS ; read PARM03 : OK'
842 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
843 & SQUEEZE_RIGHT, myThid )
844 ENDIF
845 C Check for retired parameters still being used
846 IF ( tauThetaClimRelax3Dim .NE. UNSET_RL ) THEN
847 nRetired = nRetired+1
848 WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "tauThetaClimRelax3Dim" ',
849 & 'is no longer allowed in file "data"'
850 CALL PRINT_ERROR( msgBuf, myThid )
851 WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: 3-dim. relaxation code',
852 & ' has moved to separate pkg/rbcs.'
853 CALL PRINT_ERROR( msgBuf, myThid )
854 ENDIF
855 IF ( tauSaltClimRelax3Dim .NE. UNSET_RL ) THEN
856 nRetired = nRetired+1
857 WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "tauSaltClimRelax3Dim" ',
858 & 'is no longer allowed in file "data"'
859 CALL PRINT_ERROR( msgBuf, myThid )
860 WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: 3-dim. relaxation code',
861 & ' has moved to separate pkg/rbcs.'
862 CALL PRINT_ERROR( msgBuf, myThid )
863 ENDIF
864 IF ( calendarDumps ) THEN
865 nRetired = nRetired+1
866 WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "calendarDumps" ',
867 & 'is no longer allowed in file "data"'
868 CALL PRINT_ERROR( msgBuf, myThid )
869 WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: calendarDumps',
870 & ' has moved to "data.cal"'
871 CALL PRINT_ERROR( msgBuf, myThid )
872 ENDIF
873
874 C Process "timestepping" params
875 C o Time step size
876 IF ( deltaTtracer .NE. dTtracerLev(1) .AND.
877 & deltaTtracer .NE. 0. .AND. dTtracerLev(1) .NE. 0. ) THEN
878 WRITE(msgBuf,'(A)')
879 & 'S/R INI_PARMS: deltaTtracer & dTtracerLev(1) not equal'
880 CALL PRINT_ERROR( msgBuf, myThid )
881 errCount = errCount + 1
882 ELSEIF ( dTtracerLev(1) .NE. 0. ) THEN
883 deltaTtracer = dTtracerLev(1)
884 ENDIF
885 IF ( deltaT .EQ. 0. ) deltaT = deltaTClock
886 IF ( deltaT .EQ. 0. ) deltaT = deltaTtracer
887 IF ( deltaT .EQ. 0. ) deltaT = deltaTMom
888 IF ( deltaT .EQ. 0. ) deltaT = deltaTFreeSurf
889 IF ( deltaT .EQ. 0. ) THEN
890 WRITE(msgBuf,'(2A)') 'S/R INI_PARMS: ',
891 & 'need to specify in file "data", namelist "PARM03"'
892 CALL PRINT_ERROR( msgBuf, myThid )
893 WRITE(msgBuf,'(2A)') 'S/R INI_PARMS: ',
894 & ' a model timestep (in s) deltaT or deltaTClock= ?'
895 CALL PRINT_ERROR( msgBuf, myThid )
896 errCount = errCount + 1
897 deltaT = 1.
898 ENDIF
899 IF ( deltaTMom .EQ. 0. ) deltaTMom = deltaT
900 IF ( deltaTtracer .EQ. 0. ) deltaTtracer = deltaT
901 IF ( deltaTClock .EQ. 0. ) deltaTClock = deltaT
902 DO k=1,Nr
903 IF (dTtracerLev(k).EQ.0.) dTtracerLev(k)= deltaTtracer
904 ENDDO
905 C Note that this line should set deltaFreesurf=deltaTtracer
906 C but this would change a lot of existing set-ups so we are
907 C obliged to set the default inappropriately.
908 C Be advised that when using asynchronous time stepping
909 C it is better to set deltaTreesurf=deltaTtracer
910 IF ( deltaTFreeSurf .EQ. 0. ) deltaTFreeSurf = deltaTMom
911 IF ( periodicExternalForcing ) THEN
912 IF ( externForcingCycle*externForcingPeriod .EQ. 0. ) THEN
913 WRITE(msgBuf,'(A)')
914 & 'S/R INI_PARMS: externForcingCycle,externForcingPeriod =0'
915 CALL PRINT_ERROR( msgBuf, myThid )
916 errCount = errCount + 1
917 ELSEIF ( INT(externForcingCycle/externForcingPeriod) .NE.
918 & externForcingCycle/externForcingPeriod ) THEN
919 WRITE(msgBuf,'(A)')
920 & 'S/R INI_PARMS: externForcingCycle <> N*externForcingPeriod'
921 CALL PRINT_ERROR( msgBuf, myThid )
922 errCount = errCount + 1
923 ELSEIF ( externForcingCycle.LT.externForcingPeriod ) THEN
924 WRITE(msgBuf,'(A)')
925 & 'S/R INI_PARMS: externForcingCycle < externForcingPeriod'
926 CALL PRINT_ERROR( msgBuf, myThid )
927 errCount = errCount + 1
928 ENDIF
929 IF ( externForcingPeriod.LT.deltaTClock ) THEN
930 WRITE(msgBuf,'(A)')
931 & 'S/R INI_PARMS: externForcingPeriod < deltaTClock'
932 CALL PRINT_ERROR( msgBuf, myThid )
933 errCount = errCount + 1
934 ENDIF
935 ENDIF
936 C o Adams-Bashforth time stepping:
937 IF ( momForcingOutAB .EQ. UNSET_I ) THEN
938 momForcingOutAB = 1
939 IF ( forcing_In_AB ) momForcingOutAB = 0
940 ENDIF
941 IF ( tracForcingOutAB .EQ. UNSET_I ) THEN
942 tracForcingOutAB = 1
943 IF ( forcing_In_AB ) tracForcingOutAB = 0
944 ENDIF
945 C o Convection frequency
946 IF ( cAdjFreq .LT. 0. ) THEN
947 cAdjFreq = deltaTClock
948 ENDIF
949 IF ( ivdc_kappa .NE. 0. .AND. cAdjFreq .NE. 0. ) THEN
950 WRITE(msgBuf,'(A,A)')
951 & 'S/R INI_PARMS: You have enabled both ivdc_kappa and',
952 & ' convective adjustment.'
953 CALL PRINT_ERROR( msgBuf, myThid )
954 errCount = errCount + 1
955 ENDIF
956 IF (useCDscheme) THEN
957 C o CD coupling (CD scheme):
958 IF ( tauCD .EQ. 0. _d 0 ) tauCD = deltaTMom
959 IF ( rCD .LT. 0. ) rCD = 1. _d 0 - deltaTMom/tauCD
960 IF ( epsAB_CD .EQ. UNSET_RL ) epsAB_CD = abEps
961 ENDIF
962
963 IF ( startTime.EQ.UNSET_RL .AND. nIter0.EQ.-1 ) THEN
964 C o set default value for start time & nIter0
965 startTime = baseTime
966 nIter0 = 0
967 ELSEIF ( startTime.EQ.UNSET_RL ) THEN
968 C o set start time from nIter0
969 startTime = baseTime + deltaTClock*DFLOAT(nIter0)
970 ELSEIF ( nIter0.EQ.-1 ) THEN
971 C o set nIter0 from start time
972 nIter0 = NINT( (startTime-baseTime)/deltaTClock )
973 ELSEIF ( baseTime.EQ.0. ) THEN
974 C o set base time from the 2 others
975 baseTime = startTime - deltaTClock*DFLOAT(nIter0)
976 ENDIF
977
978 nTimeSteps_l2 = 4
979 C o nTimeSteps 1
980 IF ( nTimeSteps .EQ. 0 .AND. nEndIter .NE. 0 )
981 & nTimeSteps = nEndIter-nIter0
982 C o nTimeSteps 2
983 IF ( nTimeSteps .EQ. 0 .AND. endTime .NE. 0. )
984 & nTimeSteps = NINT((endTime-startTime)/deltaTClock)
985 C o nEndIter 1
986 IF ( nEndIter .EQ. 0 .AND. nTimeSteps .NE. 0 )
987 & nEndIter = nIter0+nTimeSteps
988 C o nEndIter 2
989 IF ( nEndIter .EQ. 0 .AND. endTime .NE. 0. )
990 & nEndIter = NINT((endTime-baseTime)/deltaTClock)
991 C o End Time 1
992 IF ( endTime .EQ. 0. .AND. nTimeSteps .NE. 0 )
993 & endTime = startTime + deltaTClock*DFLOAT(nTimeSteps)
994 C o End Time 2
995 IF ( endTime .EQ. 0. .AND. nEndIter .NE. 0 )
996 & endTime = baseTime + deltaTClock*DFLOAT(nEndIter)
997
998 C o Consistent?
999 IF ( startTime .NE. baseTime+deltaTClock*DFLOAT(nIter0) ) THEN
1000 WRITE(msgBuf,'(A)')
1001 & 'S/R INI_PARMS: startTime, baseTime and nIter0 are inconsistent'
1002 CALL PRINT_ERROR( msgBuf, myThid )
1003 WRITE(msgBuf,'(A)')
1004 & 'S/R INI_PARMS: Perhaps more than two were set at once'
1005 CALL PRINT_ERROR( msgBuf, myThid )
1006 errCount = errCount + 1
1007 ENDIF
1008 IF ( nEndIter .NE. nIter0+nTimeSteps ) THEN
1009 WRITE(msgBuf,'(A)')
1010 & 'S/R INI_PARMS: nIter0, nTimeSteps and nEndIter are inconsistent'
1011 CALL PRINT_ERROR( msgBuf, myThid )
1012 WRITE(msgBuf,'(A)')
1013 & 'S/R INI_PARMS: Perhaps more than two were set at once'
1014 CALL PRINT_ERROR( msgBuf, myThid )
1015 errCount = errCount + 1
1016 ENDIF
1017 IF ( nTimeSteps .NE. NINT((endTime-startTime)/deltaTClock)
1018 & ) THEN
1019 WRITE(msgBuf,'(A)')
1020 & 'S/R INI_PARMS: both endTime and nTimeSteps have been set'
1021 CALL PRINT_ERROR( msgBuf, myThid )
1022 WRITE(msgBuf,'(A)')
1023 & 'S/R INI_PARMS: but are inconsistent'
1024 CALL PRINT_ERROR( msgBuf, myThid )
1025 errCount = errCount + 1
1026 ENDIF
1027
1028 C o Monitor (should also add CPP flag for monitor?)
1029 IF (monitorFreq.LT.0.) THEN
1030 monitorFreq=0.
1031 IF (dumpFreq.NE.0.) monitorFreq=dumpFreq
1032 IF (diagFreq.NE.0..AND.diagFreq.LT.monitorFreq)
1033 & monitorFreq=diagFreq
1034 IF (taveFreq.NE.0..AND.taveFreq.LT.monitorFreq)
1035 & monitorFreq=taveFreq
1036 IF (chkPtFreq.NE.0..AND.chkPtFreq.LT.monitorFreq)
1037 & monitorFreq=chkPtFreq
1038 IF (pChkPtFreq.NE.0..AND.pChkPtFreq.LT.monitorFreq)
1039 & monitorFreq=pChkPtFreq
1040 IF (monitorFreq.EQ.0.) monitorFreq=deltaTClock
1041 ENDIF
1042 IF ( monitorSelect.EQ.UNSET_I ) THEN
1043 monitorSelect = 2
1044 IF ( fluidIsWater ) monitorSelect = 3
1045 ENDIF
1046
1047 C-- Grid parameters
1048 C In cartesian coords distances are in metres
1049 DO k =1,Nr
1050 delZ(k) = UNSET_RL
1051 delP(k) = UNSET_RL
1052 delR(k) = UNSET_RL
1053 ENDDO
1054 C In spherical polar distances are in degrees
1055 dxSpacing = UNSET_RL
1056 dySpacing = UNSET_RL
1057 WRITE(msgBuf,'(A)') ' INI_PARMS ; starts to read PARM04'
1058 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
1059 & SQUEEZE_RIGHT, myThid )
1060 READ(UNIT=iUnit,NML=PARM04,IOSTAT=errIO)
1061 IF ( errIO .LT. 0 ) THEN
1062 WRITE(msgBuf,'(A)')
1063 & 'S/R INI_PARMS: Error reading model parameter file "data"'
1064 CALL PRINT_ERROR( msgBuf, myThid )
1065 WRITE(msgBuf,'(A)') 'S/R INI_PARMS: Problem in namelist PARM04'
1066 CALL PRINT_ERROR( msgBuf, myThid )
1067 STOP 'ABNORMAL END: S/R INI_PARMS'
1068 ELSE
1069 WRITE(msgBuf,'(A)') ' INI_PARMS ; read PARM04 : OK'
1070 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
1071 & SQUEEZE_RIGHT, myThid )
1072 ENDIF
1073
1074 C Check for retired parameters still being used
1075 IF ( Ro_SeaLevel .NE. UNSET_RL ) THEN
1076 c nRetired = nRetired+1
1077 IF ( usingPCoords ) THEN
1078 WRITE(msgBuf,'(2A)') '** WARNING ** INI_PARMS: ',
1079 & '"Ro_SeaLevel" (P @ bottom) depreciated (backward compat'
1080 ELSEIF ( usingZCoords ) THEN
1081 WRITE(msgBuf,'(2A)') '** WARNING ** INI_PARMS: ',
1082 & '"Ro_SeaLevel" (Z @ top) depreciated (backward compat'
1083 ENDIF
1084 CALL PRINT_MESSAGE( msgBuf, errorMessageUnit,
1085 & SQUEEZE_RIGHT, myThid )
1086 IF ( usingPCoords ) THEN
1087 WRITE(msgBuf,'(2A)') '** WARNING ** INI_PARMS: ',
1088 & ' only). To set vert. axis, use instead "top_Pres".'
1089 ELSEIF ( usingZCoords ) THEN
1090 WRITE(msgBuf,'(2A)') '** WARNING ** INI_PARMS: ',
1091 & ' only). To set vert. axis, use instead "seaLev_Z".'
1092 ENDIF
1093 CALL PRINT_MESSAGE( msgBuf, errorMessageUnit,
1094 & SQUEEZE_RIGHT, myThid )
1095 ENDIF
1096 IF ( rkFac .NE. UNSET_RL ) THEN
1097 nRetired = nRetired+1
1098 WRITE(msgBuf,'(A,A)')
1099 & 'S/R INI_PARMS: "rkFac" has been replaced by -rkSign ',
1100 & ' and is no longer allowed in file "data".'
1101 CALL PRINT_ERROR( msgBuf, myThid )
1102 ENDIF
1103 IF ( groundAtK1 ) THEN
1104 c nRetired = nRetired+1
1105 WRITE(msgBuf,'(A,A)')
1106 & 'S/R INI_PARMS: "groundAtK1" is set according to vertical ',
1107 & ' coordinate and is no longer allowed in file "data".'
1108 CALL PRINT_ERROR( msgBuf, myThid )
1109 ENDIF
1110 IF ( thetaMin .NE. UNSET_RL ) THEN
1111 nRetired = nRetired+1
1112 WRITE(msgBuf,'(A,A)')
1113 & 'S/R INI_PARMS: "thetaMin" no longer allowed,',
1114 & ' has been replaced by "xgOrigin"'
1115 CALL PRINT_ERROR( msgBuf, myThid )
1116 ENDIF
1117 IF ( phiMin .NE. UNSET_RL ) THEN
1118 nRetired = nRetired+1
1119 WRITE(msgBuf,'(A,A)')
1120 & 'S/R INI_PARMS: "phiMin" no longer allowed,',
1121 & ' has been replaced by "ygOrigin"'
1122 CALL PRINT_ERROR( msgBuf, myThid )
1123 ENDIF
1124
1125 C X coordinate : Check for multiple definitions
1126 goptCount = 0
1127 IF ( delX(1) .NE. UNSET_RL ) goptCount = goptCount + 1
1128 IF ( dxSpacing .NE. UNSET_RL ) goptCount = goptCount + 1
1129 IF ( delXFile .NE. ' ' ) goptCount = goptCount + 1
1130 IF ( goptCount.GT.1 ) THEN
1131 WRITE(msgBuf,'(A,A)') 'Too many specifications for delX:',
1132 & 'Specify only one of delX, dxSpacing or delXfile'
1133 CALL PRINT_ERROR( msgBuf, myThid )
1134 errCount = errCount + 1
1135 ENDIF
1136 IF ( dxSpacing .NE. UNSET_RL ) THEN
1137 DO i=1,gridNx
1138 delX(i) = dxSpacing
1139 ENDDO
1140 ENDIF
1141 C Y coordinate : Check for multiple definitions
1142 goptCount = 0
1143 IF ( delY(1) .NE. UNSET_RL ) goptCount = goptCount + 1
1144 IF ( dySpacing .NE. UNSET_RL ) goptCount = goptCount + 1
1145 IF ( delYFile .NE. ' ' ) goptCount = goptCount + 1
1146 IF ( goptCount.GT.1 ) THEN
1147 WRITE(msgBuf,'(A,A)') 'Too many specifications for delY:',
1148 & 'Specify only one of delY, dySpacing or delYfile'
1149 CALL PRINT_ERROR( msgBuf, myThid )
1150 errCount = errCount + 1
1151 ENDIF
1152 IF ( dySpacing .NE. UNSET_RL ) THEN
1153 DO j=1,gridNy
1154 delY(j) = dySpacing
1155 ENDDO
1156 ENDIF
1157
1158 C-- Check for conflicting grid definitions.
1159 goptCount = 0
1160 IF ( usingCartesianGrid ) goptCount = goptCount+1
1161 IF ( usingSphericalPolarGrid ) goptCount = goptCount+1
1162 IF ( usingCurvilinearGrid ) goptCount = goptCount+1
1163 IF ( usingCylindricalGrid ) goptCount = goptCount+1
1164 IF ( goptCount .GT. 1 ) THEN
1165 WRITE(msgBuf,'(A)')
1166 & 'S/R INI_PARMS: More than one coordinate system requested'
1167 CALL PRINT_ERROR( msgBuf, myThid )
1168 errCount = errCount + 1
1169 ENDIF
1170 IF ( goptCount .LT. 1 ) THEN
1171 C- No horizontal grid is specified => use Cartesian grid as default:
1172 WRITE(msgBuf,'(A)')
1173 & 'S/R INI_PARMS: No horizontal grid requested'
1174 CALL PRINT_MESSAGE( msgBuf, errorMessageUnit,
1175 & SQUEEZE_RIGHT, myThid )
1176 WRITE(msgBuf,'(A)')
1177 & 'S/R INI_PARMS: => Use Cartesian Grid as default'
1178 CALL PRINT_MESSAGE( msgBuf, errorMessageUnit,
1179 & SQUEEZE_RIGHT, myThid )
1180 usingCartesianGrid = .TRUE.
1181 ENDIF
1182 C- Radius of the Planet:
1183 IF ( rSphere.EQ.UNSET_RL ) THEN
1184 IF ( usingCurvilinearGrid .AND.
1185 & radius_fromHorizGrid.NE.UNSET_RL ) THEN
1186 rSphere = radius_fromHorizGrid
1187 ELSE
1188 rSphere = 6370. _d 3
1189 ENDIF
1190 ENDIF
1191 IF ( radius_fromHorizGrid.EQ.UNSET_RL ) THEN
1192 radius_fromHorizGrid = rSphere
1193 ENDIF
1194 IF ( rSphere .NE. 0. ) THEN
1195 recip_rSphere = 1. _d 0/rSphere
1196 ELSE
1197 recip_rSphere = 0.
1198 ENDIF
1199 C-- Default vertical axis origin
1200 IF ( Ro_SeaLevel .NE. UNSET_RL ) THEN
1201 IF ( usingPCoords .AND. top_Pres.NE.UNSET_RL ) THEN
1202 WRITE(msgBuf,'(2A)') 'S/R INI_PARMS: ',
1203 & 'Cannot set both "Ro_SeaLevel" and "top_Pres"'
1204 CALL PRINT_ERROR( msgBuf, myThid )
1205 errCount = errCount + 1
1206 ENDIF
1207 IF ( usingZCoords .AND. seaLev_Z.NE.UNSET_RL ) THEN
1208 WRITE(msgBuf,'(2A)') 'S/R INI_PARMS: ',
1209 & 'Cannot set both "Ro_SeaLevel" and "seaLev_Z"'
1210 CALL PRINT_ERROR( msgBuf, myThid )
1211 errCount = errCount + 1
1212 ENDIF
1213 rF(1) = Ro_SeaLevel
1214 ELSE
1215 rF(1) = UNSET_RS
1216 ENDIF
1217 IF ( top_Pres.EQ.UNSET_RL ) top_Pres = 0.
1218 IF ( seaLev_Z.EQ.UNSET_RL ) seaLev_Z = 0.
1219 C-- Default origin for X & Y axis (longitude & latitude origin):
1220 IF ( xgOrigin .EQ. UNSET_RL ) xgOrigin = 0.
1221 IF ( ygOrigin .EQ. UNSET_RL ) THEN
1222 IF ( usingSphericalPolarGrid ) THEN
1223 ygOrigin = 0.
1224 ELSEIF ( usingCartesianGrid ) THEN
1225 ygOrigin = 0.
1226 ELSEIF ( usingCylindricalGrid ) THEN
1227 ygOrigin = 0.
1228 ELSEIF ( usingCurvilinearGrid ) THEN
1229 ygOrigin = 0.
1230 ELSE
1231 WRITE(msgBuf,'(A)')
1232 & 'S/R INI_PARMS: found no coordinate system to set ygOrigin'
1233 CALL PRINT_ERROR( msgBuf, myThid )
1234 errCount = errCount + 1
1235 ENDIF
1236 ENDIF
1237 C-- Make metric term & Coriolis settings consistent with underlying grid.
1238 IF ( usingCartesianGrid ) THEN
1239 metricTerms = .FALSE.
1240 useNHMTerms = .FALSE.
1241 ENDIF
1242 IF ( usingCylindricalGrid ) THEN
1243 useNHMTerms = .FALSE.
1244 WRITE(msgBuf,'(A)') 'S/R INI_PARMS ; Cylinder OK'
1245 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
1246 & SQUEEZE_RIGHT, myThid )
1247 ENDIF
1248 IF ( usingCurvilinearGrid ) THEN
1249 metricTerms = .FALSE.
1250 ENDIF
1251 IF ( selectCoriMap.EQ.-1 ) THEN
1252 IF ( usingCartesianGrid.OR.usingCylindricalGrid ) THEN
1253 C default is to use Beta-Plane Coriolis
1254 selectCoriMap = 1
1255 ELSE
1256 C default for other grids is to use Spherical Coriolis
1257 selectCoriMap = 2
1258 ENDIF
1259 ENDIF
1260 IF ( .NOT.(nonHydrostatic.OR.quasiHydrostatic) )
1261 & use3dCoriolis = .FALSE.
1262 IF ( (selectCoriMap.EQ.0 .OR.selectCoriMap.EQ.1)
1263 & .AND. fPrime.EQ.0. ) use3dCoriolis = .FALSE.
1264
1265 C-- Grid rotation
1266 IF ( phiEuler .NE. 0. _d 0 .OR. thetaEuler .NE. 0. _d 0
1267 & .OR. psiEuler .NE. 0. _d 0 ) rotateGrid = .TRUE.
1268
1269 C-- Set default for latitude-band where relaxation to climatology applies
1270 C note: done later (once domain size is known) if using CartesianGrid
1271 IF ( latBandClimRelax .EQ. UNSET_RL) THEN
1272 IF ( usingSphericalPolarGrid ) latBandClimRelax= 180. _d 0
1273 IF ( usingCurvilinearGrid ) latBandClimRelax= 180. _d 0
1274 ENDIF
1275
1276 C-- set cell Center depth and put Interface at the middle between 2 C
1277 setCenterDr = .FALSE.
1278 DO k=1,Nr+1
1279 IF ( delRc(k).EQ.UNSET_RL ) THEN
1280 IF ( setCenterDr ) THEN
1281 WRITE(msgBuf,'(A,I4)')
1282 & 'S/R INI_PARMS: No value for delRc at k =', k
1283 CALL PRINT_ERROR( msgBuf, myThid )
1284 errCount = errCount + 1
1285 ENDIF
1286 ELSE
1287 IF ( k.EQ.1 ) setCenterDr = .TRUE.
1288 IF ( .NOT.setCenterDr ) THEN
1289 WRITE(msgBuf,'(A,I4)')
1290 & 'S/R INI_PARMS: No value for delRc at k <', k
1291 CALL PRINT_ERROR( msgBuf, myThid )
1292 errCount = errCount + 1
1293 ENDIF
1294 ENDIF
1295 ENDDO
1296 IF ( setCenterDr ) rCoordInputData = .TRUE.
1297 C-- p, z, r coord parameters
1298 setInterFDr = .FALSE.
1299 DO k = 1, Nr
1300 IF ( delZ(k) .NE. UNSET_RL ) zCoordInputData = .TRUE.
1301 IF ( delP(k) .NE. UNSET_RL ) pCoordInputData = .TRUE.
1302 IF ( delR(k) .NE. UNSET_RL ) rCoordInputData = .TRUE.
1303 IF ( delR(k) .EQ. UNSET_RL ) delR(k) = delZ(k)
1304 IF ( delR(k) .EQ. UNSET_RL ) delR(k) = delP(k)
1305 IF ( delR(k) .EQ. UNSET_RL ) THEN
1306 IF ( setInterFDr ) THEN
1307 WRITE(msgBuf,'(A,I4)')
1308 & 'S/R INI_PARMS: No value for delZ/delP/delR at k =', k
1309 CALL PRINT_ERROR( msgBuf, myThid )
1310 errCount = errCount + 1
1311 ENDIF
1312 ELSE
1313 IF ( k.EQ.1 ) setInterFDr = .TRUE.
1314 IF ( .NOT.setInterFDr ) THEN
1315 WRITE(msgBuf,'(A,I4)')
1316 & 'S/R INI_PARMS: No value for delZ/delP/delR at k <', k
1317 CALL PRINT_ERROR( msgBuf, myThid )
1318 errCount = errCount + 1
1319 ENDIF
1320 ENDIF
1321 ENDDO
1322 C Check for multiple coordinate systems
1323 coordsSet = 0
1324 IF ( zCoordInputData ) coordsSet = coordsSet + 1
1325 IF ( pCoordInputData ) coordsSet = coordsSet + 1
1326 IF ( rCoordInputData ) coordsSet = coordsSet + 1
1327 IF ( coordsSet .GT. 1 ) THEN
1328 WRITE(msgBuf,'(A)')
1329 & 'S/R INI_PARMS: Cannot mix z, p and r in the input data.'
1330 CALL PRINT_ERROR( msgBuf, myThid )
1331 errCount = errCount + 1
1332 ENDIF
1333 C- Check for double definition (file & namelist)
1334 IF ( delRcFile.NE.' ' ) THEN
1335 IF ( setCenterDr ) THEN
1336 WRITE(msgBuf,'(A)')
1337 & 'S/R INI_PARMS: Cannot set both delRc and delRcFile'
1338 CALL PRINT_ERROR( msgBuf, myThid )
1339 errCount = errCount + 1
1340 ENDIF
1341 setCenterDr = .TRUE.
1342 ENDIF
1343 IF ( delRFile.NE.' ' ) THEN
1344 IF ( setInterFDr ) THEN
1345 WRITE(msgBuf,'(A)')
1346 & 'S/R INI_PARMS: Cannot set both delR and delRFile'
1347 CALL PRINT_ERROR( msgBuf, myThid )
1348 errCount = errCount + 1
1349 ENDIF
1350 setInterFDr = .TRUE.
1351 ENDIF
1352
1353 C-- Input files
1354 WRITE(msgBuf,'(A)') ' INI_PARMS ; starts to read PARM05'
1355 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
1356 & SQUEEZE_RIGHT, myThid )
1357 READ(UNIT=iUnit,NML=PARM05) !,IOSTAT=errIO)
1358 IF ( errIO .LT. 0 ) THEN
1359 WRITE(msgBuf,'(A)')
1360 & 'S/R INI_PARMS: Error reading model parameter file "data"'
1361 CALL PRINT_ERROR( msgBuf, myThid )
1362 WRITE(msgBuf,'(A)') 'S/R INI_PARMS: Problem in namelist PARM05'
1363 CALL PRINT_ERROR( msgBuf, myThid )
1364 STOP 'ABNORMAL END: S/R INI_PARMS'
1365 ELSE
1366 WRITE(msgBuf,'(A)') ' INI_PARMS ; read PARM05 : OK'
1367 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
1368 & SQUEEZE_RIGHT, myThid )
1369 ENDIF
1370 C Check for retired parameters still being used
1371 IF ( shelfIceFile .NE. ' ' ) THEN
1372 nRetired = nRetired+1
1373 WRITE(msgBuf,'(A,A)')
1374 & 'S/R INI_PARMS: "shelfIceFile" is not allowed in "data", ',
1375 & 'substitute "SHELFICEtopoFile" in data.shelfice'
1376 CALL PRINT_ERROR( msgBuf, myThid )
1377 ENDIF
1378 IF ( dQdTFile .NE. ' ' ) THEN
1379 nRetired = nRetired+1
1380 WRITE(msgBuf,'(A,A)')
1381 & 'S/R INI_PARMS: "dQdTFile" has been retired from file "data"'
1382 CALL PRINT_ERROR( msgBuf, myThid )
1383 ENDIF
1384 C Check if two conflicting I/O directories have been specified
1385 IF (mdsioLocalDir .NE. ' ' .AND. adTapeDir .NE. ' ') THEN
1386 WRITE(msgBuf,'(A)')
1387 & 'S/R INI_PARMS: mdsioLocalDir and adTapeDir cannot be'
1388 CALL PRINT_ERROR( msgBuf, myThid )
1389 WRITE(msgBuf,'(A)')
1390 & 'S/R INI_PARMS: specified at the same time'
1391 CALL PRINT_ERROR( msgBuf, myThid )
1392 errCount = errCount + 1
1393 ENDIF
1394
1395 C Check for relaxation term:
1396 IF ( (tauThetaClimRelax.GT.0.).AND.
1397 & (thetaClimFile.EQ.' ') ) THEN
1398 WRITE(msgBuf,'(A)')
1399 & 'S/R INI_PARMS: tauThetaClimRelax > 0 but'
1400 CALL PRINT_ERROR( msgBuf, myThid )
1401 WRITE(msgBuf,'(A)')
1402 & 'S/R INI_PARMS: thetaClimFile is undefined'
1403 CALL PRINT_ERROR( msgBuf, myThid )
1404 errCount = errCount + 1
1405 ENDIF
1406 IF ( (tauSaltClimRelax.GT.0.).AND.
1407 & (saltClimFile.EQ.' ') ) THEN
1408 WRITE(msgBuf,'(A)')
1409 & 'S/R INI_PARMS: tauSaltClimRelax > 0 but'
1410 CALL PRINT_ERROR( msgBuf, myThid )
1411 WRITE(msgBuf,'(A)')
1412 & 'S/R INI_PARMS: saltClimFile is undefined'
1413 CALL PRINT_ERROR( msgBuf, myThid )
1414 errCount = errCount + 1
1415 ENDIF
1416 C Check vertical diffusivity setting:
1417 #ifdef ALLOW_3D_DIFFKR
1418 IF ( specifiedDiffKrT ) THEN
1419 WRITE(msgBuf,'(2A)') '** WARNING ** INI_PARMS: Ignores diffKr',
1420 & 'T (or Kp,Kz) setting in file "data" with ALLOW_3D_DIFFKR'
1421 CALL PRINT_MESSAGE( msgBuf, errorMessageUnit,
1422 & SQUEEZE_RIGHT, myThid )
1423 ENDIF
1424 IF ( specifiedDiffKrS .AND. diffKrFile.NE.' ' ) THEN
1425 WRITE(msgBuf,'(2A)') '** WARNING ** INI_PARMS: Ignores diffKr',
1426 & 'S (or Kp,Kz) setting in file "data" and uses diffKrFile'
1427 CALL PRINT_MESSAGE( msgBuf, errorMessageUnit,
1428 & SQUEEZE_RIGHT, myThid )
1429 ENDIF
1430 #endif
1431
1432 C-- Set Units conversion factor required to incorporate
1433 C surface forcing into z-p isomorphic equations:
1434 C mass2rUnit: from mass per unit area [kg/m2] to r-coordinate (z:=1/rho;p:=g)
1435 C rUnit2mass: from r-coordinate to mass per unit area [kg/m2] (z:=rho;p:=1/g)
1436 IF ( usingPCoords ) THEN
1437 mass2rUnit = gravity
1438 rUnit2mass = recip_gravity
1439 ELSE
1440 mass2rUnit = recip_rhoConst
1441 rUnit2mass = rhoConst
1442 ENDIF
1443
1444 C-- For backward compatibility, set temp_addMass and salt_addMass
1445 C to temp_EvPrRn and salt_EvPrRn if not set in parameter file "data"
1446 IF (temp_addMass .EQ. UNSET_RL) temp_addMass = temp_EvPrRn
1447 IF (salt_addMass .EQ. UNSET_RL) salt_addMass = salt_EvPrRn
1448
1449 C-- Make a choice (if unset) regarding using CG2D solver min.-residual sol.
1450 C for simple set-up (cartesian grid + flat bottom), default is to
1451 C use the solver minimum residual solution (cg2dUseMinResSol=1).
1452 IF ( cg2dUseMinResSol.EQ.UNSET_I ) THEN
1453 cg2dUseMinResSol = 0
1454 IF ( topoFile.EQ.' ' .AND. bathyFile.EQ.' '
1455 & .AND. usingCartesianGrid ) cg2dUseMinResSol = 1
1456 ENDIF
1457
1458 C-- Close the open data file
1459 CLOSE(iUnit)
1460
1461 WRITE(msgBuf,'(A)') ' INI_PARMS: finished reading file "data"'
1462 CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
1463 & SQUEEZE_RIGHT , myThid )
1464
1465 C-- Check whether any retired parameters were found.
1466 IF ( nRetired .GT. 0 ) THEN
1467 WRITE(msgBuf,'(A)')
1468 & 'S/R INI_PARMS: Error reading parameter file "data":'
1469 CALL PRINT_ERROR( msgBuf, myThid )
1470 WRITE(msgBuf,'(I4,A)') nRetired,
1471 & ' out of date parameters were found in the namelist'
1472 CALL PRINT_ERROR( msgBuf, myThid )
1473 errCount = errCount + 1
1474 ENDIF
1475 C-- Stop if any error was found (including retired params):
1476 IF ( errCount .GE. 1 ) THEN
1477 WRITE(msgBuf,'(A,I3,A)')
1478 & 'S/R INI_PARMS: detected', errCount,' fatal error(s)'
1479 CALL PRINT_ERROR( msgBuf, myThid )
1480 CALL ALL_PROC_DIE( 0 )
1481 STOP 'ABNORMAL END: S/R INI_PARMS'
1482 ENDIF
1483
1484 _END_MASTER(myThid)
1485
1486 C-- Everyone else must wait for the parameters to be loaded
1487 _BARRIER
1488
1489 RETURN
1490 END

  ViewVC Help
Powered by ViewVC 1.1.22