Back to home page

MITgcm

 
 

    


File indexing completed on 2025-12-03 06:08:09 UTC

view on githubraw file Latest commit 47ecd102 on 2025-12-02 23:05:14 UTC
b9b591469d Jean*0001 #include "PACKAGES_CONFIG.h"
b9377e56be Alis*0002 #include "CPP_OPTIONS.h"
b9b591469d Jean*0003 #ifdef ALLOW_EXCH2
                0004 # include "W2_OPTIONS.h"
                0005 #endif /* ALLOW_EXCH2 */
924557e60a Chri*0006 
a30418b6b9 Ed H*0007 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
9366854e02 Chri*0008 CBOP
427e24e121 Jean*0009 C !ROUTINE: INI_PARMS
                0010 C !INTERFACE:
924557e60a Chri*0011       SUBROUTINE INI_PARMS( myThid )
                0012 
427e24e121 Jean*0013 C !DESCRIPTION:
809bdccbfc Jean*0014 C     Routine to load model "parameters" from parameter file "data"
                0015 
9366854e02 Chri*0016 C     !USES:
                0017       IMPLICIT NONE
924557e60a Chri*0018 #include "SIZE.h"
                0019 #include "EEPARAMS.h"
                0020 #include "PARAMS.h"
b9b591469d Jean*0021 #ifdef ALLOW_EXCH2
                0022 # include "W2_EXCH2_SIZE.h"
                0023 # include "W2_EXCH2_TOPOLOGY.h"
                0024 #endif /* ALLOW_EXCH2 */
a37a13034c Mart*0025 #include "EOS.h"
97168926d1 Jean*0026 C- need GRID.h to save, in rF(1), half retired param Ro_SeaLevel value:
2ad79bdf32 Jean*0027 #include "GRID.h"
97168926d1 Jean*0028 #include "SET_GRID.h"
924557e60a Chri*0029 
427e24e121 Jean*0030 C !INPUT/OUTPUT PARAMETERS:
                0031 C   myThid :: Number of this instance of INI_PARMS
924557e60a Chri*0032       INTEGER myThid
                0033 
427e24e121 Jean*0034 C !LOCAL VARIABLES:
                0035 C   dxSpacing, dySpacing :: Default spacing in X and Y.
                0036 C                        :: Units are that of coordinate system
                0037 C                        :: i.e. cartesian => metres ; s. polar => degrees
                0038 C   deltaTtracer  :: Timestep for tracer equations ( s )
                0039 C   forcing_In_AB :: flag to put all forcings (Temp,Salt,Tracers,Momentum)
                0040 C                 :: contribution in (or out of) Adams-Bashforth time stepping.
                0041 C   goptCount  :: Used to count the nuber of grid options (only one is allowed!)
                0042 C   msgBuf     :: Informational/error message buffer
                0043 C   errIO      :: IO error flag
                0044 C   errCount   :: error counter (inconsitent params and other errors)
                0045 C   iUnit      :: Work variable for IO unit number
                0046 C   k, i, j    :: Loop counters
                0047 C   xxxDefault :: Default value for variable xxx
910f05e765 Chri*0048       _RL  dxSpacing
                0049       _RL  dySpacing
062a876ce5 Jean*0050       _RL deltaTtracer
924557e60a Chri*0051       CHARACTER*(MAX_LEN_MBUF) msgBuf
c2b6ed6bfd Jean*0052       LOGICAL forcing_In_AB
924557e60a Chri*0053       INTEGER goptCount
b9b591469d Jean*0054       INTEGER gridNx, gridNy
aa8f58c33d Jean*0055       INTEGER k, i, j, iUnit
90aed06079 Jean*0056       INTEGER errIO, errCount
427e24e121 Jean*0057 C   Default values for variables which have vertical coordinate system
                0058 C   dependency.
910f05e765 Chri*0059       _RL viscArDefault
                0060       _RL diffKrTDefault
                0061       _RL diffKrSDefault
                0062       _RL hFacMinDrDefault
0127add478 Alis*0063       _RL delRDefault(Nr)
427e24e121 Jean*0064 C-- Variables used to select between different coordinate systems:
                0065 C   zCoordInputData :: The vertical coordinate system in the rest of the  model
e6e223b277 Jean*0066 C   pCoordInputData :: is written in terms of r. In the model "data" file,
427e24e121 Jean*0067 C   rCoordInputData :: input params can be in terms of z, p or r.
                0068 C   coordsSet       :: e.g. delZ or delP or delR for the vertical grid spacing.
                0069 C                   :: The following rules apply:
                0070 C                   :: All parameters must use the same vertical coordinate
                0071 C                   :: system.  e.g. delZ and viscAz is legal
                0072 C                   ::           but delZ and viscAr is an error.
                0073 C                   :: Similarly specifying delZ and delP is an error.
                0074 C                   :: zCoord..., pCoord..., rCoord... are used to flag when
                0075 C                   :: z, p or r are used.
                0076 C                   :: coordsSet counts how many vertical coordinates have been
                0077 C                   :: used to specify variables. coordsSet > 1 is an error.
                0078 C   vertSetCount    :: to count number of vertical array elements which are set.
                0079 C   specifiedDiffKrT :: internal flag, true if any diffK[r,z,p,Nr]T is specified
                0080 C   specifiedDiffKrS :: internal flag, true if any diffK[r,z,p,Nr]S is specified
53f30fd0a1 Jean*0081 
910f05e765 Chri*0082       LOGICAL zCoordInputData
                0083       LOGICAL pCoordInputData
                0084       LOGICAL rCoordInputData
                0085       INTEGER coordsSet
7cd7754949 Jean*0086       INTEGER vertSetCount
13b4d13c98 Jean*0087       LOGICAL specifiedDiffKrT, specifiedDiffKrS
ede2cfe569 Chri*0088 
427e24e121 Jean*0089 C-- Variables which have vertical coordinate system dependency.
                0090 C   delZ    :: Vertical grid spacing ( m  ).
                0091 C   delP    :: Vertical grid spacing ( Pa ).
                0092 C   viscAz  :: Eddy viscosity coeff. for mixing of momentum vertically (m^2/s)
                0093 C   viscAp  :: Eddy viscosity coeff. for mixing of momentum vertically (Pa^2/s)
                0094 C   diffKzT :: Laplacian diffusion coeff. for mixing of heat vertically (m^2/s)
                0095 C   diffKpT :: Laplacian diffusion coeff. for mixing of heat vertically (Pa^2/s)
                0096 C   diffKzS :: Laplacian diffusion coeff. for mixing of salt vertically (m^2/s)
                0097 C   diffKpS :: Laplacian diffusion coeff. for mixing of salt vertically (Pa^2/s)
53f30fd0a1 Jean*0098       _RL delZ(Nr)
                0099       _RL delP(Nr)
                0100       _RL viscAz
                0101       _RL viscAp
7cd7754949 Jean*0102       _RL viscAr
53f30fd0a1 Jean*0103       _RL diffKzT
                0104       _RL diffKpT
8fc5af62b8 Jean*0105       _RL diffKrT
53f30fd0a1 Jean*0106       _RL diffKzS
                0107       _RL diffKpS
8fc5af62b8 Jean*0108       _RL diffKrS
53f30fd0a1 Jean*0109 
427e24e121 Jean*0110 C-- Retired main data file parameters. Kept here to trap use of old data files.
                0111 C   nRetired :: Counter used to trap gracefully namelists containing "retired"
                0112 C            :: parameters. These are parameters that are either no-longer used
                0113 C            :: or that have moved to a different input file and/or namelist.
742cf4499c Jean*0114 C Namelist PARM01:
427e24e121 Jean*0115 C   useConstantF  :: Coriolis coeff set to f0     (replaced by selectCoriMap=0)
                0116 C   useBetaPlaneF :: Coriolis coeff = f0 + beta.y (replaced by selectCoriMap=1)
                0117 C   useSphereF    :: Coriolis = 2.omega.sin(phi)  (replaced by selectCoriMap=2)
                0118 C   useJamartWetPoints :: for backward compat. (replaced by selectCoriScheme=1)
                0119 C                      :: Use wet-point method for Coriolis (Jamart & Ozer 1986)
                0120 C   SadournyCoriolis :: for backward compat. (replaced by selectVortScheme=2)
                0121 C   use3dCoriolis :: for backward compat. (=F,T same as select3dCoriScheme=0,1)
                0122 C   metricTerms   :: for backward compat. (=F,T same as selectMetricTerms=0,1)
                0123 C   tracerAdvScheme :: tracer advection scheme (old passive tracer code)
                0124 C   trac_EvPrRn :: tracer conc. in Rain & Evap (old passive tracer code)
                0125 C   saltDiffusion :: diffusion of salinity    on/off (flag not used)
                0126 C   tempDiffusion :: diffusion of temperature on/off (flag not used)
                0127 C   zonal_filt_lat :: Moved to package "zonal_filt"
                0128 C   gravitySign  :: direction of gravity relative to R direction
                0129 C                :: (removed from namelist and set according to z/p coordinate)
                0130 C   viscAstrain  :: replaced by standard viscosity coeff & useStrainTensionVisc
                0131 C   viscAtension :: replaced by standard viscosity coeff & useStrainTensionVisc
                0132 C   useAnisotropicViscAgridMax :: Changed to be default behavior.  Can
                0133 C                   use old method by setting useAreaViscLength=.true.
                0134 C   usePickupBeforeC35 :: to restart from old-pickup files (generated with code
                0135 C                 from before checkpoint-35, Feb 08, 2001): disabled (Jan 2007)
                0136 C   debugMode    :: to print debug msg. now read from parameter file eedata
                0137 C   allowInteriorFreezing :: Allow water at depth to freeze
                0138 C                            and rise to the surface (replaced by pkg/frazil)
                0139 C   useOldFreezing :: use the old version (before checkpoint52a_pre, 2003-11-12)
                0140 C   balanceEmPmR   :: for backward compat. (replaced by selectBalanceEmPmR=1),
                0141 C                  :: substract global mean of EmPmR at every time step
                0142 C   writeStatePrec :: Precision used for writing model state (not used)
aecc8b0f47 Mart*0143 C Namelist PARM02:
427e24e121 Jean*0144 C   cg2dChkResFreq :: Frequency with which to check 2-D con. grad solver
                0145 C                       residual (was never coded)
                0146 C   cg3dChkResFreq :: Frequency with which to check 3-D con. grad solver
                0147 C                       residual (was never coded)
742cf4499c Jean*0148 C Namelist PARM03:
427e24e121 Jean*0149 C   tauThetaClimRelax3Dim :: replaced by pkg/rbcs (3.D Relaxation B.Cs)
                0150 C   tauSaltClimRelax3Dim  :: replaced by pkg/rbcs (3.D Relaxation B.Cs)
b7411f1a84 Jean*0151 C   taveFreq       :: Frequency output for time-averaged state vars ;
                0152 C                     disabled as pkg/timeave has been removed
                0153 C   tave_lastIter  :: last (within taveFreq period) time-step fraction in time-
                0154 C                     averaded output ; disabled as pkg/timeave has been removed
427e24e121 Jean*0155 C   calendarDumps  :: moved to package "cal" (calendar)
742cf4499c Jean*0156 C Namelist PARM04:
427e24e121 Jean*0157 C   Ro_SeaLevel :: origin of the vertical R-coords axis ;
                0158 C               :: replaced by top_Pres or seaLev_Z setting
                0159 C   groundAtK1 :: put the surface(k=1) at the ground (replaced by usingPCoords)
                0160 C   rkFac      :: removed from namelist ; replaced by -rkSign
                0161 C   thetaMin   :: unfortunate variable name ; replaced by xgOrigin
                0162 C   phiMin     :: unfortunate variable name ; replaced by ygOrigin
8194bc4a99 Mart*0163 C Namelist PARM05:
427e24e121 Jean*0164 C   shelfIceFile :: File containing the topography of the shelfice draught
                0165 C                   (replaced by SHELFICEtopoFile in SHELFICE.h)
                0166 C   dQdTfile     :: File containing thermal relaxation coefficient
                0167 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
8194bc4a99 Mart*0168 
9780090eaa Jean*0169       INTEGER nRetired
c07cd3bfa8 Jean*0170       LOGICAL useConstantF, useBetaPlaneF, useSphereF
427e24e121 Jean*0171       LOGICAL useJamartWetPoints, useEnergyConservingCoriolis
7e00d7e8f9 Jean*0172       LOGICAL SadournyCoriolis
427e24e121 Jean*0173       LOGICAL use3dCoriolis, metricTerms
7a4ec66481 Jean*0174       LOGICAL tempDiffusion, saltDiffusion
                0175       INTEGER tracerAdvScheme
                0176       _RL trac_EvPrRn
2ad79bdf32 Jean*0177       _RL zonal_filt_lat
                0178 c     _RL gravitySign
8f934af0f4 Jean*0179       _RL viscAstrain, viscAtension
6d61b52b09 Jean*0180       LOGICAL useAnisotropicViscAgridMax
                0181       LOGICAL usePickupBeforeC35
83513ae28b Jean*0182       LOGICAL saveDebugMode
e793395f17 Jean*0183       LOGICAL allowInteriorFreezing, useOldFreezing
7e00d7e8f9 Jean*0184       LOGICAL balanceEmPmR
427e24e121 Jean*0185       INTEGER writeStatePrec
742cf4499c Jean*0186 C-
aecc8b0f47 Mart*0187       INTEGER cg2dChkResFreq, cg3dChkResFreq
                0188 C-
46979d29da Jean*0189       _RL tauThetaClimRelax3Dim, tauSaltClimRelax3Dim
b7411f1a84 Jean*0190       _RL taveFreq, tave_lastIter
742cf4499c Jean*0191       LOGICAL calendarDumps
                0192 C-
                0193       LOGICAL groundAtK1
