/[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.37 - (hide annotations) (download)
Mon Mar 22 15:54:03 1999 UTC (25 years, 2 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint20
Changes since 1.36: +10 -3 lines
File MIME type: text/plain
Modifications for non-hydrostatic ability + updates for open-boundaries.

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

  ViewVC Help
Powered by ViewVC 1.1.22