qugar.mesh

Mesh module for QUGaR

Functions

DOLFINx_to_lexicg_edges(dim)

Creates the permutation array to map from DOLFINx to lexicographical edges ordering.

DOLFINx_to_lexicg_faces(dim)

Creates the permutation array to map from DOLFINx to lexicographical faces ordering.

DOLFINx_to_lexicg_nodes(dim, degree)

Creates the permutation array to map from DOLFINx to lexicographical nodes ordering.

DOLFINx_to_VTK_nodes(dim, degree)

Creates the permutation array to map from DOLFINx to VTK nodes ordering.

VTK_to_DOLFINx_nodes(dim, degree)

Creates the permutation array to map from VTK to DOLFINx nodes ordering.

VTK_to_lexicg_nodes(dim, degree)

Creates the permutation array to map from VTK to lexicographical nodes ordering.

VTK_to_lexicg_faces(dim)

Creates the permutation array to map from VTK to lexicographical faces ordering.

lexicg_to_DOLFINx_edges(dim)

Creates the permutation array to map from lexicographical to DOLFINx edges ordering.

lexicg_to_DOLFINx_faces(dim)

Creates the permutation array to map from lexicographical to DOLFINx faces ordering.

lexicg_to_VTK_faces(dim)

Creates the permutation array to map from lexicographical to VTK faces ordering.

lexicg_to_DOLFINx_nodes(dim, degree)

Creates the permutation array to map from lexicographical to DOLFINx nodes ordering.

lexicg_to_VTK_nodes(dim, degree)

Creates the permutation array to map from lexicographical to VTK nodes ordering.

create_cells_to_facets_map(mesh)

Creates a map that allows to find the facets ids in a mesh from their cell ids and the local facet ids referred to that cells.

map_cells_and_local_facets_to_facets(mesh, ...)

Given cell ids and local face ids referred to those cells, returns the facet ids.

map_facets_to_cells_and_local_facets(mesh, ...)

Given the ids of facets, returns the cells and local facet ids corresponding to those facets.

create_Cartesian_mesh(comm, n_cells[, xmin, ...])

Creates a Cartesian mesh from a bounding box and the number of cells per direction.

Classes

CartesianMesh(comm, grid_cpp, degree, ...)

Cartesian mesh data structure.

TensorProductMesh(comm, nodes_coords, n_cells)

Tensor-product mesh data structure.

TensorProdIndex()

Helper class for accessing indices in arbitrary dimensions tensor-product structures.

class qugar.mesh.CartesianMesh(comm: ~mpi4py.MPI.Intracomm, grid_cpp: ~qugar.cpp.CartGridTP_2D | ~qugar.cpp.CartGridTP_3D, degree: int = 1, ghost_mode: ~dolfinx.cpp.mesh.GhostMode = GhostMode.shared_facet, dtype: type[~numpy.float32 | ~numpy.float64] = <class 'numpy.float64'>)[source]

Bases: TensorProductMesh

Cartesian mesh data structure. This class is a helper class to ease the construction and management of tensor-product (hypercube) meshes aligned with the Cartesian axes.

It is only implemented for 2D and 3D.

I.e., tensor-product meshes (see TensorProductMesh) whose points are distributed according to the Cartesian axes.

Constructor.

Note

This constructor is not intended to be called directly, but rather use the function create_Cartesian_mesh.

Parameters:
  • comm (MPI.Intracomm) – MPI communicator to be used for distributing the mesh.

  • grid_cpp (qugar.cpp.CartGridTP_2D | qugar.cpp.CartGridTP_3D) – C++ Cartesian grid object.

  • degree (int, optional) – Degree of the mesh. Defaults to 1. It must be greater than zero.

  • ghost_mode (GhostMode, optional) – Ghost mode used for mesh partitioning. Defaults to shared_facet.

  • dtype (type[np.float32 | np.float64], optional) – Type to be used in the grid. Defaults to np.float64.

property cell_breaks: list[ndarray[Any, dtype[float32 | float64]]]

Gets the coordinates of the Cartesian grid cells along each direction

Returns:

Cell coordinates along the parametric directions.

Return type:

list[npt.NDArray[np.float32 | np.float64]]

property cpp_object: CartGridTP_2D

Returns the stored (C++) Cartesian grid object.

Returns:

Stored Cartesian grid object.

Return type:

qugar.cpp.CartGridTP_2D | qugar.cpp.CartGridTP_2D

property dtype: type[float32 | float64]

Gets the type associated to the breaks.

Returns:

Type associated to the breaks.

Return type:

type[np.float32 | np.float64]

get_cell_bbox(cell_id: int32, lexicg: bool = False) ndarray[Any, dtype[float32 | float64]][source]

