ogstools.meshlib package#

class ogstools.meshlib.Boundary[source]#

Bases: ABC

Abstract base class representing a boundary within a mesh.

A boundary refers to the set of edges or faces that defines the delineation between the interior region and exterior regions of a mesh. In a 2D mesh, it is formed by a closed collection of line segments (1D). In a 3D mesh, it is formed by a closed collection of faces (2D).

abstract dim()[source]#

Get the dimension of the boundary.

Returns:

The dimension of the boundary. For example, the dimension of a boundary of a cube (3D) is 2.

Return type:

int

ogstools.meshlib.Gaussian2D(bound2D, amplitude, spread, height_offset, n)[source]#

Generate a 2D Gaussian-like surface using the provided parameters.

This method computes a 2D Gaussian-like surface by sampling the given bound and parameters.

Parameters:
  • bound2D (tuple) – Tuple of boundary coordinates (x_min, x_max, y_min, y_max).

  • amplitude (float) – Amplitude or peak value of the Gaussian curve.

  • spread (float) – Scaling factor that controls the spread or width of the Gaussian curve.

  • height_offset (float) – Constant offset added to elevate the entire surface.

  • n (int) – Number of points in each dimension for sampling.

Returns:

pyvista.PolyData: A PyVista PolyData object representing the generated surface.

Return type:

DataSet

note:
  • The larger amplitude, the taller the peak of the surface.

  • The larger spread, the wider and flatter the surface.

  • height_offset shifts the entire surface vertically.

example:

Generating a 2D Gaussian-like surface: ` bound = (-1.0, 1.0, -1.0, 1.0) amplitude = 1.0 spread = 0.5 height_offset = 0.0 n = 100 surface = MyClass.Gaussian2D(bound, amplitude, spread, height_offset, n) `

class ogstools.meshlib.Layer[source]#

Bases: Boundary

Layer(top: ogstools.meshlib.boundary_subset.Surface, bottom: ogstools.meshlib.boundary_subset.Surface, material_id: int = 0, num_subdivisions: int = 0)

top: Surface#
bottom: Surface#
material_id: int = 0#
num_subdivisions: int = 0#

Class representing a geological layer with top and bottom surfaces.

A geological layer is a distinct unit of rock or sediment that has unique properties and characteristics, associated by the material_id. It is often bounded by two surfaces: the top surface and the bottom surface. These surfaces delineate the spatial extent of the layer in the GIS system.

create_raster(resolution)[source]#

Create raster representations for the layer.

For each surface, including intermediate surfaces (num_of_subdivisions > 0), this method generates .asc files.

Parameters:

resolution (float) – The resolution for raster creation.

Returns:

A list of filenames to .asc raster files.

Return type:

list[Path]

dim()[source]#

Get the dimension of the boundary.

Returns:

The dimension of the boundary. For example, the dimension of a boundary of a cube (3D) is 2.

Return type:

int

__init__(top, bottom, material_id=0, num_subdivisions=0)#
class ogstools.meshlib.LayerSet[source]#

Bases: BoundarySet

Collection of geological layers stacked to represent subsurface arrangements.

In a geological information system, multiple layers can be stacked vertically to represent the subsurface arrangement. This class provides methods to manage and work with layered geological data.

Initializes a LayerSet. It checks if the list of provided layers are given in a top to bottom order. In neighboring layers, layers share the same surface (upper bottom == low top).

__init__(layers)[source]#

Initializes a LayerSet. It checks if the list of provided layers are given in a top to bottom order. In neighboring layers, layers share the same surface (upper bottom == low top).

bounds()[source]#
Return type:

list

filenames()[source]#
Return type:

list[Path]

classmethod from_pandas(df)[source]#

Create a LayerSet from a Pandas DataFrame.

Return type:

LayerSet

create_raster(resolution)[source]#

Create raster representations for the LayerSet.

This method generates raster files at a specified resolution for each layer’s top and bottom boundaries and returns paths to the raster files.

Return type:

tuple[Path, Path]

create_rasters(resolution)[source]#

For each surface a (temporary) raster file with given resolution is created.

Return type:

list[Path]

class ogstools.meshlib.LocationFrame[source]#

Bases: object

LocationFrame(xmin: float, xmax: float, ymin: float, ymax: float)

xmin: float#
xmax: float#
ymin: float#
ymax: float#
as_gml(filename)[source]#

Generate GML representation of the location frame.

Parameters:

filename (Path) – The filename to save the GML representation to.

Returns:

None

Return type:

None

