source: branches/numpy/anuga/config.py @ 7176

Last change on this file since 7176 was 7176, checked in by rwilson, 15 years ago

Back-merge from Numeric to numpy.

File size: 9.1 KB
Line 
1"""Module where global ANUGA model parameters and default values are set
2"""
3
4import os
5import sys
6import numpy as num
7
8
9################################################################################
10# Numerical constants
11################################################################################
12
13epsilon = 1.0e-12                    # Smallest number - used for safe division
14max_float = 1.0e36                   # Largest number - used to initialise
15                                     # (max, min) ranges
16default_smoothing_parameter = 0.001  # Default alpha for penalised
17                                     # least squares fitting
18single_precision = 1.0e-6            # Smallest single precision number
19velocity_protection = 1.0e-6         # Used to compute velocity from momentum
20                                     # See section 7.4 on Flux limiting
21                                     # in the user manual
22                           
23
24################################################################################
25# Standard filenames, directories and system parameters used by ANUGA
26################################################################################
27
28pmesh_filename = '.\\pmesh'
29version_filename = 'stored_version_info.py'
30default_datadir = '.'
31time_format = '%d/%m/%y %H:%M:%S'    # Used with timefile2netcdf
32umask = 002  # Controls file and directory permission created by anuga (UNIX)
33default_boundary_tag = 'exterior' 
34
35# Major revision number for use with create_distribution
36# and update_anuga_user_guide
37major_revision = '1.0beta'
38
39################################################################################
40# Physical constants
41################################################################################
42
43manning = 0.03  # Manning's friction coefficient
44#g = 9.80665    # Gravity - FIXME reinstate this and fix unit tests.
45g = 9.8
46#g(phi) = 9780313 * (1 + 0.0053024 sin(phi)**2 - 0.000 0059 sin(2*phi)**2) micro m/s**2, where phi is the latitude
47#The 'official' average is 9.80665
48
49eta_w = 3.0e-3 # Wind stress coefficient
50rho_a = 1.2e-3 # Atmospheric density
51rho_w = 1023   # Fluid density [kg/m^3] (rho_w = 1023 for salt water)
52
53################################################################################
54# Limiters - used with linear reconstruction of vertex
55# values from centroid values
56################################################################################
57
58# Betas [0;1] control the allowed steepness of gradient for second order
59# extrapolations. Values of 1 allow the steepes gradients while
60# lower values are more conservative. Values of 0 correspond to
61# 1'st order extrapolations.
62#
63
64# There are separate betas for the w, uh, and vh limiters
65# I think these are better SR but they conflict with the unit tests!
66beta_w      = 1.0
67beta_w_dry  = 0.2
68beta_uh     = 1.0
69beta_uh_dry = 0.2
70beta_vh     = 1.0
71beta_vh_dry = 0.2
72
73# Alpha_balance controls how limiters are balanced between deep and shallow.
74# A large value will favour the deep water limiters, allowing the a closer hug
75# to the coastline.  This will minimise 'creep' but at the same time cause
76# smaller time steps
77# Range:
78alpha_balance = 2.0 
79
80# Flag use of new limiters.
81# tight_slope_limiters = 0 means use old limiters (e.g. for some tests)
82# tight_slope_limiters = 1 means use new limiters that hug the bathymetry closer
83tight_slope_limiters = True
84
85# Use centroid velocities to reconstruct momentum at vertices in
86# very shallow water
87# This option has a first order flavour to it, but we still have second order
88# reconstruction of stage and this option only applies in
89# balance_deep_and_shallow when
90# alpha < 1 so in deeper water the full second order scheme is used.
91#
92# This option is good with tight_slope_limiters, especially for large domains.
93use_centroid_velocities = True
94       
95# FIXME (Ole) Maybe get rid of order altogether and use beta_w
96default_order = 2
97
98################################################################################
99# Timestepping
100################################################################################
101
102CFL = 1.0  # CFL condition assigned to domain.CFL - controls timestep size
103     
104# Choose type of timestepping,
105#timestepping_method = 'rk2'   # 2nd Order TVD scheme
106timestepping_method = 'euler' # 1st order euler
107
108# rk2 is a little more stable than euler, so rk2 timestepping
109# can deal with a larger beta when slope limiting the reconstructed
110# solution. The large beta is needed if solving problems sensitive
111# to numerical diffusion, like a small forced wave in an ocean
112beta_euler = 1.0
113beta_rk2   = 1.6
114
115# Option to search for signatures where isolated triangles are
116# responsible for a small global timestep.
117# Treating these by limiting their momenta may help speed up the
118# overall computation.
119# This facility is experimental.
120# protect_against_isolated_degenerate_timesteps = False
121protect_against_isolated_degenerate_timesteps = False
122
123min_timestep = 1.0e-6 # Minimal timestep accepted in ANUGA
124max_timestep = 1.0e+3
125max_smallsteps = 50   # Max number of degenerate steps allowed b4
126                      # trying first order
127
128# Perhaps minimal timestep could be based on the geometry as follows:
129# Define maximal possible speed in open water v_max, e.g. 500m/s (soundspeed?)
130# Then work out minimal internal distance in mesh r_min and set
131# min_timestep = r_min/v_max
132#
133# Max speeds are calculated in the flux function as
134#
135# lambda = v +/- sqrt(gh)
136#
137# so with 500 m/s, h ~ 500^2/g = 2500 m well out of the domain of the
138# shallow water wave equation
139#
140# The actual soundspeed can be as high as 1530m/s
141# (see http://staff.washington.edu/aganse/public.projects/clustering/clustering.html),
142# but that would only happen with h>225000m in this equation. Why ?
143# The maximal speed we specify is really related to the max speed
144# of surface pertubation
145#
146# v_max = 100 #For use in domain_ext.c
147# sound_speed = 500
148
149################################################################################
150# Ranges specific to the shallow water wave equation
151# These control maximal and minimal values of quantities
152################################################################################
153
154# Water depth below which it is considered to be 0 in the model
155minimum_allowed_height = 1.0e-3 
156
157# Water depth below which it is *stored* as 0
158minimum_storable_height = 1.0e-5
159
160# FIXME (Ole): Redefine this parameter to control maximal speeds in general
161# and associate it with protect_against_isolated_degenerate_timesteps = True
162maximum_allowed_speed = 0.0 # Maximal particle speed of water
163#maximum_allowed_speed = 1.0 # Maximal particle speed of water
164                             # Too large (100) creates 'flopping' water
165                             # Too small (0) creates 'creep'
166
167maximum_froude_number = 100.0 # To be used in limiters.
168
169################################################################################
170# Performance parameters used to invoke various optimisations
171################################################################################
172
173use_extensions = True # Use C-extensions
174use_psyco = True      # Use psyco optimisations
175
176optimise_dry_cells = True # Exclude dry and still cells from flux computation
177optimised_gradient_limiter = True # Use hardwired gradient limiter
178use_edge_limiter = False  # The edge limiter is better, but most runs have been
179                          # using vertex limiting. Validations passed with this
180                          # one True 9th May 2008, but many unit tests need
181                          # backward compatibility flag set FIXME(Ole).
182
183points_file_block_line_size = 500 # Number of lines read in from a points file
184                                  # when blocking
185
186################################################################################
187# NetCDF-specific type constants.  Used when defining NetCDF file variables.
188################################################################################
189
190netcdf_char = 'c'
191netcdf_byte = 'b'
192netcdf_int = 'i'
193netcdf_float = 'd'
194netcdf_float64 = 'd'
195netcdf_float32 = 'f'
196
197################################################################################
198# Dynamically-defined constants.
199################################################################################
200
201# Determine if we can read/write large NetCDF files
202netcdf_mode_w = 'w'
203netcdf_mode_a = 'a'
204netcdf_mode_r = 'r'
205
206# Code to set the write mode depending on
207# whether Scientific.IO supports large NetCDF files
208s = """
209import os
210from Scientific.IO.NetCDF import NetCDFFile
211fid = NetCDFFile('tmpfilenamexx', 'wl')
212fid.close()
213os.remove('tmpfilenamexx')
214"""
215
216# Need to run in a separate process due an
217# error with older versions of Scientific.IO
218if sys.platform == 'win32':
219    null = 'NUL'
220else:
221    null = '/dev/null'
222cmd = 'python -c "%s" 2> %s' % (s, null)
223err = os.system(cmd)
224
225if err != 0:
226    # The Python script s failed e.g. with a segfault
227    # which means that large file support is
228    # definitely not supported
229    pass
230else:   
231    # Try the import within this process
232    try:
233        exec(s)
234    except IOError:
235        # NetCDFFile does not segfault but it does not
236        # support large file support   
237        pass
238    else:
239        # Set the default mode to large file support
240        netcdf_mode_w = 'wl'
Note: See TracBrowser for help on using the repository browser.