/[MITgcm]/MITgcm/model/inc/PARAMS.h
ViewVC logotype

Annotation of /MITgcm/model/inc/PARAMS.h

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.49.2.4 - (hide annotations) (download)
Mon May 7 17:09:42 2001 UTC (22 years, 11 months ago) by jmc
Branch: pre38
CVS Tags: pre38-close
Changes since 1.49.2.3: +4 -10 lines
File MIME type: text/plain
Further packaging of Shapiro Filters.

1 jmc 1.49.2.4 C $Header: /u/gcmpack/models/MITgcmUV/model/inc/PARAMS.h,v 1.49.2.3 2001/04/06 12:59:45 jmc Exp $
2 jmc 1.48 C $Name: $
3 cnh 1.1 C
4     C /==========================================================\
5     C | PARAMS.h |
6     C | o Header file defining model "parameters". |
7     C |==========================================================|
8 cnh 1.33 C | The values from the model standard input file are |
9 cnh 1.1 C | stored into the variables held here. Notes describing |
10     C | the parameters can also be found here. |
11     C \==========================================================/
12 cnh 1.13
13     C Macros for special grid options
14     #include "PARAMS_MACROS.h"
15 cnh 1.1
16     C-- Contants
17     C Useful physical values
18     Real*8 PI
19     PARAMETER ( PI = 3.14159265358979323844D0 )
20     Real*8 deg2rad
21     PARAMETER ( deg2rad = 2.D0*PI/360.D0 )
22    
23 cnh 1.7 C Symbolic values
24     C precXXXX - Used to indicate what precision to use for
25     C dumping model state.
26     INTEGER precFloat32
27 adcroft 1.36 PARAMETER ( precFloat32 = 32 )
28 cnh 1.7 INTEGER precFloat64
29 adcroft 1.36 PARAMETER ( precFloat64 = 64 )
30 cnh 1.28 C UNSET_xxx - Used to indicate variables that have not been given a value
31     Real*8 UNSET_FLOAT8
32     PARAMETER ( UNSET_FLOAT8 = 1.234567D5 )
33     Real*4 UNSET_FLOAT4
34     PARAMETER ( UNSET_FLOAT4 = 1.234567E5 )
35     _RL UNSET_RL
36     PARAMETER ( UNSET_RL = 1.234567D5 )
37     _RS UNSET_RS
38     PARAMETER ( UNSET_RS = 1.234567E5 )
39     INTEGER UNSET_I
40     PARAMETER ( UNSET_I = 123456789 )
41 cnh 1.7
42     C Checkpoint data
43     INTEGER maxNoChkptLev
44     PARAMETER ( maxNoChkptLev = 2 )
45    
46 cnh 1.1 C-- COMMON /PARM_C/ Character valued parameters used by the model.
47 cnh 1.7 C checkPtSuff - List of checkpoint file suffices
48 cnh 1.16 C bathyFile - File containing bathymetry. If not defined bathymetry
49     C is taken from inline function.
50     C hydrogThetaFile - File containing initial hydrographic data for potential
51     C temperature.
52     C hydrogSaltFile - File containing initial hydrographic data for salinity.
53     C zonalWindFile - File containing zonal wind data
54     C meridWindFile - File containing meridional wind data
55 cnh 1.18 C thetaClimFile - File containing theta climataology used
56     C in relaxation term -lambda(theta-theta*)
57     C saltClimFile - File containing salt climataology used
58     C in relaxation term -lambda(salt-salt*)
59 adcroft 1.37 C surfQfile - File containing surface heat flux
60 heimbach 1.45 C surfQswfile - File containing surface shortwave radiation
61 adcroft 1.37 C EmPmRfile - File containing surface fresh water flux
62 cnh 1.26 C buoyancyRelation - Flag used to indicate which relation to use to
63     C get buoyancy.
64 cnh 1.16 COMMON /PARM_C/ checkPtSuff,
65     & bathyFile, hydrogThetaFile, hydrogSaltFile,
66 cnh 1.18 & zonalWindFile, meridWindFile, thetaClimFile,
67 adcroft 1.37 & saltClimFile, buoyancyRelation,
68 heimbach 1.45 & EmPmRfile, surfQfile, surfQswfile,
69 adcroft 1.43 & uVelInitFile, vVelInitFile, pSurfInitFile
70 cnh 1.7 CHARACTER*(5) checkPtSuff(maxNoChkptLev)
71 cnh 1.16 CHARACTER*(MAX_LEN_FNAM) bathyFile
72     CHARACTER*(MAX_LEN_FNAM) hydrogThetaFile
73     CHARACTER*(MAX_LEN_FNAM) hydrogSaltFile
74     CHARACTER*(MAX_LEN_FNAM) zonalWindFile
75     CHARACTER*(MAX_LEN_FNAM) meridWindFile
76 cnh 1.18 CHARACTER*(MAX_LEN_FNAM) thetaClimFile
77     CHARACTER*(MAX_LEN_FNAM) saltClimFile
78 adcroft 1.37 CHARACTER*(MAX_LEN_FNAM) surfQfile
79 heimbach 1.45 CHARACTER*(MAX_LEN_FNAM) surfQswfile
80 adcroft 1.37 CHARACTER*(MAX_LEN_FNAM) EmPmRfile
81 cnh 1.26 CHARACTER*(MAX_LEN_FNAM) buoyancyRelation
82 adcroft 1.42 CHARACTER*(MAX_LEN_FNAM) uVelInitFile
83     CHARACTER*(MAX_LEN_FNAM) vVelInitFile
84 adcroft 1.43 CHARACTER*(MAX_LEN_FNAM) pSurfInitFile
85 cnh 1.7
86 cnh 1.1 C-- COMMON /PARM_I/ Integer valued parameters used by the model.
87     C cg2dMaxIters - Maximum number of iterations in the
88     C two-dimensional con. grad solver.
89     C cg2dChkResFreq - Frequency with which to check residual
90     C in con. grad solver.
91 cnh 1.33 C cg3dMaxIters - Maximum number of iterations in the
92     C three-dimensional con. grad solver.
93     C cg3dChkResFreq - Frequency with which to check residual
94     C in con. grad solver.
95 cnh 1.1 C nIter0 - Start time-step number of for this run
96     C nTimeSteps - Number of timesteps to execute
97     C numStepsPerPickup - For offline setup. Frequency of pickup
98     C of flow fields.
99 cnh 1.7 C writeStatePrec - Precision used for writing model state.
100     C writeBinaryPrec - Precision used for writing binary files
101     C readBinaryPrec - Precision used for reading binary files
102     C nCheckLev - Holds current checkpoint level
103 cnh 1.34
104 cnh 1.1 COMMON /PARM_I/
105     & cg2dMaxIters,
106     & cg2dChkResFreq,
107 cnh 1.33 & cg3dMaxIters,
108     & cg3dChkResFreq,
109 adcroft 1.41 & nIter0, nTimeSteps, nEndIter,
110 cnh 1.7 & numStepsPerPickup,
111     & writeStatePrec, nCheckLev,
112 cnh 1.34 & writeBinaryPrec, readBinaryPrec,
113 jmc 1.49.2.4 & zonal_filt_sinpow, zonal_filt_cospow
114 cnh 1.1 INTEGER cg2dMaxIters
115     INTEGER cg2dChkResFreq
116 cnh 1.33 INTEGER cg3dMaxIters
117     INTEGER cg3dChkResFreq
118 cnh 1.1 INTEGER nIter0
119     INTEGER nTimeSteps
120 adcroft 1.41 INTEGER nEndIter
121 cnh 1.1 INTEGER numStepsPerPickup
122 cnh 1.7 INTEGER writeStatePrec
123     INTEGER writeBinaryPrec
124     INTEGER readBinaryPrec
125     INTEGER nCheckLev
126 adcroft 1.46 INTEGER zonal_filt_sinpow
127     INTEGER zonal_filt_cospow
128 cnh 1.1
129     C-- COMMON /PARM_L/ Logical valued parameters used by the model.
130     C usingCartesianGrid - If TRUE grid generation will be in a cartesian
131     C coordinate frame.
132     C usingSphericalPolarGrid - If TRUE grid generation will be in a
133     C spherical polar frame.
134 adcroft 1.36 C no_slip_sides - Impose "no-slip" at lateral boundaries.
135     C no_slip_bottom- Impose "no-slip" at bottom boundary.
136 adcroft 1.46 C staggerTimeStep - enable a Stagger time stepping T,S Rho then U,V
137 cnh 1.1 C momViscosity - Flag which turns momentum friction terms on and off.
138     C momAdvection - Flag which turns advection of momentum on and off.
139     C momForcing - Flag which turns external forcing of momentum on
140     C and off.
141 cnh 1.9 C momPressureForcing - Flag which turns pressure term in momentum equation
142     C on and off.
143 cnh 1.15 C metricTerms - Flag which turns metric terms on or off.
144     C usingSphericalPolarMTerms - If TRUE use spherical polar metric terms.
145 cnh 1.1 C useCoriolis - Flag which turns the coriolis terms on and off.
146     C tempDiffusion - Flag which turns diffusion of temperature on
147     C and off.
148     C tempAdvection - Flag which turns advection of temperature on
149     C and off.
150     C tempForcing - Flag which turns external forcing of temperature on
151     C and off.
152     C saltDiffusion - Flag which turns diffusion of salinit on
153     C and off.
154     C saltAdvection - Flag which turns advection of salinit on
155     C and off.
156     C saltForcing - Flag which turns external forcing of salinit on
157     C and off.
158 cnh 1.8 C rigidLid - Set to true to use rigid lid
159 jmc 1.49.2.3 C implicitFreeSurface - Set to true to use implcit free surface
160     C exactConserv - Set to true to conserve exactly the total Volume
161     C nonlinFreeSurf - Set to true to use non-linear free surface
162 cnh 1.10 C momStepping - Turns momentum equation time-stepping off
163 cnh 1.15 C tempStepping - Turns temperature equation time-stepping off
164 adcroft 1.17 C saltStepping - Turns salinity equation time-stepping off
165 cnh 1.15 C useConstantF - Coriolis parameter set to f0
166     C useBetaPlaneF - Coriolis parameter set to f0 + beta.y
167     C useSphereF - Coriolis parameter set to 2.omega.sin(phi)
168 adcroft 1.14 C implicitDiffusion - Turns implicit vertical diffusion on
169 adcroft 1.40 C implicitViscosity - Turns implicit vertical viscosity on
170 cnh 1.18 C doThetaClimRelax - Set true if relaxation to temperature
171     C climatology is required.
172     C doSaltClimRelax - Set true if relaxation to salinity
173     C climatology is required.
174 adcroft 1.19 C periodicExternalForcing - Set true if forcing is time-dependant
175 cnh 1.29 C usingPCoords - Set to indicate that we are working in pressure
176     C coords.
177     C usingZCoords - Set to indicate that we are working in height
178     C coords.
179 adcroft 1.37 C nonHydrostatic - Using non-hydrostatic terms
180 adcroft 1.38 C globalFiles - Selects between "global" and "tiled" files
181     C allowFreezing - Allows water to freeze and form ice
182 adcroft 1.46 C groundAtK1 - put the surface(k=1) at the Lower Boundary (=ground)
183 cnh 1.1 COMMON /PARM_L/ usingCartesianGrid, usingSphericalPolarGrid,
184 adcroft 1.49.2.2 & usingCurvilinearGrid,
185 adcroft 1.36 & no_slip_sides,no_slip_bottom,
186 adcroft 1.46 & staggerTimeStep,
187 cnh 1.33 & momViscosity, momAdvection, momForcing, useCoriolis,
188     & momPressureForcing,tempDiffusion, tempAdvection, tempForcing,
189 cnh 1.8 & saltDiffusion, saltAdvection, saltForcing,
190 jmc 1.49.2.3 & rigidLid, implicitFreeSurface, exactConserv, nonlinFreeSurf,
191 adcroft 1.17 & momStepping, tempStepping, saltStepping,
192 cnh 1.15 & metricTerms, usingSphericalPolarMTerms,
193     & useConstantF, useBetaPlaneF, useSphereF,
194 adcroft 1.40 & implicitDiffusion, implicitViscosity,
195     & doThetaClimRelax, doSaltClimRelax,
196 adcroft 1.35 & periodicExternalForcing, usingPCoords, usingZCoords,
197 adcroft 1.46 & nonHydrostatic, globalFiles,
198     & allowFreezing, groundAtK1
199 cnh 1.1 LOGICAL usingCartesianGrid
200     LOGICAL usingSphericalPolarGrid
201 adcroft 1.49.2.2 LOGICAL usingCurvilinearGrid
202 cnh 1.15 LOGICAL usingSphericalPolarMTerms
203 adcroft 1.36 LOGICAL no_slip_sides
204     LOGICAL no_slip_bottom
205 adcroft 1.46 LOGICAL staggerTimeStep
206 cnh 1.1 LOGICAL momViscosity
207     LOGICAL momAdvection
208     LOGICAL momForcing
209 cnh 1.9 LOGICAL momPressureForcing
210 cnh 1.1 LOGICAL useCoriolis
211     LOGICAL tempDiffusion
212     LOGICAL tempAdvection
213     LOGICAL tempForcing
214     LOGICAL saltDiffusion
215     LOGICAL saltAdvection
216     LOGICAL saltForcing
217 cnh 1.8 LOGICAL rigidLid
218 jmc 1.49.2.3 LOGICAL implicitFreeSurface
219     LOGICAL exactConserv
220     LOGICAL nonlinFreeSurf
221 cnh 1.10 LOGICAL momStepping
222     LOGICAL tempStepping
223 adcroft 1.17 LOGICAL saltStepping
224 cnh 1.15 LOGICAL metricTerms
225     LOGICAL useConstantF
226     LOGICAL useBetaPlaneF
227     LOGICAL useSphereF
228 adcroft 1.14 LOGICAL implicitDiffusion
229 adcroft 1.40 LOGICAL implicitViscosity
230 cnh 1.18 LOGICAL doThetaClimRelax
231     LOGICAL doSaltClimRelax
232 adcroft 1.19 LOGICAL periodicExternalForcing
233 cnh 1.29 LOGICAL usingPCoords
234     LOGICAL usingZCoords
235 adcroft 1.37 LOGICAL nonHydrostatic
236 adcroft 1.38 LOGICAL globalFiles
237     LOGICAL allowFreezing
238 adcroft 1.46 LOGICAL groundAtK1
239 cnh 1.1
240     C-- COMMON /PARM_R/ "Real" valued parameters used by the model.
241 jmc 1.49.2.3 C gg2dTargetResidual
242     C - Target residual for cg2d solver; no unit (RHS normalisation)
243     C cg2dTargetResWunit
244     C - Target residual for cg2d solver; W unit (No RHS normalisation)
245 cnh 1.33 C cg3dTargetResidual
246     C - Target residual for cg3d solver.
247 cnh 1.7 C cg2dpcOffDFac - Averaging weight for preconditioner off-diagonal.
248     C Note. 20th May 1998
249     C I made a weird discovery! In the model paper we argue
250     C for the form of the preconditioner used here ( see
251     C A Finite-volume, Incompressible Navier-Stokes Model
252     C ...., Marshall et. al ). The algebra gives a simple
253     C 0.5 factor for the averaging of ac and aCw to get a
254     C symmettric pre-conditioner. By using a factor of 0.51
255     C i.e. scaling the off-diagonal terms in the
256     C preconditioner down slightly I managed to get the
257     C number of iterations for convergence in a test case to
258     C drop form 192 -> 134! Need to investigate this further!
259     C For now I have introduced a parameter cg2dpcOffDFac which
260     C defaults to 0.51 but can be set at runtime.
261 cnh 1.24 C delP - Vertical grid spacing ( Pa ).
262     C delZ - Vertical grid spacing ( m ).
263 cnh 1.23 C delR - Vertical grid spacing ( units of r ).
264 cnh 1.1 C delX - Separation between cell faces (m) or (deg), depending
265     C delY on input flags.
266     C gravity - Accel. due to gravity ( m/s^2 )
267 cnh 1.32 C recip_gravity and its inverse
268 cnh 1.8 C gBaro - Accel. due to gravity used in barotropic equation ( m/s^2 )
269 cnh 1.1 C ronil - Reference density
270 cnh 1.24 C rhoConst - Vertically constant reference density
271 cnh 1.3 C startTime - Start time for model ( s )
272 cnh 1.1 C phiMin - Latitude of southern most cell face.
273     C thetaMin - Longitude of western most cell face (this
274     C is an "inert" parameter but it is included
275     C to make geographical references simple.)
276 cnh 1.3 C rSphere - Radius of sphere for a spherical polar grid ( m ).
277 cnh 1.25 C recip_RSphere - Reciprocal radius of sphere ( m ).
278 cnh 1.1 C f0 - Reference coriolis parameter ( 1/s )
279     C ( Southern edge f for beta plane )
280     C beta - df/dy ( s^-1.m^-1 )
281 cnh 1.15 C omega - Angular velocity ( rad/s )
282 cnh 1.1 C viscAh - Eddy viscosity coeff. for mixing of
283     C momentum laterally ( m^2/s )
284     C viscAz - Eddy viscosity coeff. for mixing of
285     C momentum vertically ( m^2/s )
286 cnh 1.27 C viscAp - Eddy viscosity coeff. for mixing of
287     C momentum vertically ( Pa^2/s )
288 cnh 1.24 C viscAr - Eddy viscosity coeff. for mixing of
289     C momentum vertically ( units of r^2/s )
290 cnh 1.1 C viscA4 - Biharmonic viscosity coeff. for mixing of
291     C momentum laterally ( m^4/s )
292     C diffKhT - Laplacian diffusion coeff. for mixing of
293     C heat laterally ( m^2/s )
294     C diffKzT - Laplacian diffusion coeff. for mixing of
295     C heat vertically ( m^2/s )
296 cnh 1.27 C diffKpT - Laplacian diffusion coeff. for mixing of
297     C heat vertically ( Pa^2/s )
298 cnh 1.24 C diffKrT - Laplacian diffusion coeff. for mixing of
299     C heat vertically ( units of r^2/s )
300 cnh 1.1 C diffK4T - Biharmonic diffusion coeff. for mixing of
301     C heat laterally ( m^4/s )
302     C diffKhS - Laplacian diffusion coeff. for mixing of
303     C salt laterally ( m^2/s )
304     C diffKzS - Laplacian diffusion coeff. for mixing of
305     C salt vertically ( m^2/s )
306 cnh 1.27 C diffKpS - Laplacian diffusion coeff. for mixing of
307     C salt vertically ( Pa^2/s )
308 cnh 1.24 C diffKrS - Laplacian diffusion coeff. for mixing of
309     C salt vertically ( units of r^2/s )
310 cnh 1.1 C diffK4S - Biharmonic diffusion coeff. for mixing of
311     C salt laterally ( m^4/s )
312 cnh 1.3 C deltaT - Default timestep ( s )
313 cnh 1.7 C deltaTClock - Timestep used as model "clock". This determines the
314     C IO frequencies and is used in tagging output. It can
315     C be totally different to the dynamical time. Typically
316     C it will be the deep-water timestep for accelerated runs.
317     C Frequency of checkpointing and dumping of the model state
318     C are referenced to this clock. ( s )
319 cnh 1.3 C deltaTMom - Timestep for momemtum equations ( s )
320     C deltaTtracer - Timestep for tracer equations ( s )
321 cnh 1.8 C freesurfFac - Parameter to turn implicit free surface term on or off
322     C freesurfac = 1. uses implicit free surface
323     C freesurfac = 0. uses rigid lid
324 jmc 1.48 C implicSurfPress - parameter of the Crank-Nickelson time stepping :
325     C Implicit part of Surface Pressure Gradient ( 0-1 )
326     C implicDiv2Dflow - parameter of the Crank-Nickelson time stepping :
327     C Implicit part of barotropic flow Divergence ( 0-1 )
328 cnh 1.24 C hFacMin - Minimum fraction size of a cell (affects hFacC etc...)
329     C hFacMinDz - Minimum dimesional size of a cell (affects hFacC etc..., m)
330 cnh 1.28 C hFacMinDp - Minimum dimesional size of a cell (affects hFacC etc..., Pa)
331 cnh 1.24 C hFacMinDr - Minimum dimesional size of a cell (affects hFacC etc..., units of r)
332 jmc 1.49.2.3 C hFacInf - Threshold (inf and sup) for fraction size of surface cell
333     C hFacSup that control vanishing and creating levels
334 cnh 1.24 C tauCD - CD scheme coupling timescale ( 1/s )
335     C rCD - CD scheme normalised coupling parameter ( 0-1 )
336     C startTime - Starting time for this integration ( s ).
337     C endTime - Ending time for this integration ( s ).
338     C chkPtFreq - Frequency of rolling check pointing ( s ).
339     C pChkPtFreq - Frequency of permanent check pointing ( s ).
340     C dumpFreq - Frequency with which model state is written to
341     C post-processing files ( s ).
342     C afFacMom - Advection of momentum term tracer parameter
343     C vfFacMom - Momentum viscosity tracer parameter
344     C pfFacMom - Momentum pressure forcing tracer parameter
345     C cfFacMom - Coriolis term tracer parameter
346     C foFacMom - Momentum forcing tracer parameter
347     C mtFacMom - Metric terms tracer parameter
348 adcroft 1.39 C cosPower - Power of cosine of latitude to multiply viscosity
349 cnh 1.24 C cAdjFreq - Frequency of convective adjustment
350     C
351 cnh 1.27 C taveFreq - Frequency with which time-averaged model state is written to
352     C post-processing files ( s ).
353 cnh 1.18 C tauThetaClimRelax - Relaxation to climatology time scale ( s ).
354     C lambdaThetaClimRelax - Inverse time scale for relaxation ( 1/s ).
355     C tauSaltClimRelax - Relaxation to climatology time scale ( s ).
356     C lambdaSaltClimRelax - Inverse time scale for relaxation ( 1/s ).
357 adcroft 1.19 C externForcingPeriod - Is the period of which forcing varies (eg. 1 month)
358     C externForcingCycle - Is the repeat time of the forcing (eg. 1 year)
359     C (note: externForcingCycle must be an integer
360     C number times externForcingPeriod)
361 cnh 1.30 C horiVertRatio - Ratio on units in vertical to units in horizontal.
362 cnh 1.31 C recip_horiVertRatio ( 1 if horiz in m and vertical in m ).
363 cnh 1.30 C ( g*rho if horiz in m and vertical in Pa ).
364 cnh 1.34 C latFFTFiltLo - Low latitude for FFT filtering of latitude
365     C circles ( see filter*.F )
366 adcroft 1.46 C Ro_SeaLevel - standard position of Sea-Level in "R" coordinate, used as
367     C starting value (k=1) for vertical coordinate (rf(1)=Ro_SeaLevel)
368     C bottomDragLinear - Drag coefficient built in to core dynamics
369     C " Quadratic ( linear: 1/s, quadratic: 1/m )
370 jmc 1.49.2.3 COMMON /PARM_R/ cg2dTargetResidual, cg2dTargetResWunit,
371     & cg2dpcOffDFac, cg3dTargetResidual,
372 cnh 1.24 & delP, delZ, delR, delX, delY,
373 cnh 1.33 & deltaT,deltaTmom, deltaTtracer, deltaTClock,abeps, startTime,
374 adcroft 1.49.2.1 & phiMin, thetaMin, rSphere, recip_RSphere, f0, beta,
375     & fCori, fCoriG,
376 cnh 1.24 & viscAh, viscAz, viscA4, viscAr,
377     & diffKhT, diffKzT, diffK4T, diffKrT,
378 cnh 1.25 & diffKhS, diffKzS, diffK4S, diffKrS,
379 jmc 1.48 & delT, tauCD, rCD, freeSurfFac, implicSurfPress, implicDiv2Dflow,
380 jmc 1.49.2.3 & hFacMin, hFacMinDz, hFacInf, hFacSup,
381 cnh 1.32 & gravity, recip_Gravity, gBaro, rhonil, recip_rhonil,
382 cnh 1.27 & recip_rhoConst, rhoConst, tRef, sRef,
383 adcroft 1.20 & endTime, chkPtFreq, pchkPtFreq, dumpFreq, taveFreq,
384 cnh 1.15 & afFacMom, vfFacMom, pfFacMom, cfFacMom, foFacMom, mtFacMom,
385 adcroft 1.39 & cosPower,
386 cnh 1.18 & cAdjFreq, omega, tauThetaClimRelax, lambdaThetaClimRelax,
387 adcroft 1.19 & tauSaltClimRelax, lambdaSaltClimRelax,
388 cnh 1.27 & externForcingCycle, externForcingPeriod,
389 cnh 1.28 & viscAp, diffKpT, diffKpS, hFacMinDr, hFacMinDp,
390 cnh 1.34 & theta_S, specVol_S, horiVertRatio, recip_horiVertRatio,
391 adcroft 1.46 & latFFTFiltLo, ivdc_kappa, Ro_SeaLevel, zonal_filt_lat,
392     & bottomDragLinear,bottomDragQuadratic
393 cnh 1.27
394 adcroft 1.6 _RL cg2dTargetResidual
395 jmc 1.49.2.3 _RL cg2dTargetResWunit
396 cnh 1.33 _RL cg3dTargetResidual
397 cnh 1.7 _RL cg2dpcOffDFac
398 cnh 1.25 _RL delZ(Nr)
399     _RL delP(Nr)
400     _RL delR(Nr)
401 adcroft 1.6 _RL delX(Nx)
402     _RL delY(Ny)
403     _RL deltaT
404 cnh 1.7 _RL deltaTClock
405 adcroft 1.6 _RL deltaTmom
406     _RL deltaTtracer
407     _RL abeps
408     _RL phiMin
409     _RL thetaMin
410     _RL rSphere
411 cnh 1.25 _RL recip_RSphere
412 adcroft 1.6 _RL f0
413 cnh 1.8 _RL freeSurfFac
414 jmc 1.48 _RL implicSurfPress
415     _RL implicDiv2Dflow
416 adcroft 1.22 _RL hFacMin
417     _RL hFacMinDz
418 cnh 1.28 _RL hFacMinDp
419 cnh 1.27 _RL hFacMinDr
420 jmc 1.49.2.3 _RL hFacInf
421     _RL hFacSup
422 adcroft 1.6 _RL beta
423     _RL viscAh
424     _RL viscAz
425 cnh 1.27 _RL viscAp
426 cnh 1.25 _RL viscAr
427 adcroft 1.6 _RL viscA4
428     _RL diffKhT
429 cnh 1.25 _RL diffKrT
430 adcroft 1.6 _RL diffKzT
431 cnh 1.27 _RL diffKpT
432 adcroft 1.6 _RL diffK4T
433     _RL diffKhS
434 cnh 1.25 _RL diffKrS
435 adcroft 1.6 _RL diffKzS
436 cnh 1.27 _RL diffKpS
437 adcroft 1.6 _RL diffK4S
438     _RL delt
439     _RL tauCD
440     _RL rCD
441     _RL gravity
442 cnh 1.32 _RL recip_gravity
443 cnh 1.8 _RL gBaro
444 adcroft 1.6 _RL rhonil
445 cnh 1.25 _RL recip_rhonil
446     _RL rhoConst
447 cnh 1.26 _RL recip_rhoConst
448     _RL specVol_S(Nr)
449 cnh 1.25 _RL tRef(Nr)
450 cnh 1.26 _RL theta_S(Nr)
451 cnh 1.25 _RL sRef(Nr)
452 adcroft 1.49.2.1 _RS fCori(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
453     _RS fCoriG(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy)
454 adcroft 1.6 _RL startTime
455     _RL endTime
456     _RL chkPtFreq
457 cnh 1.7 _RL pChkPtFreq
458 adcroft 1.6 _RL dumpFreq
459 adcroft 1.20 _RL taveFreq
460 cnh 1.9 _RL afFacMom
461     _RL vfFacMom
462     _RL pfFacMom
463     _RL cfFacMom
464     _RL foFacMom
465 cnh 1.15 _RL mTFacMom
466 adcroft 1.39 _RL cosPower
467 cnh 1.9 _RL cAdjFreq
468 cnh 1.15 _RL omega
469 cnh 1.18 _RL tauThetaClimRelax
470     _RL lambdaThetaClimRelax
471     _RL tauSaltClimRelax
472     _RL lambdaSaltClimRelax
473 adcroft 1.19 _RL externForcingCycle
474     _RL externForcingPeriod
475 cnh 1.30 _RL horiVertRatio
476 cnh 1.31 _RL recip_horiVertRatio
477 cnh 1.34 _RL latFFTFiltLo
478 adcroft 1.41 _RL ivdc_kappa
479 adcroft 1.46 _RL Ro_SeaLevel
480     _RL zonal_filt_lat
481     _RL bottomDragLinear
482     _RL bottomDragQuadratic
483 adcroft 1.6
484 adcroft 1.36 COMMON /PARM_A/ HeatCapacity_Cp,recip_Cp,
485 cnh 1.1 & Lamba_theta
486 adcroft 1.6 _RL HeatCapacity_Cp
487     _RL Lamba_theta
488 adcroft 1.36 _RL recip_Cp
489 adcroft 1.5
490     C Equation of State (polynomial coeffients)
491     COMMON /PARM_EOS_NL/ eosC,eosSig0,eosRefT,eosRefS
492 cnh 1.25 _RL eosC(9,Nr+1),eosSig0(Nr+1),eosRefT(Nr+1),eosRefS(Nr+1)
493 adcroft 1.5 C Linear equation of state
494     C tAlpha - Linear EOS thermal expansion coefficient ( 1/degree ).
495     C sBeta - Linear EOS haline contraction coefficient.
496 adcroft 1.12 COMMON /PARM_EOS_LIN/ tAlpha,sBeta,eosType
497 adcroft 1.6 _RL tAlpha
498     _RL sBeta
499 adcroft 1.12 character*(6) eosType
500 cnh 1.7
501 adcroft 1.44 C Logical flags for selecting packages
502 heimbach 1.45 LOGICAL useKPP
503     LOGICAL useGMRedi
504 adcroft 1.46 LOGICAL useOBCS
505 jmc 1.49 LOGICAL useAIM
506 heimbach 1.45 LOGICAL useECCO
507 jmc 1.49.2.4 LOGICAL useSHAP_FILT
508 adcroft 1.44 COMMON /PARM_PACKAGES/
509 jmc 1.49.2.4 & useKPP, useGMRedi, useOBCS, useAIM, useECCO, useSHAP_FILT
510 jmc 1.49

  ViewVC Help
Powered by ViewVC 1.1.22