/[MITgcm]/MITgcm/model/inc/PARAMS.h
ViewVC logotype

Contents of /MITgcm/model/inc/PARAMS.h

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


Revision 1.54 - (show annotations) (download)
Wed Jun 6 14:55:45 2001 UTC (22 years, 11 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint40pre1
Changes since 1.53: +3 -2 lines
File MIME type: text/plain
Added a debugMode that uses same statistics stuff as monitor.F
Can be disabled with -DEXCLUDE_DEBUGMODE. Turn on at run-time
with debugMode=.true.  Default is enabled but off.

1 C $Header: /u/gcmpack/models/MITgcmUV/model/inc/PARAMS.h,v 1.53 2001/05/30 19:33:18 adcroft Exp $
2 C $Name: $
3 C
4 C /==========================================================\
5 C | PARAMS.h |
6 C | o Header file defining model "parameters". |
7 C |==========================================================|
8 C | The values from the model standard input file are |
9 C | stored into the variables held here. Notes describing |
10 C | the parameters can also be found here. |
11 C \==========================================================/
12
13 C Macros for special grid options
14 #include "PARAMS_MACROS.h"
15
16 C-- Contants
17 C Useful physical values
18 Real*8 PI
19 PARAMETER ( PI = 3.14159265358979323844D0 )
20 Real*8 deg2rad
21 PARAMETER ( deg2rad = 2.D0*PI/360.D0 )
22
23 C Symbolic values
24 C precXXXX - Used to indicate what precision to use for
25 C dumping model state.
26 INTEGER precFloat32
27 PARAMETER ( precFloat32 = 32 )
28 INTEGER precFloat64
29 PARAMETER ( precFloat64 = 64 )
30 C UNSET_xxx - Used to indicate variables that have not been given a value
31 Real*8 UNSET_FLOAT8
32 PARAMETER ( UNSET_FLOAT8 = 1.234567D5 )
33 Real*4 UNSET_FLOAT4
34 PARAMETER ( UNSET_FLOAT4 = 1.234567E5 )
35 _RL UNSET_RL
36 PARAMETER ( UNSET_RL = 1.234567D5 )
37 _RS UNSET_RS
38 PARAMETER ( UNSET_RS = 1.234567E5 )
39 INTEGER UNSET_I
40 PARAMETER ( UNSET_I = 123456789 )
41
42 C Checkpoint data
43 INTEGER maxNoChkptLev
44 PARAMETER ( maxNoChkptLev = 2 )
45
46 C-- COMMON /PARM_C/ Character valued parameters used by the model.
47 C checkPtSuff - List of checkpoint file suffices
48 C bathyFile - File containing bathymetry. If not defined bathymetry
49 C is taken from inline function.
50 C hydrogThetaFile - File containing initial hydrographic data for potential
51 C temperature.
52 C hydrogSaltFile - File containing initial hydrographic data for salinity.
53 C zonalWindFile - File containing zonal wind data
54 C meridWindFile - File containing meridional wind data
55 C thetaClimFile - File containing theta climataology used
56 C in relaxation term -lambda(theta-theta*)
57 C saltClimFile - File containing salt climataology used
58 C in relaxation term -lambda(salt-salt*)
59 C surfQfile - File containing surface heat flux
60 C surfQswfile - File containing surface shortwave radiation
61 C dQdTfile - File containing thermal relaxation coefficient
62 C EmPmRfile - File containing surface fresh water flux
63 C buoyancyRelation - Flag used to indicate which relation to use to
64 C get buoyancy.
65 COMMON /PARM_C/ checkPtSuff,
66 & bathyFile, hydrogThetaFile, hydrogSaltFile,
67 & zonalWindFile, meridWindFile, thetaClimFile,
68 & saltClimFile, buoyancyRelation,
69 & EmPmRfile, surfQfile, surfQswfile,
70 & uVelInitFile, vVelInitFile, pSurfInitFile,
71 & dQdTfile
72 CHARACTER*(5) checkPtSuff(maxNoChkptLev)
73 CHARACTER*(MAX_LEN_FNAM) bathyFile
74 CHARACTER*(MAX_LEN_FNAM) hydrogThetaFile
75 CHARACTER*(MAX_LEN_FNAM) hydrogSaltFile
76 CHARACTER*(MAX_LEN_FNAM) zonalWindFile
77 CHARACTER*(MAX_LEN_FNAM) meridWindFile
78 CHARACTER*(MAX_LEN_FNAM) thetaClimFile
79 CHARACTER*(MAX_LEN_FNAM) saltClimFile
80 CHARACTER*(MAX_LEN_FNAM) surfQfile
81 CHARACTER*(MAX_LEN_FNAM) surfQswfile
82 CHARACTER*(MAX_LEN_FNAM) EmPmRfile
83 CHARACTER*(MAX_LEN_FNAM) buoyancyRelation
84 CHARACTER*(MAX_LEN_FNAM) uVelInitFile
85 CHARACTER*(MAX_LEN_FNAM) vVelInitFile
86 CHARACTER*(MAX_LEN_FNAM) pSurfInitFile
87 CHARACTER*(MAX_LEN_FNAM) dQdTfile
88
89 C-- COMMON /PARM_I/ Integer valued parameters used by the model.
90 C cg2dMaxIters - Maximum number of iterations in the
91 C two-dimensional con. grad solver.
92 C cg2dChkResFreq - Frequency with which to check residual
93 C in con. grad solver.
94 C cg3dMaxIters - Maximum number of iterations in the
95 C three-dimensional con. grad solver.
96 C cg3dChkResFreq - Frequency with which to check residual
97 C in con. grad solver.
98 C nIter0 - Start time-step number of for this run
99 C nTimeSteps - Number of timesteps to execute
100 C numStepsPerPickup - For offline setup. Frequency of pickup
101 C of flow fields.
102 C writeStatePrec - Precision used for writing model state.
103 C writeBinaryPrec - Precision used for writing binary files
104 C readBinaryPrec - Precision used for reading binary files
105 C nCheckLev - Holds current checkpoint level
106
107 COMMON /PARM_I/
108 & cg2dMaxIters,
109 & cg2dChkResFreq,
110 & cg3dMaxIters,
111 & cg3dChkResFreq,
112 & nIter0, nTimeSteps, nEndIter,
113 & numStepsPerPickup,
114 & writeStatePrec, nCheckLev,
115 & writeBinaryPrec, readBinaryPrec,
116 & zonal_filt_sinpow, zonal_filt_cospow
117 INTEGER cg2dMaxIters
118 INTEGER cg2dChkResFreq
119 INTEGER cg3dMaxIters
120 INTEGER cg3dChkResFreq
121 INTEGER nIter0
122 INTEGER nTimeSteps
123 INTEGER nEndIter
124 INTEGER numStepsPerPickup
125 INTEGER writeStatePrec
126 INTEGER writeBinaryPrec
127 INTEGER readBinaryPrec
128 INTEGER nCheckLev
129 INTEGER zonal_filt_sinpow
130 INTEGER zonal_filt_cospow
131
132 C-- COMMON /PARM_L/ Logical valued parameters used by the model.
133 C usingCartesianGrid - If TRUE grid generation will be in a cartesian
134 C coordinate frame.
135 C usingSphericalPolarGrid - If TRUE grid generation will be in a
136 C spherical polar frame.
137 C no_slip_sides - Impose "no-slip" at lateral boundaries.
138 C no_slip_bottom- Impose "no-slip" at bottom boundary.
139 C staggerTimeStep - enable a Stagger time stepping T,S Rho then U,V
140 C momViscosity - Flag which turns momentum friction terms on and off.
141 C momAdvection - Flag which turns advection of momentum on and off.
142 C momForcing - Flag which turns external forcing of momentum on
143 C and off.
144 C momPressureForcing - Flag which turns pressure term in momentum equation
145 C on and off.
146 C metricTerms - Flag which turns metric terms on or off.
147 C usingSphericalPolarMTerms - If TRUE use spherical polar metric terms.
148 C useCoriolis - Flag which turns the coriolis terms on and off.
149 C tempDiffusion - Flag which turns diffusion of temperature on
150 C and off.
151 C tempAdvection - Flag which turns advection of temperature on
152 C and off.
153 C tempForcing - Flag which turns external forcing of temperature on
154 C and off.
155 C saltDiffusion - Flag which turns diffusion of salinit on
156 C and off.
157 C saltAdvection - Flag which turns advection of salinit on
158 C and off.
159 C saltForcing - Flag which turns external forcing of salinit on
160 C and off.
161 C rigidLid - Set to true to use rigid lid
162 C implicitFreeSurface - Set to true to use implcit free surface
163 C exactConserv - Set to true to conserve exactly the total Volume
164 C nonlinFreeSurf - Set to true to use non-linear free surface
165 C momStepping - Turns momentum equation time-stepping off
166 C tempStepping - Turns temperature equation time-stepping off
167 C saltStepping - Turns salinity equation time-stepping off
168 C useConstantF - Coriolis parameter set to f0
169 C useBetaPlaneF - Coriolis parameter set to f0 + beta.y
170 C useSphereF - Coriolis parameter set to 2.omega.sin(phi)
171 C implicitDiffusion - Turns implicit vertical diffusion on
172 C implicitViscosity - Turns implicit vertical viscosity on
173 C doThetaClimRelax - Set true if relaxation to temperature
174 C climatology is required.
175 C doSaltClimRelax - Set true if relaxation to salinity
176 C climatology is required.
177 C periodicExternalForcing - Set true if forcing is time-dependant
178 C usingPCoords - Set to indicate that we are working in pressure
179 C coords.
180 C usingZCoords - Set to indicate that we are working in height
181 C coords.
182 C nonHydrostatic - Using non-hydrostatic terms
183 C globalFiles - Selects between "global" and "tiled" files
184 C allowFreezing - Allows water to freeze and form ice
185 C groundAtK1 - put the surface(k=1) at the Lower Boundary (=ground)
186 COMMON /PARM_L/ usingCartesianGrid, usingSphericalPolarGrid,
187 & usingCurvilinearGrid,
188 & no_slip_sides,no_slip_bottom,
189 & staggerTimeStep,
190 & momViscosity, momAdvection, momForcing, useCoriolis,
191 & momPressureForcing,tempDiffusion, tempAdvection, tempForcing,
192 & saltDiffusion, saltAdvection, saltForcing,
193 & rigidLid, implicitFreeSurface, exactConserv, nonlinFreeSurf,
194 & momStepping, tempStepping, saltStepping,
195 & metricTerms, usingSphericalPolarMTerms,
196 & useConstantF, useBetaPlaneF, useSphereF,
197 & implicitDiffusion, implicitViscosity,
198 & doThetaClimRelax, doSaltClimRelax,
199 & periodicExternalForcing, usingPCoords, usingZCoords,
200 & nonHydrostatic, globalFiles,
201 & allowFreezing, groundAtK1,
202 & usePickupBeforeC35, debugMode
203 LOGICAL usingCartesianGrid
204 LOGICAL usingSphericalPolarGrid
205 LOGICAL usingCurvilinearGrid
206 LOGICAL usingSphericalPolarMTerms
207 LOGICAL no_slip_sides
208 LOGICAL no_slip_bottom
209 LOGICAL staggerTimeStep
210 LOGICAL momViscosity
211 LOGICAL momAdvection
212 LOGICAL momForcing
213 LOGICAL momPressureForcing
214 LOGICAL useCoriolis
215 LOGICAL tempDiffusion
216 LOGICAL tempAdvection
217 LOGICAL tempForcing
218 LOGICAL saltDiffusion
219 LOGICAL saltAdvection
220 LOGICAL saltForcing
221 LOGICAL rigidLid
222 LOGICAL implicitFreeSurface
223 LOGICAL exactConserv
224 LOGICAL nonlinFreeSurf
225 LOGICAL momStepping
226 LOGICAL tempStepping
227 LOGICAL saltStepping
228 LOGICAL metricTerms
229 LOGICAL useConstantF
230 LOGICAL useBetaPlaneF
231 LOGICAL useSphereF
232 LOGICAL implicitDiffusion
233 LOGICAL implicitViscosity
234 LOGICAL doThetaClimRelax
235 LOGICAL doSaltClimRelax
236 LOGICAL periodicExternalForcing
237 LOGICAL usingPCoords
238 LOGICAL usingZCoords
239 LOGICAL nonHydrostatic
240 LOGICAL globalFiles
241 LOGICAL allowFreezing
242 LOGICAL groundAtK1
243 LOGICAL usePickupBeforeC35
244 LOGICAL debugMode
245
246 C-- COMMON /PARM_R/ "Real" valued parameters used by the model.
247 C gg2dTargetResidual
248 C - Target residual for cg2d solver; no unit (RHS normalisation)
249 C cg2dTargetResWunit
250 C - Target residual for cg2d solver; W unit (No RHS normalisation)
251 C cg3dTargetResidual
252 C - Target residual for cg3d solver.
253 C cg2dpcOffDFac - Averaging weight for preconditioner off-diagonal.
254 C Note. 20th May 1998
255 C I made a weird discovery! In the model paper we argue
256 C for the form of the preconditioner used here ( see
257 C A Finite-volume, Incompressible Navier-Stokes Model
258 C ...., Marshall et. al ). The algebra gives a simple
259 C 0.5 factor for the averaging of ac and aCw to get a
260 C symmettric pre-conditioner. By using a factor of 0.51
261 C i.e. scaling the off-diagonal terms in the
262 C preconditioner down slightly I managed to get the
263 C number of iterations for convergence in a test case to
264 C drop form 192 -> 134! Need to investigate this further!
265 C For now I have introduced a parameter cg2dpcOffDFac which
266 C defaults to 0.51 but can be set at runtime.
267 C delP - Vertical grid spacing ( Pa ).
268 C delZ - Vertical grid spacing ( m ).
269 C delR - Vertical grid spacing ( units of r ).
270 C delX - Separation between cell faces (m) or (deg), depending
271 C delY on input flags.
272 C gravity - Accel. due to gravity ( m/s^2 )
273 C recip_gravity and its inverse
274 C gBaro - Accel. due to gravity used in barotropic equation ( m/s^2 )
275 C ronil - Reference density
276 C rhoConst - Vertically constant reference density
277 C startTime - Start time for model ( s )
278 C phiMin - Latitude of southern most cell face.
279 C thetaMin - Longitude of western most cell face (this
280 C is an "inert" parameter but it is included
281 C to make geographical references simple.)
282 C rSphere - Radius of sphere for a spherical polar grid ( m ).
283 C recip_RSphere - Reciprocal radius of sphere ( m ).
284 C f0 - Reference coriolis parameter ( 1/s )
285 C ( Southern edge f for beta plane )
286 C beta - df/dy ( s^-1.m^-1 )
287 C omega - Angular velocity ( rad/s )
288 C viscAh - Eddy viscosity coeff. for mixing of
289 C momentum laterally ( m^2/s )
290 C viscAz - Eddy viscosity coeff. for mixing of
291 C momentum vertically ( m^2/s )
292 C viscAp - Eddy viscosity coeff. for mixing of
293 C momentum vertically ( Pa^2/s )
294 C viscAr - Eddy viscosity coeff. for mixing of
295 C momentum vertically ( units of r^2/s )
296 C viscA4 - Biharmonic viscosity coeff. for mixing of
297 C momentum laterally ( m^4/s )
298 C diffKhT - Laplacian diffusion coeff. for mixing of
299 C heat laterally ( m^2/s )
300 C diffKzT - Laplacian diffusion coeff. for mixing of
301 C heat vertically ( m^2/s )
302 C diffKpT - Laplacian diffusion coeff. for mixing of
303 C heat vertically ( Pa^2/s )
304 C diffKrT - Laplacian diffusion coeff. for mixing of
305 C heat vertically ( units of r^2/s )
306 C diffK4T - Biharmonic diffusion coeff. for mixing of
307 C heat laterally ( m^4/s )
308 C diffKhS - Laplacian diffusion coeff. for mixing of
309 C salt laterally ( m^2/s )
310 C diffKzS - Laplacian diffusion coeff. for mixing of
311 C salt vertically ( m^2/s )
312 C diffKpS - Laplacian diffusion coeff. for mixing of
313 C salt vertically ( Pa^2/s )
314 C diffKrS - Laplacian diffusion coeff. for mixing of
315 C salt vertically ( units of r^2/s )
316 C diffK4S - Biharmonic diffusion coeff. for mixing of
317 C salt laterally ( m^4/s )
318 C deltaT - Default timestep ( s )
319 C deltaTClock - Timestep used as model "clock". This determines the
320 C IO frequencies and is used in tagging output. It can
321 C be totally different to the dynamical time. Typically
322 C it will be the deep-water timestep for accelerated runs.
323 C Frequency of checkpointing and dumping of the model state
324 C are referenced to this clock. ( s )
325 C deltaTMom - Timestep for momemtum equations ( s )
326 C deltaTtracer - Timestep for tracer equations ( s )
327 C freesurfFac - Parameter to turn implicit free surface term on or off
328 C freesurfac = 1. uses implicit free surface
329 C freesurfac = 0. uses rigid lid
330 C implicSurfPress - parameter of the Crank-Nickelson time stepping :
331 C Implicit part of Surface Pressure Gradient ( 0-1 )
332 C implicDiv2Dflow - parameter of the Crank-Nickelson time stepping :
333 C Implicit part of barotropic flow Divergence ( 0-1 )
334 C hFacMin - Minimum fraction size of a cell (affects hFacC etc...)
335 C hFacMinDz - Minimum dimesional size of a cell (affects hFacC etc..., m)
336 C hFacMinDp - Minimum dimesional size of a cell (affects hFacC etc..., Pa)
337 C hFacMinDr - Minimum dimesional size of a cell (affects hFacC etc..., units of r)
338 C hFacInf - Threshold (inf and sup) for fraction size of surface cell
339 C hFacSup that control vanishing and creating levels
340 C tauCD - CD scheme coupling timescale ( 1/s )
341 C rCD - CD scheme normalised coupling parameter ( 0-1 )
342 C startTime - Starting time for this integration ( s ).
343 C endTime - Ending time for this integration ( s ).
344 C chkPtFreq - Frequency of rolling check pointing ( s ).
345 C pChkPtFreq - Frequency of permanent check pointing ( s ).
346 C dumpFreq - Frequency with which model state is written to
347 C post-processing files ( s ).
348 C afFacMom - Advection of momentum term tracer parameter
349 C vfFacMom - Momentum viscosity tracer parameter
350 C pfFacMom - Momentum pressure forcing tracer parameter
351 C cfFacMom - Coriolis term tracer parameter
352 C foFacMom - Momentum forcing tracer parameter
353 C mtFacMom - Metric terms tracer parameter
354 C cosPower - Power of cosine of latitude to multiply viscosity
355 C cAdjFreq - Frequency of convective adjustment
356 C
357 C taveFreq - Frequency with which time-averaged model state is written to
358 C post-processing files ( s ).
359 C tauThetaClimRelax - Relaxation to climatology time scale ( s ).
360 C lambdaThetaClimRelax - Inverse time scale for relaxation ( 1/s ).
361 C tauSaltClimRelax - Relaxation to climatology time scale ( s ).
362 C lambdaSaltClimRelax - Inverse time scale for relaxation ( 1/s ).
363 C externForcingPeriod - Is the period of which forcing varies (eg. 1 month)
364 C externForcingCycle - Is the repeat time of the forcing (eg. 1 year)
365 C (note: externForcingCycle must be an integer
366 C number times externForcingPeriod)
367 C horiVertRatio - Ratio on units in vertical to units in horizontal.
368 C recip_horiVertRatio ( 1 if horiz in m and vertical in m ).
369 C ( g*rho if horiz in m and vertical in Pa ).
370 C latFFTFiltLo - Low latitude for FFT filtering of latitude
371 C circles ( see filter*.F )
372 C Ro_SeaLevel - standard position of Sea-Level in "R" coordinate, used as
373 C starting value (k=1) for vertical coordinate (rf(1)=Ro_SeaLevel)
374 C bottomDragLinear - Drag coefficient built in to core dynamics
375 C " Quadratic ( linear: 1/s, quadratic: 1/m )
376 COMMON /PARM_R/ cg2dTargetResidual, cg2dTargetResWunit,
377 & cg2dpcOffDFac, cg3dTargetResidual,
378 & delP, delZ, delR, delX, delY,
379 & deltaT,deltaTmom, deltaTtracer, deltaTClock,abeps, startTime,
380 & phiMin, thetaMin, rSphere, recip_RSphere, f0, beta,
381 & fCori, fCoriG,
382 & viscAh, viscAz, viscA4, viscAr,
383 & diffKhT, diffKzT, diffK4T, diffKrT,
384 & diffKhS, diffKzS, diffK4S, diffKrS,
385 & delT, tauCD, rCD, freeSurfFac, implicSurfPress, implicDiv2Dflow,
386 & hFacMin, hFacMinDz, hFacInf, hFacSup,
387 & gravity, recip_Gravity, gBaro, rhonil, recip_rhonil,
388 & recip_rhoConst, rhoConst, tRef, sRef,
389 & endTime, chkPtFreq, pchkPtFreq, dumpFreq, taveFreq, monitorFreq,
390 & afFacMom, vfFacMom, pfFacMom, cfFacMom, foFacMom, mtFacMom,
391 & cosPower,
392 & cAdjFreq, omega, tauThetaClimRelax, lambdaThetaClimRelax,
393 & tauSaltClimRelax, lambdaSaltClimRelax,
394 & externForcingCycle, externForcingPeriod,
395 & viscAp, diffKpT, diffKpS, hFacMinDr, hFacMinDp,
396 & theta_S, specVol_S, horiVertRatio, recip_horiVertRatio,
397 & latFFTFiltLo, ivdc_kappa, Ro_SeaLevel, zonal_filt_lat,
398 & bottomDragLinear,bottomDragQuadratic
399
400 _RL cg2dTargetResidual
401 _RL cg2dTargetResWunit
402 _RL cg3dTargetResidual
403 _RL cg2dpcOffDFac
404 _RL delZ(Nr)
405 _RL delP(Nr)
406 _RL delR(Nr)
407 _RL delX(Nx)
408 _RL delY(Ny)
409 _RL deltaT
410 _RL deltaTClock
411 _RL deltaTmom
412 _RL deltaTtracer
413 _RL abeps
414 _RL phiMin
415 _RL thetaMin
416 _RL rSphere
417 _RL recip_RSphere
418 _RL f0
419 _RL freeSurfFac
420 _RL implicSurfPress
421 _RL implicDiv2Dflow
422 _RL hFacMin
423 _RL hFacMinDz
424 _RL hFacMinDp
425 _RL hFacMinDr
426 _RL hFacInf
427 _RL hFacSup
428 _RL beta
429 _RL viscAh
430 _RL viscAz
431 _RL viscAp
432 _RL viscAr
433 _RL viscA4
434 _RL diffKhT
435 _RL diffKrT
436 _RL diffKzT
437 _RL diffKpT
438 _RL diffK4T
439 _RL diffKhS
440 _RL diffKrS
441 _RL diffKzS
442 _RL diffKpS
443 _RL diffK4S
444 _RL delt
445 _RL tauCD
446 _RL rCD
447 _RL gravity
448 _RL recip_gravity
449 _RL gBaro
450 _RL rhonil
451 _RL recip_rhonil
452 _RL rhoConst
453 _RL recip_rhoConst
454 _RL specVol_S(Nr)
455 _RL tRef(Nr)
456 _RL theta_S(Nr)
457 _RL sRef(Nr)
458 _RS fCori(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
459 _RS fCoriG(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
460 _RL startTime
461 _RL endTime
462 _RL chkPtFreq
463 _RL pChkPtFreq
464 _RL dumpFreq
465 _RL taveFreq
466 _RL monitorFreq
467 _RL afFacMom
468 _RL vfFacMom
469 _RL pfFacMom
470 _RL cfFacMom
471 _RL foFacMom
472 _RL mTFacMom
473 _RL cosPower
474 _RL cAdjFreq
475 _RL omega
476 _RL tauThetaClimRelax
477 _RL lambdaThetaClimRelax
478 _RL tauSaltClimRelax
479 _RL lambdaSaltClimRelax
480 _RL externForcingCycle
481 _RL externForcingPeriod
482 _RL horiVertRatio
483 _RL recip_horiVertRatio
484 _RL latFFTFiltLo
485 _RL ivdc_kappa
486 _RL Ro_SeaLevel
487 _RL zonal_filt_lat
488 _RL bottomDragLinear
489 _RL bottomDragQuadratic
490
491 COMMON /PARM_A/ HeatCapacity_Cp,recip_Cp,
492 & Lamba_theta
493 _RL HeatCapacity_Cp
494 _RL Lamba_theta
495 _RL recip_Cp
496
497 C Equation of State (polynomial coeffients)
498 COMMON /PARM_EOS_NL/ eosC,eosSig0,eosRefT,eosRefS
499 _RL eosC(9,Nr+1),eosSig0(Nr+1),eosRefT(Nr+1),eosRefS(Nr+1)
500 C Linear equation of state
501 C tAlpha - Linear EOS thermal expansion coefficient ( 1/degree ).
502 C sBeta - Linear EOS haline contraction coefficient.
503 COMMON /PARM_EOS_LIN/ tAlpha,sBeta,eosType
504 _RL tAlpha
505 _RL sBeta
506 character*(6) eosType
507
508 C Logical flags for selecting packages
509 LOGICAL useKPP
510 LOGICAL useGMRedi
511 LOGICAL useOBCS
512 LOGICAL useAIM
513 LOGICAL useECCO
514 LOGICAL useSHAP_FILT
515 COMMON /PARM_PACKAGES/
516 & useKPP, useGMRedi, useOBCS, useAIM, useECCO, useSHAP_FILT
517

  ViewVC Help
Powered by ViewVC 1.1.22