2ad79bdf32 Jean*0194       _RL Ro_SeaLevel
742cf4499c Jean*0195       _RL rkFac
9780090eaa Jean*0196       _RL thetaMin, phiMin
7e00d7e8f9 Jean*0197 C-
8194bc4a99 Mart*0198       CHARACTER*(MAX_LEN_FNAM) shelfIceFile
275ad78b78 Jean*0199       CHARACTER*(MAX_LEN_FNAM) dQdTfile
924557e60a Chri*0200 
                0201 C--   Continuous equation parameters
                0202       NAMELIST /PARM01/
fe4f35d666 Alis*0203      & gravitySign, nh_Am2,
7fbd8d8c1c Jean*0204      & gravity, gBaro, gravityFile, rhonil, tAlpha, sBeta,
c07cd3bfa8 Jean*0205      & selectCoriMap, f0, beta, fPrime, omega, rotationPeriod,
8f106aecfd Bayl*0206      & viscAh, viscAhW, viscAhMax,
                0207      & viscAhGrid, viscAhGridMax, viscAhGridMin,
dc3adfb09b Jean*0208      & viscC2leith, viscC4leith, smag3D_coeff, useSmag3D,
6d61b52b09 Jean*0209      & useFullLeith, useAnisotropicViscAgridMax, useStrainTensionVisc,
f59d76b0dd Ed D*0210      & useAreaViscLength, viscC2leithD, viscC4leithD, viscC2LeithQG,
                0211      & viscC2smag, viscC4smag, viscAhD, viscAhZ, viscA4D, viscA4Z,
0b6cbae535 Mart*0212      & viscA4, viscA4W,
                0213      & viscA4Max, viscA4Grid, viscA4GridMax, viscA4GridMin,
1e9e8678dd Bayl*0214      & viscA4ReMax, viscAhReMax,
7cd7754949 Jean*0215      & cosPower, viscAstrain, viscAtension,
427e24e121 Jean*0216      & viscAz, diffKzT, diffKzS, viscAp, diffKpT, diffKpS,
                0217      & viscAr, diffKrT, diffKrS, viscArNr, diffKrNrT, diffKrNrS,
2d5bb917cc Jean*0218      & diffKhT, diffK4T, diffKhS, diffK4S, smag3D_diffCoeff,
427e24e121 Jean*0219      & diffKr4T, diffKr4S, BL79LatVary,
                0220      & diffKrBL79surf, diffKrBL79deep, diffKrBL79scl, diffKrBL79Ho,
                0221      & diffKrBLEQsurf, diffKrBLEQdeep, diffKrBLEQscl, diffKrBLEQHo,
9330bd8273 Jean*0222      & surf_pRef, tRef, sRef, tRefFile, sRefFile, rhoRefFile,
7fbd8d8c1c Jean*0223      & eosType, selectP_inEOS_Zc, integr_GeoPot, selectFindRoSurf,
aad87ebb4c Jean*0224      & HeatCapacity_Cp, celsius2K, atm_Cp, atm_Rd, atm_Rq, atm_Po,
59b35dd864 Jean*0225      & no_slip_sides, sideDragFactor, no_slip_bottom, bottomVisc_pCell,
ab47de63dc Mart*0226      & bottomDragLinear, bottomDragQuadratic,  zRoughBot,
                0227      & selectBotDragQuadr,
427e24e121 Jean*0228      & momPressureForcing, momForcing, momTidalForcing,
                0229      & momViscosity, momAdvection, vectorInvariantMomentum,
                0230      & useConstantF, useBetaPlaneF, useSphereF, useCoriolis,
                0231      & use3dCoriolis, select3dCoriScheme, selectCoriScheme,
                0232      & useCDscheme, useJamartWetPoints, useEnergyConservingCoriolis,
                0233      & useAbsVorticity, selectVortScheme, SadournyCoriolis,
                0234      & useJamartMomAdv, upwindVorticity, highOrderVorticity,
                0235      & upwindShear, selectKEscheme,
                0236      & selectMetricTerms, metricTerms, useNHMTerms, addFrictionHeating,
46918f1b26 Jean*0237      & tempDiffusion, tempAdvection, tempForcing, temp_stayPositive,
                0238      & saltDiffusion, saltAdvection, saltForcing, salt_stayPositive,
9e44938df7 Jean*0239      & implicSurfPress, implicDiv2DFlow, implicitNHPress,
7cd7754949 Jean*0240      & implicitFreeSurface, rigidLid, freeSurfFac,
                0241      & hFacMin, hFacMinDz, hFacMinDp, hFacMinDr,
1e6181f584 Davi*0242      & exactConserv, linFSConserveTr, uniformLin_PhiSurf,
                0243      & nonlinFreeSurf, hFacInf, hFacSup, select_rStar,
3fcd8a21e5 Jean*0244      & nonHydrostatic, selectNHfreeSurf, quasiHydrostatic,
627a21a418 Jean*0245      & implicitIntGravWave, staggerTimeStep, doResetHFactors,
8a28092f34 Patr*0246      & tempStepping, saltStepping, momStepping,
9e44938df7 Jean*0247      & implicitDiffusion, implicitViscosity, selectImplicitDrag,
b9d351b225 Jean*0248      & tempImplVertAdv, saltImplVertAdv, momImplVertAdv,
aad87ebb4c Jean*0249      & rhoConst, thetaConst, rhoConstFresh, buoyancyRelation,
7c642b2a2b Dimi*0250      & allowFreezing, allowInteriorFreezing, useOldFreezing, ivdc_kappa,
d1b0368d70 Davi*0251      & hMixCriteria, dRhoSmall, hMixSmooth,
427e24e121 Jean*0252      & tempAdvScheme, tempVertAdvScheme, tracerAdvScheme,
                0253      & saltAdvScheme, saltVertAdvScheme, multiDimAdvection,
d18df35fee Jean*0254      & selectAddFluid, useRealFreshWaterFlux, convertFW2Salt,
745a0098ab Jean*0255      & temp_EvPrRn, salt_EvPrRn, trac_EvPrRn,
80d98e0151 Dimi*0256      & temp_addMass, salt_addMass, zonal_filt_lat,
00c7090dc0 Mart*0257      & smoothAbsFuncRange, sIceLoadFac, selectPenetratingSW,
7e00d7e8f9 Jean*0258      & selectBalanceEmPmR, balanceEmPmR, balanceQnet, balancePrintMean,
427e24e121 Jean*0259      & balanceThetaClimRelax, balanceSaltClimRelax,
                0260      & readBinaryPrec, writeBinaryPrec, writeStatePrec, globalFiles,
                0261      & useSingleCpuIO, useSingleCpuInput, usePickupBeforeC54,
                0262      & usePickupBeforeC35, debugMode, debugLevel, plotLevel
924557e60a Chri*0263 
                0264 C--   Elliptic solver parameters
                0265       NAMELIST /PARM02/
aecc8b0f47 Mart*0266      & cg2dMaxIters, cg2dMinItersNSA, cg2dChkResFreq, cg2dUseMinResSol,
b46f9da188 Jean*0267      & cg2dTargetResidual, cg2dTargetResWunit,
                0268      & cg2dpcOffDFac, cg2dPreCondFreq,
e922d59b96 Jean*0269      & cg3dMaxIters, cg3dChkResFreq,
e6e223b277 Jean*0270      & cg3dTargetResidual, cg3dTargetResWunit,
aecc8b0f47 Mart*0271      & useNSACGSolver, useSRCGSolver, printResidualFreq
924557e60a Chri*0272 
                0273 C--   Time stepping parammeters
                0274       NAMELIST /PARM03/
6e6f314aa7 Patr*0275      & nIter0, nTimeSteps, nTimeSteps_l2, nEndIter,
f804abbd25 Jean*0276      & baseTime, startTime, endTime,
db8d49beca Jean*0277      & deltaT, deltaTClock, deltaTMom,
                0278      & deltaTtracer, dTtracerLev, deltaTFreeSurf,
b9306711eb Jean*0279      & forcing_In_AB, momForcingOutAB, tracForcingOutAB,
c2b6ed6bfd Jean*0280      & momDissip_In_AB, doAB_onGtGs,
2a618ca8c5 Jean*0281      & abEps, alph_AB, beta_AB, startFromPickupAB2, applyExchUV_early,
8039e9b985 Jean*0282      & tauCD, rCD, epsAB_CD, cAdjFreq,
f804abbd25 Jean*0283      & chkPtFreq, pChkPtFreq, pickupSuff, pickupStrictlyMatch,
8c73a5b228 Mart*0284      & writePickupAtEnd,
82d0948361 Jean*0285      & dumpFreq, dumpInitAndLast, adjDumpFreq, taveFreq, tave_lastIter,
f804abbd25 Jean*0286      & diagFreq, monitorFreq, adjMonitorFreq, monitorSelect,
2e3729af6b Jean*0287      & outputTypesInclusive, rwSuffixType,
b9306711eb Jean*0288      & tauThetaClimRelax, tauSaltClimRelax, latBandClimRelax,
a31e157718 Jean*0289      & tauThetaClimRelax3Dim, tauSaltClimRelax3Dim,
fc64655a31 Dimi*0290      & periodicExternalForcing, externForcingPeriod, externForcingCycle,
82d0948361 Jean*0291      & calendarDumps
924557e60a Chri*0292 
                0293 C--   Gridding parameters
                0294       NAMELIST /PARM04/
b9306711eb Jean*0295      & usingCartesianGrid, usingCylindricalGrid,
9780090eaa Jean*0296      & usingSphericalPolarGrid, usingCurvilinearGrid,
                0297      & xgOrigin, ygOrigin, dxSpacing, dySpacing,
                0298      & delX, delY, delXFile, delYFile, horizGridFile,
                0299      & phiEuler, thetaEuler, psiEuler,
2ad79bdf32 Jean*0300      & rSphere, radius_fromHorizGrid, deepAtmosphere, seaLev_Z,
                0301      & top_Pres, delZ, delP, delR, delRc, delRFile, delRcFile,
bb2fd3f1ad Jean*0302      & useMin4hFacEdges, interViscAr_pCell, interDiffKr_pCell,
                0303      & pCellMix_select, pCellMix_maxFac, pCellMix_delR,
                0304      & pCellMix_viscAr, pCellMix_diffKr,
                0305      & selectSigmaCoord, rSigmaBnd, hybSigmFile,
2ad79bdf32 Jean*0306      & Ro_SeaLevel, rkFac, groundAtK1, thetaMin, phiMin
924557e60a Chri*0307 
81bc00c2f0 Chri*0308 C--   Input files
                0309       NAMELIST /PARM05/
26565fca80 Jean*0310      & bathyFile, topoFile, addWwallFile, addSwallFile, shelfIceFile,
                0311      & diffKrFile, viscAhDfile, viscAhZfile, viscA4Dfile, viscA4Zfile,
4c5bb1c88e Jean*0312      & hydrogThetaFile, hydrogSaltFile,
                0313      & maskIniTemp, maskIniSalt, checkIniTemp, checkIniSalt,
88830be691 Alis*0314      & zonalWindFile, meridWindFile,
                0315      & thetaClimFile, saltClimFile,
1e273d1bf5 Jean*0316      & surfQfile, surfQnetFile, surfQswFile, EmPmRfile, saltFluxFile,
2dcaa8b9a5 Patr*0317      & uVelInitFile, vVelInitFile, pSurfInitFile,
0320e25227 Mart*0318      & dQdTFile, ploadFile, geoPotAnomFile, addMassFile,
                0319      & tCylIn, tCylOut,
90929f8806 Patr*0320      & eddyPsiXFile, eddyPsiYFile, geothermalFile,
7e00d7e8f9 Jean*0321      & lambdaThetaFile, lambdaSaltFile, wghtBalanceFile,
cbc417d429 Mart*0322      & mdsioLocalDir, adTapeDir,
9a263a84a8 Ed H*0323      & the_run_name
a30418b6b9 Ed H*0324 CEOP
81bc00c2f0 Chri*0325 
b9b591469d Jean*0326 #ifdef ALLOW_EXCH2
                0327       gridNx = exch2_mydNx(1)
                0328       gridNy = exch2_mydNy(1)
                0329 #else /* ALLOW_EXCH2 */
                0330       gridNx = Nx
                0331       gridNy = Ny
                0332 #endif /* ALLOW_EXCH2 */
                0333 
924557e60a Chri*0334       _BEGIN_MASTER(myThid)
                0335 
0127add478 Alis*0336 C     Defaults values for input parameters
                0337       CALL SET_DEFAULTS(
                0338      O   viscArDefault, diffKrTDefault, diffKrSDefault,
f15994caab Jean*0339      O   hFacMinDrDefault, delRDefault,
0127add478 Alis*0340      I   myThid )
7843dde2de jm-c 0341 
427e24e121 Jean*0342       useJamartWetPoints = .FALSE.
7843dde2de jm-c 0343       useEnergyConservingCoriolis = .FALSE.
e35a804907 Jean*0344       SadournyCoriolis = .FALSE.
427e24e121 Jean*0345       use3dCoriolis   = .TRUE.
                0346       metricTerms     = .TRUE.
7e00d7e8f9 Jean*0347       balanceEmPmR    = .FALSE.
0127add478 Alis*0348 
910f05e765 Chri*0349 C--   Initialise "which vertical coordinate system used" flags.
                0350       zCoordInputData = .FALSE.
                0351       pCoordInputData = .FALSE.
                0352       rCoordInputData = .FALSE.
                0353       coordsSet       = 0
                0354 
8fc5af62b8 Jean*0355 C--   Initialise retired parameters to unlikely value
ede2cfe569 Chri*0356       nRetired = 0
c07cd3bfa8 Jean*0357       useConstantF    = .FALSE.
                0358       useBetaPlaneF   = .FALSE.
                0359       useSphereF      = .TRUE.
7a4ec66481 Jean*0360       tempDiffusion   = .TRUE.
                0361       saltDiffusion   = .TRUE.
                0362       tracerAdvScheme = UNSET_I
                0363       trac_EvPrRn     = UNSET_RL
ede2cfe569 Chri*0364       zonal_filt_lat  = UNSET_RL
36b12bb7ff Jean*0365       gravitySign     = UNSET_RL
8f934af0f4 Jean*0366       viscAstrain     = UNSET_RL
                0367       viscAtension    = UNSET_RL
