/[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.41 - (hide annotations) (download)
Thu Aug 26 17:47:37 1999 UTC (24 years, 8 months ago) by adcroft
Branch: MAIN
Changes since 1.40: +5 -3 lines
File MIME type: text/plain
Added IVDC (Implicit Vertical Diffusion Convection).
Also facilitated a "convection counter" that is output through "diags".

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

  ViewVC Help
Powered by ViewVC 1.1.22