wolfgpu.simple_simulation

Author: HECE - University of Liege, Stéphane Champailler, Pierre Archambeau Date: 2024

Copyright (c) 2024 University of Liege. All rights reserved.

This script and its content are protected by copyright law. Unauthorized copying or distribution of this file, via any medium, is strictly prohibited.

Module Contents

class wolfgpu.simple_simulation.BoundaryConditionIndices[source]

Bases: enum.Enum

The column indices of each supported boundary conditions in the boundary conditions table. This also covers cell parameters.

IMPORTANT: These values are tied to data representation in the shader. Don’t change them!

These are separate from the BoundaryConditionsTypes because they don’t serve the same purpose. These are for data representation in wolfgpu whereas BoundaryConditionsTypes are for CPU Wolf compatibility.

BoundaryConditionsTypes and CellParameterType must all be converted to this enum when loading the sim in the gpu.

BC_TABLE_INDEX_FOR_H_ON_LEFT = 1[source]
BC_TABLE_INDEX_FOR_HMOD_ON_LEFT = 2[source]
BC_TABLE_INDEX_FOR_QX_ON_LEFT = 0[source]
BC_TABLE_INDEX_FOR_QY_ON_LEFT = 3[source]
BC_TABLE_INDEX_FOR_FROUDE_NORMAL_ON_LEFT = 8[source]
BC_TABLE_INDEX_FOR_H_ON_BOTTOM = 4[source]
BC_TABLE_INDEX_FOR_HMOD_ON_BOTTOM = 5[source]
BC_TABLE_INDEX_FOR_QX_ON_BOTTOM = 6[source]
BC_TABLE_INDEX_FOR_QY_ON_BOTTOM = 7[source]
BC_TABLE_INDEX_FOR_FROUDE_NORMAL_ON_BOTTOM = 9[source]
BC_TABLE_INDEX_FOR_BRIDGE_DECK_ELEVATION = 10[source]
BC_TABLE_INDEX_FOR_BRIDGE_ROOF_ELEVATION = 11[source]
wolfgpu.simple_simulation.BC_BLANK_VALUE[source]
wolfgpu.simple_simulation.NB_BC_TYPES[source]
class wolfgpu.simple_simulation.CellParameterType[source]

Bases: enum.Enum

Represents the type of a cell parameter.

BRIDGE_DECK_ELEVATION = 1[source]
BRIDGE_ROOF_ELEVATION = 2[source]
class wolfgpu.simple_simulation.BoundaryConditionsTypes[source]

Bases: enum.Enum

Inheritance diagram of wolfgpu.simple_simulation.BoundaryConditionsTypes

The boundary conditions types.

Note: The values match the numbers in Wolf’s simulations parameters.

H = 1[source]

Weak bondary condition over \(h+b\), the water free surface height (that is, height of water plus bathymetry level).

QX = 2[source]

Weak bondary condition over \(q_x\), the discharge over the x-axis.

QY = 3[source]

Weak bondary condition over \(q_y\), the discharge over the x-axis.

NONE = 4[source]

This is used to represent the absence of boundary condition.

QBX = 5[source]

Not used.

QBY = 6[source]

Not used.

HMOD = 7[source]

Weak bondary condition over \(h\), the water free surface height (that is, height of water, measured from the top of the bathymetry).

FROUDE_NORMAL = 8[source]

Weak boundary condition for Froude limited height. Water height becomes \(\min(h, \frac{u^2}{g \times \text{Fr}^2})\) where \(\text{Fr}\) is given.

class wolfgpu.simple_simulation.Direction[source]

Bases: enum.Enum

Inheritance diagram of wolfgpu.simple_simulation.Direction

Generic enumeration.

Derive from this class to define new enumerations.

LEFT = 1[source]
BOTTOM = 2[source]
X = 1[source]
Y = 2[source]
class wolfgpu.simple_simulation.boundary_condition_2D(i: int, j: int, ntype: BoundaryConditionsTypes, val: float, direction: Direction)[source]