__init__(xmin, xmax, ymin, ymax)#
class ogstools.meshlib.MeshSeries[source]#

Bases: Sequence[UnstructuredGrid]

A wrapper around pyvista and meshio for reading of pvd and xdmf timeseries.

Initialize a MeshSeries object

param filepath:

Path to the PVD or XDMF file.

param spatial_unit:

Unit/s of the mesh points. See note.

param time_unit:

Unit/s of the timevalues. See note.

returns:

A MeshSeries object

Note:

If given as a single string, the data is read in SI units i.e. in seconds and meters and converted to the given units. If given as a tuple, the first str corresponds to the data_unit, the second to the output_unit. E.g.: ot.MeshSeries(filepath, "km", ("a", "d")) would read in the spatial data in meters and convert to kilometers and read in the timevalues in years and convert to days.

__init__(filepath=None, spatial_unit='m', time_unit='s')[source]#

Initialize a MeshSeries object

param filepath:

Path to the PVD or XDMF file.

param spatial_unit:

Unit/s of the mesh points. See note.

param time_unit:

Unit/s of the timevalues. See note.

returns:

A MeshSeries object

Note:

If given as a single string, the data is read in SI units i.e. in seconds and meters and converted to the given units. If given as a tuple, the first str corresponds to the data_unit, the second to the output_unit. E.g.: ot.MeshSeries(filepath, "km", ("a", "d")) would read in the spatial data in meters and convert to kilometers and read in the timevalues in years and convert to days.

classmethod from_data(meshes, timevalues, spatial_unit='m', time_unit='s')[source]#

Create a MeshSeries from a list of meshes and timevalues.

Return type:

MeshSeries

extend(mesh_series)[source]#

Extends self with mesh_series. If the last element of the mesh series is within epsilon to the first element of mesh_series to extend, the duplicate element is removed

resample_temporal(timevalues)[source]#

Return a new MeshSeries interpolated to the given timevalues.

Return type:

MeshSeries

probe(points, data_name=None, interp_method='linear')[source]#

Create a new MeshSeries by probing points on an existing MeshSeries.

Parameters:
  • points (ndarray) – The points at which to probe.

  • data_name (str | Variable | list[str | Variable] | None) – Data to extract. If None, use all point data.

  • interp_method (Literal['nearest', 'linear']) – The interpolation method to use.

Returns:

A MeshSeries (Pointcloud) containing the probed data.

Return type:

MeshSeries

interpolate(mesh, data_name=None)[source]#

Create a new MeshSeries by spatial interpolation.

Parameters:
  • mesh (DataSet) – The mesh on which to interpolate.

  • data_name (str | Variable | list[str | Variable] | None) – Data to extract. If None, use all point data.

Returns:

A spatially interpolated MeshSeries.

Return type:

MeshSeries

copy(deep=True)[source]#

Create a copy of MeshSeries object. Deep copy is the default.

Parameters:

deep (bool) – switch to choose between deep (default) and shallow (self.copy(deep=False)) copy.

Returns:

Copy of self.

Return type:

MeshSeries

__getitem__(index: int) UnstructuredGrid[source]#
__getitem__(index: slice | Sequence) MeshSeries
__getitem__(index: str) ndarray
items()[source]#

Returns zipped tuples of timevalues and meshes.

Return type:

Sequence[tuple[float, UnstructuredGrid]]

aggregate_temporal(variable, func)[source]#

Aggregate data over all timesteps using a specified function.

Parameters:
  • variable (Variable | str) – The mesh variable to be aggregated.

  • func (Callable) – The aggregation function to apply. E.g. np.min, np.max, np.mean, np.median, np.sum, np.std, np.var

Returns:

A mesh with aggregated data according to the given function.

Return type:

UnstructuredGrid

clear_cache()[source]#
closest_timestep(timevalue)[source]#

Return the corresponding timestep from a timevalue.

Return type:

int

closest_timevalue(timevalue)[source]#

Return the closest timevalue to a timevalue.

Return type:

float

ip_tesselated()[source]#

Create a new MeshSeries from integration point tessellation.

Return type:

MeshSeries

mesh(timestep, lazy_eval=True)[source]#

Returns the mesh at the given timestep.

Return type:

UnstructuredGrid

rawdata_file()[source]#

Checks, if working with the raw data is possible. For example, OGS Simulation results with XDMF support efficient raw data access via h5py

Returns:

The location of the file containing the raw data. If it does not support efficient read (e.g., no efficient slicing), it returns None.

Return type:

Path | None

mesh_interp(timevalue, lazy_eval=True)[source]#

