/[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.36 - (hide annotations) (download)
Tue Dec 15 00:20:34 1998 UTC (25 years, 5 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint19
Changes since 1.35: +10 -4 lines
File MIME type: text/plain
 o Added "natural BCs" as alternative to "virtual salt flux"
 o Re-difined precFloat32 and precFloat64 to be 32 and 64
   so that their values can be meaningfuly set in the data file
 o Modified read_write.F to create an exception if readBinaryPrec
   is not set
 o Replaced CPP control of viscous BCs with run-time control
 o Tidied up input-data precision (ie. ini_depths cnh_dbg...)
 o ini_forcing.F now initialises *all* forcing arrays to zero
 o Definitively tested verification experiments 0,1,2 and 4
   (3 is atmospheric set-up which is in a state of flux)

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

  ViewVC Help
Powered by ViewVC 1.1.22