742cf4499c Jean*0368       useAnisotropicViscAgridMax=.TRUE.
6d61b52b09 Jean*0369       usePickupBeforeC35    = .FALSE.
83513ae28b Jean*0370       saveDebugMode   = debugMode
e0b3e1bdd8 Dimi*0371       allowInteriorFreezing = .FALSE.
e793395f17 Jean*0372       useOldFreezing  = .FALSE.
427e24e121 Jean*0373       writeStatePrec  = UNSET_I
7e00d7e8f9 Jean*0374       cg2dChkResFreq  = UNSET_I
aecc8b0f47 Mart*0375       cg3dChkResFreq  = UNSET_I
46979d29da Jean*0376       tauThetaClimRelax3Dim = UNSET_RL
                0377       tauSaltClimRelax3Dim  = UNSET_RL
b7411f1a84 Jean*0378       taveFreq        = UNSET_RL
                0379       tave_lastIter   = UNSET_RL
742cf4499c Jean*0380       calendarDumps   = .FALSE.
2ad79bdf32 Jean*0381       Ro_SeaLevel     = UNSET_RL
742cf4499c Jean*0382       rkFac           = UNSET_RL
                0383       groundAtK1      = .FALSE.
9780090eaa Jean*0384       thetaMin        = UNSET_RL
                0385       phiMin          = UNSET_RL
8194bc4a99 Mart*0386       shelfIceFile    = ' '
275ad78b78 Jean*0387       dQdTFile        = ' '
ede2cfe569 Chri*0388 
924557e60a Chri*0389 C--   Open the parameter file
b9306711eb Jean*0390       WRITE(msgBuf,'(A)')
aa8f58c33d Jean*0391      &     ' INI_PARMS: opening model parameter file "data"'
924557e60a Chri*0392       CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
b9306711eb Jean*0393      &                    SQUEEZE_RIGHT, myThid )
924557e60a Chri*0394 
aa8f58c33d Jean*0395       CALL OPEN_COPY_DATA_FILE( 'data', 'INI_PARMS',
                0396      O                          iUnit, myThid )
924557e60a Chri*0397 
aa8f58c33d Jean*0398 C--   Read settings from iUnit (= a copy of model parameter file "data").
                0399       errIO = 0
90aed06079 Jean*0400       errCount = 0
924557e60a Chri*0401 
                0402 C--   Set default "physical" parameters
6deb815953 Mart*0403       viscAhW  = UNSET_RL
0b6cbae535 Mart*0404       viscA4W  = UNSET_RL
915a82d1f5 Jean*0405       viscAhD  = UNSET_RL
                0406       viscAhZ  = UNSET_RL
                0407       viscA4D  = UNSET_RL
                0408       viscA4Z  = UNSET_RL
809bdccbfc Jean*0409       viscAz   = UNSET_RL
910f05e765 Chri*0410       viscAp   = UNSET_RL
7cd7754949 Jean*0411       viscAr   = UNSET_RL
910f05e765 Chri*0412       diffKzT  = UNSET_RL
                0413       diffKpT  = UNSET_RL
                0414       diffKrT  = UNSET_RL
                0415       diffKzS  = UNSET_RL
                0416       diffKpS  = UNSET_RL
                0417       diffKrS  = UNSET_RL
aad87ebb4c Jean*0418       hFacMinDr         = UNSET_RL
                0419       hFacMinDz         = UNSET_RL
                0420       hFacMinDp         = UNSET_RL
                0421       tAlpha            = UNSET_RL
                0422       sBeta             = UNSET_RL
                0423       implicitNHPress   = UNSET_RL
49efc6c1e0 Jean*0424       tempVertAdvScheme = 0
                0425       saltVertAdvScheme = 0
337b46d524 Jean*0426       plotLevel         = UNSET_I
49efc6c1e0 Jean*0427 C--   z,p,r coord input switching.
aa8f58c33d Jean*0428       WRITE(msgBuf,'(A)') ' INI_PARMS ; starts to read PARM01'
6b6f69d6b0 Jean*0429       CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
b9306711eb Jean*0430      &                    SQUEEZE_RIGHT, myThid )
88830be691 Alis*0431       READ(UNIT=iUnit,NML=PARM01) !,IOSTAT=errIO)
c206b07c1d Chri*0432       IF ( errIO .LT. 0 ) THEN
924557e60a Chri*0433        WRITE(msgBuf,'(A)')
aa8f58c33d Jean*0434      &  'S/R INI_PARMS: Error reading model parameter file "data"'
b9306711eb Jean*0435        CALL PRINT_ERROR( msgBuf, myThid )
aa8f58c33d Jean*0436        WRITE(msgBuf,'(A)') 'S/R INI_PARMS: Problem in namelist PARM01'
b9306711eb Jean*0437        CALL PRINT_ERROR( msgBuf, myThid )
924557e60a Chri*0438        STOP 'ABNORMAL END: S/R INI_PARMS'
41f1ac6de1 Jean*0439       ELSE
aa8f58c33d Jean*0440        WRITE(msgBuf,'(A)') ' INI_PARMS ; read PARM01 : OK'
41f1ac6de1 Jean*0441        CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
b9306711eb Jean*0442      &                     SQUEEZE_RIGHT, myThid )
c206b07c1d Chri*0443       ENDIF
463053c692 Jean*0444 
599be48d9a Jean*0445 C-    set the type of vertical coordinate and type of fluid
e99e7fe4e6 Jean*0446 C     according to buoyancyRelation
                0447       usingPCoords    = .FALSE.
                0448       usingZCoords    = .FALSE.
                0449       fluidIsAir      = .FALSE.
                0450       fluidIsWater    = .FALSE.
                0451       IF ( buoyancyRelation.EQ.'ATMOSPHERIC' ) THEN
                0452         usingPCoords = .TRUE.
                0453         fluidIsAir   = .TRUE.
                0454       ELSEIF ( buoyancyRelation.EQ.'OCEANICP') THEN
                0455         usingPCoords  = .TRUE.
                0456         fluidIsWater  = .TRUE.
                0457       ELSEIF ( buoyancyRelation.EQ.'OCEANIC' ) THEN
                0458         usingZCoords  = .TRUE.
                0459         fluidIsWater  = .TRUE.
                0460       ELSE
                0461         WRITE(msgBuf,'(2A)') 'S/R INI_PARMS:',
                0462      &      ' Bad value of buoyancyRelation '
b9306711eb Jean*0463         CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*0464         errCount = errCount + 1
e99e7fe4e6 Jean*0465       ENDIF
                0466 
f615684c22 Jean*0467       IF ( .NOT.rigidLid .AND.
                0468      &     .NOT.implicitFreeSurface ) THEN
                0469 C-     No barotropic solver selected => use implicitFreeSurface as default
                0470        WRITE(msgBuf,'(A)')
f76021331f Jean*0471      &  'S/R INI_PARMS: No request for barotropic solver'
f615684c22 Jean*0472        CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
b9306711eb Jean*0473      &                     SQUEEZE_RIGHT, myThid )
f615684c22 Jean*0474        WRITE(msgBuf,'(A)')
                0475      &  'S/R INI_PARMS: => Use implicitFreeSurface as default'
                0476        CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
b9306711eb Jean*0477      &                     SQUEEZE_RIGHT, myThid )
f615684c22 Jean*0478         implicitFreeSurface = .TRUE.
                0479       ENDIF
8039e9b985 Jean*0480       IF ( implicitFreeSurface ) freeSurfFac = 1. _d 0
                0481       IF ( rigidLid            ) freeSurfFac = 0. _d 0
0127add478 Alis*0482       IF ( gBaro .EQ. UNSET_RL ) gBaro=gravity
                0483       IF ( rhoConst .EQ. UNSET_RL ) rhoConst=rhoNil
01d019192b Jean*0484       IF ( rhoConstFresh .EQ. UNSET_RL ) rhoConstFresh=rhoConst
16f5093311 Jean*0485       IF ( implicitNHPress.EQ.UNSET_RL )
                0486      &     implicitNHPress = implicSurfPress
81a2602dac Jean*0487       IF ( omega .EQ. UNSET_RL ) THEN
f0bc3639a4 Jean*0488         omega = 0. _d 0
b9306711eb Jean*0489         IF ( rotationPeriod .NE. 0. _d 0 )
8039e9b985 Jean*0490      &  omega = 2. _d 0 * PI / rotationPeriod
f0bc3639a4 Jean*0491       ELSEIF ( omega .EQ. 0. _d 0 ) THEN
                0492         rotationPeriod = 0. _d 0
81a2602dac Jean*0493       ELSE
8039e9b985 Jean*0494         rotationPeriod = 2. _d 0 * PI / omega
81a2602dac Jean*0495       ENDIF
aad87ebb4c Jean*0496       IF ( atm_Rd .EQ. UNSET_RL ) THEN
463053c692 Jean*0497         atm_Rd = atm_Cp * atm_kappa
                0498       ELSE
                0499         atm_kappa = atm_Rd / atm_Cp
                0500       ENDIF
775ad00e06 Alis*0501 C--   Non-hydrostatic/quasi-hydrostatic
                0502       IF (nonHydrostatic.AND.quasiHydrostatic) THEN
                0503        WRITE(msgBuf,'(A)')
                0504      &   'Illegal: both nonHydrostatic = quasiHydrostatic = TRUE'
b9306711eb Jean*0505        CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*0506        errCount = errCount + 1
775ad00e06 Alis*0507       ENDIF
809bdccbfc Jean*0508 C--  Advection and Forcing for Temp and salt
49efc6c1e0 Jean*0509       IF (tempVertAdvScheme.EQ.0) tempVertAdvScheme = tempAdvScheme
                0510       IF (saltVertAdvScheme.EQ.0) saltVertAdvScheme = saltAdvScheme
915a82d1f5 Jean*0511 C--  horizontal viscosity (acting on Divergence or Vorticity)
                0512       IF ( viscAhD .EQ. UNSET_RL ) viscAhD = viscAh
                0513       IF ( viscAhZ .EQ. UNSET_RL ) viscAhZ = viscAh
                0514       IF ( viscA4D .EQ. UNSET_RL ) viscA4D = viscA4
                0515       IF ( viscA4Z .EQ. UNSET_RL ) viscA4Z = viscA4
2d5bb917cc Jean*0516 C--  horizontal viscosity for vertical momentum
92da82d27f Jean*0517       IF ( viscAhW .EQ. UNSET_RL ) viscAhW = viscAhD
                0518       IF ( viscA4W .EQ. UNSET_RL ) viscA4W = viscA4D
910f05e765 Chri*0519 C--   z,p,r coord input switching.
                0520       IF ( viscAz .NE. UNSET_RL ) zCoordInputData = .TRUE.
                0521       IF ( viscAp .NE. UNSET_RL ) pCoordInputData = .TRUE.
                0522       IF ( viscAr .NE. UNSET_RL ) rCoordInputData = .TRUE.
                0523       IF ( viscAr .EQ. UNSET_RL )          viscAr = viscAz
                0524       IF ( viscAr .EQ. UNSET_RL )          viscAr = viscAp
7cd7754949 Jean*0525       vertSetCount = 0
                0526       DO k=1,Nr
13b4d13c98 Jean*0527         IF ( viscArNr(k).NE.UNSET_RL ) vertSetCount = vertSetCount + 1
7cd7754949 Jean*0528       ENDDO
13b4d13c98 Jean*0529       IF ( vertSetCount.GT.0 .AND. vertSetCount.LT.Nr ) THEN
                0530         WRITE(msgBuf,'(A,2(I5,A))') 'S/R INI_PARMS: Partial setting (',
                0531      &        vertSetCount, ' /', Nr, ') of viscArNr is not allowed'
                0532         CALL PRINT_ERROR( msgBuf, myThid )
                0533         errCount = errCount + 1
                0534       ENDIF
                0535       IF ( viscAr .EQ. UNSET_RL ) THEN
                0536         viscAr = viscArDefault
                0537       ELSEIF ( vertSetCount.GT.0 ) THEN
                0538         WRITE(msgBuf,'(2A)') 'S/R INI_PARMS: Cannot set both ',
                0539      &     'viscArNr and viscAr (or Ap,Az) in param file data'
                0540         CALL PRINT_ERROR( msgBuf, myThid )
                0541         errCount = errCount + 1
                0542       ENDIF
                0543       IF ( vertSetCount.EQ.0 ) THEN
                0544         DO k=1,Nr
                0545           viscArNr(k) = viscAr
                0546         ENDDO
7cd7754949 Jean*0547       ENDIF
f9ef55fa0c Jean*0548 #ifdef ALLOW_MOM_COMMON
                0549 C-    set default scheme for quadratic bottom-drag
                0550       IF ( selectBotDragQuadr.EQ.-1 .AND. bottomDragQuadratic.NE.0. )
                0551      &  selectBotDragQuadr = 0
ab47de63dc Mart*0552       IF ( selectBotDragQuadr.EQ.-1 .AND. zRoughBot.NE.0. )
                0553      &  selectBotDragQuadr = 0
f9ef55fa0c Jean*0554 #endif /* ALLOW_MOM_COMMON */
910f05e765 Chri*0555 
2d5bb917cc Jean*0556       IF ( smag3D_diffCoeff.NE.zeroRL .AND. .NOT.useSmag3D ) THEN
                0557         WRITE(msgBuf,'(2A)') '** WARNING ** INI_PARMS: ',
                0558      &    'will not use "smag3D_diffCoeff" without useSmag3D'
                0559         CALL PRINT_MESSAGE( msgBuf, errorMessageUnit,
                0560      &                      SQUEEZE_RIGHT, myThid )
                0561         WRITE(msgBuf,'(2A)') '** WARNING ** INI_PARMS: ',
                0562      &    '==> reset "smag3D_diffCoeff" to zero'
                0563         CALL PRINT_MESSAGE( msgBuf, errorMessageUnit,
                0564      &                      SQUEEZE_RIGHT, myThid )
                0565        smag3D_diffCoeff = zeroRL
                0566       ENDIF