Type des CL générales, faibles et fortes @author Pierre Archambeau

to_dict()[source]
classmethod from_dict(js)[source]
class wolfgpu.simple_simulation.CellParameter(i: int, j: int, ntype: CellParameterType, val: float)[source]

Cell parameter

to_dict()[source]
classmethod from_dict(js)[source]
class wolfgpu.simple_simulation.SimulationDurationType[source]

Bases: enum.Enum

Inheritance diagram of wolfgpu.simple_simulation.SimulationDurationType

A type of duration.

This was introduced to abstract the duration of a simulation.

SECONDS = 0[source]

A duration expressed in seconds.

STEPS = 1[source]

A duration expressed in steps.

class wolfgpu.simple_simulation.ReportFrequencyType[source]

Bases: enum.Enum

Inheritance diagram of wolfgpu.simple_simulation.ReportFrequencyType

When specifying the report period, one can use two different units.

ITERATION = 1[source]

Report period in iterations

SIMULATION_TIME = 0[source]

Report period in simulation time (seconds)

class wolfgpu.simple_simulation.SimulationDuration(_type: SimulationDurationType, duration: float | int)[source]

In the GPU simulator, we represent durations about the simulation either as number of simulation steps either as simulation time (in seconds). This class can represent both.

to_dict()[source]
from_dict(d)[source]
classmethod from_seconds(seconds) SimulationDuration[source]
classmethod from_steps(steps) SimulationDuration[source]
classmethod from_str(d: str) SimulationDuration[source]
classmethod zero_like(d) SimulationDuration[source]
is_zero()[source]
class wolfgpu.simple_simulation.TimeStepStrategy[source]

Bases: enum.Enum

Inheritance diagram of wolfgpu.simple_simulation.TimeStepStrategy

Describes the way the time step is computed on each simulation step.

FIXED_TIME_STEP = 1[source]

The time step is constant over all simulation steps.

OPTIMIZED_TIME_STEP = 2[source]

The time step is the maximum of the minimum time steps authorized by the Courant formula.

class wolfgpu.simple_simulation.InfiltrationInterpolation[source]

Bases: enum.Enum

Inheritance diagram of wolfgpu.simple_simulation.InfiltrationInterpolation

Describes the way the infiltration are computed across an infiltration time interval.

NONE = 0[source]

The infiltration is constant over the infiltration time interval.

LINEAR = 1[source]

The infiltration is interpolated over the infiltration time interval. The initial value is the one associated with interval, the final value is the one associated with the next interval.

class wolfgpu.simple_simulation.InfiltrationChronology[source]

This is optimized for access where inserts and reads are not frequently mixed togehter.

property nb_zones: int[source]
property nb_entries: int[source]
_sorted_timeline()[source]
classmethod from_array(npa)[source]
classmethod from_list(chronology: list[float, list[float]])[source]
to_array()[source]
get_row_at_time(time: float) Tuple[float, List[float]][source]

Get infiltration chronology row active at time time.

Parameters:

time – a positive number of seconds.

get_active_entry_start_time(time: float) float[source]
set(time: float, quantities: list[float])[source]

Sets or replaces an entry in the chronology.

Parameters:
  • time – the moment at which the water must be injected. It is expressed in seconds. If there’s already something at that moment, it will be discarded.

  • quantities – A list giving the quantities of water to inject at the given time. Each quantity in the list corresponds to a infiltration zone of the infiltration map.

_get_closest_row_time_by_time(time: float)[source]

Can return None if no row is active. This happens when time is before all the chronology entries.

class wolfgpu.simple_simulation.SimpleSimulation(nx: int = 0, ny: int = 0)[source]

A simple simulation scenario.

IMPORTANT: althoug we use numpy arrays to represent various matrices, be careful: we use shapes as (nb of columns, nb of rows). This is transposed comprared to regular numpy !

property param_dx[source]

dx is the x-dimension of a mesh (expressed in meters).

