/[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.191 - (show annotations) (download)
Tue Nov 28 22:44:44 2006 UTC (17 years, 5 months ago) by jmc
Branch: MAIN
Changes since 1.190: +36 -18 lines
File MIME type: text/plain
allow to read vertical vectors: tRef,sRef,delR & delRc from (binary) file ;
 + start adding anelastic-code parameters

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

  ViewVC Help
Powered by ViewVC 1.1.22