910f05e765 Chri*0567       IF ( diffKzT .NE. UNSET_RL ) zCoordInputData  = .TRUE.
                0568       IF ( diffKpT .NE. UNSET_RL ) pCoordInputData  = .TRUE.
                0569       IF ( diffKrT .NE. UNSET_RL ) rCoordInputData  = .TRUE.
                0570       IF ( diffKrT .EQ. UNSET_RL )          diffKrT = diffKzT
                0571       IF ( diffKrT .EQ. UNSET_RL )          diffKrT = diffKpT
7cd7754949 Jean*0572       vertSetCount = 0
8fc5af62b8 Jean*0573       DO k=1,Nr
13b4d13c98 Jean*0574         IF ( diffKrNrT(k).NE.UNSET_RL ) vertSetCount = vertSetCount + 1
8fc5af62b8 Jean*0575       ENDDO
13b4d13c98 Jean*0576       IF ( vertSetCount.GT.0 .AND. vertSetCount.LT.Nr ) THEN
                0577         WRITE(msgBuf,'(A,2(I5,A))') 'S/R INI_PARMS: Partial setting (',
                0578      &        vertSetCount, ' /', Nr, ') of diffKrNrT is not allowed'
                0579         CALL PRINT_ERROR( msgBuf, myThid )
                0580         errCount = errCount + 1
                0581       ENDIF
                0582       specifiedDiffKrT = vertSetCount.EQ.Nr
                0583       IF ( diffKrT .EQ. UNSET_RL ) THEN
                0584         diffKrT = diffKrTDefault
                0585       ELSEIF ( vertSetCount.GT.0 ) THEN
                0586         WRITE(msgBuf,'(2A)') 'S/R INI_PARMS: Cannot set both ',
                0587      &     'diffKrNrT and diffKrT (or Kp,Kz) in param file data'
                0588         CALL PRINT_ERROR( msgBuf, myThid )
                0589         errCount = errCount + 1
                0590       ELSE
                0591         specifiedDiffKrT = .TRUE.
                0592       ENDIF
                0593       IF ( vertSetCount.EQ.0 ) THEN
                0594         DO k=1,Nr
                0595           diffKrNrT(k) = diffKrT
                0596         ENDDO
8fc5af62b8 Jean*0597       ENDIF
910f05e765 Chri*0598 
                0599       IF ( diffKzS .NE. UNSET_RL ) zCoordInputData  = .TRUE.
                0600       IF ( diffKpS .NE. UNSET_RL ) pCoordInputData  = .TRUE.
                0601       IF ( diffKrS .NE. UNSET_RL ) rCoordInputData  = .TRUE.
                0602       IF ( diffKrS .EQ. UNSET_RL )          diffKrS = diffKzS
                0603       IF ( diffKrS .EQ. UNSET_RL )          diffKrS = diffKpS
7cd7754949 Jean*0604       vertSetCount = 0
8fc5af62b8 Jean*0605       DO k=1,Nr
13b4d13c98 Jean*0606         IF ( diffKrNrS(k).NE.UNSET_RL ) vertSetCount = vertSetCount + 1
8fc5af62b8 Jean*0607       ENDDO
13b4d13c98 Jean*0608       IF ( vertSetCount.GT.0 .AND. vertSetCount.LT.Nr ) THEN
                0609         WRITE(msgBuf,'(A,2(I5,A))') 'S/R INI_PARMS: Partial setting (',
                0610      &        vertSetCount, ' /', Nr, ') of diffKrNrS is not allowed'
                0611         CALL PRINT_ERROR( msgBuf, myThid )
                0612         errCount = errCount + 1
                0613       ENDIF
a631d4754d Jean*0614       IF ( vertSetCount.EQ.Nr ) THEN
                0615         specifiedDiffKrS = .TRUE.
                0616         IF ( diffKrS.NE.UNSET_RL ) THEN
                0617          WRITE(msgBuf,'(2A)') 'S/R INI_PARMS: Cannot set both ',
13b4d13c98 Jean*0618      &     'diffKrNrS and diffKrS (or Kp,Kz) in param file data'
a631d4754d Jean*0619          CALL PRINT_ERROR( msgBuf, myThid )
                0620          errCount = errCount + 1
                0621         ENDIF
                0622       ELSEIF ( diffKrS.NE.UNSET_RL ) THEN
13b4d13c98 Jean*0623         specifiedDiffKrS = .TRUE.
                0624         DO k=1,Nr
                0625           diffKrNrS(k) = diffKrS
                0626         ENDDO
a631d4754d Jean*0627       ELSE
                0628         specifiedDiffKrS = .FALSE.
                0629         diffKrS = diffKrSDefault
                0630 C-    use temp diffusivity as default salt diffusivity:
                0631         DO k=1,Nr
                0632           diffKrNrS(k) = diffKrNrT(k)
                0633         ENDDO
8fc5af62b8 Jean*0634       ENDIF
910f05e765 Chri*0635 
e40c34e398 Dimi*0636       IF (diffKrBLEQsurf .EQ. UNSET_RL) diffKrBLEQsurf = diffKrBL79surf
                0637       IF (diffKrBLEQdeep .EQ. UNSET_RL) diffKrBLEQdeep = diffKrBL79deep
                0638       IF (diffKrBLEQscl  .EQ. UNSET_RL) diffKrBLEQscl  = diffKrBL79scl
                0639       IF (diffKrBLEQHo   .EQ. UNSET_RL) diffKrBLEQHo   = diffKrBL79Ho
                0640 
910f05e765 Chri*0641       IF ( hFacMinDz .NE. UNSET_RL ) zCoordInputData = .TRUE.
                0642       IF ( hFacMinDp .NE. UNSET_RL ) pCoordInputData = .TRUE.
                0643       IF ( hFacMinDr .NE. UNSET_RL ) rCoordInputData = .TRUE.
fe0a6d6615 Alis*0644       IF ( hFacMinDr .EQ. UNSET_RL ) hFacMinDr       = hFacMinDz
                0645       IF ( hFacMinDr .EQ. UNSET_RL ) hFacMinDr       = hFacMinDp
910f05e765 Chri*0646       IF ( hFacMinDr .EQ. UNSET_RL ) hFacMinDr       = hFacMinDrDefault
c0a4efc370 Chri*0647 
745a0098ab Jean*0648       IF (convertFW2Salt.EQ.UNSET_RL) THEN
                0649         convertFW2Salt = 35.
                0650         IF (useRealFreshWaterFlux) convertFW2Salt=-1
d18df35fee Jean*0651         IF ( selectAddFluid.GE.1 ) convertFW2Salt=-1
745a0098ab Jean*0652       ENDIF
47ecd1025a Jean*0653       IF ( selectPenetratingSW.EQ.UNSET_I ) THEN
                0654         selectPenetratingSW = 0
                0655 #ifdef SHORTWAVE_HEATING
                0656         IF ( fluidIsWater ) selectPenetratingSW = 1
                0657 #endif
                0658       ENDIF
745a0098ab Jean*0659 
7843dde2de jm-c 0660 C--   for backward compatibility :
                0661       IF ( selectCoriScheme.EQ.UNSET_I ) THEN
                0662         selectCoriScheme = 0
                0663         IF ( useJamartWetPoints ) selectCoriScheme = 1
                0664         IF ( useEnergyConservingCoriolis .AND.
                0665      &       .NOT.vectorInvariantMomentum )
                0666      &    selectCoriScheme = selectCoriScheme + 2
                0667       ENDIF
                0668       IF ( vectorInvariantMomentum ) THEN
                0669         IF ( useJamartWetPoints .AND. selectCoriScheme.NE.1 ) THEN
                0670           WRITE(msgBuf,'(A,I5,A)')
                0671      &     'S/R INI_PARMS: selectCoriScheme=', selectCoriScheme,
                0672      &     ' conflicts with "useJamartWetPoints"'
                0673           CALL PRINT_ERROR( msgBuf, myThid )
                0674           errCount = errCount + 1
                0675         ENDIF
                0676       ELSE
                0677         IF ( useEnergyConservingCoriolis
                0678      &       .AND. selectCoriScheme.LT.2 ) THEN
                0679           WRITE(msgBuf,'(A,I5,A)')
                0680      &     'S/R INI_PARMS: selectCoriScheme=', selectCoriScheme,
                0681      &     ' conflicts with "useEnergyConservingCoriolis"'
                0682           CALL PRINT_ERROR( msgBuf, myThid )
                0683           errCount = errCount + 1
                0684         ENDIF
                0685         IF ( useJamartWetPoints .AND. selectCoriScheme.NE.1
                0686      &                          .AND. selectCoriScheme.NE.3 ) THEN
                0687           WRITE(msgBuf,'(A,I5,A)')
                0688      &     'S/R INI_PARMS: selectCoriScheme=', selectCoriScheme,
                0689      &     ' conflicts with "useJamartWetPoints"'
                0690           CALL PRINT_ERROR( msgBuf, myThid )
                0691           errCount = errCount + 1
                0692         ENDIF
                0693       ENDIF
e35a804907 Jean*0694       IF ( SadournyCoriolis ) THEN
                0695 C--   for backward compatibility :
                0696         IF ( selectVortScheme.EQ.UNSET_I ) selectVortScheme = 2
                0697         IF ( selectVortScheme.NE.2 ) THEN
7fe6343684 Jean*0698           WRITE(msgBuf,'(A,I5,A)')
                0699      &     'S/R INI_PARMS: selectVortScheme=', selectVortScheme,
                0700      &     ' conflicts with "SadournyCoriolis"'
e35a804907 Jean*0701           CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*0702           errCount = errCount + 1
e35a804907 Jean*0703         ENDIF
                0704       ENDIF
427e24e121 Jean*0705       IF ( select3dCoriScheme.EQ.UNSET_I ) THEN
                0706 C--   for backward compatibility :
                0707         select3dCoriScheme = 0
                0708         IF ( use3dCoriolis ) select3dCoriScheme = 1
                0709       ELSEIF ( select3dCoriScheme.NE.0 .AND. .NOT.use3dCoriolis ) THEN
                0710         WRITE(msgBuf,'(A,I5,A)')
                0711      &     'S/R INI_PARMS: select3dCoriScheme=', select3dCoriScheme,
                0712      &     ' conflicts with "use3dCoriolis=F"'
                0713         CALL PRINT_ERROR( msgBuf, myThid )
                0714         errCount = errCount + 1
                0715       ENDIF
                0716       IF ( selectMetricTerms.EQ.UNSET_I ) THEN
                0717 C--   for backward compatibility :
                0718         selectMetricTerms = 0
                0719         IF ( metricTerms ) selectMetricTerms = 1
                0720       ELSEIF ( selectMetricTerms.NE.0 .AND. .NOT.metricTerms ) THEN
                0721         WRITE(msgBuf,'(A,I5,A)')
                0722      &     'S/R INI_PARMS: selectMetricTerms=', selectMetricTerms,
                0723      &     ' conflicts with "metricTerms=F"'
                0724         CALL PRINT_ERROR( msgBuf, myThid )
                0725         errCount = errCount + 1
                0726       ENDIF
7e00d7e8f9 Jean*0727       IF ( selectBalanceEmPmR.EQ.UNSET_I ) THEN
                0728         selectBalanceEmPmR = 0
                0729         IF ( balanceEmPmR ) selectBalanceEmPmR = 1
                0730       ELSEIF ( selectBalanceEmPmR.NE.1 .AND. balanceEmPmR ) THEN
                0731         WRITE(msgBuf,'(A,I5,A)')
                0732      &     'S/R INI_PARMS: selectBalanceEmPmR=', selectBalanceEmPmR,
                0733      &     ' conflicts with "balanceEmPmR"'
                0734         CALL PRINT_ERROR( msgBuf, myThid )
                0735         errCount = errCount + 1
                0736       ENDIF
e35a804907 Jean*0737 
aad87ebb4c Jean*0738       IF ( ivdc_kappa.NE.zeroRL .AND. .NOT.implicitDiffusion ) THEN
e1c6dcc4ea Jean*0739         WRITE(msgBuf,'(A,A)')
d4701cb6da Alis*0740      &  'S/R INI_PARMS: To use ivdc_kappa you must enable implicit',
                0741      &  ' vertical diffusion.'
b9306711eb Jean*0742        CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*0743        errCount = errCount + 1
d4701cb6da Alis*0744       ENDIF
                0745 
910f05e765 Chri*0746       coordsSet = 0
                0747       IF ( zCoordInputData ) coordsSet = coordsSet + 1
                0748       IF ( pCoordInputData ) coordsSet = coordsSet + 1
                0749       IF ( rCoordInputData ) coordsSet = coordsSet + 1
                0750       IF ( coordsSet .GT. 1 ) THEN
                0751        WRITE(msgBuf,'(A)')
                0752      &  'S/R INI_PARMS: Cannot mix z, p and r in the input data.'
b9306711eb Jean*0753        CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*0754        errCount = errCount + 1
c0a4efc370 Chri*0755       ENDIF
910f05e765 Chri*0756       IF ( rhoConst .LE. 0. ) THEN
                0757        WRITE(msgBuf,'(A)')
                0758      &  'S/R INI_PARMS: rhoConst must be greater than 0.'
b9306711eb Jean*0759        CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*0760        errCount = errCount + 1
                0761        recip_rhoConst = 1. _d 0
910f05e765 Chri*0762       ELSE
8039e9b985 Jean*0763        recip_rhoConst = 1. _d 0 / rhoConst
910f05e765 Chri*0764       ENDIF
53a1bebe0a Jean*0765       IF ( eosType.EQ.'LINEAR' .AND. rhoNil.LE.0. ) THEN
9bef66796c Alis*0766        WRITE(msgBuf,'(A)')
                0767      &  'S/R INI_PARMS: rhoNil must be greater than 0.'
b9306711eb Jean*0768        CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*0769        errCount = errCount + 1
9bef66796c Alis*0770       ENDIF
0127add478 Alis*0771       IF ( HeatCapacity_Cp .LE. 0. ) THEN
                0772        WRITE(msgBuf,'(A)')
                0773      &  'S/R INI_PARMS: HeatCapacity_Cp must be greater than 0.'