property param_dy[source]

dy is the y-dimension of a mesh (expressed in meters)

property h: numpy.typing.NDArray[numpy.float32][source]

An array representing the water height (expressed in meters above the bathymetry).

property qx: numpy.typing.NDArray[numpy.float32][source]

An array representing the water discharge along x-axis (expressed in m/s²).

property qy: numpy.typing.NDArray[numpy.float32][source]

An array representing the water discharge along y-axis (expressed in m/s²).

property manning: numpy.typing.NDArray[numpy.float32][source]

An array representing the Manning coefficient for each mesh.

property infiltration_zones: numpy.typing.NDArray[numpy.int32][source]

An array representing the infiltration zones for each mesh.

property nap: numpy.typing.NDArray[numpy.uint8][source]

An array representing the meshes that are part of the computation domain.

property bathymetry: numpy.typing.NDArray[numpy.float32][source]

An array representing the bathymetry each mesh (expressed in meters).

property param_nx: int[source]

Dimension of the computation domain on the x axis.

property param_ny: int[source]

Dimension of the computation domain on the y axis.

property param_froude_bc_limit_tolerance[source]

When Froude BC applies and needs to push up the water height, it will push it up to current_water_height * tolerance.

property shape: tuple[int, int][source]

The shape of the various arrays representing data in the computation domain. This is equal to (param_nx, param_ny).

property param_duration: SimulationDuration[source]

Duration of the simulation, either in steps or in simulation time (seconds). The duration is represented by an object of class SimulationDuration.

property param_report_period: SimulationDuration[source]

Duration of the period for reporting, either in steps or in simulation time (seconds). The duration is represented by an object of class SimulationDuration.

property param_timestep: float[source]

In case of a constant timestep strategy, this is the duration of a step, expressed in seconds.

property param_timestep_strategy: TimeStepStrategy[source]

Strategy used to compute the time step at each iteration. The strategy is represented by an object of class TimeStepStrategy.

property param_infiltration_lerp: InfiltrationInterpolation[source]

constant values, linear interpolation,… The interpolation type is represented by an object of class InfiltrationInterpolation.

Type:

Type of infiltration interpolation

property boundary_condition: list[boundary_condition_2D][source]
property infiltrations_chronology: list[float, list[float]][source]

Returns the chronology of infiltration.

This is a list. Each entry of the list is: [start_time, [list of values for each infiltration zone]]

The start_time is when the infiltration should start. The list of value is the quantity of water falling on each mesh of each infitlation zones (the n-th entry of the list corresponds to the n-th infiltration zone).

property bridge_roof: numpy.ndarray[source]

Return the bridge roof.

path: pathlib.Path[source]

Directory where the simulation was loaded or will be saved.

creation_date: datetime.datetime[source]

Creation or last save date.

comment: str[source]

A free text comment, no formatting supported.

param_froude_max: float[source]

Froude limit

_param_froude_bc_limit_tolerance: float[source]

When Froude BC applies and needs to push up the water height, it will push it up to current_water_height * tolerance.

param_runge_kutta: float[source]

Runge-Kutta ponderation. The value is the ponderation of the estimator. Thus 1-ponderation is applied to the corrector. None means no RK but simple Euler step.

_param_dx: float[source]

Size of a single mesh along x axis, in meters.

_param_dy: float[source]

Size of a single mesh along y axis, in meters.

_param_nx: int[source]

Size of the computation domain along the x axis, in number of meshes.

_param_ny: int[source]

Size of the computation domain along the y axis, in number of meshes.

param_base_coord_ll_x: float[source]

The computation domain coordinates are relative to “base” coordinates. This parameter is the x component of the base coordinates. The basis of the coordinates system is located, by convention, on the lower left corner of the domain (i.e. at (0,0) in the local coordinates system).

param_base_coord_ll_y: float[source]

