source: trunk/anuga_core/source/anuga/config.py @ 8418

Last change on this file since 8418 was 8418, checked in by steve, 12 years ago

Strange proble with parallel_safe

File size: 10.1 KB
Line 
1"""Module where global ANUGA model parameters and default values are set
2"""
3
4import os
5import sys
6
7
8################################################################################
9# Numerical constants
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
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                           
22
23################################################################################
24# Standard filenames, directories and system parameters used by ANUGA
25################################################################################
26
27pmesh_filename = '.\\pmesh'
28version_filename = 'stored_version_info.py'
29default_datadir = '.'
30time_format = '%d/%m/%y %H:%M:%S'    # Used with timefile2netcdf
31umask = 002  # Controls file and directory permission created by anuga (UNIX)
32default_boundary_tag = 'exterior' 
33
34# Major revision number for use with create_distribution
35# and update_anuga_user_guide
36major_revision = '1.3.0_beta'
37
38################################################################################
39# Physical constants
40################################################################################
41
42manning = 0.03  # Manning's friction coefficient
43#g = 9.80665    # Gravity - FIXME reinstate this and fix unit tests.
44g = 9.8
45#g(phi) = 9780313 * (1 + 0.0053024 sin(phi)**2 - 0.000 0059 sin(2*phi)**2)
46# 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# Note the individual beta values are set in domain.set_flow_method which also sets
58# the timestepping method
59
60beta_w = 1.0
61
62# Alpha_balance controls how limiters are balanced between deep and shallow.
63# A large value will favour the deep water limiters, allowing the a closer hug
64# to the coastline.  This will minimise 'creep' but at the same time cause
65# smaller time steps
66# Range:
67alpha_balance = 2.0 
68
69# Flag use of new limiters.
70# tight_slope_limiters = 0 means use old limiters (e.g. for some tests)
71# tight_slope_limiters = 1 means use new limiters that hug the bathymetry closer
72tight_slope_limiters = True
73
74use_edge_limiter = False  # The edge limiter is better, but most runs have been
75                          # using vertex limiting. Validations passed with this
76                          # one True 9th May 2008, but many unit tests need
77                          # backward compatibility flag set FIXME(Ole).
78
79# Use centroid velocities to reconstruct momentum at vertices in
80# very shallow water
81# This option has a first order flavour to it, but we still have second order
82# reconstruction of stage and this option only applies in
83# balance_deep_and_shallow when
84# alpha < 1 so in deeper water the full second order scheme is used.
85#
86# This option is good with tight_slope_limiters, especially for large domains.
87use_centroid_velocities = True
88       
89# FIXME (Ole) Maybe get rid of order altogether and use beta_w
90default_order = 2
91
92# Option to use velocity extrapolation instead of momentum extrapolation in the
93# routine domain.extrapolate_second_order_sw
94extrapolate_velocity_second_order=True
95
96# Option to setup compute_fluxes_method
97# Currently "original' and 'wb_1' to 'wb_3'
98compute_fluxes_method = 'wb_2'
99
100################################################################################
101# Friction Method
102################################################################################
103
104sloped_mannings_function = False
105
106################################################################################
107# Timestepping
108################################################################################
109
110CFL = 1.0  # CFL condition assigned to domain.CFL - controls timestep size
111     
112# Choose type of timestepping and spatial reconstruction method
113
114timestepping_method = 1
115
116# For shallow water we have a method that sets both timestepping and spatial reconstruction and
117# beta values. In this case the settings for timestepping_method will be overriden
118
119#flow_algorithm = '1_0'    # 1st order euler and conservative piecewise constant spatial reconstruction
120flow_algorithm = '1_5'  # 1st order euler and conservative piecewise linear spatial reconstruction
121#flow_algorithm = '1_75' # 1st order euler and more aggressive piecewise linear spatial reconstruction
122#flow_algorithm = '2_0'    # 2nd order TVD scheme and more aggressive piecewise linear spatial reconstruction
123#flow_algorithm = '2.5'  # 3rd order TVD scheme and more aggressive piecewise linear spatial reconstruction
124
125
126
127
128# rk2 is a little more stable than euler, so rk2 timestepping
129# can deal with a larger beta when slope limiting the reconstructed
130# solution. The large beta is needed if solving problems sensitive
131# to numerical diffusion, like a small forced wave in an ocean
132beta_euler = 1.0
133beta_rk2   = 1.6
134
135# Option to search for signatures where isolated triangles are
136# responsible for a small global timestep.
137# Treating these by limiting their momenta may help speed up the
138# overall computation.
139# This facility is experimental.
140# protect_against_isolated_degenerate_timesteps = False
141protect_against_isolated_degenerate_timesteps = False
142
143min_timestep = 1.0e-6 # Minimal timestep accepted in ANUGA
144max_timestep = 1.0e+3
145max_smallsteps = 50   # Max number of degenerate steps allowed b4
146                      # trying first order
147
148# Perhaps minimal timestep could be based on the geometry as follows:
149# Define maximal possible speed in open water v_max, e.g. 500m/s (soundspeed?)
150# Then work out minimal internal distance in mesh r_min and set
151# min_timestep = r_min/v_max
152#
153# Max speeds are calculated in the flux function as
154#
155# lambda = v +/- sqrt(gh)
156#
157# so with 500 m/s, h ~ 500^2/g = 2500 m well out of the domain of the
158# shallow water wave equation
159#
160# The actual soundspeed can be as high as 1530m/s
161# (see http://staff.washington.edu/aganse/public.projects/
162#            clustering/clustering.html),
163# but that would only happen with h>225000m in this equation. Why ?
164# The maximal speed we specify is really related to the max speed
165# of surface pertubation
166#
167# v_max = 100 #For use in domain_ext.c
168# sound_speed = 500
169
170################################################################################
171# Ranges specific to the shallow water wave equation
172# These control maximal and minimal values of quantities
173################################################################################
174
175# Water depth below which it is considered to be 0 in the model
176minimum_allowed_height = 1.0e-3 
177
178# Water depth below which it is *stored* as 0
179minimum_storable_height = 1.0e-05
180
181# FIXME (Ole): Redefine this parameter to control maximal speeds in general
182# and associate it with protect_against_isolated_degenerate_timesteps = True
183maximum_allowed_speed = 0.0 # Maximal particle speed of water
184#maximum_allowed_speed = 1.0 # Maximal particle speed of water
185                             # Too large (100) creates 'flopping' water
186                             # Too small (0) creates 'creep'
187
188maximum_froude_number = 100.0 # To be used in limiters.
189
190################################################################################
191# Performance parameters used to invoke various optimisations
192################################################################################
193
194use_psyco = False      # Use psyco optimisations
195
196optimise_dry_cells = True # Exclude dry and still cells from flux computation
197optimised_gradient_limiter = True # Use hardwired gradient limiter
198
199points_file_block_line_size = 500 # Number of lines read in from a points file
200                                  # when blocking
201
202################################################################################
203# NetCDF-specific type constants.  Used when defining NetCDF file variables.
204################################################################################
205
206netcdf_char = 'c'
207netcdf_byte = 'b'
208netcdf_int = 'i'
209netcdf_float = 'd'
210netcdf_float64 = 'd'
211netcdf_float32 = 'f'
212
213################################################################################
214# Dynamically-defined constants.
215################################################################################
216
217# Determine if we can read/write large NetCDF files
218netcdf_mode_w = 'w'
219netcdf_mode_a = 'a'
220netcdf_mode_r = 'r'
221
222# Code to set the write mode depending on
223# whether Scientific.IO supports large NetCDF files
224s = """
225import os, tempfile
226from Scientific.IO.NetCDF import NetCDFFile
227
228filename = tempfile.mktemp('.nc')
229
230fid = NetCDFFile(filename, 'wl')
231fid.close()
232os.remove(filename)
233"""
234
235# Need to run in a separate process due an
236# error with older versions of Scientific.IO
237if sys.platform == 'win32':
238    null = 'NUL'
239else:
240    null = '/dev/null'
241cmd = 'python -c "%s" 2> %s' % (s, null)
242err = os.system(cmd)
243
244if err != 0:
245    # The Python script s failed e.g. with a segfault
246    # which means that large file support is
247    # definitely not supported
248    pass
249else:   
250    # Try the import within this process
251    try:
252        exec(s)
253    except IOError:
254        # NetCDFFile does not segfault but it does not
255        # support large file support   
256        pass
257    else:
258        # Set the default mode to large file support
259        #netcdf_mode_w = 'w4' # Future use of HDF5       
260        netcdf_mode_w = 'wl' # Large NetCDF (new default 30/6/2009)
261        #netcdf_mode_w = 'w'   # Old style NetCDF used by OSG viewer
262
263
Note: See TracBrowser for help on using the repository browser.