b9306711eb Jean*0774        CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*0775        errCount = errCount + 1
0127add478 Alis*0776       ENDIF
b0f9f620c0 Chri*0777       IF ( gravity .LE. 0. ) THEN
                0778        WRITE(msgBuf,'(A)')
                0779      &  'S/R INI_PARMS: gravity must be greater than 0.'
b9306711eb Jean*0780        CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*0781        errCount = errCount + 1
                0782        recip_gravity = 1. _d 0
b0f9f620c0 Chri*0783       ELSE
8039e9b985 Jean*0784        recip_gravity = 1. _d 0 / gravity
b0f9f620c0 Chri*0785       ENDIF
631fe75038 Dimi*0786 
83513ae28b Jean*0787 C-    set default printResidualFreq according to debugLevel
e922d59b96 Jean*0788       printResidualFreq = -1
83513ae28b Jean*0789       IF ( debugLevel.GE.debLevE ) printResidualFreq = 1
337b46d524 Jean*0790       IF ( plotLevel.EQ.UNSET_I ) plotLevel = debugLevel
924557e60a Chri*0791 
631fe75038 Dimi*0792 C-    set useSingleCpuInput=.TRUE. if useSingleCpuIO==.TRUE.
                0793       IF ( useSingleCpuIO ) useSingleCpuInput=.TRUE.
                0794 
ede2cfe569 Chri*0795 C     Check for retired parameters still being used
                0796       nRetired = 0
c07cd3bfa8 Jean*0797       IF ( useConstantF ) THEN
                0798        nRetired = nRetired+1
                0799        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "useConstantF" ',
                0800      &  'is no longer allowed in file "data"'
                0801        CALL PRINT_ERROR( msgBuf, myThid )
                0802        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: set "selectCoriMap"',
                0803      &  ' [0,1,2,3] to impose a setting over grid default'
                0804        CALL PRINT_ERROR( msgBuf, myThid )
                0805       ENDIF
                0806       IF ( useBetaPlaneF ) THEN
                0807        nRetired = nRetired+1
                0808        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "useBetaPlaneF" ',
                0809      &  'is no longer allowed in file "data"'
                0810        CALL PRINT_ERROR( msgBuf, myThid )
                0811        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: set "selectCoriMap"',
                0812      &  ' [0,1,2,3] to impose a setting over grid default'
                0813        CALL PRINT_ERROR( msgBuf, myThid )
                0814       ENDIF
                0815       IF ( .NOT. useSphereF ) THEN
                0816        nRetired = nRetired+1
                0817        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "useSphereF" ',
                0818      &  'is no longer allowed in file "data"'
                0819        CALL PRINT_ERROR( msgBuf, myThid )
                0820        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: set "selectCoriMap"',
                0821      &  ' [0,1,2,3] to impose a setting over grid default'
                0822        CALL PRINT_ERROR( msgBuf, myThid )
                0823       ENDIF
ede2cfe569 Chri*0824       IF ( zonal_filt_lat .NE. UNSET_RL ) THEN
                0825        nRetired = nRetired+1
                0826        WRITE(msgBuf,'(A,A)')
                0827      &  'S/R INI_PARMS: Paramater "zonal_filt_lat" is',
                0828      &  ' no longer allowed in file "data".'
b9306711eb Jean*0829        CALL PRINT_ERROR( msgBuf, myThid )
ede2cfe569 Chri*0830        WRITE(msgBuf,'(A,A)')
                0831      &  'S/R INI_PARMS: Paramater "zonal_filt_lat" is',
                0832      &  ' now read from file "data.zonfilt".'
b9306711eb Jean*0833        CALL PRINT_ERROR( msgBuf, myThid )
ede2cfe569 Chri*0834       ENDIF
36b12bb7ff Jean*0835       IF ( gravitySign .NE. UNSET_RL ) THEN
                0836        nRetired = nRetired+1
                0837        WRITE(msgBuf,'(A,A)')
                0838      &  'S/R INI_PARMS: "gravitySign" is set according to vertical ',
                0839      &  ' coordinate and is no longer allowed in file "data".'
b9306711eb Jean*0840        CALL PRINT_ERROR( msgBuf, myThid )
36b12bb7ff Jean*0841       ENDIF
7a4ec66481 Jean*0842       IF ( tracerAdvScheme .NE. UNSET_I ) THEN
                0843        nRetired = nRetired+1
                0844        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "tracerAdvScheme" ',
                0845      &  '(old passive tracer code) is no longer allowed in file "data"'
b9306711eb Jean*0846        CALL PRINT_ERROR( msgBuf, myThid )
7a4ec66481 Jean*0847       ENDIF
                0848       IF ( trac_EvPrRn .NE. UNSET_RL ) THEN
                0849        nRetired = nRetired+1
                0850        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "trac_EvPrRn" ',
                0851      &  '(old passive tracer code) is no longer allowed in file "data"'
b9306711eb Jean*0852        CALL PRINT_ERROR( msgBuf, myThid )
7a4ec66481 Jean*0853       ENDIF
                0854       IF ( .NOT. tempDiffusion ) THEN
                0855        nRetired = nRetired+1
                0856        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "tempDiffusion" ',
                0857      &  'is no longer allowed in file "data"'
b9306711eb Jean*0858        CALL PRINT_ERROR( msgBuf, myThid )
7a4ec66481 Jean*0859        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: to turn off diffusion',
                0860      &  ' => set diffusivity to zero'
b9306711eb Jean*0861        CALL PRINT_ERROR( msgBuf, myThid )
7a4ec66481 Jean*0862       ENDIF
                0863       IF ( .NOT. saltDiffusion ) THEN
                0864        nRetired = nRetired+1
                0865        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "saltDiffusion" ',
                0866      &  'is no longer allowed in file "data"'
b9306711eb Jean*0867        CALL PRINT_ERROR( msgBuf, myThid )
7a4ec66481 Jean*0868        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: to turn off diffusion',
                0869      &  ' => set diffusivity to zero'
b9306711eb Jean*0870        CALL PRINT_ERROR( msgBuf, myThid )
7a4ec66481 Jean*0871       ENDIF
8f934af0f4 Jean*0872       IF ( viscAstrain .NE. UNSET_RL ) THEN
                0873        nRetired = nRetired+1
                0874        WRITE(msgBuf,'(A,A)')
                0875      &  'S/R INI_PARMS: "viscAstrain" ',
                0876      &  'is no longer allowed in file "data"'
b9306711eb Jean*0877        CALL PRINT_ERROR( msgBuf, myThid )
8f934af0f4 Jean*0878        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: to use Strain & Tension',
                0879      &  ' formulation => set useStrainTensionVisc to TRUE'
b9306711eb Jean*0880        CALL PRINT_ERROR( msgBuf, myThid )
8f934af0f4 Jean*0881       ENDIF
                0882       IF ( viscAtension .NE. UNSET_RL ) THEN
                0883        nRetired = nRetired+1
                0884        WRITE(msgBuf,'(A,A)')
                0885      &  'S/R INI_PARMS: "viscAtension" ',
                0886      &  'is no longer allowed in file "data"'
b9306711eb Jean*0887        CALL PRINT_ERROR( msgBuf, myThid )
8f934af0f4 Jean*0888        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: to use Strain & Tension',
                0889      &  ' formulation => set useStrainTensionVisc to TRUE'
b9306711eb Jean*0890        CALL PRINT_ERROR( msgBuf, myThid )
8f934af0f4 Jean*0891       ENDIF
6d61b52b09 Jean*0892       IF ( .NOT.useAnisotropicViscAgridMax ) THEN
cef83a8f90 Bayl*0893        nRetired = nRetired+1
                0894        WRITE(msgBuf,'(A,A)')
                0895      &  'S/R INI_PARMS: "useAnisotropicViscAgridMax" ',
                0896      &  'is not allowed in "data" substitute useAreaViscLength=true'
b9306711eb Jean*0897        CALL PRINT_ERROR( msgBuf, myThid )
cef83a8f90 Bayl*0898       ENDIF
6d61b52b09 Jean*0899       IF ( usePickupBeforeC35 ) THEN
                0900        nRetired = nRetired+1
                0901        WRITE(msgBuf,'(A,A)')
                0902      &  'S/R INI_PARMS: "usePickupBeforeC35" ',
                0903      &  'is no longer supported & not longer allowed in file "data"'
                0904        CALL PRINT_ERROR( msgBuf, myThid )
                0905       ENDIF
83513ae28b Jean*0906       IF ( debugMode.NEQV.saveDebugMode ) THEN
                0907        nRetired = nRetired+1
                0908        WRITE(msgBuf,'(A,A)')
                0909      &  'S/R INI_PARMS: "debugMode" has been moved to "eedata"',
                0910      &  ' and is no longer allowed in file "data"'
                0911        CALL PRINT_ERROR( msgBuf, myThid )
                0912       ENDIF
e0b3e1bdd8 Dimi*0913       IF ( allowInteriorFreezing ) THEN
                0914        nRetired = nRetired+1
                0915        WRITE(msgBuf,'(A,A)')
                0916      &  'S/R INI_PARMS: "allowInteriorFreezing" has been replaced',
                0917      &  ' by pkg/frazil and is no longer allowed in file "data"'
                0918        CALL PRINT_ERROR( msgBuf, myThid )
                0919       ENDIF
e793395f17 Jean*0920       IF ( useOldFreezing ) THEN
                0921        nRetired = nRetired+1
                0922        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "useOldFreezing" ',
                0923      &  'is no longer supported & not longer allowed in file "data"'
                0924        CALL PRINT_ERROR( msgBuf, myThid )
                0925       ENDIF
427e24e121 Jean*0926       IF ( writeStatePrec .NE. UNSET_I ) THEN
                0927        nRetired = nRetired+1
                0928        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "writeStatePrec" ',
                0929      &  '(un-used) is no longer allowed in file "data"'
                0930        CALL PRINT_ERROR( msgBuf, myThid )
                0931       ENDIF
ede2cfe569 Chri*0932 
924557e60a Chri*0933 C--   Elliptic solver parameters
aa8f58c33d Jean*0934       WRITE(msgBuf,'(A)') ' INI_PARMS ; starts to read PARM02'
6b6f69d6b0 Jean*0935       CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
b9306711eb Jean*0936      &                    SQUEEZE_RIGHT, myThid )
88830be691 Alis*0937       READ(UNIT=iUnit,NML=PARM02) !,IOSTAT=errIO)
c206b07c1d Chri*0938       IF ( errIO .LT. 0 ) THEN
924557e60a Chri*0939        WRITE(msgBuf,'(A)')
aa8f58c33d Jean*0940      &  'S/R INI_PARMS: Error reading model parameter file "data"'
b9306711eb Jean*0941        CALL PRINT_ERROR( msgBuf, myThid )
aa8f58c33d Jean*0942        WRITE(msgBuf,'(A)') 'S/R INI_PARMS: Problem in namelist PARM02'
b9306711eb Jean*0943        CALL PRINT_ERROR( msgBuf, myThid )
924557e60a Chri*0944        STOP 'ABNORMAL END: S/R INI_PARMS'
41f1ac6de1 Jean*0945       ELSE
aa8f58c33d Jean*0946        WRITE(msgBuf,'(A)') ' INI_PARMS ; read PARM02 : OK'
41f1ac6de1 Jean*0947        CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
b9306711eb Jean*0948      &                     SQUEEZE_RIGHT, myThid )
809bdccbfc Jean*0949       ENDIF
aecc8b0f47 Mart*0950 C     Check for retired parameters still being used
                0951       IF ( cg2dChkResFreq .NE. UNSET_I ) THEN
                0952        nRetired = nRetired+1
                0953        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: unused "cg2dChkResFreq"',
                0954      &  ' is no longer allowed in file "data"'
                0955        CALL PRINT_ERROR( msgBuf, myThid )
                0956       ENDIF
                0957       IF ( cg3dChkResFreq .NE. UNSET_I ) THEN
                0958        nRetired = nRetired+1
                0959        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: unused "cg3dChkResFreq"',
                0960      &  ' is no longer allowed in file "data"'
                0961        CALL PRINT_ERROR( msgBuf, myThid )
                0962       ENDIF
924557e60a Chri*0963 
                0964 C--   Time stepping parameters
8039e9b985 Jean*0965       rCD               = -1. _d 0
                0966       epsAB_CD          = UNSET_RL
36b12bb7ff Jean*0967       latBandClimRelax  = UNSET_RL
062a876ce5 Jean*0968       deltaTtracer      = 0. _d 0
c2b6ed6bfd Jean*0969       forcing_In_AB     = .TRUE.
aa8f58c33d Jean*0970       WRITE(msgBuf,'(A)') ' INI_PARMS ; starts to read PARM03'
6b6f69d6b0 Jean*0971       CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
b9306711eb Jean*0972      &                    SQUEEZE_RIGHT, myThid )
88830be691 Alis*0973       READ(UNIT=iUnit,NML=PARM03) !,IOSTAT=errIO)
c206b07c1d Chri*0974       IF ( errIO .LT. 0 ) THEN
924557e60a Chri*0975        WRITE(msgBuf,'(A)')
aa8f58c33d Jean*0976      &  'S/R INI_PARMS: Error reading model parameter file "data"'
b9306711eb Jean*0977        CALL PRINT_ERROR( msgBuf, myThid )
aa8f58c33d Jean*0978        WRITE(msgBuf,'(A)') 'S/R INI_PARMS: Problem in namelist PARM03'
b9306711eb Jean*0979        CALL PRINT_ERROR( msgBuf, myThid )
924557e60a Chri*0980        STOP 'ABNORMAL END: S/R INI_PARMS'
41f1ac6de1 Jean*0981       ELSE
aa8f58c33d Jean*0982        WRITE(msgBuf,'(A)') ' INI_PARMS ; read PARM03 : OK'
41f1ac6de1 Jean*0983        CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
b9306711eb Jean*0984      &                     SQUEEZE_RIGHT, myThid )
809bdccbfc Jean*0985       ENDIF
46979d29da Jean*0986 C     Check for retired parameters still being used
                0987       IF ( tauThetaClimRelax3Dim .NE. UNSET_RL ) THEN
                0988        nRetired = nRetired+1
                0989        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "tauThetaClimRelax3Dim" ',
                0990      &  'is no longer allowed in file "data"'
