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

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

Some changes to netcdf config

File size: 10.2 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' and 'tsunami'
98compute_fluxes_method = 'wb_2'
99
100# Option to setup distribute_to_vertices_and_edges_method
101# Currently "original' and 'tsunami'
102distribute_to_vertices_and_edges_method = 'original'
103
104################################################################################
105# Friction Method
106################################################################################
107
108sloped_mannings_function = False
109
110################################################################################
111# Timestepping
112################################################################################
113
114CFL = 1.0  # CFL condition assigned to domain.CFL - controls timestep size
115     
116# Choose type of timestepping and spatial reconstruction method
117
118timestepping_method = 1
119
120# For shallow water we have a method that sets both timestepping and spatial reconstruction and
121# beta values. In this case the settings for timestepping_method will be overriden
122
123#flow_algorithm = '1_0'    # 1st order euler and conservative piecewise constant spatial reconstruction
124flow_algorithm = '1_5'  # 1st order euler and conservative piecewise linear spatial reconstruction
125#flow_algorithm = '1_75' # 1st order euler and more aggressive piecewise linear spatial reconstruction
126#flow_algorithm = '2_0'    # 2nd order TVD scheme and more aggressive piecewise linear spatial reconstruction
127#flow_algorithm = '2.5'  # 3rd order TVD scheme and more aggressive piecewise linear spatial reconstruction
128
129
130
131
132# rk2 is a little more stable than euler, so rk2 timestepping
133# can deal with a larger beta when slope limiting the reconstructed
134# solution. The large beta is needed if solving problems sensitive
135# to numerical diffusion, like a small forced wave in an ocean
136beta_euler = 1.0
137beta_rk2   = 1.6
138
139# Option to search for signatures where isolated triangles are
140# responsible for a small global timestep.
141# Treating these by limiting their momenta may help speed up the
142# overall computation.
143# This facility is experimental.
144# protect_against_isolated_degenerate_timesteps = False
145protect_against_isolated_degenerate_timesteps = False
146
147min_timestep = 1.0e-6 # Minimal timestep accepted in ANUGA
148max_timestep = 1.0e+3
149max_smallsteps = 50   # Max number of degenerate steps allowed b4
150                      # trying first order
151
152# Perhaps minimal timestep could be based on the geometry as follows:
153# Define maximal possible speed in open water v_max, e.g. 500m/s (soundspeed?)
154# Then work out minimal internal distance in mesh r_min and set
155# min_timestep = r_min/v_max
156#
157# Max speeds are calculated in the flux function as
158#
159# lambda = v +/- sqrt(gh)
160#
161# so with 500 m/s, h ~ 500^2/g = 2500 m well out of the domain of the
162# shallow water wave equation
163#
164# The actual soundspeed can be as high as 1530m/s
165# (see http://staff.washington.edu/aganse/public.projects/
166#            clustering/clustering.html),
167# but that would only happen with h>225000m in this equation. Why ?
168# The maximal speed we specify is really related to the max speed
169# of surface pertubation
170#
171# v_max = 100 #For use in domain_ext.c
172# sound_speed = 500
173
174################################################################################
175# Ranges specific to the shallow water wave equation
176# These control maximal and minimal values of quantities
177################################################################################
178
179# Water depth below which it is considered to be 0 in the model
180minimum_allowed_height = 1.0e-03 
181
182# Water depth below which it is *stored* as 0
183minimum_storable_height = 1.0e-05
184
185# FIXME (Ole): Redefine this parameter to control maximal speeds in general
186# and associate it with protect_against_isolated_degenerate_timesteps = True
187maximum_allowed_speed = 0.0 # Maximal particle speed of water
188#maximum_allowed_speed = 1.0 # Maximal particle speed of water
189                             # Too large (100) creates 'flopping' water
190                             # Too small (0) creates 'creep'
191
192maximum_froude_number = 100.0 # To be used in limiters.
193
194################################################################################
195# Performance parameters used to invoke various optimisations
196################################################################################
197
198use_psyco = False      # Use psyco optimisations
199
200optimise_dry_cells = True # Exclude dry and still cells from flux computation
201optimised_gradient_limiter = True # Use hardwired gradient limiter
202
203points_file_block_line_size = 5e7 # Number of lines read in from a points file
204                                  # when blocking
205
206################################################################################
207# NetCDF-specific type constants.  Used when defining NetCDF file variables.
208################################################################################
209
210netcdf_char = 'c'
211netcdf_byte = 'b'
212netcdf_int = 'i'
213netcdf_float = 'd'
214netcdf_float64 = 'd'
215netcdf_float32 = 'f'
216
217################################################################################
218# Dynamically-defined constants.
219################################################################################
220
221# Determine if we can read/write large NetCDF files
222netcdf_mode_w = 'wl'
223netcdf_mode_a = 'r+'
224netcdf_mode_r = 'r'
225
226# Code to set the write mode depending on
227# whether Scientific.IO supports large NetCDF files
228s = """
229import os, tempfile
230from anuga.file.netcdf import NetCDFFile
231
232filename = tempfile.mktemp('.nc')
233
234fid = NetCDFFile(filename, 'wl')
235fid.close()
236os.remove(filename)
237"""
238
239"""
240# Need to run in a separate process due an
241# error with older versions of Scientific.IO
242if sys.platform == 'win32':
243    null = 'NUL'
244else:
245    null = '/dev/null'
246cmd = 'python -c "%s" 2> %s' % (s, null)
247err = os.system(cmd)
248
249if err != 0:
250    # The Python script s failed e.g. with a segfault
251    # which means that large file support is
252    # definitely not supported
253    pass
254else:   
255    # Try the import within this process
256    try:
257        exec(s)
258    except IOError:
259        # NetCDFFile does not segfault but it does not
260        # support large file support   
261        pass
262    else:
263        # Set the default mode to large file support
264        #netcdf_mode_w = 'w4' # Future use of HDF5       
265        netcdf_mode_w = 'wl' # Large NetCDF (new default 30/6/2009)
266        #netcdf_mode_w = 'w'   # Old style NetCDF used by OSG viewer
267
268"""
Note: See TracBrowser for help on using the repository browser.