/[MITgcm]/MITgcm_contrib/snarayan/divided_adjoint/model/src/ini_parms.F
ViewVC logotype

Annotation of /MITgcm_contrib/snarayan/divided_adjoint/model/src/ini_parms.F

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


Revision 1.1 - (hide annotations) (download)
Tue May 19 21:32:23 2015 UTC (10 years, 2 months ago) by snarayan
Branch: MAIN
CVS Tags: HEAD
Divided adjoints

Files for the implementing divided adjoints using OpenAD

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

  ViewVC Help
Powered by ViewVC 1.1.22