source: anuga_core/source/anuga/config.py @ 7073

Last change on this file since 7073 was 7055, checked in by ole, 16 years ago

Added a few comments to config.py

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