Return the temporal interpolated mesh for a given timevalue.

Return type:

UnstructuredGrid

property timevalues: ndarray#

Return the timevalues.

property timesteps: list#

Return the OGS simulation timesteps of the timeseries data. Not to be confused with timevalues which returns a list of times usually given in time units.

values(variable: str | Variable) ndarray[source]#
values(variable: list[str | Variable]) list[ndarray]

Get the data in the MeshSeries for all timesteps.

Adheres to time slicing via __getitem__ and an applied pyvista filter via transform if the applied filter produced ‘vtkOriginalPointIds’ or ‘vtkOriginalCellIds’ (e.g. clip(…, crinkle=True), extract_cells(…), threshold(…).)

Parameters:

variable – Data to read/process from the MeshSeries. Can also be a list of str or Variable.

Returns:

A numpy array of shape (n_timesteps, n_points/c_cells). If given an argument of type Variable is given, its transform function is applied on the data. If a list of str or Variable is given, a list of the individual values is returned.

property point_data: DataDict#

Useful for reading or setting point_data for the entire meshseries.

property cell_data: DataDict#

Useful for reading or setting cell_data for the entire meshseries.

property field_data: DataDict#

Useful for reading or setting field_data for the entire meshseries.

time_of_min(variable)[source]#

Returns a Mesh with the time of the variable minimum as data.

Return type:

UnstructuredGrid

time_of_max(variable)[source]#

Returns a Mesh with the time of the variable maximum as data.

Return type:

UnstructuredGrid

aggregate_spatial(variable, func)[source]#

Aggregate data over domain per timestep using a specified function.

Parameters:
  • variable (Variable | str) – The mesh variable to be aggregated.

  • func (Callable) – The aggregation function to apply. E.g. np.min, np.max, np.mean, np.median, np.sum, np.std, np.var

Returns:

A numpy array with aggregated data.

Return type:

ndarray

probe_values(points, data_name, interp_method='linear')[source]#

Return the sampled data of the MeshSeries at observation points.

Similar to probe() but returns the data directly instead of creating a new MeshSeries.

Parameters:
  • points (ndarray | list) – The observation points to sample at.

  • data_name (str | Variable | list[str | Variable]) – Data to sample. If provided as a Variable, the output will transformed accordingly. Can also be a list of str or Variable.

  • interp_method (Literal['nearest', 'linear']) – Interpolation method, defaults to linear

Returns:

numpy array/s of interpolated data at observation points with the following shape:

  • multiple points: (n_timesteps, n_points, [n_components])

  • single points: (n_timesteps, [n_components])

If data_name is a list, a corresponding list of arrays is returned.

Return type:

ndarray | list[ndarray]

plot_line(var1=None, var2=None, ax=None, sort=True, outer_legend=False, **kwargs)#

Plot some data of a (1D) dataset.

You can pass “x”, “y” or “z” to either of x_var or y_var to specify which spatial dimension should be used for the corresponding axis. By passing “time” the timevalues will be use for this axis. You can also pass two data variables for a phase plot. if no value is given, automatic detection of spatial axis is tried.

>>> line(ms, ot.variables.temperature)          # temperature over time
>>> line(ms, ot.variables.temperature, "time")  # time over temperature
>>> line(ms, "pressure", "temperature")     # temperature over pressure
>>> line(mesh, ot.variables.temperature)    # temperature over x, y or z
>>> line(mesh, "y", "temperature")          # temperature over y
>>> line(mesh, ot.variables.pressure, "y")  # y over pressure
>>> line(mesh)  # z=const: y over x, y=const: z over x, x=const: z over y
Parameters:
  • var1 (str | Variable | None) – Variable for the x-axis if var2 is given else for y-axis.

  • var2 (str | Variable | None) – Variable for the y-axis if var1 is given.

  • ax (Axes | None) – The matplotlib axis to use for plotting, if None a new figure will be created.

  • sort (bool) – Automatically sort the values along the dimension of the mesh with the largest extent

Outer_legend:

Draw legend to the right next to the plot area. By default False (legend stays inside). User can pass a tuple of two floats (x, y), which will be passed to bbox_to_anchor parameter in matplotlib legend call. True will pass the default values (1.05, 1.0).

Return type:

Figure | None

Keyword Arguments:
  • figsize: figure size (default=[16, 10])

  • color: color of the line

  • linewidth: width of the line

  • linestyle: style of the line

  • label: label in the legend

  • grid: if True, show grid

  • monospace: if True, the legend uses a monospace font

  • loc: location of the legend (default=”upper right”)

  • annotate: string to be annotate at the center of the mesh

  • clip_on: If True, clip the output to stay within the Axes.

    (default=False)

  • all other kwargs get passed to matplotlib’s plot function

