/[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.165 - (show annotations) (download)
Fri Sep 16 19:33:05 2005 UTC (18 years, 8 months ago) by baylor
Branch: MAIN
CVS Tags: checkpoint57s_post
Changes since 1.164: +11 -8 lines
File MIME type: text/plain
Parameter and config_summary changes for mom_calc_visc.F

1 C $Header: /u/gcmpack/MITgcm/model/inc/PARAMS.h,v 1.164 2005/08/24 23:08:35 jmc Exp $
2 C $Name: $
3 C
4
5 CBOP
6 C !ROUTINE: PARAMS.h
7 C !INTERFACE:
8 C #include PARAMS.h
9
10 C !DESCRIPTION:
11 C Header file defining model "parameters". The values from the
12 C model standard input file are stored into the variables held
13 C here. Notes describing the parameters can also be found here.
14
15 CEOP
16
17 C Macros for special grid options
18 #include "PARAMS_MACROS.h"
19
20 C-- Contants
21 C Useful physical values
22 Real*8 PI
23 PARAMETER ( PI = 3.14159265358979323844D0 )
24 Real*8 deg2rad
25 PARAMETER ( deg2rad = 2.D0*PI/360.D0 )
26
27 C Symbolic values
28 C precXXXX :: Used to indicate what precision to use for
29 C dumping model state.
30 INTEGER precFloat32
31 PARAMETER ( precFloat32 = 32 )
32 INTEGER precFloat64
33 PARAMETER ( precFloat64 = 64 )
34 C UNSET_xxx :: Used to indicate variables that have not been given a value
35 Real*8 UNSET_FLOAT8
36 PARAMETER ( UNSET_FLOAT8 = 1.234567D5 )
37 Real*4 UNSET_FLOAT4
38 PARAMETER ( UNSET_FLOAT4 = 1.234567E5 )
39 _RL UNSET_RL
40 PARAMETER ( UNSET_RL = 1.234567D5 )
41 _RS UNSET_RS
42 PARAMETER ( UNSET_RS = 1.234567E5 )
43 INTEGER UNSET_I
44 PARAMETER ( UNSET_I = 123456789 )
45
46 C Checkpoint data
47 INTEGER maxNoChkptLev
48 PARAMETER ( maxNoChkptLev = 2 )
49
50 C-- COMMON /PARM_C/ Character valued parameters used by the model.
51 C checkPtSuff :: List of checkpoint file suffices
52 C delXFile :: File containing X-spacing grid definition (1.D array)
53 C delYFile :: File containing Y-spacing grid definition (1.D array)
54 C horizGridFile :: File containing horizontal-grid definition
55 C (only when using curvilinear_grid)
56 C bathyFile :: File containing bathymetry. If not defined bathymetry
57 C is taken from inline function.
58 C topoFile :: File containing the topography of the surface (unit=m)
59 C (mainly used for the atmosphere = ground height).
60 C hydrogThetaFile :: File containing initial hydrographic data for potential
61 C temperature.
62 C hydrogSaltFile :: File containing initial hydrographic data for salinity.
63 C zonalWindFile :: File containing zonal wind data
64 C meridWindFile :: File containing meridional wind data
65 C thetaClimFile :: File containing theta climataology used
66 C in relaxation term -lambda(theta-theta*)
67 C saltClimFile :: File containing salt climataology used
68 C in relaxation term -lambda(salt-salt*)
69 C surfQfile :: File containing surface heat flux, excluding SW
70 C (old version, kept for backward compatibility)
71 C surfQnetFile :: File containing surface net heat flux
72 C surfQswFile :: File containing surface shortwave radiation
73 C dQdTfile :: File containing thermal relaxation coefficient
74 C EmPmRfile :: File containing surface fresh water flux
75 C saltFluxFile :: File containing surface salt flux
76 C pLoadFile :: File containing pressure loading
77 C eddyTauxFile :: File containing zonal Eddy stress data
78 C eddyTauyFile :: File containing meridional Eddy stress data
79 C buoyancyRelation :: Flag used to indicate which relation to use to
80 C get buoyancy.
81 C eosType :: choose the equation of state:
82 C LINEAR, POLY3, UNESCO, JMD95Z, JMD95P, MDJWF, IDEALGAS
83 C the_run_name :: string identifying the name of the model "run"
84 COMMON /PARM_C/ checkPtSuff,
85 & delXFile, delYFile, horizGridFile,
86 & bathyFile, topoFile,
87 & hydrogThetaFile, hydrogSaltFile,
88 & zonalWindFile, meridWindFile, thetaClimFile,
89 & saltClimFile, buoyancyRelation,
90 & EmPmRfile, saltFluxFile,
91 & surfQfile, surfQnetFile, surfQswFile,
92 & lambdaThetaFile, lambdaSaltFile,
93 & uVelInitFile, vVelInitFile, pSurfInitFile,
94 & dQdTfile, ploadFile,
95 & eddyTauxFile, eddyTauyFile,
96 & eosType, pickupSuff,
97 & mdsioLocalDir,
98 & the_run_name
99 CHARACTER*(5) checkPtSuff(maxNoChkptLev)
100 CHARACTER*(MAX_LEN_FNAM) delXFile
101 CHARACTER*(MAX_LEN_FNAM) delYFile
102 CHARACTER*(MAX_LEN_FNAM) horizGridFile
103 CHARACTER*(MAX_LEN_FNAM) bathyFile, topoFile
104 CHARACTER*(MAX_LEN_FNAM) hydrogThetaFile
105 CHARACTER*(MAX_LEN_FNAM) hydrogSaltFile
106 CHARACTER*(MAX_LEN_FNAM) zonalWindFile
107 CHARACTER*(MAX_LEN_FNAM) meridWindFile
108 CHARACTER*(MAX_LEN_FNAM) thetaClimFile
109 CHARACTER*(MAX_LEN_FNAM) saltClimFile
110 CHARACTER*(MAX_LEN_FNAM) surfQfile
111 CHARACTER*(MAX_LEN_FNAM) surfQnetFile
112 CHARACTER*(MAX_LEN_FNAM) surfQswFile
113 CHARACTER*(MAX_LEN_FNAM) EmPmRfile
114 CHARACTER*(MAX_LEN_FNAM) saltFluxFile
115 CHARACTER*(MAX_LEN_FNAM) buoyancyRelation
116 CHARACTER*(MAX_LEN_FNAM) uVelInitFile
117 CHARACTER*(MAX_LEN_FNAM) vVelInitFile
118 CHARACTER*(MAX_LEN_FNAM) pSurfInitFile
119 CHARACTER*(MAX_LEN_FNAM) dQdTfile
120 CHARACTER*(MAX_LEN_FNAM) ploadFile
121 CHARACTER*(MAX_LEN_FNAM) eddyTauxFile
122 CHARACTER*(MAX_LEN_FNAM) eddyTauyFile
123 CHARACTER*(MAX_LEN_FNAM) lambdaThetaFile
124 CHARACTER*(MAX_LEN_FNAM) lambdaSaltFile
125 CHARACTER*(MAX_LEN_FNAM) mdsioLocalDir
126 CHARACTER*(MAX_LEN_FNAM) the_run_name
127 CHARACTER*(6) eosType
128 CHARACTER*(10) pickupSuff
129
130 C-- COMMON /PARM_I/ Integer valued parameters used by the model.
131 C cg2dMaxIters :: Maximum number of iterations in the
132 C two-dimensional con. grad solver.
133 C cg2dChkResFreq :: Frequency with which to check residual
134 C in con. grad solver.
135 C cg2dPreCondFreq :: Frequency for updating cg2d preconditioner
136 C (non-linear free-surf.)
137 C cg3dMaxIters :: Maximum number of iterations in the
138 C three-dimensional con. grad solver.
139 C cg3dChkResFreq :: Frequency with which to check residual
140 C in con. grad solver.
141 C nIter0 :: Start time-step number of for this run
142 C nTimeSteps :: Number of timesteps to execute
143 C numStepsPerPickup :: For offline setup. Frequency of pickup
144 C of flow fields.
145 C writeStatePrec :: Precision used for writing model state.
146 C writeBinaryPrec :: Precision used for writing binary files
147 C readBinaryPrec :: Precision used for reading binary files
148 C nCheckLev :: Holds current checkpoint level
149 C nonlinFreeSurf :: option related to non-linear free surface
150 C =0 Linear free surface ; >0 Non-linear
151 C select_rStar :: option related to r* vertical coordinate
152 C =0 (default) use r coord. ; > 0 use r*
153 C tempAdvScheme :: Temp. Horiz.Advection scheme selector
154 C tempVertAdvScheme :: Temp. Vert. Advection scheme selector
155 C saltAdvScheme :: Salt. Horiz.advection scheme selector
156 C saltVertAdvScheme :: Salt. Vert. Advection scheme selector
157 C debugLevel :: debug level selector: higher -> more writing
158
159 COMMON /PARM_I/
160 & cg2dMaxIters,
161 & cg2dChkResFreq, cg2dPreCondFreq,
162 & cg3dMaxIters,
163 & cg3dChkResFreq,
164 & nIter0, nTimeSteps, nEndIter,
165 & numStepsPerPickup,
166 & writeStatePrec, nCheckLev,
167 & writeBinaryPrec, readBinaryPrec,
168 & nonlinFreeSurf, select_rStar,
169 & tempAdvScheme, tempVertAdvScheme,
170 & saltAdvScheme, saltVertAdvScheme,
171 & debugLevel
172 INTEGER cg2dMaxIters
173 INTEGER cg2dChkResFreq
174 INTEGER cg2dPreCondFreq
175 INTEGER cg3dMaxIters
176 INTEGER cg3dChkResFreq
177 INTEGER nIter0
178 INTEGER nTimeSteps
179 INTEGER nEndIter
180 INTEGER numStepsPerPickup
181 INTEGER writeStatePrec
182 INTEGER writeBinaryPrec
183 INTEGER readBinaryPrec
184 INTEGER nCheckLev
185 INTEGER nonlinFreeSurf
186 INTEGER select_rStar
187 INTEGER tempAdvScheme, tempVertAdvScheme
188 INTEGER saltAdvScheme, saltVertAdvScheme
189 INTEGER debugLevel
190
191 C
192 INTEGER debLevZero
193 PARAMETER(debLevZero=0)
194 INTEGER debLevA
195 PARAMETER(debLevA=1)
196 INTEGER debLevB
197 PARAMETER(debLevB=2)
198
199 C-- COMMON /PARM_L/ Logical valued parameters used by the model.
200 C usingCartesianGrid :: If TRUE grid generation will be in a cartesian
201 C coordinate frame.
202 C usingSphericalPolarGrid :: If TRUE grid generation will be in a
203 C spherical polar frame.
204 C usingCylindricalGrid :: If TRUE grid generation will be Cylindrical
205 C no_slip_sides :: Impose "no-slip" at lateral boundaries.
206 C no_slip_bottom :: Impose "no-slip" at bottom boundary.
207 C staggerTimeStep :: enable a Stagger time stepping T,S Rho then U,V
208 C momViscosity :: Flag which turns momentum friction terms on and off.
209 C momAdvection :: Flag which turns advection of momentum on and off.
210 C momForcing :: Flag which turns external forcing of momentum on
211 C and off.
212 C momPressureForcing :: Flag which turns pressure term in momentum equation
213 C on and off.
214 C metricTerms :: Flag which turns metric terms on or off.
215 C usingSphericalPolarMTerms :: If TRUE use spherical polar metric terms.
216 C useNHMTerms :: If TRUE use non-hydrostatic metric terms.
217 C useCoriolis :: Flag which turns the coriolis terms on and off.
218 C tempAdvection :: Flag which turns advection of temperature on
219 C and off.
220 C tempForcing :: Flag which turns external forcing of temperature on
221 C and off.
222 C saltAdvection :: Flag which turns advection of salinity on
223 C and off.
224 C saltForcing :: Flag which turns external forcing of salinity on
225 C and off.
226 C useRealFreshWaterFlux :: if True (=Natural BCS), treats P+R-E flux
227 C as a real Fresh Water (=> changes the Sea Level)
228 C if F, converts P+R-E to salt flux (no SL effect)
229 C useFullLeith :: Set to true to use full Leith viscosity (may be unstable
230 C on irregular grids)
231 C useAnisotropicViscAGridMax :: Set to true to use Alistair's latest
232 C anisotropic length scale. It is used only for maximum viscosity
233 C calculations. Alistair recommends a value of viscA*GridMax=.25
234 C useStrainTensionVisc:: Set to true to use Strain-Tension viscous terms
235 C rigidLid :: Set to true to use rigid lid
236 C implicitFreeSurface :: Set to true to use implcit free surface
237 C exactConserv :: Set to true to conserve exactly the total Volume
238 C uniformLin_PhiSurf :: Set to true to use a uniform Bo_surf in the
239 C linear relation Phi_surf = Bo_surf*eta
240 C momStepping :: Turns momentum equation time-stepping off
241 C tempStepping :: Turns temperature equation time-stepping off
242 C saltStepping :: Turns salinity equation time-stepping off
243 C useConstantF :: Coriolis parameter set to f0
244 C useBetaPlaneF :: Coriolis parameter set to f0 + beta.y
245 C useSphereF :: Coriolis parameter set to 2.omega.sin(phi)
246 C useCDscheme :: use CD-scheme to calculate Coriolis terms.
247 C useJamartWetPoints :: Use wet-point method for Coriolis (Jamart and Ozer, 1986)
248 C useJamartMomAdv :: Use wet-point method for V.I. non-linear term
249 C SadournyCoriolis :: use the enstrophy conserving scheme by Sadourny
250 C upwindVorticity :: bias interpolation of vorticity in the Coriolis term
251 C highOrderVorticity :: use 3rd/4th order interp. of vorticity (V.I., advection)
252 C upwindShear :: use 1rst order upwind interp. (V.I., vertical advection)
253 C useAbsVorticity :: work with f+zeta in Coriolis terms
254 C implicitDiffusion :: Turns implicit vertical diffusion on
255 C implicitViscosity :: Turns implicit vertical viscosity on
256 C tempImplVertAdv :: Turns on implicit vertical advection for Temperature
257 C saltImplVertAdv :: Turns on implicit vertical advection for Salinity
258 C momImplVertAdv :: Turns on implicit vertical advection for Momentum
259 C multiDimAdvection :: Flag that enable multi-dimension advection
260 C useMultiDimAdvec :: True if multi-dim advection is used at least once
261 C forcing_In_AB :: if False, put forcing (Temp,Salt,Tracers) contribution
262 C out off Adams-Bashforth time stepping.
263 C startFromPickupAB2 :: with AB-3 code, start from an AB-2 pickup
264 C doThetaClimRelax :: Set true if relaxation to temperature
265 C climatology is required.
266 C doSaltClimRelax :: Set true if relaxation to salinity
267 C climatology is required.
268 C periodicExternalForcing :: Set true if forcing is time-dependant
269 C usingPCoords :: Set to indicate that we are working in a pressure
270 C type coordinate (p or p*).
271 C usingZCoords :: Set to indicate that we are working in a height
272 C type coordinate (z or z*)
273 C fluidIsAir :: Set to indicate that the fluid major constituent
274 C is air
275 C fluidIsWater :: Set to indicate that the fluid major constituent
276 C is water
277 C useDynP_inEos_Zc :: use the dynamical pressure in EOS (with Z-coord.)
278 C this requires specific code for restart & exchange
279 C setCenterDr :: set cell Center depth and put Interface at the middle
280 C nonHydrostatic :: Using non-hydrostatic terms
281 C quasiHydrostatic :: Using non-hydrostatic terms in hydrostatic algorithm
282 C globalFiles :: Selects between "global" and "tiled" files
283 C useSingleCpuIO :: On SGI platforms, option globalFiles is either
284 C slow (f77) or does not work (f90). When
285 C useSingleCpuIO is set, mdsio_writefield.F
286 C outputs from master mpi process only.
287 C allowFreezing :: Allows surface water to freeze and form ice
288 C useOldFreezing :: use the old version (before checkpoint52a_pre, 2003-11-12)
289 C pickup_write_mdsio :: use mdsio to write pickups
290 C pickup_read_mdsio :: use mdsio to read pickups
291 C pickup_write_immed :: echo the pickup immediately (for conversion)
292 C timeave_mdsio :: use mdsio for timeave output
293 C snapshot_mdsio :: use mdsio for "snapshot" (dumpfreq/diagfreq) output
294 C monitor_stdio :: use stdio for monitor output
295 C calendarDumps :: When set, approximate months (30-31 days) and years (360-372 days)
296 C for parameters chkPtFreq, pChkPtFreq, taveFreq, SEAICE_taveFreq,
297 C KPP_taveFreq, and freq in pkg/diagnostics are converted to exact
298 C calendar months and years. Requires pkg/cal.
299 C dumpInitAndLast :: dumps model state to files at Initial (nIter0)
300 C & Last iteration, in addition multiple of dumpFreq iter.
301 COMMON /PARM_L/ usingCartesianGrid, usingSphericalPolarGrid,
302 & usingCurvilinearGrid, usingCylindricalGrid,
303 & no_slip_sides,no_slip_bottom,
304 & staggerTimeStep,
305 & momViscosity, momAdvection, momForcing, useCoriolis,
306 & momPressureForcing, vectorInvariantMomentum,
307 & tempAdvection, tempForcing,
308 & saltAdvection, saltForcing,
309 & useRealFreshWaterFlux,
310 & useFullLeith, useAnisotropicViscAGridMax, useStrainTensionVisc,
311 & rigidLid, implicitFreeSurface, exactConserv, uniformLin_PhiSurf,
312 & momStepping, tempStepping, saltStepping,
313 & metricTerms, usingSphericalPolarMTerms, useNHMTerms,
314 & useConstantF, useBetaPlaneF, useSphereF,
315 & useCDscheme,
316 & useEnergyConservingCoriolis, useJamartWetPoints, useJamartMomAdv,
317 & SadournyCoriolis, upwindVorticity, highOrderVorticity,
318 & useAbsVorticity, upwindShear,
319 & implicitDiffusion, implicitViscosity,
320 & tempImplVertAdv, saltImplVertAdv, momImplVertAdv,
321 & multiDimAdvection, useMultiDimAdvec, forcing_In_AB,
322 & doThetaClimRelax, doSaltClimRelax, doTr1ClimRelax,
323 & periodicExternalForcing,
324 & fluidIsAir, fluidIsWater,
325 & usingPCoords, usingZCoords, useDynP_inEos_Zc, setCenterDr,
326 & nonHydrostatic, quasiHydrostatic, globalFiles, useSingleCpuIO,
327 & allowFreezing, useOldFreezing,
328 & usePickupBeforeC35, usePickupBeforeC54, startFromPickupAB2,
329 & pickup_read_mdsio, pickup_write_mdsio, pickup_write_immed,
330 & timeave_mdsio, snapshot_mdsio, monitor_stdio,
331 & outputTypesInclusive, dumpInitAndLast, debugMode,
332 & inAdMode, inAdTrue, inAdFalse, inAdExact,
333 & calendarDumps
334
335 LOGICAL usingCartesianGrid
336 LOGICAL usingSphericalPolarGrid
337 LOGICAL usingCylindricalGrid
338 LOGICAL usingCurvilinearGrid
339 LOGICAL usingSphericalPolarMTerms
340 LOGICAL useNHMTerms
341 LOGICAL no_slip_sides
342 LOGICAL no_slip_bottom
343 LOGICAL staggerTimeStep
344 LOGICAL momViscosity
345 LOGICAL momAdvection
346 LOGICAL momForcing
347 LOGICAL momPressureForcing
348 LOGICAL useCoriolis
349 LOGICAL vectorInvariantMomentum
350 LOGICAL tempAdvection
351 LOGICAL tempForcing
352 LOGICAL saltAdvection
353 LOGICAL saltForcing
354 LOGICAL useRealFreshWaterFlux
355 LOGICAL useFullLeith
356 LOGICAL useAnisotropicViscAGridMax
357 LOGICAL useStrainTensionVisc
358 LOGICAL rigidLid
359 LOGICAL implicitFreeSurface
360 LOGICAL exactConserv
361 LOGICAL uniformLin_PhiSurf
362 LOGICAL momStepping
363 LOGICAL tempStepping
364 LOGICAL saltStepping
365 LOGICAL metricTerms
366 LOGICAL useConstantF
367 LOGICAL useBetaPlaneF
368 LOGICAL useSphereF
369 LOGICAL useCDscheme
370 LOGICAL useEnergyConservingCoriolis
371 LOGICAL useJamartWetPoints
372 LOGICAL useJamartMomAdv
373 LOGICAL SadournyCoriolis
374 LOGICAL upwindVorticity
375 LOGICAL highOrderVorticity
376 LOGICAL useAbsVorticity
377 LOGICAL upwindShear
378 LOGICAL implicitDiffusion
379 LOGICAL implicitViscosity
380 LOGICAL tempImplVertAdv
381 LOGICAL saltImplVertAdv
382 LOGICAL momImplVertAdv
383 LOGICAL multiDimAdvection
384 LOGICAL useMultiDimAdvec
385 LOGICAL forcing_In_AB
386 LOGICAL doThetaClimRelax
387 LOGICAL doSaltClimRelax
388 LOGICAL doTr1ClimRelax
389 LOGICAL periodicExternalForcing
390 LOGICAL fluidIsAir
391 LOGICAL fluidIsWater
392 LOGICAL usingPCoords
393 LOGICAL usingZCoords
394 LOGICAL useDynP_inEos_Zc
395 LOGICAL setCenterDr
396 LOGICAL nonHydrostatic
397 LOGICAL quasiHydrostatic
398 LOGICAL globalFiles
399 LOGICAL useSingleCpuIO
400 LOGICAL allowFreezing
401 LOGICAL useOldFreezing
402 LOGICAL usePickupBeforeC35
403 LOGICAL usePickupBeforeC54
404 LOGICAL startFromPickupAB2
405 LOGICAL dumpInitAndLast
406 LOGICAL debugMode
407 LOGICAL pickup_read_mdsio, pickup_write_mdsio
408 LOGICAL pickup_write_immed
409 LOGICAL timeave_mdsio, snapshot_mdsio, monitor_stdio
410 LOGICAL outputTypesInclusive
411 LOGICAL inAdMode, inAdTrue, inAdFalse, inAdExact
412 LOGICAL calendarDumps
413
414 C-- COMMON /PARM_R/ "Real" valued parameters used by the model.
415 C cg2dTargetResidual
416 C :: Target residual for cg2d solver; no unit (RHS normalisation)
417 C cg2dTargetResWunit
418 C :: Target residual for cg2d solver; W unit (No RHS normalisation)
419 C cg3dTargetResidual
420 C :: Target residual for cg3d solver.
421 C cg2dpcOffDFac :: Averaging weight for preconditioner off-diagonal.
422 C Note. 20th May 1998
423 C I made a weird discovery! In the model paper we argue
424 C for the form of the preconditioner used here ( see
425 C A Finite-volume, Incompressible Navier-Stokes Model
426 C ...., Marshall et. al ). The algebra gives a simple
427 C 0.5 factor for the averaging of ac and aCw to get a
428 C symmettric pre-conditioner. By using a factor of 0.51
429 C i.e. scaling the off-diagonal terms in the
430 C preconditioner down slightly I managed to get the
431 C number of iterations for convergence in a test case to
432 C drop form 192 -> 134! Need to investigate this further!
433 C For now I have introduced a parameter cg2dpcOffDFac which
434 C defaults to 0.51 but can be set at runtime.
435 C delR :: Vertical grid spacing ( units of r ).
436 C delRc :: Vertical grid spacing between cell centers (r unit).
437 C delX :: Separation between cell faces (m) or (deg), depending
438 C delY on input flags.
439 C gravity :: Accel. due to gravity ( m/s^2 )
440 C recip_gravity and its inverse
441 C gBaro :: Accel. due to gravity used in barotropic equation ( m/s^2 )
442 C rhoNil :: Reference density for the linear equation of state
443 C rhoConst :: Vertically constant reference density
444 C rhoConstFresh :: Constant reference density for fresh water (rain)
445 C tRef :: reference vertical profile for potential temperature
446 C sRef :: reference vertical profile for salinity/specific humidity
447 C phiMin :: Latitude of southern most cell face.
448 C thetaMin :: Longitude of western most cell face (this
449 C is an "inert" parameter but it is included
450 C to make geographical references simple.)
451 C rSphere :: Radius of sphere for a spherical polar grid ( m ).
452 C recip_RSphere :: Reciprocal radius of sphere ( m ).
453 C f0 :: Reference coriolis parameter ( 1/s )
454 C ( Southern edge f for beta plane )
455 C beta :: df/dy ( s^-1.m^-1 )
456 C omega :: Angular velocity ( rad/s )
457 C rotationPeriod :: Rotation period (s) (= 2.pi/omega)
458 C viscAh :: Eddy viscosity coeff. for mixing of
459 C momentum laterally ( m^2/s )
460 C viscAhW :: Eddy viscosity coeff. for mixing of vertical
461 C momentum laterally, no effect for hydrostatic
462 C model, defaults to viscAh if unset ( m^2/s )
463 C viscAr :: Eddy viscosity coeff. for mixing of
464 C momentum vertically ( units of r^2/s )
465 C viscA4 :: Biharmonic viscosity coeff. for mixing of
466 C momentum laterally ( m^4/s )
467 C viscA4W :: Biharmonic viscosity coeff. for mixing of vertical
468 C momentum laterally, no effect for hydrostatic
469 C model, defaults to viscA4 if unset ( m^2/s )
470 C viscAhD :: Eddy viscosity coeff. for mixing of momentum laterally
471 C (act on Divergence part) ( m^2/s )
472 C viscAhZ :: Eddy viscosity coeff. for mixing of momentum laterally
473 C (act on Vorticity part) ( m^2/s )
474 C viscA4D :: Biharmonic viscosity coeff. for mixing of momentum laterally
475 C (act on Divergence part) ( m^4/s )
476 C viscA4Z :: Biharmonic viscosity coeff. for mixing of momentum laterally
477 C (act on Vorticity part) ( m^4/s )
478 C viscC2leith :: Leith non-dimensional viscosity factor (grad(vort))
479 C viscC2leithD :: Modified Leith non-dimensional viscosity factor (grad(div))
480 C viscC2smag :: Smagorinsky non-dimensional viscosity factor (harmonic)
481 C viscC4smag :: Smagorinsky non-dimensional viscosity factor (biharmonic)
482 C viscAhMax :: Maximum eddy viscosity coeff. for mixing of
483 C momentum laterally ( m^2/s )
484 C viscAhGridMax:: maximum and minimum harmonic viscosity coefficients ...
485 C viscAhGridMin:: in terms of non-dimensional grid-size dependent viscosity
486 C viscA4Max :: Maximum biharmonic viscosity coeff. for mixing of
487 C momentum laterally ( m^4/s )
488 C viscAhGrid:: non-dimensional grid-size dependent viscosity
489 C viscA4Grid:: non-dimensional grid-size dependent bi-harmonic viscosity
490 C viscA4GridMax:: maximum and minimum biharmonic viscosity coefficients ...
491 C viscA4GridMin:: in terms of non-dimensional grid-size dependent viscosity
492 C viscC4leith :: Leith non-dimensional viscosity factor (grad(vort))
493 C viscC4leithD :: Modified Leith non-dimensional viscosity factor (grad(div))
494 C diffKhT :: Laplacian diffusion coeff. for mixing of
495 C heat laterally ( m^2/s )
496 C diffKrNrT :: vertical profile of Laplacian diffusion coeff.
497 C for mixing of heat vertically ( units of r^2/s )
498 C diffK4T :: Biharmonic diffusion coeff. for mixing of
499 C heat laterally ( m^4/s )
500 C diffKhS :: Laplacian diffusion coeff. for mixing of
501 C salt laterally ( m^2/s )
502 C diffKrNrS :: vertical profile of Laplacian diffusion coeff.
503 C for mixing of salt vertically ( units of r^2/s ),
504 C diffK4S :: Biharmonic diffusion coeff. for mixing of
505 C salt laterally ( m^4/s )
506 C diffKrBL79surf :: T/S surface diffusivity (m^2/s) Bryan and Lewis, 1979
507 C diffKrBL79deep :: T/S deep diffusivity (m^2/s) Bryan and Lewis, 1979
508 C diffKrBL79scl :: depth scale for arctan fn (m) Bryan and Lewis, 1979
509 C diffKrBL79Ho :: depth offset for arctan fn (m) Bryan and Lewis, 1979
510 C deltaT :: Default timestep ( s )
511 C deltaTClock :: Timestep used as model "clock". This determines the
512 C IO frequencies and is used in tagging output. It can
513 C be totally different to the dynamical time. Typically
514 C it will be the deep-water timestep for accelerated runs.
515 C Frequency of checkpointing and dumping of the model state
516 C are referenced to this clock. ( s )
517 C deltaTMom :: Timestep for momemtum equations ( s )
518 C dTtracerLev :: Timestep for tracer equations ( s ), function of level k
519 C deltaTfreesurf :: Timestep for free-surface equation ( s )
520 C freesurfFac :: Parameter to turn implicit free surface term on or off
521 C freesurfac = 1. uses implicit free surface
522 C freesurfac = 0. uses rigid lid
523 C abEps :: Adams-Bashforth-2 stabilizing weight
524 C alph_AB :: Adams-Bashforth-3 primary factor
525 C beta_AB :: Adams-Bashforth-3 secondary factor
526 C implicSurfPress :: parameter of the Crank-Nickelson time stepping :
527 C Implicit part of Surface Pressure Gradient ( 0-1 )
528 C implicDiv2Dflow :: parameter of the Crank-Nickelson time stepping :
529 C Implicit part of barotropic flow Divergence ( 0-1 )
530 C hFacMin :: Minimum fraction size of a cell (affects hFacC etc...)
531 C hFacMinDz :: Minimum dimesional size of a cell (affects hFacC etc..., m)
532 C hFacMinDp :: Minimum dimesional size of a cell (affects hFacC etc..., Pa)
533 C hFacMinDr :: Minimum dimesional size of a cell (affects hFacC etc..., units of r)
534 C hFacInf :: Threshold (inf and sup) for fraction size of surface cell
535 C hFacSup that control vanishing and creating levels
536 C tauCD :: CD scheme coupling timescale ( 1/s )
537 C rCD :: CD scheme normalised coupling parameter ( 0-1 )
538 C baseTime :: model base time (time origin) = time @ iteration zero
539 C startTime :: Starting time for this integration ( s ).
540 C endTime :: Ending time for this integration ( s ).
541 C chkPtFreq :: Frequency of rolling check pointing ( s ).
542 C pChkPtFreq :: Frequency of permanent check pointing ( s ).
543 C dumpFreq :: Frequency with which model state is written to
544 C post-processing files ( s ).
545 C diagFreq :: Frequency with which model writes diagnostic output
546 C of intermediate quantities.
547 C afFacMom :: Advection of momentum term tracer parameter
548 C vfFacMom :: Momentum viscosity tracer parameter
549 C pfFacMom :: Momentum pressure forcing tracer parameter
550 C cfFacMom :: Coriolis term tracer parameter
551 C foFacMom :: Momentum forcing tracer parameter
552 C mtFacMom :: Metric terms tracer parameter
553 C cosPower :: Power of cosine of latitude to multiply viscosity
554 C cAdjFreq :: Frequency of convective adjustment
555 C
556 C taveFreq :: Frequency with which time-averaged model state
557 C is written to post-processing files ( s ).
558 C tave_lastIter :: (for state variable only) fraction of the last time
559 C step (of each taveFreq period) put in the time average.
560 C (fraction for 1rst iter = 1 - tave_lastIter)
561 C tauThetaClimRelax :: Relaxation to climatology time scale ( s ).
562 C tauSaltClimRelax :: Relaxation to climatology time scale ( s ).
563 C latBandClimRelax :: latitude band where Relaxation to Clim. is applied,
564 C i.e. where |yC| <= latBandClimRelax
565 C externForcingPeriod :: Is the period of which forcing varies (eg. 1 month)
566 C externForcingCycle :: Is the repeat time of the forcing (eg. 1 year)
567 C (note: externForcingCycle must be an integer
568 C number times externForcingPeriod)
569 C convertFW2Salt :: salinity, used to convert Fresh-Water Flux to Salt Flux
570 C (use model surface (local) value if set to -1)
571 C temp_EvPrRn :: temperature of Rain & Evap.
572 C salt_EvPrRn :: salinity of Rain & Evap.
573 C (notes: a) tracer content of Rain/Evap only used if both
574 C NonLin_FrSurf & useRealFreshWater are set.
575 C b) use model surface (local) value if set to UNSET_RL)
576 C horiVertRatio :: Ratio on units in vertical to units in horizontal.
577 C recip_horiVertRatio ( 1 if horiz in m and vertical in m ).
578 C ( g*rho if horiz in m and vertical in Pa ).
579 C Ro_SeaLevel :: standard position of Sea-Level in "R" coordinate, used as
580 C starting value (k=1) for vertical coordinate (rf(1)=Ro_SeaLevel)
581 C bottomDragLinear :: Drag coefficient built in to core dynamics
582 C --"-"-- Quadratic ( linear: 1/s, quadratic: 1/m )
583 COMMON /PARM_R/ cg2dTargetResidual, cg2dTargetResWunit,
584 & cg2dpcOffDFac, cg3dTargetResidual,
585 & delR, delRc, delX, delY,
586 & deltaT, deltaTmom, dTtracerLev, deltaTfreesurf, deltaTClock,
587 & abEps, alph_AB, beta_AB,
588 & phiMin, thetaMin, rSphere, recip_RSphere, f0, beta,
589 & viscAh, viscAhW, viscAhMax,
590 & viscAhGrid, viscAhGridMax, viscAhGridMin,
591 & viscC2leith, viscC2leithD,
592 & viscC2smag, viscC4smag,
593 & viscAhD, viscAhZ, viscA4D, viscA4Z,
594 & viscA4, viscA4W,
595 & viscA4Max, viscA4Grid, viscA4GridMax, viscA4GridMin,
596 & viscC4leith, viscC4leithD, viscAr,
597 & diffKhT, diffK4T, diffKrNrT,
598 & diffKhS, diffK4S, diffKrNrS,
599 & diffKrBL79surf, diffKrBL79deep, diffKrBL79scl, diffKrBL79Ho,
600 & delT, tauCD, rCD, freeSurfFac, implicSurfPress, implicDiv2Dflow,
601 & hFacMin, hFacMinDz, hFacInf, hFacSup,
602 & gravity, recip_Gravity, gBaro, rhonil, recip_rhonil,
603 & recip_rhoConst, rhoConst,
604 & rhoConstFresh, convertEmP2rUnit, tRef, sRef,
605 & baseTime, startTime, endTime,
606 & chkPtFreq, pchkPtFreq, dumpFreq, adjDumpFreq,
607 & diagFreq, taveFreq, tave_lastIter, monitorFreq, adjMonitorFreq,
608 & afFacMom, vfFacMom, pfFacMom, cfFacMom, foFacMom, mtFacMom,
609 & cosPower, cAdjFreq, omega, rotationPeriod,
610 & tauThetaClimRelax,
611 & tauSaltClimRelax,
612 & tauTr1ClimRelax, lambdaTr1ClimRelax, latBandClimRelax,
613 & externForcingCycle, externForcingPeriod,
614 & convertFW2Salt, temp_EvPrRn, salt_EvPrRn,
615 & hFacMinDr, hFacMinDp,
616 & horiVertRatio, recip_horiVertRatio,
617 & ivdc_kappa, Ro_SeaLevel,
618 & bottomDragLinear,bottomDragQuadratic,nh_Am2,
619 & tCylIn, tCylOut
620
621 _RL nh_Am2
622 _RL cg2dTargetResidual
623 _RL cg2dTargetResWunit
624 _RL cg3dTargetResidual
625 _RL cg2dpcOffDFac
626 _RL delR(Nr)
627 _RL delRc(Nr+1)
628 _RL delX(Nx)
629 _RL delY(Ny)
630 _RL deltaT
631 _RL deltaTClock
632 _RL deltaTmom
633 _RL dTtracerLev(Nr)
634 _RL deltaTfreesurf
635 _RL abEps, alph_AB, beta_AB
636 _RL phiMin
637 _RL thetaMin
638 _RL rSphere
639 _RL recip_RSphere
640 _RL f0
641 _RL freeSurfFac
642 _RL implicSurfPress
643 _RL implicDiv2Dflow
644 _RL hFacMin
645 _RL hFacMinDz
646 _RL hFacMinDp
647 _RL hFacMinDr
648 _RL hFacInf
649 _RL hFacSup
650 _RL beta
651 _RL viscAh
652 _RL viscAhW
653 _RL viscAhD
654 _RL viscAhZ
655 _RL viscAhMax
656 _RL viscAhGrid
657 _RL viscAhGridMax
658 _RL viscAhGridMin
659 _RL viscC2leith
660 _RL viscC2leithD
661 _RL viscC2smag
662 _RL viscC4smag
663 _RL viscAr
664 _RL viscA4
665 _RL viscA4W
666 _RL viscA4D
667 _RL viscA4Z
668 _RL viscA4Max
669 _RL viscA4Grid, viscA4GridMax, viscA4GridMin
670 _RL viscC4leith
671 _RL viscC4leithD
672 _RL diffKhT
673 _RL diffKrNrT(Nr)
674 _RL diffK4T
675 _RL diffKhS
676 _RL diffKrNrS(Nr)
677 _RL diffK4S
678 _RL diffKrBL79surf
679 _RL diffKrBL79deep
680 _RL diffKrBL79scl
681 _RL diffKrBL79Ho
682 _RL delt
683 _RL tauCD
684 _RL rCD
685 _RL gravity
686 _RL recip_gravity
687 _RL gBaro
688 _RL rhonil
689 _RL recip_rhonil
690 _RL rhoConst
691 _RL recip_rhoConst
692 _RL rhoConstFresh
693 _RL convertEmP2rUnit
694 _RL tRef(Nr)
695 _RL sRef(Nr)
696 _RL baseTime
697 _RL startTime
698 _RL endTime
699 _RL chkPtFreq
700 _RL pChkPtFreq
701 _RL dumpFreq
702 _RL adjDumpFreq
703 _RL diagFreq
704 _RL taveFreq
705 _RL tave_lastIter
706 _RL monitorFreq
707 _RL adjMonitorFreq
708 _RL afFacMom
709 _RL vfFacMom
710 _RL pfFacMom
711 _RL cfFacMom
712 _RL foFacMom
713 _RL mTFacMom
714 _RL cosPower
715 _RL cAdjFreq
716 _RL omega
717 _RL rotationPeriod
718 _RL tauThetaClimRelax
719 _RL tauSaltClimRelax
720 _RL tauTr1ClimRelax
721 _RL lambdaTr1ClimRelax
722 _RL latBandClimRelax
723 _RL externForcingCycle
724 _RL externForcingPeriod
725 _RL convertFW2Salt
726 _RL temp_EvPrRn
727 _RL salt_EvPrRn
728 _RL horiVertRatio
729 _RL recip_horiVertRatio
730 _RL ivdc_kappa
731 _RL Ro_SeaLevel
732 _RL bottomDragLinear
733 _RL bottomDragQuadratic
734 _RL tCylIn
735 _RL tCylOut
736
737 C-- COMMON /PARM_A/ Thermodynamics constants ?
738 COMMON /PARM_A/ HeatCapacity_Cp,recip_Cp
739 _RL HeatCapacity_Cp
740 _RL recip_Cp
741
742 C-- COMMON /PARM_ATM/ Atmospheric physical parameters (Ideal Gas EOS, ...)
743 C celsius2K :: convert centigrade (Celsius) degree to Kelvin
744 C atm_Po :: standard reference pressure
745 C atm_Cp :: specific heat (Cp) of the (dry) air at constant pressure
746 C atm_Rd :: gas constant for dry air
747 C atm_kappa :: kappa = R/Cp (R: constant of Ideal Gas EOS)
748 C atm_Rq :: water vapour specific volume anomaly relative to dry air
749 C (e.g. typical value = (29/18 -1) 10^-3 with q [g/kg])
750 C integr_GeoPot :: option to select the way we integrate the geopotential
751 C (still a subject of discussions ...)
752 C selectFindRoSurf :: select the way surf. ref. pressure (=Ro_surf) is
753 C derived from the orography. Implemented: 0,1 (see INI_P_GROUND)
754 COMMON /PARM_ATM/
755 & celsius2K,
756 & atm_Cp, atm_Rd, atm_kappa, atm_Rq, atm_Po,
757 & integr_GeoPot, selectFindRoSurf
758 _RL celsius2K
759 _RL atm_Po, atm_Cp, atm_Rd, atm_kappa, atm_Rq
760 INTEGER integr_GeoPot, selectFindRoSurf
761
762 C Logical flags for selecting packages
763 LOGICAL useOPPS
764 LOGICAL usePP81
765 LOGICAL useMY82
766 LOGICAL useGGL90
767 LOGICAL useKPP
768 LOGICAL useGMRedi
769 LOGICAL useOBCS
770 LOGICAL useAIM
771 LOGICAL useLand
772 LOGICAL useGrdchk
773 LOGICAL useECCO
774 LOGICAL useSHAP_FILT
775 LOGICAL useZONAL_FILT
776 LOGICAL useFLT
777 LOGICAL usePTRACERS
778 LOGICAL useGCHEM
779 LOGICAL useSBO
780 LOGICAL useSEAICE
781 LOGICAL useBulkForce
782 LOGICAL useThSIce
783 LOGICAL usefizhi
784 LOGICAL usegridalt
785 LOGICAL usediagnostics
786 LOGICAL useEBM
787 LOGICAL useMNC
788 LOGICAL useMATRIX
789 LOGICAL useRunClock
790 COMMON /PARM_PACKAGES/
791 & useKPP, useGMRedi, useOBCS, useAIM, useLand, useECCO,
792 & useSHAP_FILT, useZONAL_FILT, useGrdchk, useFLT,
793 & usePTRACERS, useGCHEM,
794 & useSBO, useSEAICE, useThSIce, useBulkForce,
795 & usefizhi, usegridalt, usediagnostics, useEBM, useMNC,
796 & usePP81, useMY82, useOPPS, useGGL90, useMATRIX,
797 & useRunClock
798
799 CEH3 ;;; Local Variables: ***
800 CEH3 ;;; mode:fortran ***
801 CEH3 ;;; End: ***

  ViewVC Help
Powered by ViewVC 1.1.22