ogstools.meshes package#

class ogstools.meshes.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, domain_key='domain')[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.

  • domain_key (str) – String which is only in the domain filepath

Return type:

Self

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:

Self

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

Self

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

Create Meshes by extract boundaries of domain mesh.

The provided mesh must be already conforming to OGS standards.

Parameters:
  • mesh (UnstructuredGrid) – The domain mesh

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

Self

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(include_domain=False)[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’}

modify_names(prefix='', suffix='')[source]#

Add prefix and/or suffix to names of meshes. Separators (underscore, etc.) have to be provided by user as part of prefix and/or suffix parameters.

Parameters:
  • prefix (str) – Prefix to be added at the beginning of the mesh name

  • suffix (str) – Suffix to be added after the mesh name

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 (Path | str | 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, **kwargs)[source]#

Save all meshes.

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

Parameters:
  • meshes_path (Path | str | 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

remove_material(mat_id, tolerance=1e-12)[source]#

Remove material from meshes and update integration point data.

Parameters:
  • mat_id (int | Sequence[int]) – MaterialID/s to be removed from domain, subdomain elements, which only belonged to this material are also removed. If given as a sequence, then it must be of length 2 and all ids in between are removed.

  • tolerance (float) – Absolute distance threshold to check if subdomain nodes still have a corresponding domain node after removal of the designated material.

Submodules#