/[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.143 - (show annotations) (download)
Thu Feb 10 05:25:37 2005 UTC (19 years, 3 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57d_post
Changes since 1.142: +6 -3 lines
File MIME type: text/plain
Add flag(s) inAdExact to help D. grow seaice in the Sahara.

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

  ViewVC Help
Powered by ViewVC 1.1.22