| |
- builtins.object
-
- Age
- Domain
- Loop
- ModifiedShallowWater2d
- ShallowWater2d
- ShallowWater2dWD
- ShallowWaterTracer2d
class Age(builtins.object) |
|
Create a system so solve the water age or partial age. To solve it, use add_age from class Loop |
|
Methods defined here:
- __init__(self, domain, temporal_scheme, hydro_sol, name='Water', offline=False, initial_time=None, final_time=None)
- keyword arguments:
* domain
object of class Domain
* temporal_scheme
scheme for the temporal integration
* "explicit"
explicit Euler scheme
* "implicit"
implicit order 2 Runge-Kutta
* hydro_sol
hydrodynamics solution
Format:
- equation : ShallowWater2d object (online mode)
- string : path to the folder containing the hydrodynamics folder 'offline_sw2d' (offline mode).
* name
name of the tracer studied (default: "Water")
* offline
resolve shallow water equation (False) or use file to provide hydrodynamics (True) (default: False)
* initial_time
initial time of the simulation (default: 0)
Format:
- float [in seconds]
- date string in format "year-month-day hours:minutes:seconds" (e.g. "2000-06-20 10:05:00")
* final_time
final time of the simulation (default: 0)
Format:
- float [in seconds]
- date string in format "year-month-day hours:minutes:seconds" (e.g. "2000-06-20 22:42:17")
- set_age(self)
- Sets the aging term over the whole domain (standard age)
- set_boundary_coast(self, physical_tag)
- Closed boundary condition: no flux condition at the boundary
keyword argument:
* physical_tag
tag of the part of the boundary on which this boundary condition is applied.
- set_boundary_open(self, physical_tag, inflow=True)
- Open boundary condition: inflow of water (of interest or not)
keyword arguments:
* physical_tag
tag of the part of the boundary on which this boundary condition is applied.
* inflow
boolean stating if the inflowing water is of interest (C=1) or not (C=0) (default: True)
- set_diffusivity(self, mode, okubo_coefficient=0.03, constant_diffusivity=1e-06)
- Diffusivity
keyword arguments:
* mode
type of diffusivity
* "okubo"
Okubo scheme
* "constant"
constant value
* okubo_coefficient
coefficient for Okubo scheme [in m^0.85/s] (default: 0.03 m^0.85/s)
* constant_diffusivity
constant value [in m^2/s] (default: 1e-6 m^2/s)
- set_initial_concentration(self, initial_concentration)
- Initial concentration
keyword argument:
* initial_concentration
concentration initial condition [-] (.msh or .nc format) (default: 0 )
- set_partial_age(self, subdomain, name=None)
- Sets the aging term over a subregion of the domain (partial age). Warning : must be called before setting boundary conditions
keyword arguments:
* subdomain
netcdf or .msh file defining the subregion of interest: 1 inside the subdomain, 0 outside
* name
name of the age concentration (default: 'Age name'_age_concentration)
Data descriptors defined here:
- __dict__
- dictionary for instance variables (if defined)
- __weakref__
- list of weak references to the object (if defined)
|
class Domain(builtins.object) |
|
Create the numerical domain |
|
Methods defined here:
- __init__(self, mesh, bathy, g=9.81, density=1025, solve_on_sphere=False, order=1, periodic_mesh='')
- keyword arguments:
* mesh
path to the mesh file (.msh format). The partitioned mesh will automatically be loaded in multiprocessing
* bathy
bathymetric file [in meters, positive] (.msh, .idx or .nc format)
if bathy is None, the z coordinates of the mesh nodes will be used as bathymetry
* g
map of the mean gravitational acceleration (.msh, .idx or .nc format) (default: 9.81 m/s^2)
* density
density of the liquid (default: 1025 kg/m^3)
* solve_on_sphere
boolean stating if you want to solve the shallow water equation on a sphere or on any 2d manifold instead of solving it on a plane (default: False)
* order
Polynomial order of the discretization
* periodic_mesh
path to the file defining the mapping between the edges of the boundaries to reconnet periodically
Data descriptors defined here:
- __dict__
- dictionary for instance variables (if defined)
- __weakref__
- list of weak references to the object (if defined)
|
class Loop(builtins.object) |
|
Temporal solver |
|
Methods defined here:
- __init__(self, maximum_time_step=3600, export_time=3600, path='./output', evaluateEveryTimeStep=False, export_format='msh')
- keyword arguments:
* maximum_time_step
maximum value for the time step allowed by the user in order to represent the physical processes. [in seconds] (default: 3600 s)
The real time step is computed from numerical restriction of the temporal scheme (if implicit, this value is used as time step).
* export_time
time step for results exportation [in seconds] (default: 3600 s)
* path
path to the output directory [string] (default: ./output)
* evaluateEveryTimeStep
flag to define if quantities (export_value_at_points and compute_mass) are evaluated every time step (True) or only at export time (False) (default: False)
* export_format
format of the output files (Default: 'msh', others: 'vtk' and 'xdmf')
- add_age(self, age, options=None, absolute_tolerance=1e-06, relative_tolerance=0.001)
- Add age system to the temporal solver
keywords arguments:
* age
Age object to add to the temporal solver
* options
options for petSC solver (done by SLIM by default)
* absolute_tolerance
absolute tolerance of the non linear solver for the implicit temporal scheme (default: 1e-6)
* relative_tolerance
relative tolerance of the non linear solver for the implicit temporal scheme (default: 1e-3)
- add_equation(self, equation, options=None, absolute_tolerance=1e-06, relative_tolerance=0.001)
- Add equation to the temporal solver
keywords arguments:
* equation
equation to add to the temporal solver
* options
options for petSC solver (done by SLIM by default)
* absolute_tolerance
absolute tolerance of the non linear solver for the implicit temporal scheme (default: 1e-6)
* relative_tolerance
relative tolerance of the non linear solver for the implicit temporal scheme (default: 1e-3)
- archive_results(self, name, delete=True)
- archive the results in a .tar file
keywords arguments:
* name
name of the archive file
* delete
delete the directory after archiving it successfully (default: True)
- compress_results(self, name, delete=True)
- compress and archive the results in a .tar.gz file
keywords arguments:
* name
name of the final compressed file
* delete
delete the directory after archiving and compressing it successfully (default: True)
- compute_dt(self)
- compute the optimal time step
- evaluatePoints(self)
- export_checkpoint(self)
- Export the solutions at checkpoints without loss of precision enabling an accurate restart
- export_on_structured_grid(self, equation_list, x_min, x_max, y_min, y_max, d_x, d_y, file_name='output.nc')
- Export the solution of equation_name on a structured grid in (longitude, latitude) in netcdf format
keywords arguments:
* equation_list
list of equation(s) whose the solution will be exported on a structured grid
* x_min
minimum x
* x_max
maximum x
* y_min
minimum y
* y_max
maximum y
* d_x
step between x
* d_y
step between y
* file_name
name of the output file (default: output.nc)
- export_solutions(self)
- Export the solutions
- get_time_step(self)
- Get the time step of the equation
- iterate(self, dt=None, allow_failure=False)
- Time integration loop with exports
keywords arguments:
* dt
you can set a time step different from the automatically computed time step
* allow_failure
return False if newton does not converge instead of a Fatal Error
- restart(self, index, time, use_safe_iterate=False, reload_dir=None)
- Restart run from a time-step already resolved
keywords arguments:
* index
index used for restart
* time
time corresponding to this index
* use_safe_iterate
recursively split the iteration into sub time-step if the non-linear solver does not converge (Default: False)
* reload_dir
path to the reload directory, if different from the output directory (Default: None, same as output directory and erase previous files if existing !)
- run(self, use_safe_iterate=False)
- Time integration loop with exports
keywords arguments:
* use_safe_iterate
recursively split the iteration into sub time-step if the non-linear solver does not converge
- safe_iterate(self, dt)
- set_export_options(self, float_size=8, continuous=False, checkpoint_ratio=-1)
- keyword arguments:
* float_size
Precision of the float (Default: 8, double precision; other : 4, simple precision)
* continuous
Boolean to state if the output is continuous or not (Default: False)
* checkpoint_ratio
ratio between output time exports and checkpoint exports in order to restart the simulation (default: -1, no checkpoint)
- set_periodic_offline(self, index_start=1, n_index_per_period=-1)
- Set the periodicity information of the already computed hydrodynamics for the tracer equation computations
keywords arguments:
* index_start
index of the solution exported at each time step starting the period (default: 1)
* n_index_per_period
number of file per period for offline method (default: -1)
Data descriptors defined here:
- __dict__
- dictionary for instance variables (if defined)
- __weakref__
- list of weak references to the object (if defined)
|
class ModifiedShallowWater2d(builtins.object) |
|
Create the modifed shallow water equation with various options |
|
Methods defined here:
- __init__(self, domain, temporal_scheme, linear_equation=False, export_every_sub_time_step=False, initial_time=None, final_time=None)
- keyword arguments:
* domain
object of class Domain
* temporal_scheme
scheme for the temporal integration
* "explicit"
explicit Euler scheme
* "implicit"
implicit order 2 Runge-Kutta
* "multirate"
multirate scheme (see Seny et al. 2012, 2014)
* "steady"
solve the steady equations
Warning: This temporal scheme is not yet implemented for tracers
* linear_equation
boolean stating if you want tos solve the linear version of the shallow water equation (Default: False)
* export_every_sub_time_step
boolean stating if you want to export the hydrodynamics every sub time step for further use of offline equation (Default: False)
* initial_time
initial time of the simulation (specific to each equation) (Default: 0)
Format:
- float [in seconds]
- date string in format "year-month-day hours:minutes:seconds" (e.g. "2000-02-15 10:05:00")
* final_time
final time of the simulation (default: 0)
Format:
- float [in seconds]
- date string in format "year-month-day hours:minutes:seconds" (e.g. "2007-06-25 22:42:17")
- compute_mass(self, output_name)
- compute the total mass of the domain
keyword arguments:
* output_name
name of the output file
- export_value_at_point(self, output_name, x, y, z=None)
- set the coordinates where a time serie of the solution is exported
keyword arguments:
* output_name
name of the output file
* coord_x
vector containing the x coordinates where the solution is evaluated
* coord_y
vector containing the y coordinates where the solution is evaluated
* coord_z
vector containing the z coordinates where the solution is evaluated (Default: None, replaced by vector of 0)
- set_boundary_coast(self, physical_tag, slip=False)
- Boundary condition: normal speed set to zero and tangential speed set to zero or free (no slip or free slip condition)
keyword arguments:
* physical_tag
tag of the boundary on which this condition is applied.
* slip
flag whether applying slip (True) or no slip (False) condition at the boundary (Default: False)
- set_boundary_open(self, physical_tag, p=None, ux=None, uy=None, uz=None, transport_flux=False, ramp_period=0)
- Open boundary condition: one can impose the discharge or the sea surface elevation and/or the velocity (ux and uy have to be defined at the same time but they can be zero).
keyword arguments:
* physical_tag
tag of the part of the boundary on which this boundary condition is applied.
* p
netcdf or .msh file containing the sea surface pressure which will be applied at the boundary [in Pa] (default: None, i.e. the sea surface pressure computed by the model is set to this boundary)
* ux
netcdf or .msh file containing the velocity along the x axis of the local basis or the velocity along the x axis of the global basis if you solve the equation on a sphere. This velocity will be applied at the boundary [in m/s] (default: None, i.e. the velocity computed by the model is set to this boundary)
* uy
netcdf or .msh file containing the velocity along the y axis of the local basis or the velocity along the y axis of the global basis if you solve the equation on a sphere. This velocity will be applied at the boundary [in m/s] (default: None, i.e. the velocity computed by the model is set to this boundary)
* uz
netcdf or .msh file containing the velocity along the z axis of the global basis if you solve the equation on a sphere. This velocity will be applied at the boundary [in m/s] (default: None, i.e. the velocity computed by the model is set to this boundary)
* transport_flux
set whether ux and uy (and uz) are velocities (False) or transport (True) (Default: False)
* ramp_period
period for the linear ramp before applying the nudging [in s] (Default: 0s, no ramp)
- set_coriolis(self, coriolis)
- Add a coriolis term (2*Omega*cos(latitude)) in the shallow water equation
keyword argument:
* coriolis
netcdf or .msh file containing the coriolis term for the whole domain.
- set_dissipation(self, mode, coefficient=None)
- Dissipation due to bottom friction: linear, basic quadratic scheme or Chezy-Manning scheme
keyword arguments:
* mode
scheme used for the bottom friction
* "linear"
linear dissipation tau_b/(rho*H) = gamma*u (default: gamma = 1e-6 s^-1)
* "quadratic"
basic quadratic dissipation tau_b/(rho*H) = (C_d/H)*|u|*u (default: C_d = 2.5e-3)
* "manning"
Chezy-Manning-Strickler formulation tau_b/(rho*H) = (n*n*g/(H^4/3))*|u|*u (default: n = 0.03 sm^(-1/3))
* "slope"
Bottom slope formulation tau_b/(rho*H) = C_s*||u.grad(h)||*u/H (default: C_s = 10)
* coefficient
netcdf or .msh file containing the coefficient value (gamma, C_d, n or C_s depending on the mode) over the whole domain (default: None, the default value of the mode will be used)
- set_forcing(self, forcing_x, forcing_y, forcing_z=None)
- Add a momentum source term derived from a gravitational potential in the shallow water equation.
keyword arguments:
* forcing_x
netcdf or .msh file containing the forcing term along the x-axis in the local basis
(except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain [in ms^-2].
* forcing_y
netcdf or .msh file containing the forcing term along the y-axis in the local basis
(except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain [in ms^-2].
* forcing_z
netcdf or .msh file containing the forcing term along the z-axis in the local basis
(except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain [in ms^-2].
- set_initial_condition(self, p, ux, uy, uz=None)
- Initial conditions
keyword arguments:
* p
surface pressure [in Pa] (.msh or .nc format)
* ux
velocity along the x-axis in the local basis or in the global basis (x,y,z) if you want to solve on a sphere [in m/s] (.idx, .msh or .nc format)
* uy
velocity along the y-axis in the local basis or in the global basis (x,y,z) if you want to solve on a sphere [in m/s] (.idx, .msh or .nc format)
* uz
velocity along the z-axis in the global basis (x,y,z) if you want to solve on a sphere [in m/s] (.idx, .msh or .nc format) (default: None)
- set_relaxation_factor(self, epsilonS, epsilonT)
- Set the relaxation term for the surface pressure in the continuity equation
keyword arguments:
* epsilonS
coefficient to control the spatial relaxation term
* epsilonT
coefficient to control the temporal relaxation term
- set_sea_surface_elevation(self, sse, sse_temporal_variation)
- Add the sea surface elevation temporal variation as a source in the continuity equation and set the sea surface elevation.
keyword arguments:
* sse
netcdf or .msh file containing the sea surface elevation in [m]
* sse_temporal_variation
netcdf or .msh file containing the temporal variation of the sea surface elevation
- set_viscosity(self, mode, smagorinsky_coefficient=0.1, constant_viscosity=1e-06)
- Turbulent viscosity : Smagorinsky scheme or constant value
keyword arguments:
* mode
type of turbulent viscosity
* "smagorinsky"
Smagorinsky scheme (see Smagorinsky 1962)
* "constant"
constant value
* smagorinsky_coefficient
coefficient value for Smagorinsky scheme (default: 0.1)
* constant_viscosity
constant value [in m^2/s] (default: 1e-6 m^2/s)
- set_wind_stress(self, mode, wind_x, wind_y, density_air=-1, wind_z=None)
- Add a wind stress term in the shallow water equation.
keyword arguments:
* mode
type of wind forcing
* "speed"
wind speed given in m s^-1. wind_stress will be computed with Smith-Banke formula.
stress = C_D(speed)*density_air*speed
Setting the air density is mandatory
* "moon"
wind speed given in m s^-1. wind_stress will be computed with Moon et al. parametrization
stress = C_D(speed)*density_air*speed with C_D = k*k*(log(10/z0))**(-2)
Setting the air density is mandatory
* "stress"
wind stress given in kg m^-1 s^-1. (air density will not be used)
* wind_x
netcdf or .msh file containing the wind term along the x-axis in the local basis
(except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain.
* wind_y
netcdf or .msh file containing the wind term along the y-axis in the local basis
(except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain.
* density_air
density of the ambiant air (only necessary for "speed" mode) (default: -1)
* wind_z
netcdf or .msh file containing the wind term along the z-axis in the global basis (x,y,z) for the whole domain (only for equation on a sphere). (default: None)
- update_time(self)
- update the time depending variables which are not given in a netcdf file
Data descriptors defined here:
- __dict__
- dictionary for instance variables (if defined)
- __weakref__
- list of weak references to the object (if defined)
|
class ShallowWater2d(builtins.object) |
|
Create the shallow water equation with various options |
|
Methods defined here:
- __init__(self, domain, temporal_scheme, linear_equation=False, export_every_sub_time_step=False, initial_time=None, final_time=None)
- keyword arguments:
* domain
object of class Domain
* temporal_scheme
scheme for the temporal integration
* "explicit"
explicit Euler scheme
* "implicit"
implicit order 2 Runge-Kutta
* "multirate"
multirate scheme (see Seny et al. 2012, 2014)
* "steady"
solve the steady equations
Warning: This temporal scheme is not yet implemented for tracers
* linear_equation
boolean stating if you want tos solve the linear version of the shallow water equation (Default: False)
* export_every_sub_time_step
boolean stating if you want to export the hydrodynamics every sub time step for further use of offline equation (Default: False)
* initial_time
initial time of the simulation (specific to each equation) (Default: 0)
Format:
- float [in seconds]
- date string in format "year-month-day hours:minutes:seconds" (e.g. "2000-02-15 10:05:00")
* final_time
final time of the simulation (default: 0)
Format:
- float [in seconds]
- date string in format "year-month-day hours:minutes:seconds" (e.g. "2007-06-25 22:42:17")
- compute_mass(self, output_name)
- compute the total mass of the domain
keyword arguments:
* output_name
name of the output file
- export_value_at_point(self, output_name, x, y, z=None)
- set the coordinates where a time serie of the solution is exported
keyword arguments:
* output_name
name of the output file
* coord_x
vector containing the x coordinates where the solution is evaluated
* coord_y
vector containing the y coordinates where the solution is evaluated
* coord_z
vector containing the z coordinates where the solution is evaluated (Default: None, replaced by vector of 0)
- set_boundary_coast(self, physical_tag, slip=False)
- Boundary condition: normal speed set to zero and tangential speed set to zero or free (no slip or free slip condition)
keyword arguments:
* physical_tag
tag of the part of the boundary on which this boundary condition is applied.
* slip
flag wheter applying slip (True) or no slip (False) condition at the boundary (Default: False)
- set_boundary_open(self, physical_tag, discharge=None, sse=None, ux=None, uy=None, uz=None, transport_flux=False, ramp_period=0)
- Open boundary condition: one can impose the discharge or the sea surface elevation and/or the velocity (ux and uy have to be defined at the same time but they can be zero).
keyword arguments:
* physical_tag
tag of the part of the boundary on which this boundary condition is applied.
* discharge
netcdf or .msh file containing the discharge which will be applied at the boundary [in m^3s^-1] (default: None)
* sse
netcdf or .msh file containing the sea surface elevation which will be applied at the boundary [in m] (default: None, i.e. the sea surface elevation computed by the model is set to this boundary)
* ux
netcdf or .msh file containing the velocity along the x axis of the local basis or the velocity along the x axis of the global basis if you solve the equation on a sphere. This velocity will be applied at the boundary [in m/s] (default: None, i.e. the velocity computed by the model is set to this boundary)
* uy
netcdf or .msh file containing the velocity along the y axis of the local basis or the velocity along the y axis of the global basis if you solve the equation on a sphere. This velocity will be applied at the boundary [in m/s] (default: None, i.e. the velocity computed by the model is set to this boundary)
* uz
netcdf or .msh file containing the velocity along the z axis of the global basis if you solve the equation on a sphere. This velocity will be applied at the boundary [in m/s] (default: None, i.e. the velocity computed by the model is set to this boundary)
* transport_flux
set whether ux and uy (and uz) are velocities (False) or transport (True) (Default: False)
* ramp_period
period for the linear ramp before applying the nudging [in s] (Default: 0s, no ramp)
- set_coriolis(self, coriolis)
- Add a coriolis term (2*Omega*cos(latitude)) in the shallow water equation
keyword argument:
* coriolis
netcdf or .msh file containing the coriolis term for the whole domain.
- set_dissipation(self, mode, coefficient=None)
- Dissipation due to bottom friction: linear, basic quadratic scheme or Chezy-Manning scheme
keyword arguments:
* mode
scheme used for the bottom friction
* "linear"
linear dissipation tau_b/(rho*H) = gamma*u (default: gamma = 1e-6 s^-1)
* "quadratic"
basic quadratic dissipation tau_b/(rho*H) = (C_d/H)*|u|*u (default: C_d = 2.5e-3)
* "manning"
Chezy-Manning-Strickler formulation tau_b/(rho*H) = (n*n*g/(H^4/3))*|u|*u (default: n = 0.03 sm^(-1/3))
* "slope"
Bottom slope formulation tau_b/(rho*H) = C_s*||u.grad(h)||*u/H (default: C_s = 10)
* coefficient
netcdf or .msh file containing the coefficient value (gamma, C_d, n or C_s depending on the mode) over the whole domain (default: None, the default value of the mode will be used)
- set_forcing(self, forcing_x, forcing_y, forcing_z=None)
- Add a momentum source term derived from a gravitational potential in the shallow water equation.
keyword arguments:
* forcing_x
netcdf or .msh file containing the forcing term along the x-axis in the local basis
(except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain [in ms^-2].
* forcing_y
netcdf or .msh file containing the forcing term along the y-axis in the local basis
(except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain [in ms^-2].
* forcing_z
netcdf or .msh file containing the forcing term along the z-axis in the local basis
(except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain [in ms^-2].
- set_initial_condition(self, sse, ux, uy, uz=None)
- Initial conditions
keyword arguments:
* sse
sea surface elevation [in m] (.msh or .nc format)
* ux
velocity along the x-axis in the local basis or in the global basis (x,y,z) if you want to solve on a sphere [in m/s] (.idx, .msh or .nc format)
* uy
velocity along the y-axis in the local basis or in the global basis (x,y,z) if you want to solve on a sphere [in m/s] (.idx, .msh or .nc format)
* uz
velocity along the z-axis in the global basis (x,y,z) if you want to solve on a sphere [in m/s] (.idx, .msh or .nc format) (default: None)
- set_nudging(self, coefficient, external_ux, external_uy, external_uz=None, ramp_period=0, transport_flux=False)
- Add a nudging term as a source term in the shallow water equation.
keyword arguments:
* coefficient
netcdf or .msh file containing the nudging coefficient for the whole domain [in s^-1]
* external_ux
netcdf or .msh file containing the external velocity along the x-axis expressed in the local basis (except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain [in m/s]
* external_uy
netcdf or .msh file containing the external velocity along the y-axis expressed in the local basis (except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain [in m/s]
* external_uz
netcdf or .msh file containing the external velocity along the z-axis expressed in the global basis (x,y,z) for the whole domain [in m/s]
* ramp_period
period for the linear ramp before applying the nudging [in s] (Default: 0s, no ramp)
* transport_flux
set whether ux and uy (and uz) are velocities (False) or transport (True) (Default: False)
- set_viscosity(self, mode, smagorinsky_coefficient=0.1, constant_viscosity=1e-06)
- Turbulent viscosity : Smagorinsky scheme or constant value
keyword arguments:
* mode
type of turbulent viscosity
* "smagorinsky"
Smagorinsky scheme (see Smagorinsky 1962)
* "constant"
constant value
* smagorinsky_coefficient
coefficient value for Smagorinsky scheme (default: 0.1)
* constant_viscosity
constant value [in m^2/s] (default: 1e-6 m^2/s)
- set_wind_stress(self, mode, wind_x, wind_y, density_air=-1, wind_z=None)
- Add a wind stress term in the shallow water equation.
keyword arguments:
* mode
type of wind forcing
* "speed"
wind speed given in m s^-1. wind_stress will be computed with Smith-Banke formula.
stress = C_D(speed)*density_air*speed
Setting the air density is mandatory
* "moon"
wind speed given in m s^-1. wind_stress will be computed with Moon et al. parametrization
stress = C_D(speed)*density_air*speed with C_D = k*k*(log(10/z0))**(-2)
Setting the air density is mandatory
* "stress"
wind stress given in kg m^-1 s^-1. (air density will not be used)
* wind_x
netcdf or .msh file containing the wind term along the x-axis in the local basis
(except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain.
* wind_y
netcdf or .msh file containing the wind term along the y-axis in the local basis
(except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain.
* density_air
density of the ambiant air (only necessary for "speed" mode) (default: -1)
* wind_z
netcdf or .msh file containing the wind term along the z-axis in the global basis (x,y,z) for the whole domain (only for equation on a sphere). (default: None)
- update_time(self)
- update the time depending variables which are not given in a netcdf file
Data descriptors defined here:
- __dict__
- dictionary for instance variables (if defined)
- __weakref__
- list of weak references to the object (if defined)
|
class ShallowWater2dWD(builtins.object) |
|
Create the shallow water equation with various options |
|
Methods defined here:
- __init__(self, domain, temporal_scheme, export_every_sub_time_step=False, initial_time=None, final_time=None, hlimDry=0.03)
- keyword arguments:
* domain
object of class Domain
* temporal_scheme
scheme for the temporal integration
* "explicit"
explicit Euler scheme
* "implicit"
implicit Euler scheme
* "steady"
solve the steady equations
Warning: This temporal scheme is not yet implemented for tracers
* export_every_sub_time_step
boolean stating if you want to export the hydrodynamics every sub time step for further use of offline equation (Default: False)
* initial_time
initial time of the simulation (specific to each equation) (Default: 0)
Format:
- float [in seconds]
- date string in format "year-month-day hours:minutes:seconds" (e.g. "2000-02-15 10:05:00")
* final_time
final time of the simulation (default: 0)
Format:
- float [in seconds]
- date string in format "year-month-day hours:minutes:seconds" (e.g. "2007-06-25 22:42:17")
* hlimDry
threshold to define the limit before dry area
- add_dissipation(self, mode, coefficient=None)
- Dissipation due to bottom friction: linear, basic quadratic scheme, Chezy-Manning scheme or slope dependent formulation
keyword arguments:
* mode
scheme used for the bottom friction
* "linear"
linear dissipation tau_b/(rho*H) = gamma*u (default: gamma = 1e-6 s^-1)
* "quadratic"
basic quadratic dissipation tau_b/(rho*H) = (C_d/H)*|u|*u (default: C_d = 2.5e-3)
* "manning"
Chezy-Manning-Strickler formulation tau_b/(rho*H) = (n*n*g/(H^4/3))*|u|*u (default: n = 0.03 sm^(-1/3))
* "slope"
Bottom slope formulation tau_b/(rho*H) = C_s*||u.grad(h)||*u/H (default: C_s = 10)
* coefficient
netcdf or .msh file containing the coefficient value (gamma, C_d, n or C_s depending on the mode) over the whole domain (default: None, the default value of the mode will be used)
- add_viscosity(self, mode, smagorinsky_coefficient=0.1, constant_viscosity=1e-06)
- Turbulent viscosity : Smagorinsky scheme or constant value
keyword arguments:
* mode
type of turbulent viscosity
* "smagorinsky"
Smagorinsky scheme (see Smagorinsky 1962)
* "constant"
constant value
* smagorinsky_coefficient
coefficient value for Smagorinsky scheme (default: 0.1)
* constant_viscosity
constant value [in m^2/s] (default: 1e-6 m^2/s)
- compute_mass(self, output_name)
- compute the total mass of the domain
keyword arguments:
* output_name
name of the output file
- export_value_at_point(self, output_name, x, y, z=None, conservative=True)
- set the coordinates where a time serie of the solution is exported
keyword arguments:
* output_name
name of the output file
* x
vector containing the x coordinates where the solution is evaluated
* y
vector containing the y coordinates where the solution is evaluated
* z
vector containing the z coordinates where the solution is evaluated (Default: None, replaced by vector of 0)
* conservative
Boolean to output conservative variables, H and Hu (True), or non conservative variables, eta and u (False) (Default: True)
- set_atmospheric_pressure(self, atmospheric_pressure)
- Add a momentum source term accounting for the variation in the atmospheric pressure in the shallow water equations.
keyword argument:
* atmospheric_pressure
netcdf or .msh file containing the atmospheric pressure [in Pa]
- set_boundary_coast(self, physical_tag, slip=False)
- Boundary condition: normal speed set to zero and tangential speed set to zero or free (no slip or free slip condition)
keyword arguments:
* physical_tag
tag of the part of the boundary on which this boundary condition is applied.
* slip
flag wheter applying slip (True) or no slip (False) condition at the boundary (Default: False)
- set_boundary_open(self, physical_tag, discharge=None, sse=None, ux=None, uy=None, uz=None, transport_flux=False, ramp_period=0)
- Open boundary condition: one can impose the discharge or the sea surface elevation and/or the velocity (ux and uy have to be defined at the same time but they can be zero).
keyword arguments:
* physical_tag
tag of the part of the boundary on which this boundary condition is applied.
* discharge
netcdf or .msh file containing the discharge which will be applied at the boundary [in m^3s^-1] (default: None)
* sse
netcdf or .msh file containing the sea surface elevation which will be applied at the boundary [in m] (default: None, i.e. the sea surface elevation computed by the model is set to this boundary)
* ux
netcdf or .msh file containing the velocity along the x axis of the local basis or the velocity along the x axis of the global basis if you solve the equation on a sphere. This velocity will be applied at the boundary [in m/s] (default: None, i.e. the velocity computed by the model is set to this boundary)
* uy
netcdf or .msh file containing the velocity along the y axis of the local basis or the velocity along the y axis of the global basis if you solve the equation on a sphere. This velocity will be applied at the boundary [in m/s] (default: None, i.e. the velocity computed by the model is set to this boundary)
* uz
netcdf or .msh file containing the velocity along the z axis of the global basis if you solve the equation on a sphere. This velocity will be applied at the boundary [in m/s] (default: None, i.e. the velocity computed by the model is set to this boundary)
* transport_flux
set whether ux and uy (and uz) are velocities (False) or transport (True) (Default: False)
* ramp_period
period for the linear ramp before applying the nudging [in s] (Default: 0s, no ramp)
- set_coriolis(self, coriolis)
- Add a coriolis term (2*Omega*cos(latitude)) in the shallow water equation
keyword argument:
* coriolis
netcdf or .msh file containing the coriolis term for the whole domain.
- set_forcing(self, forcing_x, forcing_y, forcing_z=None)
- Add a momentum source term derived from a gravitational potential in the shallow water equation.
keyword arguments:
* forcing_x
netcdf or .msh file containing the forcing term along the x-axis in the local basis
(except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain [in ms^-2].
* forcing_y
netcdf or .msh file containing the forcing term along the y-axis in the local basis
(except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain [in ms^-2].
* forcing_z
netcdf or .msh file containing the forcing term along the z-axis in the local basis
(except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain [in ms^-2].
- set_initial_condition(self, sse, ux, uy, uz=None)
- Initial conditions
keyword arguments:
* sse
sea surface elevation [in m] (.msh or .nc format)
* ux
velocity along the x-axis in the local basis or in the global basis (x,y,z) if you want to solve on a sphere [in m/s] (.idx, .msh or .nc format)
* uy
velocity along the y-axis in the local basis or in the global basis (x,y,z) if you want to solve on a sphere [in m/s] (.idx, .msh or .nc format)
* uz
velocity along the z-axis in the global basis (x,y,z) if you want to solve on a sphere [in m/s] (.idx, .msh or .nc format) (default: None)
- set_nudging(self, coefficient, external_ux, external_uy, ramp_period=0, transport_flux=False)
- Add a nudging term as a source term in the conservative shallow water equation.
keyword arguments:
* coefficient
netcdf or .msh file containing the nudging coefficient for the whole domain [in s^-1]
* external_ux
netcdf or .msh file containing the external velocity along the x-axis [in m/s]
* external_uy
netcdf or .msh file containing the external velocity along the y-axis [in m/s]
* ramp_period
period for the linear ramp before applying the nudging [in s] (Default: 0s, no ramp)
* transport_flux
set whether ux and uy (and uz) are velocities (False) or transport (True) (Default: False)
- set_wind_stress(self, mode, wind_x, wind_y, density_air=-1, wind_z=None)
- Add a wind stress term in the shallow water equation.
keyword arguments:
* mode
type of wind forcing
* "speed"
wind speed given in m s^-1. wind_stress will be computed with Smith-Banke formula.
stress = C_D(speed)*density_air*speed
Setting the air density is mandatory
* "moon"
wind speed given in m s^-1. wind_stress will be computed with Moon et al. parametrization
stress = C_D(speed)*density_air*speed with C_D = k*k*(log(10/z0))**(-2)
Setting the air density is mandatory
* "stress"
wind stress given in kg m^-1 s^-1. (air density will not be used)
* wind_x
netcdf or .msh file containing the wind term along the x-axis in the local basis
(except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain.
* wind_y
netcdf or .msh file containing the wind term along the y-axis in the local basis
(except for equation on a sphere for which it has to be expressed in the global basis (x,y,z)) for the whole domain.
* density_air
density of the ambiant air (only necessary for "speed" mode) (default: -1)
* wind_z
netcdf or .msh file containing the wind term along the z-axis in the global basis (x,y,z) for the whole domain (only for equation on a sphere). (default: None)
- update_time(self)
- update the time depending variables which are not given in a netcdf file
Data descriptors defined here:
- __dict__
- dictionary for instance variables (if defined)
- __weakref__
- list of weak references to the object (if defined)
|
class ShallowWaterTracer2d(builtins.object) |
|
Create the shallow water equation for a tracer |
|
Methods defined here:
- __init__(self, domain, temporal_scheme, hydro_sol, linear_equation=False, name='Tracer', offline=False, index_start=0, initial_time=None, final_time=None)
- keyword arguments:
* domain
object of class Domain
* temporal_scheme
scheme for the temporal integration
* "explicit"
explicit Euler scheme
* "implicit"
implicit order 2 Runge-Kutta
* "steady"
solve the steady equations
* hydro_sol
hydrodynamics solution
Format:
- equation : ShallowWater2d object (online mode)
- string : path to the folder containing the hydrodynamics folder 'offline_sw2d' (offline mode).
* linear_equation
boolean stating linear shallow water tracer equation (default: False)
* name
name of the Tracer (default: "Tracer")
* offline
resolve shallow water equation (False) or use file to provide hydrodynamics (True) (default: False)
* index_start
index of the first file for the offline mode (default: 0)
* initial_time
initial time of the simulation (default: 0)
Format:
- float [in seconds]
- date string in format "year-month-day hours:minutes:seconds" (e.g. "2000-06-20 10:05:00")
* final_time
final time of the simulation (default: 0)
Format:
- float [in seconds]
- date string in format "year-month-day hours:minutes:seconds" (e.g. "2000-06-20 22:42:17")
- compute_mass(self, output_name)
- Compute the integral of the tracer over the whole domain
keyword arguments:
* output_name
name of the output file
- set_boundary_coast(self, physical_tag)
- Boundary condition: no flux condition at the boundary
keyword argument:
* physical_tag
tag of the part of the boundary on which this boundary condition is applied.
- set_boundary_flux(self, physical_tag, flux)
- Boundary condition for a tracer flux: one can impose the flux throw the boundary condition, in [c s^-1], (with the tracer units [c]).
keyword arguments:
* physical_tag
tag of the part of the boundary on which this boundary condition is applied.
* flux
netcdf or .msh file containing the flux which will be applied at the boundary in [c s^-1]
- set_boundary_open(self, physical_tag, concentration)
- Open boundary condition: one can impose the external concentration
keyword arguments:
* physical_tag
tag of the part of the boundary on which this boundary condition is applied.
* concentration
netcdf or .msh file containing the concentration which will be applied at the boundary [in c]
- set_diffusivity(self, mode, okubo_coefficient=0.03, constant_diffusivity=1e-06)
- Diffusivity
keyword arguments:
* mode
type of diffusivity
* "okubo"
Okubo scheme
* "constant"
constant value
* okubo_coefficient
coefficient for Okubo scheme [in m^0.85/s] (default: 0.03 m^0.85/s)
* constant_diffusivity
constant value [in m^2/s] (default: 1e-6 m^2/s)
- set_initial_condition(self, initial_condition)
- Initial conditions
keyword argument:
* tracer
tracer initial condition [-] (.msh or .nc format) (default: 0 )
- set_sediment(self, initial_bottom_concentration, windU, windV, factorW=0.01, wMax=0.0001, u0=0.3, v0=0.3, w0=10, a01=0.028, a02=0.144, a1=1e-06, n=4.0, omega2=2.4538, eros0Fact=0.245, hsig=-1)
- Define wind velocity and all parameters for the sediment module.
keyword arguments:
* initial_bottom_concentration
netcdf or .msh file containing the initial sediment concentration at the bottom for the whole domain [in c/m^2].
* windU
netcdf or .msh file containing the surface wind velocity along the x-axis in the local basis [in ms^-1].
* windV
netcdf or .msh file containing the surface wind velocity along the y-axis in the local basis [in ms^-1].
* factorW
dependence factor between the settling viscosity and the suspended sediment concentration (default: 0.01 m^4/kg/s)
* wMax
maximum settling viscosity (default: 1e-4 m/s)
* u0
threshold velocity for erosion (default: 0.3 m/s, see Lambrechts et al., 2010, eq 1 and 6)
* v0
threshold velocity for settling (default: 0.3 m/s, see Lambrechts et al., 2010, eq 1 and 6)
* w0
wind speed threshold (default: 10 m/s)
* a01
empirical constant that depends only on the characteristics of the fine sediment on the seafloor (default: 2.8e-2).
It is used to compute A1 in eq 8 of Lambrechts et al., 2010: A1= a01 + (a02-a01) * (0.5 - atan(10*alpha)/(2*atan(10))) where alpha is the wind orientation factor
* a02
empirical constant that depends only on the characteristics of the fine sediment on the seafloor (default: 1.44e-1)
It is used to compute A1 in eq 8 of Lambrechts et al., 2010: A1= a01 + (a02-a01) * (0.5 - atan(10*alpha)/(2*atan(10))) where alpha is the wind orientation factor
* a1
empirical constant that depends only on the characteristics of the fine sediment on the seafloor (default: 1e-6, A2 in eq 9 of Lambrechts et al., 2010)
* n
constant in erosion flux parametrisation (default: 4, see Lambrechts et al., 2010, eq 1)
* omega2
square of the wave frequency (default: 2.4538 Hz, see Lambrechts et al., 2010, eq 8-10)
* eros0Fact
erosion factor multiplying F (defualt: 0.245, see Lambrechts et al., 2010, eq 10)
* hsig
significant wave height (default: -1 and significant wave height computed from wind velocity)
- set_source(self, tracer_source)
- Add a tracer source term in the tracer equation.
keyword argument:
* tracer_source
netcdf or .msh file containing the source term for the whole domain [in c*s^-1].
- update_time(self)
- update the time depending variables which are not given in a netcdf file
Data descriptors defined here:
- __dict__
- dictionary for instance variables (if defined)
- __weakref__
- list of weak references to the object (if defined)
| |