/[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.62 - (hide annotations) (download)
Thu Sep 13 17:43:55 2001 UTC (22 years, 8 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint40
Changes since 1.61: +3 -2 lines
File MIME type: text/plain
Added package "flt".
 o pkg/flt
 o verification/flt_example
 o visualization of trajectories supplied
 o works but output not available to testscript

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

  ViewVC Help
Powered by ViewVC 1.1.22