The computation domain coordinates are relative to “base” coordinates. This parameter is the y component of the base coordinates. The basis of the coordinates system is located, by convention, on the lower left corner of the domain (i.e. at (0,0) in the local coordinates system).

_param_timestep_strategy: TimeStepStrategy[source]

Strategy used to compute the time step at each iteration.

_param_timestep: float[source]

In case of a constant time step strategy, this is the duration of the time step. It will be ignored in other cases.

_param_infiltration_lerp: InfiltrationInterpolation[source]

constant values, linear interpolation,…

Type:

Type of infiltration interpolation

_param_duration: SimulationDuration[source]

Duration of the simulation, either in steps or in simulation time (seconds)

_param_report_period: SimulationDuration[source]

Duration of the period for reporting, either in steps or in simulation time (seconds)

_h: numpy.typing.NDArray[numpy.float32][source]

Matrix of water height, in meters.

_qx: numpy.typing.NDArray[numpy.float32][source]

Matrix of quantity of movement along x axis, in m²/s.

_qy: numpy.typing.NDArray[numpy.float32][source]

Matrix of quantity of movement along y axis, in m²/s.

_bathymetry: numpy.typing.NDArray[numpy.float32][source]

Bathymetry, in meters.

_manning: numpy.typing.NDArray[numpy.float32][source]

Manning coefficients.

_infiltration_zones: numpy.typing.NDArray[numpy.int32][source]

Infiltration zones indices. Indices are numbered from one. 0 means no infiltration in the mesh.

_nap: numpy.typing.NDArray[numpy.uint8][source]

Matrix representing which mesh are part of the computation domain or not. 1 means a mesh is part of the computation domain, 0 means it’s not.

_bridge_roof_elevation: numpy.typing.NDArray[numpy.float32][source]

Elevation of the bridge roof.

_wp: wolfhece.PyParams.Wolf_Param[source]

Attributes of the simulation as Wolf_Param instance.

_check_np_array_type(a: numpy.array, dtype=np.float32, check_size=True)[source]
force_Euler_scheme()[source]

Force the simulation to use a simple Euler scheme.

add_boundary_condition(i: int, j: int, bc_type: BoundaryConditionsTypes, bc_value: float, border: Direction)[source]

Add a boundary condition to the model.

Parameters:
  • j (i,) – position of the boundary condition (1-based)

  • bc_type – type of boundary condition

  • bc_value – value/parameter of the boundary condition

  • border – mesh’s border on which the boundary condtion applies

number_of_bc_along_X()[source]
number_of_bc_along_Y()[source]
add_cell_parameter(i: int, j: int, bc_type: CellParameterType, bc_value: float)[source]

Add a cell parameter to the model.

Parameters:
  • i – i position of the boundary condition (1-based)

  • j – j position of the boundary condition (1-based)

  • bc_type – type of cell parameter

  • bc_value – value/parameter of the boundary condition

  • border – mesh’s border on which the boundary condtion applies

get_cell_parameter(i: int, j: int, bc_type: CellParameterType) float[source]

Get the value of a cell parameter at a given position.

Parameters:
  • j (i,) – position of the boundary condition (1-based)

  • bc_type – type of cell parameter

get_all_cell_parameters(bc_type: CellParameterType) List[CellParameter][source]

Get all the cell parameters of a given type.

Parameters:

bc_type – type of cell parameter

clear_cell_parameters(bc_type: CellParameterType = None)[source]

Clears all the cell parameters. Useful if one wants to build a new simulation starting with an existing one.

Parameters:

bc_type – if given, only the cell parameters of this type will be cleared.

clear_boundary_conditions()[source]

Clears all the boundary conditions. Useful if one wants to build a new simulation starting with an existing one.

add_infiltration(time_start: float, zones_values: list[float])[source]

Add an infiltration point in the infiltration chronology.

Parameters:
  • time_start – start time of the infiltration (in seconds)

  • zones_values – an array representing the quantity of water per second to add to each infiltration zones, starting at time_start. The first item corrsesponds to the zone 1, second item corresponds to zone 2 and so on.