b9306711eb Jean*0991        CALL PRINT_ERROR( msgBuf, myThid )
46979d29da Jean*0992        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: 3-dim. relaxation code',
                0993      & ' has moved to separate pkg/rbcs.'
b9306711eb Jean*0994        CALL PRINT_ERROR( msgBuf, myThid )
46979d29da Jean*0995       ENDIF
                0996       IF ( tauSaltClimRelax3Dim .NE. UNSET_RL ) THEN
                0997        nRetired = nRetired+1
                0998        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "tauSaltClimRelax3Dim" ',
                0999      &  'is no longer allowed in file "data"'
b9306711eb Jean*1000        CALL PRINT_ERROR( msgBuf, myThid )
46979d29da Jean*1001        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: 3-dim. relaxation code',
                1002      & ' has moved to separate pkg/rbcs.'
b9306711eb Jean*1003        CALL PRINT_ERROR( msgBuf, myThid )
46979d29da Jean*1004       ENDIF
b7411f1a84 Jean*1005       IF ( taveFreq .NE. UNSET_RL ) THEN
                1006        nRetired = nRetired+1
                1007        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "taveFreq" ',
                1008      &  'is no longer allowed in file "data"'
                1009        CALL PRINT_ERROR( msgBuf, myThid )
                1010       ENDIF
                1011       IF ( tave_lastIter .NE. UNSET_RL ) THEN
                1012        nRetired = nRetired+1
                1013        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "tave_lastIter" ',
                1014      &  'is no longer allowed in file "data"'
                1015        CALL PRINT_ERROR( msgBuf, myThid )
                1016       ENDIF
                1017       IF ( taveFreq.NE.UNSET_RL .OR. tave_lastIter.NE.UNSET_RL ) THEN
                1018        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: ',
                1019      & ' since "pkg/timeave" has been removed.'
                1020        CALL PRINT_ERROR( msgBuf, myThid )
                1021       ENDIF
742cf4499c Jean*1022       IF ( calendarDumps ) THEN
                1023        nRetired = nRetired+1
                1024        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: "calendarDumps" ',
                1025      &  'is no longer allowed in file "data"'
b9306711eb Jean*1026        CALL PRINT_ERROR( msgBuf, myThid )
742cf4499c Jean*1027        WRITE(msgBuf,'(A,A)') 'S/R INI_PARMS: calendarDumps',
                1028      & ' has moved to "data.cal"'
b9306711eb Jean*1029        CALL PRINT_ERROR( msgBuf, myThid )
742cf4499c Jean*1030       ENDIF
46979d29da Jean*1031 
7a7a4899b4 Chri*1032 C     Process "timestepping" params
                1033 C     o Time step size
062a876ce5 Jean*1034       IF ( deltaTtracer .NE. dTtracerLev(1) .AND.
                1035      &     deltaTtracer .NE. 0. .AND. dTtracerLev(1) .NE. 0. ) THEN
                1036         WRITE(msgBuf,'(A)')
                1037      &   'S/R INI_PARMS: deltaTtracer & dTtracerLev(1) not equal'
b9306711eb Jean*1038         CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1039         errCount = errCount + 1
062a876ce5 Jean*1040       ELSEIF ( dTtracerLev(1) .NE. 0. ) THEN
                1041         deltaTtracer = dTtracerLev(1)
                1042       ENDIF
ad330e45ef Jean*1043       IF ( deltaT       .EQ. 0. ) deltaT       = deltaTClock
7a7a4899b4 Chri*1044       IF ( deltaT       .EQ. 0. ) deltaT       = deltaTtracer
db8d49beca Jean*1045       IF ( deltaT       .EQ. 0. ) deltaT       = deltaTMom
                1046       IF ( deltaT       .EQ. 0. ) deltaT       = deltaTFreeSurf
90aed06079 Jean*1047       IF ( deltaT       .EQ. 0. ) THEN
                1048         WRITE(msgBuf,'(2A)') 'S/R INI_PARMS: ',
                1049      &   'need to specify in file "data", namelist "PARM03"'
                1050         CALL PRINT_ERROR( msgBuf, myThid )
                1051         WRITE(msgBuf,'(2A)') 'S/R INI_PARMS: ',
                1052      &   ' a model timestep (in s) deltaT or deltaTClock= ?'
                1053         CALL PRINT_ERROR( msgBuf, myThid )
                1054         errCount = errCount + 1
                1055         deltaT = 1.
                1056       ENDIF
db8d49beca Jean*1057       IF ( deltaTMom    .EQ. 0. ) deltaTMom    = deltaT
7a7a4899b4 Chri*1058       IF ( deltaTtracer .EQ. 0. ) deltaTtracer = deltaT
66dc79a095 Chri*1059       IF ( deltaTClock  .EQ. 0. ) deltaTClock  = deltaT
062a876ce5 Jean*1060       DO k=1,Nr
                1061        IF (dTtracerLev(k).EQ.0.) dTtracerLev(k)= deltaTtracer
                1062       ENDDO
ad330e45ef Jean*1063 C Note that this line should set deltaFreesurf=deltaTtracer
d48ed2b20a Alis*1064 C but this would change a lot of existing set-ups so we are
                1065 C obliged to set the default inappropriately.
                1066 C Be advised that when using asynchronous time stepping
                1067 C it is better to set deltaTreesurf=deltaTtracer
db8d49beca Jean*1068       IF ( deltaTFreeSurf .EQ. 0. ) deltaTFreeSurf  = deltaTMom
4edf45584c Alis*1069       IF ( periodicExternalForcing ) THEN
                1070        IF ( externForcingCycle*externForcingPeriod .EQ. 0. ) THEN
                1071         WRITE(msgBuf,'(A)')
                1072      &   'S/R INI_PARMS: externForcingCycle,externForcingPeriod =0'
b9306711eb Jean*1073         CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1074         errCount = errCount + 1
                1075        ELSEIF ( INT(externForcingCycle/externForcingPeriod) .NE.
                1076      &              externForcingCycle/externForcingPeriod ) THEN
4edf45584c Alis*1077         WRITE(msgBuf,'(A)')
                1078      &   'S/R INI_PARMS: externForcingCycle <> N*externForcingPeriod'
b9306711eb Jean*1079         CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1080         errCount = errCount + 1
                1081        ELSEIF ( externForcingCycle.LT.externForcingPeriod ) THEN
4edf45584c Alis*1082         WRITE(msgBuf,'(A)')
                1083      &   'S/R INI_PARMS: externForcingCycle < externForcingPeriod'
b9306711eb Jean*1084         CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1085         errCount = errCount + 1
4edf45584c Alis*1086        ENDIF
db8d49beca Jean*1087        IF ( externForcingPeriod.LT.deltaTClock ) THEN
4edf45584c Alis*1088         WRITE(msgBuf,'(A)')
db8d49beca Jean*1089      &   'S/R INI_PARMS: externForcingPeriod < deltaTClock'
b9306711eb Jean*1090         CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1091         errCount = errCount + 1
4edf45584c Alis*1092        ENDIF
                1093       ENDIF
c2b6ed6bfd Jean*1094 C     o Adams-Bashforth time stepping:
b9306711eb Jean*1095       IF ( momForcingOutAB  .EQ. UNSET_I ) THEN
c2b6ed6bfd Jean*1096         momForcingOutAB  = 1
                1097         IF ( forcing_In_AB ) momForcingOutAB  = 0
                1098       ENDIF
b9306711eb Jean*1099       IF ( tracForcingOutAB .EQ. UNSET_I ) THEN
c2b6ed6bfd Jean*1100         tracForcingOutAB = 1
                1101         IF ( forcing_In_AB ) tracForcingOutAB = 0
                1102       ENDIF
cf8488c0fd Chri*1103 C     o Convection frequency
                1104       IF ( cAdjFreq .LT. 0. ) THEN
                1105        cAdjFreq = deltaTClock
                1106       ENDIF
d4701cb6da Alis*1107       IF ( ivdc_kappa .NE. 0. .AND. cAdjFreq .NE. 0. ) THEN
                1108        WRITE(msgBuf,'(A,A)')
                1109      &  'S/R INI_PARMS: You have enabled both ivdc_kappa and',
                1110      &  ' convective adjustment.'
b9306711eb Jean*1111        CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1112        errCount = errCount + 1
d4701cb6da Alis*1113       ENDIF
a97467b673 Jean*1114       IF (useCDscheme) THEN
                1115 C     o CD coupling (CD scheme):
db8d49beca Jean*1116         IF ( tauCD .EQ. 0. _d 0 ) tauCD = deltaTMom
a97467b673 Jean*1117         IF ( rCD .LT. 0. ) rCD = 1. _d 0 - deltaTMom/tauCD
8039e9b985 Jean*1118         IF ( epsAB_CD .EQ. UNSET_RL ) epsAB_CD = abEps
b05b067368 Chri*1119       ENDIF
88830be691 Alis*1120 
19eda8b250 Jean*1121       IF ( startTime.EQ.UNSET_RL .AND. nIter0.EQ.-1 ) THEN
                1122 C     o set default value for start time & nIter0
                1123         startTime = baseTime
                1124         nIter0 = 0
                1125       ELSEIF ( startTime.EQ.UNSET_RL ) THEN
                1126 C     o set start time from nIter0
d987e9fc35 Jean*1127          startTime = baseTime + deltaTClock*DFLOAT(nIter0)
19eda8b250 Jean*1128       ELSEIF ( nIter0.EQ.-1 ) THEN
                1129 C     o set nIter0 from start time
                1130          nIter0 = NINT( (startTime-baseTime)/deltaTClock )
                1131       ELSEIF ( baseTime.EQ.0. ) THEN
                1132 C     o set base time from the 2 others
d987e9fc35 Jean*1133          baseTime = startTime - deltaTClock*DFLOAT(nIter0)
19eda8b250 Jean*1134       ENDIF
d4701cb6da Alis*1135 
6e6f314aa7 Patr*1136       nTimeSteps_l2 = 4
d4701cb6da Alis*1137 C     o nTimeSteps 1
                1138       IF ( nTimeSteps .EQ. 0 .AND. nEndIter .NE. 0 )
                1139      &     nTimeSteps = nEndIter-nIter0
                1140 C     o nTimeSteps 2
88830be691 Alis*1141       IF ( nTimeSteps .EQ. 0 .AND. endTime .NE. 0. )
db8d49beca Jean*1142      &     nTimeSteps = NINT((endTime-startTime)/deltaTClock)
d4701cb6da Alis*1143 C     o nEndIter 1
                1144       IF ( nEndIter .EQ. 0 .AND. nTimeSteps .NE. 0 )
                1145      &     nEndIter = nIter0+nTimeSteps
                1146 C     o nEndIter 2
                1147       IF ( nEndIter .EQ. 0 .AND. endTime .NE. 0. )
db8d49beca Jean*1148      &     nEndIter = NINT((endTime-baseTime)/deltaTClock)
d4701cb6da Alis*1149 C     o End Time 1
                1150       IF ( endTime .EQ. 0. .AND. nTimeSteps .NE. 0 )
d987e9fc35 Jean*1151      &     endTime = startTime + deltaTClock*DFLOAT(nTimeSteps)
d4701cb6da Alis*1152 C     o End Time 2
                1153       IF ( endTime .EQ. 0. .AND. nEndIter .NE. 0 )
d987e9fc35 Jean*1154      &     endTime = baseTime + deltaTClock*DFLOAT(nEndIter)
d4701cb6da Alis*1155 
88830be691 Alis*1156 C     o Consistent?
d987e9fc35 Jean*1157       IF ( startTime .NE. baseTime+deltaTClock*DFLOAT(nIter0) ) THEN
16c445e700 Jean*1158        WRITE(msgBuf,'(A)')
                1159      & 'S/R INI_PARMS: startTime, baseTime and nIter0 are inconsistent'
b9306711eb Jean*1160        CALL PRINT_ERROR( msgBuf, myThid )
16c445e700 Jean*1161        WRITE(msgBuf,'(A)')
                1162      & 'S/R INI_PARMS: Perhaps more than two were set at once'
b9306711eb Jean*1163        CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1164        errCount = errCount + 1
16c445e700 Jean*1165       ENDIF
d4701cb6da Alis*1166       IF ( nEndIter .NE. nIter0+nTimeSteps ) THEN
                1167        WRITE(msgBuf,'(A)')
                1168      & 'S/R INI_PARMS: nIter0, nTimeSteps and nEndIter are inconsistent'
b9306711eb Jean*1169        CALL PRINT_ERROR( msgBuf, myThid )
d4701cb6da Alis*1170        WRITE(msgBuf,'(A)')
                1171      & 'S/R INI_PARMS: Perhaps more than two were set at once'
b9306711eb Jean*1172        CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1173        errCount = errCount + 1
d4701cb6da Alis*1174       ENDIF
b9306711eb Jean*1175       IF ( nTimeSteps .NE. NINT((endTime-startTime)/deltaTClock)
                1176      &   ) THEN
                1177        WRITE(msgBuf,'(A)')
                1178      & 'S/R INI_PARMS: both endTime and nTimeSteps have been set'
                1179        CALL PRINT_ERROR( msgBuf, myThid )
                1180        WRITE(msgBuf,'(A)')
                1181      & 'S/R INI_PARMS: but are inconsistent'
                1182        CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1183        errCount = errCount + 1
7a7a4899b4 Chri*1184       ENDIF
ea212a0263 Alis*1185 
dc684458c1 Alis*1186 C     o Monitor (should also add CPP flag for monitor?)
                1187       IF (monitorFreq.LT.0.) THEN
                1188        monitorFreq=0.
7abc7f9caf Alis*1189        IF (dumpFreq.NE.0.) monitorFreq=dumpFreq
3a279374db Alis*1190        IF (diagFreq.NE.0..AND.diagFreq.LT.monitorFreq)
                1191      &         monitorFreq=diagFreq
