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

Last change on this file since 7340 was 7302, checked in by ole, 16 years ago

Get ready for Numpy release v 1.1beta_7302

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