wolfgpu.simple_simulation

Module Contents

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.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]
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.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]
property param_dy[source]
property h: numpy.typing.NDArray[numpy.float32][source]
property qx: numpy.typing.NDArray[numpy.float32][source]
property qy: numpy.typing.NDArray[numpy.float32][source]
property manning: numpy.typing.NDArray[numpy.float32][source]
property bathymetry: numpy.typing.NDArray[numpy.float32][source]
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 shape: tuple[int, int][source]
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 infitlration zone]]

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_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.

_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]
add_boundary_condition(i: int, j: int, bc_type: BoundaryConditionsTypes, bc_value: float, border: Direction)[source]

Add a boundary condition to the model.

i, j: 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

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.

Param:

time_start start time of the infiltration (in seconds)

Param:

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]
_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_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)[source]

Load a simulation.

Parameters:

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

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.

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.