7abc7f9caf Alis*1192        IF (chkPtFreq.NE.0..AND.chkPtFreq.LT.monitorFreq)
                1193      &         monitorFreq=chkPtFreq
                1194        IF (pChkPtFreq.NE.0..AND.pChkPtFreq.LT.monitorFreq)
                1195      &         monitorFreq=pChkPtFreq
db8d49beca Jean*1196        IF (monitorFreq.EQ.0.) monitorFreq=deltaTClock
dc684458c1 Alis*1197       ENDIF
f804abbd25 Jean*1198       IF ( monitorSelect.EQ.UNSET_I ) THEN
                1199         monitorSelect = 2
                1200         IF ( fluidIsWater ) monitorSelect = 3
                1201       ENDIF
dc684458c1 Alis*1202 
924557e60a Chri*1203 C--   Grid parameters
                1204 C     In cartesian coords distances are in metres
7cd7754949 Jean*1205       DO k =1,Nr
                1206        delZ(k) = UNSET_RL
                1207        delP(k) = UNSET_RL
                1208        delR(k) = UNSET_RL
924557e60a Chri*1209       ENDDO
                1210 C     In spherical polar distances are in degrees
0127add478 Alis*1211       dxSpacing = UNSET_RL
                1212       dySpacing = UNSET_RL
bb2fd3f1ad Jean*1213 C-    pCell-Mix  parameters:
                1214       interViscAr_pCell = .FALSE.
                1215       interDiffKr_pCell = .FALSE.
                1216       pCellMix_select = 0
                1217       pCellMix_maxFac = 1. _d 4
                1218       pCellMix_delR = 0.
                1219       DO k=1,Nr
                1220         pCellMix_viscAr(k) = viscArNr(k)
                1221         pCellMix_diffKr(k) = diffKrNrT(k)
                1222       ENDDO
                1223 
aa8f58c33d Jean*1224       WRITE(msgBuf,'(A)') ' INI_PARMS ; starts to read PARM04'
6b6f69d6b0 Jean*1225       CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
b9306711eb Jean*1226      &                    SQUEEZE_RIGHT, myThid )
b6e98c1b91 Jean*1227       READ(UNIT=iUnit,NML=PARM04) !,IOSTAT=errIO)
c206b07c1d Chri*1228       IF ( errIO .LT. 0 ) THEN
924557e60a Chri*1229        WRITE(msgBuf,'(A)')
aa8f58c33d Jean*1230      &  'S/R INI_PARMS: Error reading model parameter file "data"'
b9306711eb Jean*1231        CALL PRINT_ERROR( msgBuf, myThid )
aa8f58c33d Jean*1232        WRITE(msgBuf,'(A)') 'S/R INI_PARMS: Problem in namelist PARM04'
b9306711eb Jean*1233        CALL PRINT_ERROR( msgBuf, myThid )
924557e60a Chri*1234        STOP 'ABNORMAL END: S/R INI_PARMS'
41f1ac6de1 Jean*1235       ELSE
aa8f58c33d Jean*1236        WRITE(msgBuf,'(A)') ' INI_PARMS ; read PARM04 : OK'
41f1ac6de1 Jean*1237        CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
b9306711eb Jean*1238      &                     SQUEEZE_RIGHT, myThid )
809bdccbfc Jean*1239       ENDIF
9bec0e0933 Alis*1240 
da40681a8e Jean*1241 C     Check for retired parameters still being used
2ad79bdf32 Jean*1242       IF ( Ro_SeaLevel .NE. UNSET_RL ) THEN
                1243 c      nRetired = nRetired+1
                1244        IF ( usingPCoords ) THEN
                1245         WRITE(msgBuf,'(2A)') '** WARNING ** INI_PARMS: ',
                1246      &   '"Ro_SeaLevel" (P @ bottom) depreciated (backward compat'
                1247        ELSEIF ( usingZCoords ) THEN
                1248         WRITE(msgBuf,'(2A)') '** WARNING ** INI_PARMS: ',
                1249      &   '"Ro_SeaLevel" (Z @ top) depreciated (backward compat'
                1250        ENDIF
                1251        CALL PRINT_MESSAGE( msgBuf, errorMessageUnit,
                1252      &                     SQUEEZE_RIGHT, myThid )
                1253        IF ( usingPCoords ) THEN
                1254         WRITE(msgBuf,'(2A)') '** WARNING ** INI_PARMS: ',
                1255      &   ' only). To set vert. axis, use instead "top_Pres".'
                1256        ELSEIF ( usingZCoords ) THEN
                1257         WRITE(msgBuf,'(2A)') '** WARNING ** INI_PARMS: ',
                1258      &   ' only). To set vert. axis, use instead "seaLev_Z".'
                1259        ENDIF
                1260        CALL PRINT_MESSAGE( msgBuf, errorMessageUnit,
                1261      &                     SQUEEZE_RIGHT, myThid )
                1262       ENDIF
da40681a8e Jean*1263       IF ( rkFac .NE. UNSET_RL ) THEN
                1264        nRetired = nRetired+1
                1265        WRITE(msgBuf,'(A,A)')
                1266      &  'S/R INI_PARMS: "rkFac" has been replaced by -rkSign ',
                1267      &  ' and is no longer allowed in file "data".'
b9306711eb Jean*1268        CALL PRINT_ERROR( msgBuf, myThid )
da40681a8e Jean*1269       ENDIF
                1270       IF ( groundAtK1 ) THEN
                1271 c      nRetired = nRetired+1
                1272        WRITE(msgBuf,'(A,A)')
                1273      &  'S/R INI_PARMS: "groundAtK1" is set according to vertical ',
                1274      &  ' coordinate and is no longer allowed in file "data".'
b9306711eb Jean*1275        CALL PRINT_ERROR( msgBuf, myThid )
da40681a8e Jean*1276       ENDIF
9780090eaa Jean*1277       IF ( thetaMin .NE. UNSET_RL ) THEN
                1278        nRetired = nRetired+1
                1279        WRITE(msgBuf,'(A,A)')
                1280      &  'S/R INI_PARMS: "thetaMin" no longer allowed,',
                1281      &  ' has been replaced by "xgOrigin"'
                1282        CALL PRINT_ERROR( msgBuf, myThid )
                1283       ENDIF
                1284       IF ( phiMin .NE. UNSET_RL ) THEN
                1285        nRetired = nRetired+1
                1286        WRITE(msgBuf,'(A,A)')
                1287      &  'S/R INI_PARMS: "phiMin" no longer allowed,',
                1288      &  ' has been replaced by "ygOrigin"'
                1289        CALL PRINT_ERROR( msgBuf, myThid )
                1290       ENDIF
da40681a8e Jean*1291 
f76021331f Jean*1292 C     X coordinate : Check for multiple definitions
                1293       goptCount = 0
                1294       IF ( delX(1)   .NE. UNSET_RL ) goptCount = goptCount + 1
                1295       IF ( dxSpacing .NE. UNSET_RL ) goptCount = goptCount + 1
                1296       IF ( delXFile  .NE. ' ' )      goptCount = goptCount + 1
                1297       IF ( goptCount.GT.1 ) THEN
b9306711eb Jean*1298        WRITE(msgBuf,'(A,A)') 'Too many specifications for delX:',
9bec0e0933 Alis*1299      &   'Specify only one of delX, dxSpacing or delXfile'
b9306711eb Jean*1300        CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1301        errCount = errCount + 1
9bec0e0933 Alis*1302       ENDIF
0127add478 Alis*1303       IF ( dxSpacing .NE. UNSET_RL ) THEN
b9b591469d Jean*1304        DO i=1,gridNx
0127add478 Alis*1305         delX(i) = dxSpacing
                1306        ENDDO
                1307       ENDIF
f76021331f Jean*1308 C     Y coordinate : Check for multiple definitions
                1309       goptCount = 0
                1310       IF ( delY(1)   .NE. UNSET_RL ) goptCount = goptCount + 1
                1311       IF ( dySpacing .NE. UNSET_RL ) goptCount = goptCount + 1
                1312       IF ( delYFile  .NE. ' ' )      goptCount = goptCount + 1
                1313       IF ( goptCount.GT.1 ) THEN
b9306711eb Jean*1314        WRITE(msgBuf,'(A,A)') 'Too many specifications for delY:',
9bec0e0933 Alis*1315      &   'Specify only one of delY, dySpacing or delYfile'
b9306711eb Jean*1316        CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1317        errCount = errCount + 1
9bec0e0933 Alis*1318       ENDIF
0127add478 Alis*1319       IF ( dySpacing .NE. UNSET_RL ) THEN
b9b591469d Jean*1320        DO j=1,gridNy
f76021331f Jean*1321         delY(j) = dySpacing
0127add478 Alis*1322        ENDDO
                1323       ENDIF
2c7e43108f Jean*1324 
910f05e765 Chri*1325 C--   Check for conflicting grid definitions.
924557e60a Chri*1326       goptCount = 0
                1327       IF ( usingCartesianGrid )      goptCount = goptCount+1
                1328       IF ( usingSphericalPolarGrid ) goptCount = goptCount+1
aea29c8517 Alis*1329       IF ( usingCurvilinearGrid )    goptCount = goptCount+1
0ac260a803 Andr*1330       IF ( usingCylindricalGrid )    goptCount = goptCount+1
aea29c8517 Alis*1331       IF ( goptCount .GT. 1 ) THEN
924557e60a Chri*1332        WRITE(msgBuf,'(A)')
                1333      &  'S/R INI_PARMS: More than one coordinate system requested'
b9306711eb Jean*1334        CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1335        errCount = errCount + 1
924557e60a Chri*1336       ENDIF
aea29c8517 Alis*1337       IF ( goptCount .LT. 1 ) THEN
f2459ca318 Jean*1338 C-     No horizontal grid is specified => use Cartesian grid as default:
aea29c8517 Alis*1339        WRITE(msgBuf,'(A)')
f2459ca318 Jean*1340      &  'S/R INI_PARMS: No horizontal grid requested'
                1341        CALL PRINT_MESSAGE( msgBuf, errorMessageUnit,
b9306711eb Jean*1342      &                     SQUEEZE_RIGHT, myThid )
f2459ca318 Jean*1343        WRITE(msgBuf,'(A)')
                1344      &  'S/R INI_PARMS: => Use Cartesian Grid as default'
                1345        CALL PRINT_MESSAGE( msgBuf, errorMessageUnit,
b9306711eb Jean*1346      &                     SQUEEZE_RIGHT, myThid )
f2459ca318 Jean*1347        usingCartesianGrid = .TRUE.
aea29c8517 Alis*1348       ENDIF
9744e36521 Jean*1349 C-    Radius of the Planet:
                1350       IF ( rSphere.EQ.UNSET_RL ) THEN
                1351         IF ( usingCurvilinearGrid .AND.
                1352      &       radius_fromHorizGrid.NE.UNSET_RL ) THEN
                1353            rSphere = radius_fromHorizGrid
                1354         ELSE
                1355            rSphere = 6370. _d 3
                1356         ENDIF
                1357       ENDIF
                1358       IF ( radius_fromHorizGrid.EQ.UNSET_RL ) THEN
                1359            radius_fromHorizGrid = rSphere
                1360       ENDIF
                1361       IF ( rSphere .NE. 0. ) THEN
                1362        recip_rSphere = 1. _d 0/rSphere
                1363       ELSE
                1364        recip_rSphere = 0.
                1365       ENDIF
427e24e121 Jean*1366 C--   Grid rotation
                1367       IF ( phiEuler .NE. 0. _d 0 .OR. thetaEuler .NE. 0. _d 0
                1368      &     .OR. psiEuler .NE. 0. _d 0 ) rotateGrid = .TRUE.
                1369 
2ad79bdf32 Jean*1370 C--   Default vertical axis origin
                1371       IF ( Ro_SeaLevel .NE. UNSET_RL ) THEN
                1372        IF ( usingPCoords .AND. top_Pres.NE.UNSET_RL ) THEN
ef7fa20135 Jean*1373         WRITE(msgBuf,'(2A)') 'S/R INI_PARMS: ',
                1374      &       'Cannot set both "Ro_SeaLevel" and "top_Pres"'
                1375         CALL PRINT_ERROR( msgBuf, myThid )
                1376         errCount = errCount + 1
2ad79bdf32 Jean*1377        ENDIF
                1378        IF ( usingZCoords .AND. seaLev_Z.NE.UNSET_RL ) THEN
ef7fa20135 Jean*1379         WRITE(msgBuf,'(2A)') 'S/R INI_PARMS: ',
                1380      &       'Cannot set both "Ro_SeaLevel" and "seaLev_Z"'
                1381         CALL PRINT_ERROR( msgBuf, myThid )
                1382         errCount = errCount + 1
2ad79bdf32 Jean*1383        ENDIF
                1384        rF(1) = Ro_SeaLevel
                1385       ELSE
                1386        rF(1) = UNSET_RS
                1387       ENDIF
                1388       IF ( top_Pres.EQ.UNSET_RL ) top_Pres = 0.
                1389       IF ( seaLev_Z.EQ.UNSET_RL ) seaLev_Z = 0.
9780090eaa Jean*1390 C--   Default origin for  X & Y axis (longitude & latitude origin):
                1391       IF ( xgOrigin .EQ. UNSET_RL ) xgOrigin = 0.
                1392       IF ( ygOrigin .EQ. UNSET_RL ) THEN
                1393         IF ( usingSphericalPolarGrid )  THEN
                1394           ygOrigin = 0.
                1395         ELSEIF ( usingCartesianGrid )   THEN
                1396           ygOrigin = 0.
                1397         ELSEIF ( usingCylindricalGrid ) THEN
                1398           ygOrigin = 0.
                1399         ELSEIF ( usingCurvilinearGrid ) THEN
                1400           ygOrigin = 0.
                1401         ELSE
                1402           WRITE(msgBuf,'(A)')
                1403      &    'S/R INI_PARMS: found no coordinate system to set ygOrigin'
                1404           CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1405           errCount = errCount + 1
9780090eaa Jean*1406         ENDIF
                1407       ENDIF