Note:

Using loc=”best” will take a long time, if you plot lines on top of a contourplot, as matplotlib is calculating the best position against all the underlying cells.

plot_time_slice(x, y, variable, time_logscale=False, fig=None, ax=None, cbar=True, **kwargs)[source]#

Create a heatmap for a variable over time and space.

Parameters:
  • x (Literal['x', 'y', 'z', 'time']) – What to display on the x-axis (x, y, z or time)

  • y (Literal['x', 'y', 'z', 'time']) – What to display on the y-axis (x, y, z or time)

  • variable (str | Variable) – The variable to be visualized.

  • time_logscale (bool) – Should log-scaling be applied to the time-axis?

  • fig (Figure | None) – matplotlib figure to use for plotting.

  • ax (Axes | None) – matplotlib axis to use for plotting.

  • cbar (bool) – If True, adds a colorbar.

Return type:

Figure | None

Keyword Arguments:
  • cb_labelsize: colorbar labelsize

  • cb_loc: colorbar location (‘left’ or ‘right’)

  • cb_pad: colorbar padding

  • cmap: colormap

  • vmin: minimum value for colorbar

  • vmax: maximum value for colorbar

  • num_levels: number of levels for colorbar

  • figsize: figure size

  • dpi: resolution

  • log_scaled: logarithmic scaling

property mesh_func: Callable[[UnstructuredGrid], UnstructuredGrid]#

Returns stored transformation function or identity if not given.

transform(mesh_func=lambda mesh: ...)[source]#

Apply a transformation function to the underlying mesh.

Parameters:

mesh_func (Callable[[UnstructuredGrid], UnstructuredGrid]) – A function which expects to read a mesh and return a mesh. Useful for slicing / clipping / thresholding.

Returns:

A deep copy of this MeshSeries with transformed meshes.

Return type:

MeshSeries

scale(spatial=1.0, time=1.0)[source]#

Scale the spatial coordinates and timevalues.

Useful to convert to other units, e.g. “m” to “km” or “s” to “a”. Converts from SI units (i.e. ‘m’ and ‘s’) to the given arguments.

Parameters:
  • spatial (int | float | str) – Float factor or str for target unit.

  • time (int | float | str) – Float factor or str for target unit.

Returns:

None.

Return type:

None

classmethod difference(ms_a, ms_b, variable=None)[source]#

Compute difference of variables between the two MeshSeries instances from which this method is called and a second MeshSeries instance passed as method parameter. Returns new instance of MeshSeries: ms = ms_a - ms_b

Parameters:
  • base_ms – The mesh from which data is to be subtracted.

  • subtract_mesh – The mesh whose data is to be subtracted.

  • variable (Variable | str | None) – The variable of interest. If not given, all point and cell_data will be processed raw.

Returns:

