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.Mesh[source]#

Bases: UnstructuredGrid

A wrapper around pyvista.UnstructuredGrid.

Contains additional data and functions mainly for postprocessing.

Initialize a Mesh object

param pv_mesh:

Underlying pyvista mesh.

filepath: Path | None = None#
difference(subtract_mesh, variable=None)#

Compute the difference of variables between two meshes.

Parameters:
  • 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

depth(use_coords=False)#

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

p_fluid()#

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

plot_contourf(variable, fig=None, ax=None, **kwargs)#

Plot the variable field of meshes with default settings.

The resulting figure adheres to the configurations in plot.setup. For 2D, the whole domain, for 3D a set of slices is displayed.

Parameters:
  • variable (Variable | str) – The field to be visualized on all meshes

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

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

Keyword Arguments:
  • cb_labelsize: colorbar labelsize

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

  • cb_pad: colorbar padding

  • cmap: colormap

  • dpi: resolution

  • figsize: figure size

  • fontsize size for labels and captions

  • levels: user defined levels

  • log_scaled: logarithmic scaling

  • show_edges: show element edges

  • show_max: mark the location of the maximum value

  • show_min: mark the location of the minimum value

  • show_region_bounds: show the edges of the different regions

  • vmin: minimum value

  • vmax: maximum value

Return type:

Figure | None

plot_quiver(ax, variable, projection=None, glyph_type='arrow')#

Plot arrows or lines corresponding to vectors on a matplotlib axis.

Parameters:
  • ax (Axes) – Matplotlib axis to plot onto

  • variable (Vector) – Vector variable to visualize

  • projection (int | None) – Index of flat dimension (e.g. 2 for z axis), gets automatically determined if not given

  • glyph_type (Literal['arrow', 'line']) – Whether to plot arrows or lines.

plot_streamlines(ax, variable, projection=None)#

Plot the vector streamlines on a matplotlib axis.

Parameters:
  • ax (Axes) – Matplotlib axis to plot onto

  • variable (Vector) – Vector variable to visualize

  • projection (int | None) – Index of flat dimension (e.g. 2 for z axis), gets automatically determined if not given

to_ip_mesh()[source]#

Create a mesh with cells centered around integration points.

Return type:

Mesh

to_ip_point_cloud()[source]#

Convert integration point data to a pyvista point cloud.

Return type:

Mesh

__init__(pv_mesh=None, **kwargs)[source]#

Initialize a Mesh object

param pv_mesh:

Underlying pyvista mesh.

classmethod read(filepath)[source]#

Initialize a Mesh object

param filepath:

Path to the mesh or shapefile file.

returns:

A Mesh object

Return type:

Mesh

classmethod read_shape(simplify=False, mesh_generator='triangle', cellsize=None)#

Generate a pyvista Unstructured grid from a shapefile.

Parameters:
  • simplify (bool) – With the Douglas-Peucker algorithm the geometry is simplified. The original line is split into smaller parts. All points with a distance smaller than half the cellsize are removed. Endpoints are preserved. More infos at https://geopandas.org/en/stable/docs/reference/api/geopandas.GeoSeries.simplify.html.

  • mesh_generator (str) – Choose between ‘triangle’ and ‘gmsh’ to generate the mesh.

  • cellsize (int | None) – Size of the cells in the mesh - only needed for simplify algorithm. If None - cellsize is 1/100 of larger bound (x or y).

Returns:

pv.UnstructuredGrid

Return type:

UnstructuredGrid

reindex_material_ids()[source]#
classmethod read_feflow(feflow_file)[source]#

Initialize a Mesh object read from a FEFLOW file. This mesh stores all model specific information such as boundary conditions or material parameters.

param feflow_file:

Path to the feflow file.

returns:

A Mesh object

Return type:

Mesh

class ogstools.meshlib.MeshSeries[source]#

Bases: object

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.

returns:

A MeshSeries object

__init__(filepath)[source]#

Initialize a MeshSeries object

param filepath:

Path to the PVD or XDMF file.

returns:

