/[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.112 - (show annotations) (download)
Wed Jun 9 14:03:35 2004 UTC (19 years, 10 months ago) by adcroft
Branch: MAIN
Changes since 1.111: +10 -1 lines
File MIME type: text/plain
Added vertical diffusivity profile (T/S) due to Bryan and Lewis, 1979.
New parameters:
 diffKrBL79surf - surface diffusivity
 diffKrBL79deep - deep diffusivity
 diffKrBL79Ho   - turning depth for arctan function
 diffKrBL79scl  - depth scale for arctan function
This diffusivity is added to all other diffusivities. The defaults are
set so as to give zero diffusivity.

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

  ViewVC Help
Powered by ViewVC 1.1.22