/[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.205 - (show annotations) (download)
Fri Oct 19 14:34:13 2007 UTC (16 years, 7 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint59i
Changes since 1.204: +30 -37 lines
File MIME type: text/plain
prepare for "clever pickup" implementation:
new header file: RESTART.h for internal parameters related to restart process
(move parameters from PARAMS.h & GAD.h to new header file RESTART.h)

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

  ViewVC Help
Powered by ViewVC 1.1.22