qugar.mesh
Mesh module for QUGaR
Functions
Creates the permutation array to map from DOLFINx to lexicographical edges ordering. |
|
Creates the permutation array to map from DOLFINx to lexicographical faces ordering. |
|
|
Creates the permutation array to map from DOLFINx to lexicographical nodes ordering. |
|
Creates the permutation array to map from DOLFINx to VTK nodes ordering. |
|
Creates the permutation array to map from VTK to DOLFINx nodes ordering. |
|
Creates the permutation array to map from VTK to lexicographical nodes ordering. |
|
Creates the permutation array to map from VTK to lexicographical faces ordering. |
Creates the permutation array to map from lexicographical to DOLFINx edges ordering. |
|
Creates the permutation array to map from lexicographical to DOLFINx faces ordering. |
|
|
Creates the permutation array to map from lexicographical to VTK faces ordering. |
|
Creates the permutation array to map from lexicographical to DOLFINx nodes ordering. |
|
Creates the permutation array to map from lexicographical to VTK nodes ordering. |
|
Creates a Cartesian mesh from a bounding box and the number of cells per direction. |
|
Creates an uniftted Cartesian mesh generated with a bounding box and the number of cells per direction, and an implicit function describing the domain |
Classes
|
Enriched mesh data structure. |
|
Class for storing an unfitted domains and access its cut, full, and empty cells and facets, and create custom quadratures associated to those cells and facets. |
|
Cartesian mesh data structure. |
|
Tensor-product mesh data structure. |
Helper class for accessing indices in arbitrary dimensions tensor-product structures. |
|
|
Class for storing an unfitted Cartesian mesh domain and access its cut, full, and empty cells and facets. |
- class qugar.mesh.CartesianMesh(comm: ~mpi4py.MPI.Comm, cart_grid_tp_cpp: ~qugar.cpp.CartGridTP_2D | ~qugar.cpp.CartGridTP_3D, degree: int = 1, active_cells: ~numpy.ndarray[~typing.Any, ~numpy.dtype[~numpy.int64]] | None = None, ghost_mode: ~dolfinx.cpp.mesh.GhostMode = GhostMode.none, 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.Comm) – MPI communicator to be used for distributing the mesh.
cart_grid_tp_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.
active_cells (Optional[npt.NDArray[np.int64]]) – Array storing the active cells of the Cartesian mesh cart_grid_tp_cpp. If not set, all the cells are considered active.
ghost_mode (GhostMode, optional) – Ghost mode used for mesh partitioning. Defaults to none.
dtype (type[np.float32 | np.float64], optional) – Type to be used in the grid. Defaults to np.float64.
- property cart_grid_tp_cpp_object: CartGridTP_2D
Returns the stored (C++) Cartesian grid object.
- Returns:
Stored Cartesian grid object.
- Return type:
- 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]]
- get_cell_bbox(dlf_local_cell_id: int32) ndarray[Any, dtype[float32 | float64]] [source]
Computes the bounding box of a single cell.
- Parameters:
dlf_local_cell_id (np.int32) – Local DOLFINx cell whose bounding box is created.
- 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.Mesh(comm: Comm, nodes_coords: ndarray[Any, dtype[float32 | float64]], conn: ndarray[Any, dtype[int64]], cell_type: CellType, degree: int = 1, active_cells: ndarray[Any, dtype[int64]] | None = None, ghost_mode: GhostMode = GhostMode.none, merge_nodes: bool = False, merge_tol: type[float32 | float64] | None = None)[source]
Bases:
Mesh
Enriched mesh data structure.
This class derives from dolfinx.mesh.Mesh, easing the management of of some mesh quantities.
For instance, it 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.
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.
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.get_original_node_ids or self.get_original_cell_ids. Or, from the original numbering used to create the mesh to the local DOLFINx numbering using self.get_DOLFINx_local_node_ids or self.get_DOLFINx_local_cell_ids, or the global with self.get_DOLFINx_global_node_ids or self.get_DOLFINx_global_cell_ids`,
Constructor.
- Parameters:
comm (MPI.Comm) – MPI communicator to be used for distributing the mesh.
nodes_coords (npt.NDArray[np.float32 | np.float64]) – Nodes coordinates. The rows correspond to the different nodes and the columns to the coordinates of each point.
conn (npt.NDArray[np.int64]) – Connectivity of the mesh.
cell_type (CellType) – Type of the mesh cell.
degree (int) – Degree of the mesh.
active_cells (Optional[npt.NDArray[np.int64]]) – Array of active cells. If not set, all the cells are considered active. The cells are referred to the original numbering used to create the mesh.
ghost_mode (GhostMode, optional) – Ghost mode used for mesh partitioning. Defaults to none.
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.
- property degree: int
Gets the mesh’ degree.
- Returns:
Mesh’s degree.
- Return type:
int
- property dtype: dtype[float32] | dtype[float64]
Gets the type associated to the breaks.
- Returns:
Type associated to the breaks.
- Return type:
np.dtype[np.float32] | np.dtype[np.float64]
- property gdim: int
Gets the geometrical dimension of the mesh.
- Returns:
Mesh’s geometrical dimension.
- Return type:
int
- get_DOLFINx_global_cell_ids(orig_cell_ids: ndarray[Any, dtype[int64]]) ndarray[Any, dtype[int64]] [source]
Transforms the given local orig_cell_ids into the corresponding global ids of the underlying DOLFINx mesh.
- Parameters:
orig_cell_ids (npt.NDArray[np.int64]) – Cell indices to be transformed. They are the original indices used to create the mesh.
Note
All the indices in orig_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 indices belonging to the current subdomain (process).
- Return type:
npt.NDArray[np.int64]
- get_DOLFINx_local_cell_ids(orig_cell_ids: ndarray[Any, dtype[int64]]) ndarray[Any, dtype[int32]] [source]
Transforms the given orig_cell_ids from the original numbering into the corresponding local ids of the underlying DOLFINx mesh.
- Parameters:
orig_cell_ids (npt.NDArray[np.int64]) – Cell indices to be transformed. They are the original indices used to create the mesh.
- 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 associated input original indices are not contained in the subdomain. The length of the output array is the same as the input.
- Return type:
npt.NDArray[np.int32]
- get_DOLFINx_local_node_ids(orig_node_ids: ndarray[Any, dtype[int64]]) ndarray[Any, dtype[int32]] [source]
Transforms the given orig_node_ids from the original numbering into the corresponding local ids of the underlying DOLFINx mesh.
- Parameters:
orig_node_ids (npt.NDArray[np.int64]) – Node indices to be transformed. They are referred to the original numbering used to create the mesh.
- 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 original 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]
- get_DOLFINx_local_to_global_cell_ids(dlf_local_cell_ids: ndarray[Any, dtype[int32]]) ndarray[Any, dtype[int64]] [source]
Transforms the given local dlf_local_cell_ids (DOLFINx) local cell ids into the corresponding global ones associated to the current process.
- Parameters:
dlf_local_cell_ids (npt.NDArray[np.int32]) – DOLFINx local cell indices to be transformed.
- Returns:
Global indices of the DOLFINx cells.
- Return type:
npt.NDArray[np.int64]
- get_all_original_cell_ids() ndarray[Any, dtype[int64]] [source]
Retrieves the ids of all the cells from the mesh using the original numbering.
- Returns:
An array containing the original cell IDs.
- Return type:
npt.nDArray[np.int64]
- get_cell_facets(dlf_local_cell_id: int32) ndarray[Any, dtype[int32]] [source]
Extracts the facet ids of the given cell.
- Parameters:
dlf_local_cell_id (np.int32) – Id of the queried cell. It is a DOLFINx (local) cell id.
- Returns:
Array of DOLFINx cell facets.
- Return type:
npt.NDArray[np.int32]
- get_original_cell_ids(dlf_local_cell_ids: ndarray[Any, dtype[int32]]) ndarray[Any, dtype[int64]] [source]
Transforms given dlf_local_cell_ids into the original numbering used to create the mesh.
- Parameters:
dlf_local_cell_ids (npt.NDArray[np.int32]) – Cell indices to be transformed. They arelocal indices referred to the underlying DOLFINx mesh.
- Returns:
Ids of the cells in the original numbering. The length of the output array is the same as the input.
- Return type:
npt.NDArray[np.int64]
- get_original_node_ids(dlf_local_node_ids: ndarray[Any, dtype[int32]]) ndarray[Any, dtype[int64]] [source]
Transforms given dlf_local_node_ids from the (local) DOLFINx mesh numbering to original numbering used to create the mesh.
- Parameters:
dlf_local_node_ids (npt.NDArray[np.int32]) – Node indices to be transformed. They are referred to the underlying (local) DOLFINx mesh.
- Returns:
Ids of the nodes in the original numbering.
- Return type:
npt.NDArray[np.int64]
- property has_inactive_cells: bool
Checks if the mesh has inactive cells.
- Returns:
Whether the mesh has inactive cells.
- Return type:
bool
- 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 nodes.
- Return type:
npt.NDArray[np.int64]
- Type:
Gets the map of original merged nodes
- property num_global_cells: int64
Gets the total (global) number of cells of the DOLFINx mesh.
- Returns:
Number of global cells of the DOLFINx mesh.
- 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 original_active_cells: ndarray[Any, dtype[int64]]
Gets the original active cells of the mesh.
- Returns:
Original active cells.
- Return type:
npt.NDArray[np.int64]
- property tdim: int
Gets the topological dimension of the mesh.
- Returns:
Mesh’s topological dimension (1, 2, or 3).
- Return type:
int
- 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: Comm, nodes_coords: ndarray[Any, dtype[float32 | float64]], n_cells: ndarray[Any, dtype[int32]] | list[int32] | list[int], degree: int = 1, active_cells: ndarray[Any, dtype[int64]] | None = None, ghost_mode: GhostMode = GhostMode.none, merge_nodes: bool = False, merge_tol: type[float32 | float64] | None = None)[source]
Bases:
Mesh
Tensor-product mesh data structure.
This class derives from dolfinx.mesh.Mesh, easing 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.get_original_node_ids or self.get_original_cell_ids. Or, from the lexicographical to the local DOLFINx numbering using self.get_DOLFINx_local_node_ids or self.get_DOLFINx_local_cell_ids, or the global with self.get_DOLFINx_global_node_ids or self.get_DOLFINx_global_cell_ids,
Constructor.
- Parameters:
comm (MPI.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.
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.
active_cells (Optional[npt.NDArray[np.int64]]) – Array storing the active cells of the mesh. If not set, all the cells are considered active. Defaults to None.
ghost_mode (GhostMode, optional) – Ghost mode used for mesh partitioning. Defaults to none.
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() 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.
The created tags follow the lexicographical ordering of the tensor-product mesh. I.e., umin:0, umax:1, vmin:2, vmax:3, …
- Returns:
Generated tags.
- Return type:
dolfinx.mesh.MeshTags
- get_cell_tensor_index(dlf_local_cell_id: int32) ndarray[Any, dtype[int32]] [source]
Computes the tensor index referred to the tensor-product mesh of a given (local) DOLFINx cell index.
- Parameters:
dlf_local_cell_id (np.int32) – Local DOLFINx cell index to be transformed.
- Returns:
Cell tensor index
- Return type:
npt.NDArray[np.int32]
- get_face_DOLFINx_cells(tp_face_id: int) ndarray[Any, dtype[int32]] [source]
Extracts the (local) DOLFINx indices of the cells touching the given face of the tensor-product mesh.
- Parameters:
tp_face_id (int) – Id of the face whose cells are extracted. This id follows the lexicographical ordering.
- 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_facets(face_id: int) ndarray[Any, dtype[int32]] [source]
Extracts the facets indices of the given face of the tensor-product mesh.
- Parameters:
face_id (int) – Id of the face whose facets are extracted. This id follows the lexicographical ordering.
- 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]
- 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_cells_tp: int64
Gets the total number of cells of the underlying tensor-product mesh.
- Returns:
Number of cells of the tensor-product mesh.
- 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]
- class qugar.mesh.UnfittedCartMesh(comm: ~mpi4py.MPI.Comm, cpp_unf_domain: ~qugar.cpp.UnfittedDomain_2D | ~qugar.cpp.UnfittedDomain_3D, cpp_cart_grid_tp: ~qugar.cpp.CartGridTP_2D | ~qugar.cpp.CartGridTP_3D, degree: int = 1, exclude_empty_cells: bool = True, ghost_mode: ~dolfinx.cpp.mesh.GhostMode = GhostMode.none, dtype: type[~numpy.float32 | ~numpy.float64] = <class 'numpy.float64'>)[source]
Bases:
CartesianMesh
,UnfittedDomain
Class for storing an unfitted Cartesian mesh domain and access its cut, full, and empty cells and facets.
Constructor.
Note
This constructor is not intended to be called directly, but rather use the function create_Cartesian_mesh.
- Parameters:
comm (MPI.Comm) – MPI communicator to be used for distributing the mesh.
cpp_unf_domain (UnfittedDomain_2D | UnfittedDomain_3D) – C++ unfitted domain object.
cpp_cart_grid_tp (CartGridTP_2D | CartGridTP_3D) – C++ Cartesian grid object.
degree (int, optional) – Degree of the mesh. Defaults to 1. It must be greater than zero.
exclude_empty_cells (bool, optional) – If True, empty cells are excluded from the mesh. Defaults to True.
ghost_mode (GhostMode, optional) – Ghost mode used for mesh partitioning. Defaults to none.
dtype (type[np.float32 | np.float64], optional) – Type to be used in the grid. Defaults to np.float64.
- class qugar.mesh.UnfittedDomain(mesh: Mesh, cpp_unf_domain_object: UnfittedDomain_2D | UnfittedDomain_3D)[source]
Bases:
UnfittedDomainABC
Class for storing an unfitted domains and access its cut, full, and empty cells and facets, and create custom quadratures associated to those cells and facets.
Note
This class is not intended to be used directly, but rather as a class to be inherited by other classes, such as UnfittedCartMesh.
Constructor.
- Parameters:
mesh (Mesh) – Mesh object associated to the unfitted domain.
cpp_unf_domain_object (UnfittedDomain_2D | UnfittedDomain_3D) – Already generated unfitted domain binary object.
- property cpp_unf_domain_object: UnfittedDomain_2D | UnfittedDomain_3D
Gets the internal C++ unfitted domain.
- Returns:
Internal C++ unfitted domain.
- Return type:
- create_quad_custom_cells(degree: int, cells: ndarray[Any, dtype[int32]]) CustomQuad [source]
Returns the custom quadratures for the given cells.
For empty cells, no quadrature is generated and will have 0 points associated to them. For full cells, the quadrature will be the standard one associated to the cell type. While for cut cells a custom quadrature for the cell’s interior will be generated.
Note
This call requires the generation of the quadratures on the fly, what can be potentially expensive.
- Parameters:
degree (int) – Expected degree of exactness of the quadrature to be generated.
cells (npt.NDArray[np.int32]) – Array of DOLFINx cell ids (local to current MPI process) for which quadratures will be generated.
- Returns:
Generated custom quadrature.
- Return type:
CustomQuadInterface
- create_quad_custom_facets(degree: int, facets: MeshFacets, ext_integral: bool) CustomQuadFacet [source]
Returns the custom quadratures for the given facets.
For empty facets, no quadrature is generated and will have 0 points associated to them. For full facets, the quadrature will be the standard one associated to the facet type. While for cut facets a custom quadrature for the facet’s interior will be generated.
Note
This call requires the generation of the quadratures on the fly, what can be potentially expensive.
- Parameters:
degree (int) – Expected degree of exactness of the quadrature to be generated.
facets (MeshFacets) – MeshFacets object containing the DOLFINx (local) facets for which quadratures will be generated.
ext_integral (bool) – Whether exterior integrals are to be computed using the generated quadratures. If True, the quadrature will be generated for the facets that belong to the domain’s boundary (either the mesh’s domain or the unfitted boundary). Otherwise, the quadrature will be generated only for the parts of the interior facets that are inside the domain (but not its boundary).
- Returns:
Generated custom facet quadratures.
- Return type:
CustomQuadFacetInterface
- create_quad_unf_boundaries(degree: int, cells: ndarray[Any, dtype[int32]]) CustomQuadUnfBoundary [source]
Returns the custom quadrature for unfitted boundaries for the given cells.
Note
Some unfitted boundary parts may lay over facets. The quadrature corresponding to those parts will not be included in the quadrature generated by this function, but in the one generated with the method create_quad_custom_facets.
Note
For cells not containing unfitted boundaries, no quadrature is generated and will have 0 points associated to them.
Note
This call requires the generation of the quadratures on the fly, what can be potentially expensive.
- Parameters:
degree (int) – Expected degree of exactness of the quadrature to be generated.
cells (npt.NDArray[np.int32]) – Array of DOLFINx cell ids (local to current MPI process) for which quadratures will be generated.
- Returns:
Generated custom quadrature for unfitted boundaries.
- Return type:
CustomQuadUnfBoundaryInterface
- get_cut_cells() ndarray[Any, dtype[int32]] [source]
Gets the ids of the cut cells.
- Returns:
Array of cut cells associated to the current process following the DOLFINx local numbering. The cell ids are sorted.
- Return type:
npt.NDArray[np.int32]
- get_cut_facets(ext_integral: bool = True) MeshFacets [source]
Gets the cut facets as a MeshFacets object following the DOLFINx local numbering.
The list of facets will be filtered for exterior or interior integrals according to the argument ext_integral.
Note
The selection of facets is performed differently depending on whether exterior or interior integrals are to be computed. For interior integrals we consider interior facets (shared by two cells) that are cut (i.e., they are partially inside the domain). For exterior integrals we consider either interior or exterior (that belong to a single cell) facets that partially belong to the domain’s boundary (either the mesh’s domain or the unfitted boundary). In the case of exterior integrals, if a facet is fully contained in an unfitted boundary, it is considered as full.
- Parameters:
ext_integral (bool) – Whether the list of facets is retrieved for computing exterior or interior integrals (see note above).
- Returns:
Cut facets (following DOLFINx local ordering).
- Return type:
MeshFacets
- get_empty_cells() ndarray[Any, dtype[int32]] [source]
Gets the ids of the empty cells.
- Returns:
Array of empty cells associated to the current process following the DOLFINx local numbering. The cell ids are sorted.
- Return type:
npt.NDArray[np.int32]
- get_empty_facets(ext_integral: bool = True) MeshFacets [source]
Gets the empty facets as a MeshFacets object following the DOLFINx local numbering.
The list of facets will be filtered to only exterior or interior facets according to the argument exterior.
Note
The selection of facets is performed differently depending on whether exterior or interior integrals are to be computed. For interior integrals we consider interior facets (shared by two cells) that are not contained in the domain (they may be contained (fully or partially) in the unfitted boundary). For exterior integrals we consider exterior facets (that belong to a single cell) that are not contained in the domain or its boundary.
- Parameters:
ext_integral (bool) – Whether the list of facets is retrieved for computing exterior or interior integrals (see note above).
- Returns:
Empty facets (following DOLFINx local ordering).
- Return type:
MeshFacets
- get_full_cells() ndarray[Any, dtype[int32]] [source]
Gets the ids of the full cells.
- Returns:
Array of full cells associated to the current process following the DOLFINx local numbering. The cell ids are sorted.
- Return type:
npt.NDArray[np.int32]
- get_full_facets(ext_integral: bool = True) MeshFacets [source]
Gets the full facets as a MeshFacets object following the DOLFINx local numbering.
The list of facets will be filtered for exterior or interior integrals according to the argument ext_integral.
Note
The selection of facets is performed differently depending on whether exterior or interior integrals are to be computed. For interior integrals we consider interior facets (shared by two cells) that are fully inside the domain (and not contained in any unfitted boundary). For exterior integrals we consider facets (interior or exterior) that are fully contained in the domain’s boundary (either the mesh’s domain or the unfitted boundary).
- Parameters:
ext_integral (bool) – Whether the list of facets is retrieved for computing exterior or interior integrals (see note above).
- Returns:
Full facets (following DOLFINx local ordering).
- Return type:
MeshFacets
- static get_num_Gauss_points(degree: int) int [source]
Returns the number of points per direction to be used for a given degree of exactness.
It computes the number of points assuming the Gauss-Legendre quadrature. The degree of exactness is 2 times the number of points minus one.
It does not consider the map of of the quadrature tile. A slightly higher number of points/degree may be needed. See the Appendix of the following paper for more details:
- Antolin, P., Wei, X. and Buffa, A., 2022. Robust numerical integration on curved
polyhedra based on folded decompositions. Computer Methods in Applied Mechanics and Engineering, 395, p.114948.
- Parameters:
degree (int) – Expected degree of exactness of the quadrature.
- Returns:
Number of quadrature points per direction for every quadrature tile.
- 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: Comm, 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.none) 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 none.
- Returns:
Generated Cartesian mesh.
- Return type:
- qugar.mesh.create_unfitted_impl_Cartesian_mesh(comm: Comm, impl_func: ImplicitFunc, n_cells: ndarray[Any, dtype[int32]] | list[int32] | list[int] | int, xmin: ndarray[Any, dtype[float32 | float64]] | None = None, xmax: ndarray[Any, dtype[float32 | float64]] | None = None, degree: int = 1, exclude_empty_cells: bool = True, ghost_mode: GhostMode = GhostMode.none, dtype: type[float32 | float64] | None = None) UnfittedCartMesh [source]
Creates an uniftted Cartesian mesh generated with a bounding box and the number of cells per direction, and an implicit function describing the domain
- Parameters:
comm – MPI communicator to be used for distributing the mesh.
impl_func (ImplicitFunc_2D | ImplicitFunc_3D) – Implicit function that describes the domain.
n_cells (npt.NDArray[np.int32] | list[np.int32] | list[int] | 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.
exclude_empty_cells (bool, optional) – If True, empty cells are excluded from the mesh. Defaults to True.
ghost_mode (GhostMode, optional) – Ghost mode used for mesh partitioning. Defaults to none.
dtype (Optional[type[np.float32 | np.float64]]) – Type to be used in the grid. If not provided, it will be inferred from xmin and/or xmax, and if not provided either it will be set to np.float64.
- Returns:
Generated unfitted Cartesian mesh.
- Return type:
- 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.).