A MeshSeries object

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) Mesh[source]#
__getitem__(index: slice | Sequence) MeshSeries
aggregate_over_time(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:

Mesh

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:

Mesh

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

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

Return the temporal interpolated mesh for a given timevalue.

Return type:

Mesh

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

Get the data in the MeshSeries for all timesteps.

Adheres to time slicing via __get_item__ 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 (str | Variable) – Variable to read/process from the MeshSeries.

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.

Return type:

ndarray

time_of_min(variable)[source]#

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

Return type:

Mesh

time_of_max(variable)[source]#

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

Return type:

Mesh

aggregate_over_domain(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

plot_domain_aggregate(variable, func, ax=None, **kwargs)[source]#

Plot the transient aggregated data over the domain per timestep.

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

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

  • kwargs (Any) – Keyword args passed to matplotlib’s plot function.

Returns:

A matplotlib Figure or None if plotting on existing axis.

Return type:

Figure | None

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

Probe the MeshSeries at observation points.

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

  • data_name (str) – Name of the data to sample.

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

Returns:

numpy array of interpolated data at observation points.

Return type:

ndarray

plot_probe(points, variable, variable_abscissa=None, labels=None, interp_method='linear', colors=None, linestyles=None, ax=None, fill_between=False, **kwargs)[source]#

Plot the transient variable on the observation points in the MeshSeries.

param points:

The points to sample at.

param variable:

The variable to be sampled.

param labels:

The labels for each observation point.

param interp_method:

Choose the interpolation method, defaults to linear for xdmf MeshSeries and probefilter for pvd MeshSeries.

param interp_backend:

Interpolation backend for PVD MeshSeries.

param kwargs:

Keyword arguments passed to matplotlib’s plot function.

returns:

A matplotlib Figure

Return type:

Figure | None

animate(variable, timevalues=None, plot_func=lambda *_: ..., **kwargs)[source]#

Create an animation for a variable with given timevalues.

Parameters:
  • variable (Variable) – the field to be visualized on all timevalues

  • timevalues (Sequence | None) – the timevalues to animate

  • plot_func (Callable[[Axes, float], None]) – A function which expects to read a matplotlib Axes and the time value of the current frame. Useful to customize the plot in the animation.

Keyword Arguments:

See ogstools.plot.contourf

Return type:

FuncAnimation

plot_time_slice(variable, points, y_axis='auto', interpolate=True, time_logscale=False, fig=None, ax=None, cbar=True, **kwargs)[source]#
Parameters:
  • variable (Variable | str) – The variable to be visualized.

  • points (ndarray) – The points along which the data is sampled over time.

  • y_axis (Literal['x', 'y', 'z', 'dist', 'auto']) – The component of the sampling points which labels the y-axis. By default, if only one coordinate of the points is changing, this axis is taken, otherwise the distance along the line is taken.

  • interpolate (bool) – Smoothen the result be interpolation.

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

  • cb_loc – Colorbar location. If None, omit colorbar.

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

Return type:

Figure

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

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[[Mesh], Mesh]) – 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”. If given as tuple of strings, the latter units will also be set in ot.plot.setup.spatial_unt and ot.plot.setup.time_unit for plotting.

Parameters:
  • spatial (float | tuple[str, str]) – Float factor or a tuple of str (from_unit, to_unit).

  • time (float | tuple[str, str]) – Float factor or a tuple of str (from_unit, to_unit).

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)[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.

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.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) – The first list/array of meshes to be subtracted from.

  • meshes_2 (list | ndarray) – 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.meshes_from_gmsh(filename, dim=0, reindex=True, log=True)[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 | list[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

Returns:

A dictionary of names and corresponding meshes

Return type:

dict[str, Mesh]

ogstools.meshlib.read_shape(shapefile, simplify=False, mesh_generator='triangle', cellsize=None)[source]#

Generate a pyvista Unstructured grid from a shapefile.

Parameters:
  • shapefile (str | Path) – Shapefile to be meshed.

  • simplify (bool) – With the Douglas-Peucker algorithm the geometry is simplified. The original line is split into smaller parts. All points with a distance smaller than half the cellsize are removed. Endpoints are preserved. More infos at https://geopandas.org/en/stable/docs/reference/api/geopandas.GeoSeries.simplify.html.

  • mesh_generator (str) – Choose between ‘triangle’ and ‘gmsh’ to generate the mesh.

  • cellsize (int | None) – Size of the cells in the mesh - only needed for simplify algorithm. If None - cellsize is 1/100 of larger bound (x or y).

Returns:

pv.UnstructuredGrid

Return type:

UnstructuredGrid

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]#
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#