Computes the bounding box of a single cell.

Parameters:
  • cell_id (np.int32) – Cell whose bounding box is created.

  • lexicg (bool, optional) – If True, the given cell_id is assumed to be referred to the lexicographical mesh numbering convention. Otherwise, if False, cell_id follows the numbering of the underlying DOLFINx mesh.

Returns:

Computed bounding box

of the cell. It is a 2D with two rows and as many columns as coordinates. Both rows represent the minimum and maximum coordintes of the bounding box, respectively.

Return type:

npt.NDArray[np.float32 | np.float64]

qugar.mesh.DOLFINx_to_VTK_nodes(dim: int, degree: int) ndarray[Any, dtype[int32]][source]

Creates the permutation array to map from DOLFINx to VTK nodes ordering.

Parameters:
  • dim (int) – Parametric dimension of the cells (1D, 2D, or 3D).

  • degree (int) – Cell’s degree. It must be greater than 0.

Returns:

Permutation array from DOLFINx to VTK

such that a_dlf[i] = a_vtk[perm_array[i]] or a_vtk[i] = perm_array[a_dlf[i]].

Return type:

npt.NDArray[np.int32]

Note

The VTK numbering is referred to the Arbitrary-order Lagrange Finite Elements as defined in https://www.kitware.com//modeling-arbitrary-order-lagrange-finite-elements-in-the-visualization-toolkit These are the cells with cell ids 68 (1D), 70 (2D), or 72 (3D), and not the linear ones (3, 9, 12, respec.) or quadratic (21, 28, 29, respec.).

qugar.mesh.DOLFINx_to_lexicg_edges(dim: int) ndarray[Any, dtype[int32]][source]

Creates the permutation array to map from DOLFINx to lexicographical edges ordering.

See https://github.com/FEniCS/basix/#supported-elements

Parameters:

dim (int) – Parametric dimension of the cells (1D, 2D, or 3D).

Returns:

Permutation array from DOLFINx to

lexicographical such that edges_lex[i] = edges_dlf[perm_array[i]] or edges_dlf[i] = perm_array[edges_lex[i]].

Return type:

npt.NDArray[np.int32]

qugar.mesh.DOLFINx_to_lexicg_faces(dim: int) ndarray[Any, dtype[int32]][source]

Creates the permutation array to map from DOLFINx to lexicographical faces ordering.

See https://github.com/FEniCS/basix/#supported-elements

Parameters:

dim (int) – Parametric dimension of the cells (1D, 2D, or 3D).

Returns:

Permutation array from DOLFINx to

lexicographical such that faces_lex[i] = faces_dlf[perm_array[i]] or faces_dlf[i] = perm_array[faces_lex[i]].

Return type:

npt.NDArray[np.int32]

qugar.mesh.DOLFINx_to_lexicg_nodes(dim: int, degree: int) ndarray[Any, dtype[int32]][source]

Creates the permutation array to map from DOLFINx to lexicographical nodes ordering.

Parameters:
  • dim (int) – Parametric dimension of the cells (1D, 2D, or 3D).

  • degree (int) – Cell’s degree. It must be greater than 0.

Returns:

Permutation array from DOLFINx to

lexicographical, i.e., such that a_dlf[i] = a_lex[perm_array[i]] or a_lex[i] = perm_array[a_dlf[i]].

Return type:

npt.NDArray[np.int32]

class qugar.mesh.TensorProdIndex[source]

Bases: object

Helper class for accessing indices in arbitrary dimensions tensor-product structures.

The indices in the tensor product structure follow the lexicographical order, I.e., the direction 0 iterates the fastest, i.e., is the inner-most, and the direction dim-1 iterates the slowest, i.e., the outer-most. Where dim is the number of dimensions.

I.e., in 2D this would be

^ y | 8 — 9 —10 — 11 | | | | 4 — 5 — 6 — 7 | | | | 0 — 1 — 2 — 3 –> x

When requested for specific corner, face, or edge indices, the local index of the corresponding entity (corner, face, edge) can be specified in lexicographical ordering or following DOLFINx ordering. See https://github.com/FEniCS/basix/#supported-elements

static get_all_corners(n_inds: ndarray[Any, dtype[int32]]) ndarray[Any, dtype[int64]][source]

Extracts the indices of the 2^dim corners of the hypercube.

Parameters:

n_inds (npt.NDArray[np.int32]) – Number of indices per dimension.

Returns:

Array storing the sorted indices of the corners.

Return type:

npt.NDArray[np.int64]

static get_all_edges(n_inds: ndarray[Any, dtype[int32]]) ndarray[Any, dtype[int64]][source]

Extracts the edge indices of the hypercube.

