/[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.65.4.1.2.1 - (show annotations) (download)
Wed Feb 12 06:45:54 2003 UTC (21 years, 2 months ago) by dimitri
Branch: c24_e25_ice
CVS Tags: ecco_ice2
Changes since 1.65.4.1: +4 -3 lines
File MIME type: text/plain
02/12/03
  Added pkg/seaice and modified pkg/exf as per release1_p12_pre

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

  ViewVC Help
Powered by ViewVC 1.1.22