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

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

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


Revision 1.42 - (show annotations) (download)
Mon Aug 30 18:29:21 1999 UTC (24 years, 8 months ago) by adcroft
Branch: MAIN
CVS Tags: checkpoint25, checkpoint26
Changes since 1.41: +5 -2 lines
File MIME type: text/plain
Corrected interaction between OBCs and algorithm. The
positioning of set_obcs() within the time-stepping sequence
is crucial for stable open-boundaries. Forcing the boundaries
with time-dependent flow previously led to horrible results...

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

  ViewVC Help
Powered by ViewVC 1.1.22