MeshSeries containing the difference of variable` or of all datasets between both MeshSeries.

Return type:

MeshSeries

extract(index, preference='points')[source]#

Extract a subset of the domain by point or cell indices.

Parameters:
  • index (slice | int | ndarray | list) – Indices of points or cells to extract.

  • preference (Literal['points', 'cells']) – Selected entities.

Returns:

A MeshSeries with the selected domain subset.

Return type:

MeshSeries

save(filename, deep=True, ascii=False)[source]#

Save mesh series to disk.

Parameters:
  • filename (str) – Filename to save the series to. Extension specifies the file type. Currently only PVD is supported.

  • deep (bool) – Specifies whether VTU/H5 files should be written.

  • ascii (bool) – Specifies if ascii or binary format should be used, defaults to binary (False) - True for ascii.

remove_array(name, data_type='field', skip_last=False)[source]#

Removes an array from all time slices of the mesh series.

Parameters:
  • name (str) – Array name

  • data_type (str) – Data type of the array. Could be either field, cell or point

  • skip_last (bool) – Skips the last time slice (e.g. for restart purposes).

class ogstools.meshlib.Meshes[source]#

Bases: MutableMapping

OGS input mesh. Refers to prj - file section <meshes>

Initialize a Meshes object.
param meshes:

List of Mesh objects (pyvista UnstructuredGrid) The first mesh is the domain mesh. All following meshes represent subdomains, and their points must align with points on the domain mesh. If needed, the domain mesh itself can also be added again as a subdomain.

returns:

A Meshes object

__init__(meshes)[source]#
Initialize a Meshes object.
param meshes:

List of Mesh objects (pyvista UnstructuredGrid) The first mesh is the domain mesh. All following meshes represent subdomains, and their points must align with points on the domain mesh. If needed, the domain mesh itself can also be added again as a subdomain.

returns:

A Meshes object

__getitem__(key)[source]#
Return type:

UnstructuredGrid

classmethod from_files(filepaths)[source]#

Initialize a Meshes object from a Sequence of existing files.

Parameters:

filepaths (Sequence[str | Path]) – Sequence of Mesh files (.vtu) The first mesh is the domain mesh. All following meshes represent subdomains, and their points must align with points on the domain mesh.

Return type:

Meshes

classmethod from_yaml(geometry_file)[source]#
Return type:

Meshes

classmethod from_gmsh(filename, dim=0, reindex=True, log=True, meshname='domain')[source]#

Generates pyvista unstructured grids from a gmsh mesh (.msh).

Extracts domain-, boundary- and physical group-submeshes.

Parameters:
  • filename (Path) – Gmsh mesh file (.msh) as input data

  • dim (int | Sequence[int]) – Spatial dimension (1, 2 or 3), trying automatic detection, if not given. If multiple dimensions are provided, all elements of these dimensions are embedded in the resulting domain mesh.

  • reindex (bool) – Physical groups / regions / Material IDs to be renumbered consecutively beginning with zero.

  • log (bool) – If False, silence log messages

  • meshname (str) – The name of the domain mesh and used as a prefix for subdomain meshes.

Returns:

A dictionary of names and corresponding meshes

Return type:

Meshes

classmethod from_mesh(mesh, threshold_angle=15.0, domain_name='domain')[source]#

Extract 1D boundaries of a 2D mesh.

Parameters:
  • mesh (UnstructuredGrid) – The 2D domain

  • threshold_angle (float | None) – If None, the boundary will be split by the assumption of vertical lateral boundaries. Otherwise it represents the angle (in degrees) between neighbouring elements which - if exceeded - determines the corners of the boundary mesh.

  • domain_name (str) – The name of the domain mesh.

Returns:

A Meshes object.

Return type:

Meshes

add_gml_subdomains(domain_path, gml_path, out_dir=None, tolerance=1e-12)[source]#

Add Meshes from geometry definition in the gml file to subdomains.

Parameters:
  • gml_file – Path to the gml file.

  • out_dir (Path | None) – Where to write the gml meshes (default: gml dir)

  • tolerance (float) – search length for node search algorithm

sort()[source]#

Sort the subdomains alphanumerically.

property domain: UnstructuredGrid#

Get the domain mesh.

By convention, the domain mesh is the first mesh in the dictionary of meshes when the Meshes object was constructed. The domain mesh is expected to be constant. e.g. Do not: myobject.domain = pv.Sphere()

Returns:

The domain mesh

property domain_name: str#

Get the name of the domain mesh.

By convention, the domain mesh is the first mesh in the dictionary of meshes when the Meshes object was constructed.

Returns:

The name of the domain mesh

property subdomains: dict[str, UnstructuredGrid]#

Get the subdomain meshes.

By convention, all meshes except the first one are considered subdomains. This returns a list of those subdomain meshes.

Returns:

A dictionary of {name: Mesh} for all subdomains

identify_subdomain()[source]#
rename_subdomains(rename_map)[source]#

Rename subdomain meshes according to the provided mapping.

Parameters:

rename_map (dict[str, str]) – A dictionary mapping old subdomain names -> new names. e.g. {‘left’:’subdomain_left’}

rename_subdomains_legacy()[source]#

Add to the name physical_group to restore legacy convention

static create_metis(domain_file, output_path, dry_run=False)[source]#

Creates a metis files. This file is needed to partition the OGS input mesh (for parallel OGS compution) using the OGS cmd line tool partmesh.

Parameters:
  • domain_file (Path | str) – A Path to existing domain mesh file (.vtu extension)

  • output – A Path to existing folder. Here the resulting metis file will be stored (.mesh)

  • dry_run (bool) – If True: Metis file is not written If False: Metis file is written

Returns:

Path to the generated metis file.

Return type:

Path

static create_partitioning(num_partitions, domain_file, subdomain_files, metis_file=None, dry_run=False)[source]#

Creates a subfolder in the metis_file’ folder. Puts .bin files into this folder that are needed as input files for running OGS parallel (MPI). Wrapper around command line tool partmesh, adding file checks, dry-run option, normalized behaviour for partition == 1 Only use this function directly when you want to bypass creating the Meshes object (e.g. files for domain and subdomains are already present)

Parameters:
  • num_partitions (int) – List of integers or a single integer that indicate the number of partitions similar to the OGS binary tool partmesh. The serial mesh will always be generated Example 1: num_partitions = [1,2,4,8,16] Example 2: num_partitions = 2

  • domain_file (Path | str) – A Path to existing domain mesh file (.vtu extension)

  • subdomain_files (Sequence[Path | str]) – A list of Path to existing subdomain files (.vtu extensions)

  • metis_file (str | Path | None) – A Path to existing metis partitioned file (.mesh extension).

  • dry_run (bool) – If True: Writes no files, but returns the list of files expected to be created If False: Writes files and returns the list of created files

Returns:

A list of Paths pointing to the saved mesh files, if num_partitions are given (also just [1]), then it returns A dict, with keys representing the number of partitions and values A list of Paths (like above)

Return type:

list[Path]

save(meshes_path=None, overwrite=False, num_partitions=None, dry_run=False)[source]#

Save all meshes.

This function will perform identifySubdomains, if not yet been done.

Parameters:
  • meshes_path (str | Path | None) – Optional path to the directory where meshes should be saved. It must already exist (will not be created). If None, a temporary directory will be used.

  • overwrite (bool) – If True, existing mesh files will be overwritten. If False, an error is raised if any file already exists.

  • num_partitions (int | Sequence[int] | None) – List of integers or a single integer that indicate the number of partitions similar to the OGS binary tool partmesh. The serial mesh will always be generated Example 1: num_partitions = [1,2,4,8,16] Example 2: num_partitions = 2

  • dry_run (bool) – If True: Writes no files, but returns the list of files expected to be created If False: Writes files and returns the list of created files

Returns:

A list of Paths pointing to the saved mesh files, if num_partitions are given (also just [1]), then it returns A dict, with keys representing the number of partitions and values A list of Paths (like above)

Return type:

list[Path] | dict[int, list[Path]]

plot(**kwargs)[source]#

Plot the domain mesh and the subdomains.

keyword arguments: see ogstools.plot.contourf()

Return type:

Figure

class ogstools.meshlib.Raster[source]#

Bases: object

Class representing a raster representation of a location frame.

This class provides methods to create and save a raster representation based on a specified location frame and resolution.

frame: LocationFrame#
resolution: float#
__init__(frame, resolution)#
as_vtu(outfilevtu)[source]#

Create and save a raster representation as a VTK unstructured grid.

Parameters:

outfilevtu (Path) – The path to save the VTK unstructured grid.

Returns:

The path to the saved VTK unstructured grid representation.

Return type:

Path

class ogstools.meshlib.Surface[source]#

Bases: object

A surface is a sub group of a polygon mesh (2D). A surface is not closed and therefore does not represent a volume. (Geological) layers (stratigraphic units) can be defined by an upper and lower surface. By convention, properties (material_id and resolution ), actually associated to the stratigraphic unit layer, are given together with the lower boundary (class Surface) of a stratigraphic unit (class Layer).

Initialize a surface mesh. Either from pyvista or from a file.

property material_id: int#
__init__(input, material_id)[source]#

Initialize a surface mesh. Either from pyvista or from a file.

create_raster_file(resolution)[source]#

Creates a raster file specific to resolution. If outfile is not specified, the extension is replaced by .asc

:returns the path and filename of the created file (.asc)

Return type:

Path

ogstools.meshlib.cuboid(lengths=1.0, n_edge_cells=1, n_layers=1, structured_grid=True, order=1, mixed_elements=False, out_name=Path('unit_cube.msh'), msh_version=None)[source]#
ogstools.meshlib.depth(mesh, use_coords=False)[source]#

Returns the depth values of the mesh.

For 2D, the last axis of the plane wherein the mesh is lying is used as the vertical axis (i.e. y if the mesh is in the xy-plane, z if it is in the xz-plane), for 3D, the z-axes is used. If use_coords is True, returns the negative coordinate value of the vertical axis. Otherwise, the vertical distance to the top facing edges surfaces are returned.

Return type:

ndarray

ogstools.meshlib.difference(base_mesh, subtract_mesh, variable=None)[source]#

Compute the difference of variables between two meshes.

Parameters:
  • base_mesh (Mesh) – The mesh to subtract from.

  • subtract_mesh (Mesh) – The mesh whose data is to be subtracted.

  • variable (Variable | str | None) – The variable of interest. If not given, all point and cell_data will be processed raw.

Returns:

A new mesh containing the difference of variable or of all datasets between both meshes.

Return type:

Mesh

ogstools.meshlib.difference_matrix(meshes_1, meshes_2=None, variable=None)[source]#

Compute the difference between all combinations of two meshes from one or two arrays based on a specified variable.

Parameters:
  • meshes_1 (list | ndarray) – The first list/array of meshes to be subtracted from.

  • meshes_2 (list | ndarray | None) – The second list/array of meshes, it is subtracted from the first list/array of meshes - meshes_1 (optional).

  • variable (Variable | str | None) – The variable of interest. If not given, all point and cell_data will be processed raw.

Returns:

An array of meshes containing the differences of variable or all datasets between meshes_1 and meshes_2 for all possible combinations.

Return type:

ndarray

ogstools.meshlib.difference_pairwise(meshes_1, meshes_2, variable=None)[source]#

Compute pairwise difference between meshes from two lists/arrays (they have to be of the same length).

Parameters:
  • meshes_1 (list | ndarray | MeshSeries) – The first list/array of meshes to be subtracted from.

  • meshes_2 (list | ndarray | MeshSeries) – The second list/array of meshes whose data is subtracted from the first list/array of meshes - meshes_1.

  • variable (Variable | str | None) – The variable of interest. If not given, all point and cell_data will be processed raw.

Returns:

An array of meshes containing the differences of variable or all datasets between meshes_1 and meshes_2.

Return type:

ndarray

ogstools.meshlib.mesh_from_simulator(simulation, name, node_properties=None, cell_properties=None, field_properties=None)[source]#

Constructs a pyvista mesh from a running simulation. It always contains points (geometry) and cells (topology) and optionally the given node-based or cell-based properties Properties must be added afterwards

param simulator:

Initialized and not finalized simulator object

param name:

Name of the submesh (e.g. domain, left, … )

param node_properties:

Given properties will be added to the mesh None or [] -> no properties will be added

param cell_properties:

Given properties will be added to the mesh None or [] -> no properties will be added

returns:

A Mesh (Pyvista Unstructured Grid) object

Return type:

UnstructuredGrid

ogstools.meshlib.meshes_from_gmsh(filename, dim=0, reindex=True, log=True, meshname='domain')[source]#

Deprecated since version 0.7.1.

Use Meshes. from_gmsh() instead. Generates pyvista unstructured grids from a gmsh mesh (.msh).

Extracts domain, boundary- and submeshes by definition of named physical groups within gmsh.

Parameters:
  • filename (Path) – Gmsh mesh file (.msh) as input data

  • dim (int | Sequence[int]) – Spatial dimension (1, 2 or 3), trying automatic detection, if not given. If multiple dimensions are provided, all elements of these dimensions are embedded in the resulting domain mesh.

  • reindex (bool) – Physical groups / regions / Material IDs to be renumbered consecutively beginning with zero.

  • log (bool) – If False, silence log messages

  • meshame – The name of the domain mesh.

Returns:

A dictionary of names and corresponding meshes

Return type:

dict[str, UnstructuredGrid]

ogstools.meshlib.p_fluid(mesh)[source]#

Return the fluid pressure in the mesh.

If “depth” is given in the mesh’s point _data, it is used return a hypothetical water column defined as:

\[p_{fl} = 1000 \frac{kg}{m^3} 9.81 \frac{m}{s^2} h\]

where h is the depth below surface. Otherwise, If “pressure” is given in the mesh, return the “pressure” data of the mesh. If that is also not the case, the hypothetical water column from above is returned with the depth being calculated via ogstools.meshlib.geo.depth().

Return type:

PlainQuantity

ogstools.meshlib.rect(lengths=1.0, n_edge_cells=1, n_layers=1, structured_grid=True, order=1, mixed_elements=False, jiggle=0.0, out_name=Path('rect.msh'), msh_version=None, layer_ids=None)[source]#

Generates a rectangular mesh using gmsh.

Parameters:
  • lengths (float | tuple[float, float]) – Length of the rectangle in x and y direction. Provide a tuple (x, y) or a scalar for a square. All values must be >= 1e-7 and <= 1e12.

  • n_edge_cells (int | tuple[int, int]) – Number of edge cells in x and y direction. Provide a tuple (x, y) or a scalar for a square. All values must be >= 1.

  • n_layers (int) – Number of layers in y direction. Must be >= 1.

  • structured_grid (bool) – If True, the mesh will be structured. If False, the mesh will be unstructured.

  • order (int) – Order of the mesh elements. 1 for linear, 2 for quadratic.

  • mixed_elements (bool) – If True, the mesh will be mixed elements. If False, the mesh will be structured.

  • jiggle (float) – Amount of random displacement to apply to the mesh nodes. Default is 0.0 (no displacement).

  • out_name (Path | str) – Name of the output mesh file. Default is “rect.msh”.

  • msh_version (float | None) – Version of the GMSH mesh file format. Default is None (use the default version).

  • layer_ids (list | None) – List of layer IDs for the physical groups. If None, the IDs will be generated automatically.

ogstools.meshlib.reindex_material_ids(mesh)[source]#
ogstools.meshlib.to_ip_mesh(mesh)[source]#

Create a mesh with cells centered around integration points.

Return type:

UnstructuredGrid

ogstools.meshlib.to_ip_point_cloud(mesh)[source]#

Convert integration point data to a pyvista point cloud.

Return type:

UnstructuredGrid

ogstools.meshlib.to_region_prism(layer_set, resolution)[source]#

Convert a layered geological structure into a RegionSet using prism meshing.

This function takes a boundary_set.LayerSet and converts it into a region.RegionSet object using prism or tetrahedral meshing technique. The function will use prism elements for meshing if possible; otherwise, it will use tetrahedral elements.

Parameters:
  • layer_set (LayerSet) – A boundary_set.LayerSet.

  • resolution (float) – The desired resolution in [meter] for meshing. It must greater than 0.

Returns:

A boundary_set.LayerSet object containing the meshed representation of the geological structure.

Return type:

RegionSet

raises:

ValueError: If an error occurs during the meshing process.

example:

layer_set = LayerSet(…) resolution = 0.1 region_set = to_region_prism(layer_set, resolution)

ogstools.meshlib.to_region_simplified(layer_set, xy_resolution, rank)[source]#

Convert a layered geological structure to a simplified meshed region.

This function converts a layered geological structure represented by a LayerSet into a simplified meshed region using the specified xy_resolution and rank.

Parameters:
  • (LayerSet) (layer_set) – A LayerSet object representing the layered geological structure.

  • (float) (xy_resolution) – The desired spatial resolution of the mesh in the XY plane.

  • (int) (rank) – The rank of the mesh (2 for 2D, 3 for 3D).

Returns:

A RegionSet object containing the simplified meshed representation of the geological structure.

Return type:

RegionSet

raises:

AssertionError: If the length of the bounds retrieved from the layer_set is not 6.

example:

layer_set = LayerSet(…) xy_resolution = 0.1 # Example resolution in XY plane rank = 2 # Mesh will be 2D region_set = to_region_simplified(layer_set, xy_resolution, rank)

ogstools.meshlib.to_region_tetraeder(layer_set, resolution)[source]#

Convert a layered geological structure to a tetrahedral meshed region.

This function converts a layered geological structure represented by a LayerSet into a tetrahedral meshed region using the specified resolution.

Parameters:
  • layer_set (LayerSet) – A LayerSet object representing the layered geological structure.

  • resolution (int) – The desired resolution for meshing.

Returns:

A RegionSet object containing the tetrahedral meshed representation of the geological structure.

Return type:

RegionSet

raises:

ValueError: If an error occurs during the meshing process.

notes:
  • The LayerSet object contains information about the layers in the geological structure.

  • The resolution parameter determines the desired spatial resolution of the mesh.

  • The function utilizes tetrahedral meshing using Tetgen software to create the meshed representation.

  • The resulting mesh is tetrahedral, and material IDs are assigned to mesh cells based on the geological layers.

example:

layer_set = LayerSet(…) resolution = 1 # Example resolution for meshing region_set = to_region_tetraeder(layer_set, resolution)

ogstools.meshlib.to_region_voxel(layer_set, resolution)[source]#

Convert a layered geological structure to a voxelized mesh.

This function converts a layered geological structure represented by a LayerSet into a voxelized mesh using the specified resolution.

Parameters:
  • layer_set (LayerSet) – A LayerSet object representing the layered geological structure.

  • resolution (list) – A list of [x_resolution, y_resolution, z_resolution] for voxelization.

Returns:

A Mesh object containing the voxelized mesh representation of the geological structure.

Return type:

RegionSet

raises:

ValueError: If an error occurs during the voxelization process.

example:

layer_set = LayerSet(…) resolution = [0.1, 0.1, 0.1] # Example voxelization resolutions in x, y, and z dimensions voxel_mesh = to_region_voxel(layer_set, resolution)

Subpackages#

Submodules#