/[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.2.8.2.1 - (show annotations) (download)
Thu Aug 25 16:38:10 2005 UTC (18 years, 8 months ago) by dimitri
Branch: release1_50yr
Changes since 1.65.2.8: +13 -3 lines
File MIME type: text/plain
updating model/inc/PARAMS.h model/src/ini_parms.F and
model/src/set_defaults.F for release1_50yr

1 C $Header: /u/gcmpack/MITgcm/model/inc/PARAMS.h,v 1.65.2.8 2003/05/02 14:25:29 dimitri Exp $
2 C $Name: release1_50yr $
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 & debugLevel
131 INTEGER cg2dMaxIters
132 INTEGER cg2dChkResFreq
133 INTEGER cg3dMaxIters
134 INTEGER cg3dChkResFreq
135 INTEGER nIter0
136 INTEGER nTimeSteps
137 INTEGER nEndIter
138 INTEGER numStepsPerPickup
139 INTEGER writeStatePrec
140 INTEGER writeBinaryPrec
141 INTEGER readBinaryPrec
142 INTEGER nCheckLev
143 INTEGER nonlinFreeSurf
144 INTEGER tempAdvScheme
145 INTEGER saltAdvScheme
146 INTEGER tracerAdvScheme
147 INTEGER debugLevel
148
149 C
150 INTEGER debLevZero
151 PARAMETER(debLevZero=0)
152 INTEGER debLevA
153 PARAMETER(debLevA=1)
154 INTEGER debLevB
155 PARAMETER(debLevB=2)
156
157 C-- COMMON /PARM_L/ Logical valued parameters used by the model.
158 C usingCartesianGrid :: If TRUE grid generation will be in a cartesian
159 C coordinate frame.
160 C usingSphericalPolarGrid :: If TRUE grid generation will be in a
161 C spherical polar frame.
162 C no_slip_sides :: Impose "no-slip" at lateral boundaries.
163 C no_slip_bottom :: Impose "no-slip" at bottom boundary.
164 C staggerTimeStep :: enable a Stagger time stepping T,S Rho then U,V
165 C momViscosity :: Flag which turns momentum friction terms on and off.
166 C momAdvection :: Flag which turns advection of momentum on and off.
167 C momForcing :: Flag which turns external forcing of momentum on
168 C and off.
169 C momPressureForcing :: Flag which turns pressure term in momentum equation
170 C on and off.
171 C metricTerms :: Flag which turns metric terms on or off.
172 C usingSphericalPolarMTerms :: If TRUE use spherical polar metric terms.
173 C useCoriolis :: Flag which turns the coriolis terms on and off.
174 C tempDiffusion :: Flag which turns diffusion of temperature on
175 C and off.
176 C tempAdvection :: Flag which turns advection of temperature on
177 C and off.
178 C tempForcing :: Flag which turns external forcing of temperature on
179 C and off.
180 C saltDiffusion :: Flag which turns diffusion of salinity on
181 C and off.
182 C saltAdvection :: Flag which turns advection of salinity on
183 C and off.
184 C saltForcing :: Flag which turns external forcing of salinity on
185 C and off.
186 C useRealFreshWaterFlux :: if true (=Natural BCS), treats P+R-E flux
187 C as a real Fresh Water (=> changes the seal level)
188 C if false, converts P+R-E to virtual salt flux
189 C rigidLid :: Set to true to use rigid lid
190 C implicitFreeSurface :: Set to true to use implcit free surface
191 C exactConserv :: Set to true to conserve exactly the total Volume
192 C uniformLin_PhiSurf :: Set to true to use a uniform Bo_surf in the
193 C linear relation Phi_surf = Bo_surf*eta
194 C momStepping :: Turns momentum equation time-stepping off
195 C tempStepping :: Turns temperature equation time-stepping off
196 C saltStepping :: Turns salinity equation time-stepping off
197 C tr1Stepping :: Turns passive tracer 1 time-stepping on/off
198 C useConstantF :: Coriolis parameter set to f0
199 C useBetaPlaneF :: Coriolis parameter set to f0 + beta.y
200 C useSphereF :: Coriolis parameter set to 2.omega.sin(phi)
201 C implicitDiffusion :: Turns implicit vertical diffusion on
202 C implicitViscosity :: Turns implicit vertical viscosity on
203 C doThetaClimRelax :: Set true if relaxation to temperature
204 C climatology is required.
205 C doSaltClimRelax :: Set true if relaxation to salinity
206 C climatology is required.
207 C periodicExternalForcing :: Set true if forcing is time-dependant
208 C usingPCoords :: Set to indicate that we are working in pressure
209 C coords.
210 C usingZCoords :: Set to indicate that we are working in height
211 C coords.
212 C nonHydrostatic :: Using non-hydrostatic terms
213 C globalFiles :: Selects between "global" and "tiled" files
214 C useSingleCpuIO :: On SGI platforms, option globalFiles is either
215 C slow (f77) or does not work (f90). When
216 C useSingleCpuIO is set, mdsio_writefield.F
217 C outputs from master mpi process only.
218 C allowFreezing :: Allows water to freeze and form ice
219 C groundAtK1 :: put the surface(k=1) at the Lower Boundary (=ground)
220 C useJamartWetPoints :: Use wet-point method for Coriolis (Jamart and Ozer, 1986)
221 COMMON /PARM_L/ usingCartesianGrid, usingSphericalPolarGrid,
222 & usingCurvilinearGrid,
223 & no_slip_sides,no_slip_bottom,
224 & staggerTimeStep,
225 & momViscosity, momAdvection, momForcing, useCoriolis,
226 & momPressureForcing, vectorInvariantMomentum,
227 & tempDiffusion, tempAdvection, tempForcing,
228 & saltDiffusion, saltAdvection, saltForcing,
229 & useRealFreshWaterFlux,
230 & rigidLid, implicitFreeSurface, exactConserv, uniformLin_PhiSurf,
231 & momStepping, tempStepping, saltStepping, tr1Stepping,
232 & metricTerms, usingSphericalPolarMTerms,
233 & useConstantF, useBetaPlaneF, useSphereF,
234 & implicitDiffusion, implicitViscosity,
235 & doThetaClimRelax, doSaltClimRelax, doTr1ClimRelax,
236 & periodicExternalForcing, usingPCoords, usingZCoords,
237 & nonHydrostatic, globalFiles, useSingleCpuIO,
238 & allowFreezing, groundAtK1,
239 & usePickupBeforeC35, debugMode,
240 & readPickupWithTracer, writePickupWithTracer,
241 & multiDimAdvection, useEnergyConservingCoriolis,
242 & useJamartWetPoints
243 LOGICAL usingCartesianGrid
244 LOGICAL usingSphericalPolarGrid
245 LOGICAL usingCurvilinearGrid
246 LOGICAL usingSphericalPolarMTerms
247 LOGICAL no_slip_sides
248 LOGICAL no_slip_bottom
249 LOGICAL staggerTimeStep
250 LOGICAL momViscosity
251 LOGICAL momAdvection
252 LOGICAL momForcing
253 LOGICAL momPressureForcing
254 LOGICAL useCoriolis
255 LOGICAL vectorInvariantMomentum
256 LOGICAL tempDiffusion
257 LOGICAL tempAdvection
258 LOGICAL tempForcing
259 LOGICAL saltDiffusion
260 LOGICAL saltAdvection
261 LOGICAL saltForcing
262 LOGICAL useRealFreshWaterFlux
263 LOGICAL rigidLid
264 LOGICAL implicitFreeSurface
265 LOGICAL exactConserv
266 LOGICAL uniformLin_PhiSurf
267 LOGICAL momStepping
268 LOGICAL tempStepping
269 LOGICAL saltStepping
270 LOGICAL tr1Stepping
271 LOGICAL metricTerms
272 LOGICAL useConstantF
273 LOGICAL useBetaPlaneF
274 LOGICAL useSphereF
275 LOGICAL implicitDiffusion
276 LOGICAL implicitViscosity
277 LOGICAL doThetaClimRelax
278 LOGICAL doSaltClimRelax
279 LOGICAL doTr1ClimRelax
280 LOGICAL periodicExternalForcing
281 LOGICAL usingPCoords
282 LOGICAL usingZCoords
283 LOGICAL nonHydrostatic
284 LOGICAL globalFiles
285 LOGICAL useSingleCpuIO
286 LOGICAL allowFreezing
287 LOGICAL groundAtK1
288 LOGICAL usePickupBeforeC35
289 LOGICAL debugMode
290 LOGICAL readPickupWithTracer
291 LOGICAL writePickupWithTracer
292 LOGICAL multiDimAdvection
293 LOGICAL useEnergyConservingCoriolis
294 LOGICAL useJamartWetPoints
295
296 C-- COMMON /PARM_R/ "Real" valued parameters used by the model.
297 C gg2dTargetResidual
298 C :: Target residual for cg2d solver; no unit (RHS normalisation)
299 C cg2dTargetResWunit
300 C :: Target residual for cg2d solver; W unit (No RHS normalisation)
301 C cg3dTargetResidual
302 C :: Target residual for cg3d solver.
303 C cg2dpcOffDFac :: Averaging weight for preconditioner off-diagonal.
304 C Note. 20th May 1998
305 C I made a weird discovery! In the model paper we argue
306 C for the form of the preconditioner used here ( see
307 C A Finite-volume, Incompressible Navier-Stokes Model
308 C ...., Marshall et. al ). The algebra gives a simple
309 C 0.5 factor for the averaging of ac and aCw to get a
310 C symmettric pre-conditioner. By using a factor of 0.51
311 C i.e. scaling the off-diagonal terms in the
312 C preconditioner down slightly I managed to get the
313 C number of iterations for convergence in a test case to
314 C drop form 192 -> 134! Need to investigate this further!
315 C For now I have introduced a parameter cg2dpcOffDFac which
316 C defaults to 0.51 but can be set at runtime.
317 C delP :: Vertical grid spacing ( Pa ).
318 C delZ :: Vertical grid spacing ( m ).
319 C delR :: Vertical grid spacing ( units of r ).
320 C delX :: Separation between cell faces (m) or (deg), depending
321 C delY on input flags.
322 C gravity :: Accel. due to gravity ( m/s^2 )
323 C recip_gravity and its inverse
324 C gBaro :: Accel. due to gravity used in barotropic equation ( m/s^2 )
325 C ronil :: Reference density
326 C rhoConst :: Vertically constant reference density
327 C startTime :: Start time for model ( s )
328 C phiMin :: Latitude of southern most cell face.
329 C thetaMin :: Longitude of western most cell face (this
330 C is an "inert" parameter but it is included
331 C to make geographical references simple.)
332 C rSphere :: Radius of sphere for a spherical polar grid ( m ).
333 C recip_RSphere :: Reciprocal radius of sphere ( m ).
334 C f0 :: Reference coriolis parameter ( 1/s )
335 C ( Southern edge f for beta plane )
336 C beta :: df/dy ( s^-1.m^-1 )
337 C omega :: Angular velocity ( rad/s )
338 C viscAh :: Eddy viscosity coeff. for mixing of
339 C momentum laterally ( m^2/s )
340 C viscAz :: Eddy viscosity coeff. for mixing of
341 C momentum vertically ( m^2/s )
342 C viscAp :: Eddy viscosity coeff. for mixing of
343 C momentum vertically ( Pa^2/s )
344 C viscAr :: Eddy viscosity coeff. for mixing of
345 C momentum vertically ( units of r^2/s )
346 C viscA4 :: Biharmonic viscosity coeff. for mixing of
347 C momentum laterally ( m^4/s )
348 C diffKhT :: Laplacian diffusion coeff. for mixing of
349 C heat laterally ( m^2/s )
350 C diffKzT :: Laplacian diffusion coeff. for mixing of
351 C heat vertically ( m^2/s )
352 C diffKpT :: Laplacian diffusion coeff. for mixing of
353 C heat vertically ( Pa^2/s )
354 C diffKrT :: Laplacian diffusion coeff. for mixing of
355 C heat vertically ( units of r^2/s )
356 C diffK4T :: Biharmonic diffusion coeff. for mixing of
357 C heat laterally ( m^4/s )
358 C diffKhS :: Laplacian diffusion coeff. for mixing of
359 C salt laterally ( m^2/s )
360 C diffKzS :: Laplacian diffusion coeff. for mixing of
361 C salt vertically ( m^2/s )
362 C diffKpS :: Laplacian diffusion coeff. for mixing of
363 C salt vertically ( Pa^2/s )
364 C diffKrS :: Laplacian diffusion coeff. for mixing of
365 C salt vertically ( units of r^2/s )
366 C diffK4S :: Biharmonic diffusion coeff. for mixing of
367 C salt laterally ( m^4/s )
368 C deltaT :: Default timestep ( s )
369 C deltaTClock :: Timestep used as model "clock". This determines the
370 C IO frequencies and is used in tagging output. It can
371 C be totally different to the dynamical time. Typically
372 C it will be the deep-water timestep for accelerated runs.
373 C Frequency of checkpointing and dumping of the model state
374 C are referenced to this clock. ( s )
375 C deltaTMom :: Timestep for momemtum equations ( s )
376 C deltaTtracer :: Timestep for tracer equations ( s )
377 C freesurfFac :: Parameter to turn implicit free surface term on or off
378 C freesurfac = 1. uses implicit free surface
379 C freesurfac = 0. uses rigid lid
380 C implicSurfPress :: parameter of the Crank-Nickelson time stepping :
381 C Implicit part of Surface Pressure Gradient ( 0-1 )
382 C implicDiv2Dflow :: parameter of the Crank-Nickelson time stepping :
383 C Implicit part of barotropic flow Divergence ( 0-1 )
384 C hFacMin :: Minimum fraction size of a cell (affects hFacC etc...)
385 C hFacMinDz :: Minimum dimesional size of a cell (affects hFacC etc..., m)
386 C hFacMinDp :: Minimum dimesional size of a cell (affects hFacC etc..., Pa)
387 C hFacMinDr :: Minimum dimesional size of a cell (affects hFacC etc..., units of r)
388 C hFacInf :: Threshold (inf and sup) for fraction size of surface cell
389 C hFacSup that control vanishing and creating levels
390 C tauCD :: CD scheme coupling timescale ( 1/s )
391 C rCD :: CD scheme normalised coupling parameter ( 0-1 )
392 C startTime :: Starting time for this integration ( s ).
393 C endTime :: Ending time for this integration ( s ).
394 C chkPtFreq :: Frequency of rolling check pointing ( s ).
395 C pChkPtFreq :: Frequency of permanent check pointing ( s ).
396 C dumpFreq :: Frequency with which model state is written to
397 C post-processing files ( s ).
398 C tr1dumpFreq :: dumpFreq for passive tracer field (s)
399 C diagFreq :: Frequency with which model writes diagnostic output
400 C of intermediate quantities.
401 C afFacMom :: Advection of momentum term tracer parameter
402 C vfFacMom :: Momentum viscosity tracer parameter
403 C pfFacMom :: Momentum pressure forcing tracer parameter
404 C cfFacMom :: Coriolis term tracer parameter
405 C foFacMom :: Momentum forcing tracer parameter
406 C mtFacMom :: Metric terms tracer parameter
407 C cosPower :: Power of cosine of latitude to multiply viscosity
408 C cAdjFreq :: Frequency of convective adjustment
409 C
410 C taveFreq :: Frequency with which time-averaged model state is written to
411 C post-processing files ( s ).
412 C tauThetaClimRelax :: Relaxation to climatology time scale ( s ).
413 C lambdaThetaClimRelax :: Inverse time scale for relaxation ( 1/s ).
414 C tauSaltClimRelax :: Relaxation to climatology time scale ( s ).
415 C lambdaSaltClimRelax :: Inverse time scale for relaxation ( 1/s ).
416 C externForcingPeriod :: Is the period of which forcing varies (eg. 1 month)
417 C externForcingCycle :: Is the repeat time of the forcing (eg. 1 year)
418 C (note: externForcingCycle must be an integer
419 C number times externForcingPeriod)
420 C convertFW2Salt :: salinity, used to convert Fresh-Water Flux to Salt Flux
421 C (use model surface (local) value if set to -1)
422 C temp_EvPrRn :: temperature of Rain & Evap.
423 C salt_EvPrRn :: salinity of Rain & Evap.
424 C trac_EvPrRn :: tracer concentration in Rain & Evap.
425 C (notes: a) tracer content of Rain/Evap only used if both
426 C NonLin_FrSurf & useRealFreshWater are set.
427 C b) use model surface (local) value if set to UNSET_RL)
428 C horiVertRatio :: Ratio on units in vertical to units in horizontal.
429 C recip_horiVertRatio ( 1 if horiz in m and vertical in m ).
430 C ( g*rho if horiz in m and vertical in Pa ).
431 C Ro_SeaLevel :: standard position of Sea-Level in "R" coordinate, used as
432 C starting value (k=1) for vertical coordinate (rf(1)=Ro_SeaLevel)
433 C bottomDragLinear :: Drag coefficient built in to core dynamics
434 C --"-"-- Quadratic ( linear: 1/s, quadratic: 1/m )
435 COMMON /PARM_R/ cg2dTargetResidual, cg2dTargetResWunit,
436 & cg2dpcOffDFac, cg3dTargetResidual,
437 & delP, delZ, delR, delX, delY,
438 & deltaT,deltaTmom, deltaTtracer, deltaTClock,abeps, startTime,
439 & phiMin, thetaMin, rSphere, recip_RSphere, f0, beta,
440 & fCori, fCoriG,
441 & viscAh, viscAz, viscA4, viscAr, viscAstrain, viscAtension,
442 & diffKhT, diffKzT, diffK4T, diffKrT,
443 & diffKhS, diffKzS, diffK4S, diffKrS,
444 & delT, tauCD, rCD, freeSurfFac, implicSurfPress, implicDiv2Dflow,
445 & hFacMin, hFacMinDz, hFacInf, hFacSup,
446 & gravity, recip_Gravity, gBaro, rhonil, recip_rhonil,
447 & recip_rhoConst, rhoConst, tRef, sRef,
448 & endTime, chkPtFreq, pchkPtFreq, dumpFreq, tr1dumpFreq,
449 & diagFreq, taveFreq, monitorFreq,
450 & afFacMom, vfFacMom, pfFacMom, cfFacMom, foFacMom, mtFacMom,
451 & cosPower, cAdjFreq, omega,
452 & tauThetaClimRelax, lambdaThetaClimRelax,
453 & tauSaltClimRelax, lambdaSaltClimRelax,
454 & tauTr1ClimRelax, lambdaTr1ClimRelax,
455 & externForcingCycle, externForcingPeriod,
456 & convertFW2Salt, temp_EvPrRn, salt_EvPrRn, trac_EvPrRn,
457 & viscAp, diffKpT, diffKpS, hFacMinDr, hFacMinDp,
458 & horiVertRatio, recip_horiVertRatio,
459 & ivdc_kappa, Ro_SeaLevel,
460 & bottomDragLinear,bottomDragQuadratic
461
462 _RL cg2dTargetResidual
463 _RL cg2dTargetResWunit
464 _RL cg3dTargetResidual
465 _RL cg2dpcOffDFac
466 _RL delZ(Nr)
467 _RL delP(Nr)
468 _RL delR(Nr)
469 _RL delX(Nx)
470 _RL delY(Ny)
471 _RL deltaT
472 _RL deltaTClock
473 _RL deltaTmom
474 _RL deltaTtracer
475 _RL abeps
476 _RL phiMin
477 _RL thetaMin
478 _RL rSphere
479 _RL recip_RSphere
480 _RL f0
481 _RL freeSurfFac
482 _RL implicSurfPress
483 _RL implicDiv2Dflow
484 _RL hFacMin
485 _RL hFacMinDz
486 _RL hFacMinDp
487 _RL hFacMinDr
488 _RL hFacInf
489 _RL hFacSup
490 _RL beta
491 _RL viscAh
492 _RL viscAstrain
493 _RL viscAtension
494 _RL viscAz
495 _RL viscAp
496 _RL viscAr
497 _RL viscA4
498 _RL diffKhT
499 _RL diffKrT
500 _RL diffKzT
501 _RL diffKpT
502 _RL diffK4T
503 _RL diffKhS
504 _RL diffKrS
505 _RL diffKzS
506 _RL diffKpS
507 _RL diffK4S
508 _RL delt
509 _RL tauCD
510 _RL rCD
511 _RL gravity
512 _RL recip_gravity
513 _RL gBaro
514 _RL rhonil
515 _RL recip_rhonil
516 _RL rhoConst
517 _RL recip_rhoConst
518 _RL tRef(Nr)
519 _RL sRef(Nr)
520 _RS fCori(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
521 _RS fCoriG(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
522 _RL startTime
523 _RL endTime
524 _RL chkPtFreq
525 _RL pChkPtFreq
526 _RL dumpFreq
527 _RL tr1dumpFreq
528 _RL diagFreq
529 _RL taveFreq
530 _RL monitorFreq
531 _RL afFacMom
532 _RL vfFacMom
533 _RL pfFacMom
534 _RL cfFacMom
535 _RL foFacMom
536 _RL mTFacMom
537 _RL cosPower
538 _RL cAdjFreq
539 _RL omega
540 _RL tauThetaClimRelax
541 _RL lambdaThetaClimRelax
542 _RL tauSaltClimRelax
543 _RL lambdaSaltClimRelax
544 _RL tauTr1ClimRelax
545 _RL lambdaTr1ClimRelax
546 _RL externForcingCycle
547 _RL externForcingPeriod
548 _RL convertFW2Salt
549 _RL temp_EvPrRn
550 _RL salt_EvPrRn
551 _RL trac_EvPrRn
552 _RL horiVertRatio
553 _RL recip_horiVertRatio
554 _RL ivdc_kappa
555 _RL Ro_SeaLevel
556 _RL bottomDragLinear
557 _RL bottomDragQuadratic
558
559 COMMON /PARM_A/ HeatCapacity_Cp,recip_Cp,
560 & Lamba_theta
561 _RL HeatCapacity_Cp
562 _RL Lamba_theta
563 _RL recip_Cp
564
565 C Equation of State (polynomial coeffients)
566 COMMON /PARM_EOS_NL/ eosC,eosSig0,eosRefT,eosRefS
567 _RL eosC(9,Nr+1),eosSig0(Nr+1),eosRefT(Nr+1),eosRefS(Nr+1)
568 C Linear equation of state
569 C tAlpha :: Linear EOS thermal expansion coefficient ( 1/degree ).
570 C sBeta :: Linear EOS haline contraction coefficient.
571 COMMON /PARM_EOS_LIN/ tAlpha,sBeta,eosType
572 _RL tAlpha
573 _RL sBeta
574 character*(6) eosType
575
576 C Atmospheric physical parameters (Ideal Gas EOS, ...)
577 C atm_po :: standard reference pressure
578 C atm_cp :: specific heat (Cp) of the (dry) air at constant pressure
579 C atm_kappa :: kappa = R/Cp (R: constant of Ideal Gas EOS)
580 C Integr_GeoPot :: option to select the way we integrate the geopotential
581 C (still a subject of discussions ...)
582 COMMON /PARM_ATM/ atm_cp, atm_kappa, atm_po,
583 & Integr_GeoPot
584 _RL atm_cp, atm_kappa, atm_po
585 INTEGER Integr_GeoPot
586
587 C Logical flags for selecting packages
588 LOGICAL useKPP
589 LOGICAL useGMRedi
590 LOGICAL useOBCS
591 LOGICAL useAIM
592 LOGICAL useGrdchk
593 LOGICAL useECCO
594 LOGICAL useSHAP_FILT
595 LOGICAL useZONAL_FILT
596 LOGICAL useFLT
597 LOGICAL useSBO
598 LOGICAL useSEAICE
599 LOGICAL useDirectSolver
600 COMMON /PARM_PACKAGES/
601 & useKPP, useGMRedi, useOBCS, useAIM, useECCO,
602 & useSHAP_FILT, useZONAL_FILT, useGrdchk, useFLT,
603 & useSBO, useSEAICE, useDirectSolver

  ViewVC Help
Powered by ViewVC 1.1.22