In 3D the edge indices are those that correspond to the intersection of two faces. In 1D and 2D, they are the same as the face indices.

Warning

This method is not implemented yet for 3D.

Parameters:

n_inds (npt.NDArray[np.int32]) – Number of indices per dimension.

Returns:

Array storing the sorted edge indices.

Return type:

npt.NDArray[np.int64]

static get_all_faces(n_inds: ndarray[Any, dtype[int32]]) ndarray[Any, dtype[int64]][source]

Extracts the face indices of the hypercube.

The face indices are those that are on the boundaries.

Parameters:

n_inds (npt.NDArray[np.int32]) – Number of indices per dimension.

Returns:

Array storing the sorted boundary indices.

Return type:

npt.NDArray[np.int64]

static get_all_internal(n_inds: ndarray[Any, dtype[int32]]) ndarray[Any, dtype[int64]][source]

Extracts the internal indices of the hypercube.

The internal indices are those that are not on the boundaries.

Parameters:

n_inds (npt.NDArray[np.int32]) – Number of indices per dimension.

Returns:

Array storing the sorted internal indices.

Return type:

npt.NDArray[np.int64]

static get_corner(n_inds: ndarray[Any, dtype[int32]], local_corner_id: int | int32) int64[source]

Get the index associated to the given corner.

Note

The local ordering of corners in DOLFINx coincides with the lexicographical ordering. See https://github.com/FEniCS/basix/#supported-elements