7fbd8d8c1c Jean*1408 
c28ce1627a Jean*1409 C--   set cell Center depth and put Interface at the middle between 2 C
                1410       setCenterDr = .FALSE.
7cd7754949 Jean*1411       DO k=1,Nr+1
                1412        IF ( delRc(k).EQ.UNSET_RL ) THEN
b9306711eb Jean*1413         IF ( setCenterDr ) THEN
                1414          WRITE(msgBuf,'(A,I4)')
7cd7754949 Jean*1415      &        'S/R INI_PARMS: No value for delRc at k =', k
b9306711eb Jean*1416          CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1417          errCount = errCount + 1
b9306711eb Jean*1418         ENDIF
                1419        ELSE
                1420         IF ( k.EQ.1 ) setCenterDr = .TRUE.
                1421         IF ( .NOT.setCenterDr ) THEN
                1422          WRITE(msgBuf,'(A,I4)')
7cd7754949 Jean*1423      &        'S/R INI_PARMS: No value for delRc at k <', k
b9306711eb Jean*1424          CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1425          errCount = errCount + 1
b9306711eb Jean*1426         ENDIF
                1427        ENDIF
c28ce1627a Jean*1428       ENDDO
b9306711eb Jean*1429       IF ( setCenterDr ) rCoordInputData = .TRUE.
910f05e765 Chri*1430 C--   p, z, r coord parameters
b9306711eb Jean*1431       setInterFDr = .FALSE.
7cd7754949 Jean*1432       DO k = 1, Nr
                1433        IF ( delZ(k) .NE. UNSET_RL ) zCoordInputData = .TRUE.
                1434        IF ( delP(k) .NE. UNSET_RL ) pCoordInputData = .TRUE.
                1435        IF ( delR(k) .NE. UNSET_RL ) rCoordInputData = .TRUE.
                1436        IF ( delR(k) .EQ. UNSET_RL ) delR(k) = delZ(k)
                1437        IF ( delR(k) .EQ. UNSET_RL ) delR(k) = delP(k)
                1438        IF ( delR(k) .EQ. UNSET_RL ) THEN
b9306711eb Jean*1439         IF ( setInterFDr ) THEN
88830be691 Alis*1440          WRITE(msgBuf,'(A,I4)')
7cd7754949 Jean*1441      &        'S/R INI_PARMS: No value for delZ/delP/delR at k =', k
b9306711eb Jean*1442          CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1443          errCount = errCount + 1
b9306711eb Jean*1444         ENDIF
                1445        ELSE
                1446         IF ( k.EQ.1 ) setInterFDr = .TRUE.
                1447         IF ( .NOT.setInterFDr ) THEN
                1448          WRITE(msgBuf,'(A,I4)')
7cd7754949 Jean*1449      &        'S/R INI_PARMS: No value for delZ/delP/delR at k <', k
b9306711eb Jean*1450          CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1451          errCount = errCount + 1
b9306711eb Jean*1452         ENDIF
88830be691 Alis*1453        ENDIF
910f05e765 Chri*1454       ENDDO
                1455 C     Check for multiple coordinate systems
b9306711eb Jean*1456       coordsSet = 0
910f05e765 Chri*1457       IF ( zCoordInputData ) coordsSet = coordsSet + 1
                1458       IF ( pCoordInputData ) coordsSet = coordsSet + 1
                1459       IF ( rCoordInputData ) coordsSet = coordsSet + 1
                1460       IF ( coordsSet .GT. 1 ) THEN
                1461        WRITE(msgBuf,'(A)')
                1462      &  'S/R INI_PARMS: Cannot mix z, p and r in the input data.'
b9306711eb Jean*1463        CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1464        errCount = errCount + 1
910f05e765 Chri*1465       ENDIF
b9306711eb Jean*1466 C-    Check for double definition (file & namelist)
                1467       IF ( delRcFile.NE.' ' ) THEN
                1468         IF ( setCenterDr ) THEN
                1469           WRITE(msgBuf,'(A)')
                1470      &         'S/R INI_PARMS: Cannot set both delRc and delRcFile'
                1471           CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1472           errCount = errCount + 1
b9306711eb Jean*1473         ENDIF
                1474         setCenterDr = .TRUE.
                1475       ENDIF
                1476       IF ( delRFile.NE.' ' ) THEN
                1477         IF ( setInterFDr ) THEN
                1478           WRITE(msgBuf,'(A)')
                1479      &         'S/R INI_PARMS: Cannot set both delR and delRFile'
                1480           CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1481           errCount = errCount + 1
b9306711eb Jean*1482         ENDIF
                1483         setInterFDr = .TRUE.
                1484       ENDIF
b05b067368 Chri*1485 
81bc00c2f0 Chri*1486 C--   Input files
aa8f58c33d Jean*1487       WRITE(msgBuf,'(A)') ' INI_PARMS ; starts to read PARM05'
6b6f69d6b0 Jean*1488       CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
b9306711eb Jean*1489      &                    SQUEEZE_RIGHT, myThid )
88830be691 Alis*1490       READ(UNIT=iUnit,NML=PARM05) !,IOSTAT=errIO)
809bdccbfc Jean*1491       IF ( errIO .LT. 0 ) THEN
81bc00c2f0 Chri*1492        WRITE(msgBuf,'(A)')
aa8f58c33d Jean*1493      &  'S/R INI_PARMS: Error reading model parameter file "data"'
b9306711eb Jean*1494        CALL PRINT_ERROR( msgBuf, myThid )
aa8f58c33d Jean*1495        WRITE(msgBuf,'(A)') 'S/R INI_PARMS: Problem in namelist PARM05'
b9306711eb Jean*1496        CALL PRINT_ERROR( msgBuf, myThid )
81bc00c2f0 Chri*1497        STOP 'ABNORMAL END: S/R INI_PARMS'
41f1ac6de1 Jean*1498       ELSE
aa8f58c33d Jean*1499        WRITE(msgBuf,'(A)') ' INI_PARMS ; read PARM05 : OK'
41f1ac6de1 Jean*1500        CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
b9306711eb Jean*1501      &                     SQUEEZE_RIGHT, myThid )
809bdccbfc Jean*1502       ENDIF
8194bc4a99 Mart*1503 C     Check for retired parameters still being used
                1504       IF ( shelfIceFile .NE. ' ' ) THEN
                1505        nRetired = nRetired+1
                1506        WRITE(msgBuf,'(A,A)')
                1507      &  'S/R INI_PARMS: "shelfIceFile" is not allowed in "data", ',
                1508      &  'substitute "SHELFICEtopoFile" in data.shelfice'
                1509        CALL PRINT_ERROR( msgBuf, myThid )
                1510       ENDIF
275ad78b78 Jean*1511       IF ( dQdTFile .NE. ' ' ) THEN
                1512        nRetired = nRetired+1
                1513        WRITE(msgBuf,'(A,A)')
                1514      &  'S/R INI_PARMS: "dQdTFile" has been retired from file "data"'
                1515        CALL PRINT_ERROR( msgBuf, myThid )
                1516       ENDIF
cbc417d429 Mart*1517 C     Check if two conflicting I/O directories have been specified
                1518       IF (mdsioLocalDir .NE. ' ' .AND. adTapeDir .NE. ' ') THEN
                1519          WRITE(msgBuf,'(A)')
                1520      &   'S/R INI_PARMS: mdsioLocalDir and adTapeDir cannot be'
                1521          CALL PRINT_ERROR( msgBuf, myThid )
                1522          WRITE(msgBuf,'(A)')
                1523      &   'S/R INI_PARMS: specified at the same time'
                1524          CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1525          errCount = errCount + 1
cbc417d429 Mart*1526       ENDIF
81bc00c2f0 Chri*1527 
7b50f71112 Gael*1528 C     Check for relaxation term:
                1529        IF ( (tauThetaClimRelax.GT.0.).AND.
                1530      &      (thetaClimFile.EQ.' ') ) THEN
                1531          WRITE(msgBuf,'(A)')
                1532      &   'S/R INI_PARMS: tauThetaClimRelax > 0 but'
                1533          CALL PRINT_ERROR( msgBuf, myThid )
                1534          WRITE(msgBuf,'(A)')
                1535      &   'S/R INI_PARMS: thetaClimFile is undefined'
                1536          CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1537          errCount = errCount + 1
7b50f71112 Gael*1538        ENDIF
                1539        IF ( (tauSaltClimRelax.GT.0.).AND.
                1540      &      (saltClimFile.EQ.' ') ) THEN
                1541          WRITE(msgBuf,'(A)')
                1542      &   'S/R INI_PARMS: tauSaltClimRelax > 0 but'
                1543          CALL PRINT_ERROR( msgBuf, myThid )
                1544          WRITE(msgBuf,'(A)')
                1545      &   'S/R INI_PARMS: saltClimFile is undefined'
                1546          CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1547          errCount = errCount + 1
7b50f71112 Gael*1548        ENDIF
13b4d13c98 Jean*1549 C     Check vertical diffusivity setting:
                1550 #ifdef ALLOW_3D_DIFFKR
                1551        IF ( specifiedDiffKrT ) THEN
3fcdbd6e66 Jean*1552          WRITE(msgBuf,'(2A)') '** WARNING ** INI_PARMS: Ignores diffKr',
                1553      &   'T (or Kp,Kz) setting in file "data" with ALLOW_3D_DIFFKR'
                1554          CALL PRINT_MESSAGE( msgBuf, errorMessageUnit,
                1555      &                       SQUEEZE_RIGHT, myThid )
13b4d13c98 Jean*1556        ENDIF
                1557        IF ( specifiedDiffKrS .AND. diffKrFile.NE.' ' ) THEN
3fcdbd6e66 Jean*1558          WRITE(msgBuf,'(2A)') '** WARNING ** INI_PARMS: Ignores diffKr',
                1559      &   'S (or Kp,Kz) setting in file "data" and uses diffKrFile'
                1560          CALL PRINT_MESSAGE( msgBuf, errorMessageUnit,
                1561      &                       SQUEEZE_RIGHT, myThid )
13b4d13c98 Jean*1562        ENDIF
                1563 #endif
7b50f71112 Gael*1564 
599be48d9a Jean*1565 C--   Set Units conversion factor required to incorporate
                1566 C     surface forcing into z-p isomorphic equations:
427e24e121 Jean*1567 C     mass2rUnit: from mass per unit area [kg/m2] to r-coordinate
                1568 C     rUnit2mass: from r-coordinate to mass per unit area [kg/m2]
e99e7fe4e6 Jean*1569       IF ( usingPCoords ) THEN
599be48d9a Jean*1570          mass2rUnit = gravity
                1571          rUnit2mass = recip_gravity
                1572       ELSE
                1573          mass2rUnit = recip_rhoConst
                1574          rUnit2mass = rhoConst
9800a0be87 Jean*1575       ENDIF
688d11fba8 Alis*1576 
80d98e0151 Dimi*1577 C--   For backward compatibility, set temp_addMass and salt_addMass
b46f9da188 Jean*1578 C     to temp_EvPrRn and salt_EvPrRn if not set in parameter file "data"
80d98e0151 Dimi*1579       IF (temp_addMass .EQ. UNSET_RL) temp_addMass = temp_EvPrRn
                1580       IF (salt_addMass .EQ. UNSET_RL) salt_addMass = salt_EvPrRn
                1581 
b46f9da188 Jean*1582 C--   Make a choice (if unset) regarding using CG2D solver min.-residual sol.
                1583 C      for simple set-up (cartesian grid + flat bottom), default is to
                1584 C      use the solver minimum residual solution (cg2dUseMinResSol=1).
                1585       IF ( cg2dUseMinResSol.EQ.UNSET_I ) THEN
                1586         cg2dUseMinResSol = 0
                1587         IF ( topoFile.EQ.' ' .AND. bathyFile.EQ.' '
                1588      &       .AND. usingCartesianGrid ) cg2dUseMinResSol = 1
                1589       ENDIF
                1590 
aa8f58c33d Jean*1591 C--   Close the open data file
a6485b3244 Mart*1592 #ifdef SINGLE_DISK_IO
9ba4d7e4bd Chri*1593       CLOSE(iUnit)
a6485b3244 Mart*1594 #else
                1595       CLOSE(iUnit,STATUS='DELETE')
                1596 #endif /* SINGLE_DISK_IO */
9ba4d7e4bd Chri*1597 
aa8f58c33d Jean*1598       WRITE(msgBuf,'(A)') ' INI_PARMS: finished reading file "data"'
                1599       CALL PRINT_MESSAGE( msgBuf, standardMessageUnit,
                1600      &                    SQUEEZE_RIGHT , myThid )
                1601 
ede2cfe569 Chri*1602 C--   Check whether any retired parameters were found.
46979d29da Jean*1603       IF ( nRetired .GT. 0 ) THEN
ede2cfe569 Chri*1604        WRITE(msgBuf,'(A)')
90aed06079 Jean*1605      &  'S/R INI_PARMS: Error reading parameter file "data":'
b9306711eb Jean*1606        CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1607        WRITE(msgBuf,'(I4,A)') nRetired,
                1608      &      ' out of date parameters were found in the namelist'
b9306711eb Jean*1609        CALL PRINT_ERROR( msgBuf, myThid )
90aed06079 Jean*1610        errCount = errCount + 1
                1611       ENDIF
                1612 C--   Stop if any error was found (including retired params):
                1613       IF ( errCount .GE. 1 ) THEN
                1614         WRITE(msgBuf,'(A,I3,A)')
                1615      &       'S/R INI_PARMS: detected', errCount,' fatal error(s)'
                1616         CALL PRINT_ERROR( msgBuf, myThid )
                1617         CALL ALL_PROC_DIE( 0 )
                1618         STOP 'ABNORMAL END: S/R INI_PARMS'
ede2cfe569 Chri*1619       ENDIF
                1620 
924557e60a Chri*1621       _END_MASTER(myThid)
                1622 
                1623 C--   Everyone else must wait for the parameters to be loaded
                1624       _BARRIER
2c7e43108f Jean*1625 
924557e60a Chri*1626       RETURN
                1627       END