clear_infiltration_chronology()[source]

Clear the infiltration chronology. Useful if one wants to build a new simulation starting with an existing one.

_infiltrations_timings_as_gpu_array()[source]

Compute the infiltration timing array in a form suitable for the GPU processing.

_check_infiltrations(infiltration_zones: numpy.ndarray, infiltration_timings, interp_infiltration: InfiltrationInterpolation) list[str][source]

Check if the infiltration data is consistent.

plot_infiltration(figax=None, toshow: bool = True)[source]

Plot the infiltration chronology.

Parameters:
  • figax – a tuple (fig, ax) of matplotlib objects. If None, a new figure is created.

  • toshow – if True, the figure is shown.

Returns:

the figure and the axis.

check_warnings() None | list[str][source]

Check if the model currently set in this object is valid. A number of tests are made, they are not exhaustive.

Returns None if no problem was found. Returns a list of warnings messages else.

check_errors() None | list[str][source]

Check if the model currently set in this object is valid. A number of tests are made, they are not exhaustive.

Returns None if no problem was found. Returns a list of error messages else.

_fail_if_invalid()[source]
classmethod load(simulation_dir: pathlib.Path | str)[source]

Load a simulation.

Parameters:

store_dir (Path) – where the simulation model is stored (a directory).

copy() SimpleSimulation[source]

Create a new simulation from an existing one.

save(store_dir: pathlib.Path = None)[source]

Save the model to a directory. If store_dir is None, the path used to load the model is used. If the model was not loaded from a directory, an exception is raised.

_save_json(store_dir: pathlib.Path = None)[source]

Save the parameters to a json file.

to_wolfparam() wolfhece.PyParams.Wolf_Param[source]

Return a Wolf_Param object that represents the parameters of the simulation.

add_to_wolfparam(dest: wolfhece.PyParams.Wolf_Param) wolfhece.PyParams.Wolf_Param[source]

Return a Wolf_Param object that represents the parameters of the simulation.

from_wolfparam(wp: wolfhece.PyParams.Wolf_Param)[source]

Callback called when the Wolf_Param object is modified and apply routine is called.

write_initial_condition_from_record(recpath: pathlib.Path = None, id_rec: int = None, destpath: pathlib.Path = None) int[source]

Define the initial condition of the simulation from a previous record.

Parameters:
  • recpath – the path to the records. if None, the default path is used and ‘simul_gpu_results’ as result directory.

  • id_rec – the index of the record you want to start from.

  • destpath – the path where to save the initial condition. If None, the current path is used.

Returns:

0 if everything went well, -1 if the record doesn’t exist, -2 if the record path doesn’t exist.

copy_ic_from_simulation(source: SimpleSimulation)[source]

Copy the initial condition from a simulation to the current simulation.

copy_ic_from_record(source: pathlib.Path, id_rec: int = None)[source]

Copy the initial condition from a record to the current simulation.

has_infiltration()[source]
mask_nap_contour(nb_cells: int = 1)[source]

Mask the contour of the NAP array.

make_void_arrays()[source]

Create void arrays from (nx, ny) shape.

Arrays :
  • nap

  • h

  • qx

  • qy

  • manning

  • bathymetry

make_void_infiltration_zones()[source]

Create a void infiltration zones array.

has_bridge() bool[source]

Return True if the simulation has a bridge roof.

make_void_bridge_roof()[source]

Create a void bridge roof array.

convert_roof_array2params(roof_array: numpy.ndarray = None)[source]

Convert a roof array to a list of parameters.

It will remove all the existing bridge roof parameters before add the new ones.

_params_have_duplicated_values() str[source]

Check if the parameters are not duplicated.

ensure_nullvalue_outside_nap(nodata=99999.0)[source]
Force the nodata where nap is 0 for:
  • topo/bathymetry

  • manning

Force zero value for initial conditions where nap is 0:
  • h

  • qx

  • qy