Parameters:
  • n_inds (npt.NDArray[np.int32]) – Number of indices per dimension.

  • local_corner_id (int | np.int32) – Corner whose index is extracted. It must be a value in the range [0,2^dim[ following the lexicographical ordering.

Returns:

Index of the corner.

Return type:

np.int64

static get_edge(n_inds: ndarray[Any, dtype[int32]], edge_id: int | int32, lexicg=False) ndarray[Any, dtype[int64]][source]

Get the indices associated to the given edge.

Parameters:
  • n_inds (npt.NDArray[np.int32]) – Number of indices per dimension.

  • edge_id (int | np.int32) – Local edge id whose indices are extracted.

  • lexicg – If True, the given edge_id is referred to the lexicographical ordering, otherwise, to the DOLFINx one. See https://github.com/FEniCS/basix/#supported-elements

Returns:

Indices associated to the edge.

Return type:

npt.NDArray[np.int64]

static get_face(n_inds: ndarray[Any, dtype[int32]], face_id: int | int32, lexicg=False) ndarray[Any, dtype[int64]][source]

Get the indices associated to the given face.

Parameters:
  • n_inds (npt.NDArray[np.int32]) – Number of indices per dimension.

  • face_id (int | np.int32) – Face whose indices is extracted. It must be a value in the range [0, 2 * dim[, where dim is the dimension of the domain.

  • lexicg (bool) – If True, the given face_id is referred to the lexicographical ordering, otherwise, to the DOLFINx one. See https://github.com/FEniCS/basix/#supported-elements

Returns:

Indices associated to the face.

Return type:

npt.NDArray[np.int64]

static get_flat_index(n_inds: ndarray[Any, dtype[int32]], tid: ndarray[Any, dtype[int32]]) int64[source]

Transform a given tensor index to a flat one considering n_inds indices per direction.

Parameters:
  • n_inds (npt.NDArray[np.int32]) – Number of indices per dimension.

  • tid (npt.NDArray[np.int32]) – Tensor indices to be transformed.

Returns:

Computed flat index.

Return type:

np.int64

static get_tensor_index(n_inds: ndarray[Any, dtype[int32]], id: int | int32 | int64) ndarray[Any, dtype[int32]][source]

Transform a given index to a tensor one considering n_inds indices per direction.

Parameters:
  • n_inds (npt.NDArray[np.int32]) – Number of indices per dimension.

  • id (int | np.int32 | np.int64) – Flat index to be transformed.

Returns:

Computed tensor index.

Return type:

npt.NDArray[np.int32]

class qugar.mesh.TensorProductMesh(comm: Intracomm, nodes_coords: ndarray[Any, dtype[float32 | float64]], n_cells: ndarray[Any, dtype[int32]] | list[int32] | list[int], degree: int = 1, ghost_mode: GhostMode = GhostMode.shared_facet, merge_nodes: bool = False, merge_tol: type[float32 | float64] | None = None)[source]

Bases: object

Tensor-product mesh data structure.

This class is a helper class to ease the management of structured tensor-product meshes of intervals (1D), quadrangles (2D), or hexahedra (3D). The geometric dimension can be arbitrarity according to the DOLFINx capabilities.

The class provides methods for easily quering the id of a cell in the reference mesh or its corresponding counterpart referred to a mesh created with DOLFINx.

It also provides functionalities for extracting the facets of a cell, and the facets that belong to the mesh’s boundary, among others.

The mesh can be partitioned among different processes using a MPI communicator passed to the contructor. So, when queried about indices of vertices, facets, cells, etc., the returns are referred to the indices present in the current subdomain of the mesh.

The created mesh follows the denoted as lexicographical convention for numbering cells and nodes. This ordering is such that the parametric direction 0 iterates the fastest, i.e., is the inner-most, and the direction dim-1 iterates the slowest, i.e., the outer-most. The nodes and cells numbering for a 2D case would be:

8 — 9 —10 —11 . — . — . — . | | | | | 3 | 4 | 5 | 4 — 5 — 6 — 7 . — . — . — . | | | | | 0 | 1 | 2 | 0 — 1 — 2 — 3 . — . — . — .

Note that after the mesh creation, DOLFINx renumbers and partitions the mesh, so, this numeration will be different. It is always possible to retrieve the original numbering using the maps self.dolfinx_to_lexicg_nodes or self.dolfinx_to_lexicg_cells. Or, from the lexicographical to the DOLFINx numbering using self.lexicg_to_dolfinx_nodes or self.lexicg_to_dolfinx_cells.

Constructor.

Parameters:
  • comm (MPI.Intracomm) – MPI communicator to be used for distributing the mesh.

  • n_cells (npt.NDArray[np.int32] | list[np.int32] | list[int]) – Number of cells per direction in the mesh.

  • nodes_coords (npt.NDArray[np.float32 | np.float64]) – Nodes coordinates. The rows correspond to the different nodes (following lexicographical ordering), and the columns to the coordinates of each point.

  • degree (int) – Degree of the mesh.

  • ghost_mode (GhostMode, optional) – Ghost mode used for mesh partitioning. Defaults to shared_facet.

  • merge_nodes (bool, optional) – If True, coincident nodes will be merged together into a single one. Otherwise, duplicated nodes will not be merged. Defaults to False.

  • merge_tol (Optional[type[np.float32 | np.float64]]) – Absolute tolerance to be used for seeking coincident nodes. If not set, and if merge_nodes is set to True, this tolerance will be automatically computed in the function merge_coincident_points_in_mesh.

create_facet_bdry_tags(lexicg: bool = False) MeshTags[source]

Creates face facet tags for the tensor product mesh.

It creates the tags for the facets of the underlying DOLFINx mesh, setting a different tag for each boundary face of the domain.

Parameters:

lexicg (bool, optional) – Defines how the tag values for the different boundary facets are set. It set to True, the lexicographical ordering is used (i.e., umin:0, umax:1, vmin:2, vmax:3, …). Otherwise, the faces correspond to the FEniCSx cell faces ordering, but applied to the full hypercube. See https://github.com/FEniCS/basix/#supported-elements It defaults to False.

Returns:

Generated tags.

Return type:

dolfinx.mesh.MeshTags

property degree: int

Gets the mesh’ degree.

Returns:

Mesh’s degree.

Return type:

int

property dolfinx_mesh: Mesh

Gets the underlying DOLFINx mesh.

Returns:

Underlying DOLFINx mesh.

Return type:

dolfinx.mesh.Mesh

property gdim: int

Gets the geometrical dimension of the mesh.

Returns:

Mesh’s geometrical dimension.

Return type:

int

get_DOLFINx_global_cell_ids(cell_ids: ndarray[Any, dtype[int32 | int64]] | int32 | int64, lexicg: bool = True) ndarray[Any, dtype[int32]][source]

Transforms the given local cell_ids into the corresponding global ids of the underlying DOLFINx mesh.

Parameters:
  • cell_ids (npt.NDArray[np.int32 | np.int64] | np.int32 | np.int64) – Cell indices to be transformed. They may be local indices referred to the underlying DOLFINx mesh (if lexicg is set to False) or lexicographical indices (if lexicg is True). If lexicg is set to False, all the local DOLFINx indices must belong to the current subdomain (process).

  • lexicg (bool, optional) – Whether cell_ids follow the tensor-product lexicographical numbering or the DOLFINx one. Defaults to True.

Note

All the indices in cell_ids must be contained in the current subdomain.

Returns:

Indices of the cells in the underlying DOLFINx mesh. These indices correspond to the global numbering associated to the current subdomain (process).

Return type:

npt.NDArray[np.int32]

get_DOLFINx_local_cell_ids(cell_ids: ndarray[Any, dtype[int32 | int64]] | int32 | int64, lexicg: bool = True) ndarray[Any, dtype[int32]] | int32[source]

Transforms the given cell_ids from the tensor-product mesh into the corresponding local ids of the underlying DOLFINx mesh.

Parameters:
  • cell_ids (npt.NDArray[np.int32 | np.int64] | np.int32 | np.int64) – Cell indices to be transformed. They may be local indices referred to the underlying DOLFINx mesh (if lexicg is set to False) or lexicographical indices (if lexicg is True). If lexicg is set to False, all the local DOLFINx indices must belong to the current subdomain (process).

  • lexicg (bool, optional) – Whether cell_ids follow the tensor-product lexicographical numbering or the DOLFINx one. Defaults to True.

Returns:

Indices of the cells in the underlying DOLFINx mesh. These indices correspond to the local numbering associated to the current subdomain (process). Note that some indices may be set to -1 in the case the input indices are lexicographical and some of them are not contained in the subdomain. The length of the output array is the same as the input.

Return type:

npt.NDArray[np.int32] | np.int32

get_DOLFINx_node_ids(node_ids: ndarray[Any, dtype[int32 | int64]] | int32 | int64, lexicg: bool = True) ndarray[Any, dtype[int32]] | int32[source]

Transforms the given node_ids from the tensor-product mesh into the corresponding local ids of the underlying DOLFINx mesh.

Parameters:
  • node_ids (npt.NDArray[np.int32 | np.int64] | np.int32 | np.int64) – Node indices to be transformed. They may be local indices referred to the underlying DOLFINx mesh (if lexicg is set to False) or lexicographical indices (if lexicg is True). If lexicg is set to False, all the local DOLFINx indices must belong to the current subdomain (process).

  • lexicg (bool, optional) – Whether node_ids follow the tensor-product lexicographical numbering or the DOLFINx one. Defaults to True.

Returns:

Indices of the nodes in the underlying DOLFINx mesh. These indices correspond to the local numbering associated to the current subdomain (process). Note that some indices may be set to -1 in the case the input indices are lexicographical and some of them are not contained in the subdomain. The length of the output array is the same as the input.

Return type:

npt.NDArray[np.int32] | np.int32

get_cell_facets(cell_id: int32, out_lexicg: bool = False) ndarray[Any, dtype[int32]][source]

Extracts the facet ids of the given cell.

Parameters:
  • cell_id (np.int32) – Id of the queried cell. It is a DOLFINx (local) cell id. It must be contained in the subdomain.

  • out_lexicg (bool, optional) – Whether the returned facet ids are reordered according to the lexicographical ordering of the cell’s bounding box. I.e., according to the directions umin, umax, vmin, vmax, … Otherwise, they follow the FEniCSx convention. See https://github.com/FEniCS/basix/#supported-elements Defaults to False.

Returns:

Array of cell facets.

Return type:

npt.NDArray[np.int32]

get_cell_tensor_index(cell_id: int64 | int32, lexicg: bool = False) ndarray[Any, dtype[int32]][source]

Computes the tensor index referred to the tensor-product mesh of a given cell index.

Parameters:
  • cell_id (np.int64 | np.int32) – Cell index to be transformed. It must be contained in the subdomain (process).

  • lexicg (bool, optional) – Whether cell_id follows the tensor-product lexicographical numbering or the one associated to the DOLINFx mesh. Defaults to False.

Returns:

Cell tensor index

Return type:

npt.NDArray[np.int32]

get_cells_exterior_facets(cell_ids: ndarray[Any, dtype[int32 | int64]] | int32 | int64, lexicg: bool = False) ndarray[Any, dtype[int32]][source]

Extracts the exterior facets of associated to a list of cells for the tensor-product mesh.

The exterior facets are those that belong to one single cell (they are not at the interface between two cells, but on the mesh boundary).

Parameters:
  • cell_ids (npt.NDArray[np.int32 | np.int64] | np.int32 | np.int64) – Array of cells whose exterior facets are extracted. They must be all contained in the subdomain (process).

  • lexicg (bool, optional) – Whether cell_ids follow the tensor-product lexicographical numbering or the one associated to the DOLINFx mesh. Defaults to False.

Returns:

Sorted unique array of exterior facets present in the current subdomain.

Return type:

npt.NDArray[np.int32]

get_face_DOLFINx_cells(face_id: int, lexicg: bool = False) ndarray[Any, dtype[int32]][source]

Extracts the DOLFINx indices of the cells touching the given face of the tensor-product mesh.

Parameters:
Returns:

Array storing the local cells indices referred to the underlying DOLFINx mesh numbering for the current subdomain (process). This array is unique and sorted. Note that it may be empty if the cells associated to the face do not belong to the current subdomain.

Return type:

npt.NDArray[np.int32]

get_face_DOLFINx_nodes(face_id: int, face_lexicg: bool = False) ndarray[Any, dtype[int32]][source]

Extracts the nodes of a given boundary face that belong to the subdomain of the tensor-product mesh.

Parameters:
  • face_id (int) – Id of the face whose nodes are extracted.

  • face_lexicg (bool, optional) – Whether face_id follows the tensor-product lexicographical numbering or the FEniCSx one. See https://github.com/FEniCS/basix/#supported-elements Defaults to False.

  • nodes_lexicg (bool, optional) – Whether the returned nodes indices must follow the tensor-product lexicographical numbering or the one associated to the underlying DOLINFx mesh. Defaults to False.

Returns:

Sorted unique array of local DOLFINx nodes of the boundary face that belong to the current subdomain (process).

Return type:

npt.NDArray[np.int32]

get_face_facets(face_id: int, lexicg: bool = False) ndarray[Any, dtype[int32]][source]

Extracts the facets indices of the given face of the tensor-product mesh.

Parameters:
Returns:

Array storing the facets indices referred to the underlying DOLFINx mesh numbering. This array is unique and sorted. Note that it may be empty if the facets associated to the face do not belong to the current subdomain.

Return type:

npt.NDArray[np.int32]

get_lexicg_cell_ids(cell_ids: ndarray[Any, dtype[int32 | int64]] | int32 | int64, lexicg: bool = False) ndarray[Any, dtype[int64]] | int64[source]

Transforms given cell_ids from the tensor-product mesh into the tensor-product mesh numbering, that follows a lexicographical indexing.

Parameters:
  • cell_ids (npt.NDArray[np.int32 | np.int64] | np.int32 | np.int64) – Cell indices to be transformed. They may be local indices referred to the underlying DOLFINx mesh (if lexicg is set to False) or lexicographical indices (if lexicg is True). If lexicg is set to False, all the local DOLFINx indices must belong to the current subdomain (process).

  • lexicg (bool, optional) – Whether cell_ids follow the tensor-product lexicographical numbering or the DOLFINx one. Defaults to True.

Returns:

Ids of the cells in the tensor-product lexicographical indexing. Note that some indices may be set to -1 if their associated cells do not belong to current subdomain. The length of the output array is the same as the input.

Return type:

npt.NDArray[np.int64] | np.int64

get_lexicg_node_ids(node_ids: ndarray[Any, dtype[int32]] | int32, lexicg: bool = False) ndarray[Any, dtype[int64]] | int64[source]

Transforms given node_ids from the tensor-product mesh into the tensor-product mesh numbering, that follows a lexicographical indexing.

Parameters:
  • node_ids (npt.NDArray[np.int32 | np.int64] | np.int32 | np.int64) – Node indices to be transformed. They may be local indices referred to the underlying DOLFINx mesh (if lexicg is set to False) or lexicographical indices (if lexicg is True). If lexicg is set to False, all the local DOLFINx indices must belong to the current subdomain (process).

  • lexicg (bool, optional) – Whether node_ids follow the tensor-product lexicographical numbering or the DOLFINx one. Defaults to True.

Returns:

Ids of the nodes in the tensor-product lexicographical indexing. Note that some indices may be set to -1 if their associated cells do not belong to current subdomain. The length of the output array is the same as the input.

Return type:

npt.NDArray[np.int64] | np.int64

property has_merged_nodes: bool

Checks if the mesh has merged nodes.

Returns:

Whether the mesh has merged nodes.

Return type:

bool

property merged_nodes_map: ndarray[Any, dtype[int64]]

from original ids to unique ones.

Returns:

Map from old to new lexicographical nodes.

Return type:

npt.NDArray[np.int64]

Type:

Gets the map of lexicographical merged nodes

property num_cells_dir: ndarray[Any, dtype[int32]]

Gets the number of cells per direction in the mesh.

Returns:

Number of cells per direction.

Return type:

npt.NDArray[np.int32]

property num_global_cells: int64

Gets the total (global) number of cells.

Returns:

Number of global cells.

Return type:

np.int64

property num_local_cells: int64

Gets the local number of cells (associated to the current subdomain, i.e., MPI process).

Returns:

Number of local cells.

Return type:

np.int64

property num_pts_dir: ndarray[Any, dtype[int32]]

Gets the number of points per direction in the mesh taking into account the number of cells per direction and the mesh’s degree.

Returns:

Number of points per direction.

Return type:

npt.NDArray[np.int32]

property tdim: int

Gets the topological dimension of the mesh.

Returns:

Mesh’s topological dimension (1, 2, or 3).

Return type:

int

qugar.mesh.VTK_to_DOLFINx_nodes(dim: int, degree: int) ndarray[Any, dtype[int32]][source]

Creates the permutation array to map from VTK to DOLFINx nodes ordering.

Parameters:
  • dim (int) – Parametric dimension of the cells (1D, 2D, or 3D).

  • degree (int) – Cell’s degree. It must be greater than 0.

Returns:

Permutation array from VTK to DOLFINx

such that a_vtk[i] = a_dlf[perm_array[i]] or a_dlf[i] = perm_array[a_vtk[i]].

Return type:

npt.NDArray[np.int32]

Note

The VTK numbering is referred to the Arbitrary-order Lagrange Finite Elements as defined in https://www.kitware.com//modeling-arbitrary-order-lagrange-finite-elements-in-the-visualization-toolkit These are the cells with cell ids 68 (1D), 70 (2D), or 72 (3D), and not the linear ones (3, 9, 12, respec.) or quadratic (21, 28, 29, respec.).

qugar.mesh.VTK_to_lexicg_faces(dim: int) ndarray[Any, dtype[int32]][source]

Creates the permutation array to map from VTK to lexicographical faces ordering.

Parameters:

dim (int) – Parametric dimension of the cells (1D, 2D, or 3D).

Returns:

Permutation array from VTK to

lexicographical such that faces_lex[i] = faces_vtk[perm_array[i]] or faces_vtk[i] = perm_array[faces_lex[i]].

Return type:

npt.NDArray[np.int32]

qugar.mesh.VTK_to_lexicg_nodes(dim: int, degree: int) ndarray[Any, dtype[int32]][source]

Creates the permutation array to map from VTK to lexicographical nodes ordering.

Parameters:
  • dim (int) – Parametric dimension of the cells (1D, 2D, or 3D).

  • degree (int) – Cell’s degree. It must be greater than 0.

Returns:

Permutation array from VTK to

lexicographical such that a_vtk[i] = a_lex[perm_array[i]].

Return type:

npt.NDArray[np.int32]

Note

The VTK numbering is referred to the Arbitrary-order Lagrange Finite Elements as defined in https://www.kitware.com//modeling-arbitrary-order-lagrange-finite-elements-in-the-visualization-toolkit These are the cells with cell ids 68 (1D), 70 (2D), or 72 (3D), and not the linear ones (3, 9, 12, respec.) or quadratic (21, 28, 29, respec.).

qugar.mesh.create_Cartesian_mesh(comm: Intracomm, n_cells: ndarray[Any, dtype[int32]] | list[int32] | list[int], xmin: ndarray[Any, dtype[float32 | float64]] | None = None, xmax: ndarray[Any, dtype[float32 | float64]] | None = None, degree: int = 1, ghost_mode: GhostMode = GhostMode.shared_facet) CartesianMesh[source]

Creates a Cartesian mesh from a bounding box and the number of cells per direction.

Parameters:
  • comm – MPI communicator to be used for distributing the mesh.

  • n_cells (npt.NDArray[np.int32] | list[np.int32] | list[int]) – Number of cells per direction in the mesh.

  • xmin (Optional[npt.NDArray[np.float32 | np.float64]]) – Minimum coordinates of the mesh’s bounding box. Defaults to a vector of zeros with double floating precision.

  • xmax (Optional[npt.NDArray[np.float32 | np.float64]]) – Maximum coordinates of the mesh’s bounding box. Defaults to a vector of ones with double floating precision.

  • degree (int, optional) – Degree of the mesh. Defaults to 1.

  • ghost_mode (GhostMode, optional) – Ghost mode used for mesh partitioning. Defaults to shared_facet.

Returns:

Generated Cartesian mesh.

Return type:

CartesianMesh

qugar.mesh.create_cells_to_facets_map(mesh: Mesh) ndarray[Any, dtype[int32]][source]

Creates a map that allows to find the facets ids in a mesh from their cell ids and the local facet ids referred to that cells.

Parameters:

mesh (dolfinx.mesh.Mesh) – Mesh from which facet ids are extracted.

Returns:

Map from the cells and local facet

ids to the facet ids. It is a 2D array where the first column corresponds to the cells and the second one to the local facet ids.

Thus, given the cell and local facet ids of a particular facet, the facet id can be accessed as facet_id = facets_map[cell_id, local_facet_id], where facets_map is the generated 2D array.

Return type:

npt.NDArray[np.int32]

qugar.mesh.lexicg_to_DOLFINx_edges(dim: int) ndarray[Any, dtype[int32]][source]

Creates the permutation array to map from lexicographical to DOLFINx edges ordering.

See https://github.com/FEniCS/basix/#supported-elements

Parameters:

dim (int) – Parametric dimension of the cells (1D, 2D, or 3D).

Returns:

Permutation array from lexicographical

to DOLFINx such that edges_dlf[i] = edges_lex[perm_array[i]] or edges_lex[i] = perm_array[edges_dlf[i]].

Return type:

npt.NDArray[np.int32]

qugar.mesh.lexicg_to_DOLFINx_faces(dim: int) ndarray[Any, dtype[int32]][source]

Creates the permutation array to map from lexicographical to DOLFINx faces ordering.

See https://github.com/FEniCS/basix/#supported-elements

Parameters:

dim (int) – Parametric dimension of the cells (1D, 2D, or 3D).

Returns:

Permutation array from lexicographical

to DOLFINx such that faces_dlf[i] = faces_lex[perm_array[i]] or faces_lex[i] = perm_array[faces_dlf[i]].

Return type:

npt.NDArray[np.int32]

qugar.mesh.lexicg_to_DOLFINx_nodes(dim: int, degree: int) ndarray[Any, dtype[int32]][source]

Creates the permutation array to map from lexicographical to DOLFINx nodes ordering.

Parameters:
  • dim (int) – Parametric dimension of the cells (1D, 2D, or 3D).

  • degree (int) – Cell’s degree. It must be greater than 0.

Returns:

Permutation array from lexicographical to

DOLFINx, i.e., such that a_lex[i] = a_dlf[perm_array[i]] or a_dlf[i] = perm_array[a_lex[i]].

Return type:

npt.NDArray[np.int32]

qugar.mesh.lexicg_to_VTK_faces(dim: int) ndarray[Any, dtype[int32]][source]

Creates the permutation array to map from lexicographical to VTK faces ordering.

Parameters:

dim (int) – Parametric dimension of the cells (1D, 2D, or 3D).

Returns:

Permutation array from lexicographical

to VTK such that faces_vtk[i] = faces_lex[perm_array[i]] or faces_lex[i] = perm_array[faces_vtk[i]].

Return type:

npt.NDArray[np.int32]

qugar.mesh.lexicg_to_VTK_nodes(dim: int, degree: int) ndarray[Any, dtype[int32]][source]

Creates the permutation array to map from lexicographical to VTK nodes ordering.

Parameters:
  • dim (int) – Parametric dimension of the cells (1D, 2D, or 3D).

  • degree (int) – Cell’s degree. It must be greater than 0.

Returns:

Permutation array from lexicographical to

VTK such that a_lex[i] = a_vtk[perm_array[i]] or a_vtk[i] = perm_array[a_lex[i]].

Return type:

npt.NDArray[np.int32]

Note

The VTK numbering is referred to the Arbitrary-order Lagrange Finite Elements as defined in https://www.kitware.com//modeling-arbitrary-order-lagrange-finite-elements-in-the-visualization-toolkit These are the cells with cell ids 68 (1D), 70 (2D), or 72 (3D), and not the linear ones (3, 9, 12, respec.) or quadratic (21, 28, 29, respec.).

qugar.mesh.map_cells_and_local_facets_to_facets(mesh: Mesh, cells: ndarray[Any, dtype[int32]], local_facets: ndarray[Any, dtype[int32]]) ndarray[Any, dtype[int32]][source]

Given cell ids and local face ids referred to those cells, returns the facet ids.

Parameters:
  • mesh (dolfinx.mesh.Mesh) – Mesh to which the facets belong to.

  • cells (npt.NDArray[np.int32]) – Array of DOLFINx cell ids (local to the current process) to be transformed.

  • local_facets (npt.NDArray[np.int32]) – Array of local facets referred to the cells. They follow the FEniCSx convention. See https://github.com/FEniCS/basix/#supported-elements

Returns:

Computed DOLFINx facet ids (local to the current process).

Return type:

npt.NDArray[np.int32]

qugar.mesh.map_facets_to_cells_and_local_facets(mesh: Mesh, facets: ndarray[Any, dtype[int32]]) tuple[ndarray[Any, dtype[int32]], ndarray[Any, dtype[int32]]][source]

Given the ids of facets, returns the cells and local facet ids corresponding to those facets.

Note

Interior facets belong to more whan one cell, in those cases only one cell (and local facet) is returned for that particular facet. The one chosen depends on the way in which that information is stored in the mesh connectivity.

Parameters:
  • mesh (dolfinx.mesh.Mesh) – Mesh to which the facets belong to.

  • facets (npt.NDArray[np.int32]) – Array of DOLFINx facets (local to the current process) to be transformed.

Returns:

Cells

and local facet ids the associated to the facets. The first entry of the tuple corresponds to the cells and the second to the the local facets.

The local facets indices follow the FEniCSx convention. See https://github.com/FEniCS/basix/#supported-elements

Return type:

tuple[npt.NDArray[np.int32], npt.NDArray[np.int32]]