qugar.quad

Quadrature module for QUGaR

Classes

CellState(*values)

CustomQuad(cells, n_pts_per_entity, points, ...)

Data class for storing custom quadratures for a collection of cells.

CustomQuadFacet(facets, n_pts_per_entity, ...)

Data class for storing custom quadratures for a collection of facets.

CustomQuadUnfBoundary(cells, ...)

Data class for storing custom quadratures for unfitted custom boundaries.

QuadGenerator(*args, **kwargs)

Protocol class for generating quadratures for unfitted domains.

class qugar.quad.CellState(*values)[source]

Bases: Enum

cut = 1
empty = 3
full = 2
class qugar.quad.CustomQuad(cells: ndarray[Any, dtype[int32]], n_pts_per_entity: ndarray[Any, dtype[int32]], points: ndarray[Any, dtype[float64]], weights: ndarray[Any, dtype[float64]])[source]

Bases: NamedTuple

Data class for storing custom quadratures for a collection of cells.

The quadratures are defined by custom points and weights referred to the [0,1]^tdim unit domain of every cell, where tdim is the cell’s topological dimension.

All the points (and weights) are stored contiguously for all the cells involved. Thus, the way of differentiating which points belong to which cells is through n_pts_per_entity.

Parameters:
  • cells (npt.NDArray[np.int32]) – Array of cell ids to which the quadrature points and weights are associated to. This is a unique array. The ids are local to the MPI rank.

  • n_pts_per_entity (npt.NDArray[np.int32]) – Array indicating the number of quadrature points and weights for every cell in cells.

  • points (npt.NDArray[np.float64]) – Array of custom quadrature points for all the cells. This is a 2D array that has as many rows as points and as many columns as coordinates. It stores all the points contiguously according the ordering of cells.

  • weights (npt.NDArray[np.float64]) – Array of custom quadrature weights for all the cells. This is a 1D array with one entry per point. It stores all the weights contiguously according the ordering of cells.

Create new instance of CustomQuad(cells, n_pts_per_entity, points, weights)

cells: ndarray[Any, dtype[int32]]

Alias for field number 0

n_pts_per_entity: ndarray[Any, dtype[int32]]

Alias for field number 1

points: ndarray[Any, dtype[float64]]

Alias for field number 2

weights: ndarray[Any, dtype[float64]]

Alias for field number 3

class qugar.quad.CustomQuadFacet(facets: MeshFacets, n_pts_per_entity: ndarray[Any, dtype[int32]], points: ndarray[Any, dtype[float64]], weights: ndarray[Any, dtype[float64]])[source]

Bases: NamedTuple

Data class for storing custom quadratures for a collection of facets.

It contains the same attributes as CustomQuad plus an additional array for identifying the local ids of the facets.

The quadratures are defined by custom points and weights referred to the [0,1]^(tdim-1) unit domain of every facet, where tdim is the cell’s topological dimension, and therefore (tdim-1) is the facet’s topological dimension.

Thus, the points have as many coordinates as the facets of a cell. E.g., for triangles or quadrilaterals, points will have 1 coordinate, while for for tetrahedra or hexahedra they will have 2.

Parameters:
  • facest (MeshFacets) – Mesh facets object that stores the information of the facets. The facets are local to the MPI rank.

  • n_pts_per_entity (npt.NDArray[np.int32]) – Array indicating the number of quadrature points and weights for every facet.

  • points (npt.NDArray[np.float64]) – Array of custom quadrature points for all the cells. This is a 2D array that has as many rows as points and as many columns as coordinates in the facet. It stores all the points contiguously according the ordering of cells and facets.

  • weights (npt.NDArray[np.float64]) – Array of custom quadrature weights for all the facets. This is a 1D array with one entry per point. It stores all the weights contiguously according the ordering of cells and facets.

Create new instance of CustomQuadFacet(facets, n_pts_per_entity, points, weights)

facets: MeshFacets

Alias for field number 0

n_pts_per_entity: ndarray[Any, dtype[int32]]

Alias for field number 1

points: ndarray[Any, dtype[float64]]

Alias for field number 2

weights: ndarray[Any, dtype[float64]]

Alias for field number 3

class qugar.quad.CustomQuadUnfBoundary(cells: ndarray[Any, dtype[int32]], n_pts_per_entity: ndarray[Any, dtype[int32]], points: ndarray[Any, dtype[float64]], weights: ndarray[Any, dtype[float64]], normals: ndarray[Any, dtype[float64]])[source]

Bases: NamedTuple

Data class for storing custom quadratures for unfitted custom boundaries. I.e., boundaries that are not on the exterior boundary of a domain, but inside, as those derived from trimming operations.

It contains the same attributes as CustomQuad plus an additional array for storing the unit outer normals on the boundary.

The normal vectors (and points) have as many coordinates as the cell. E.g., for triangles or quadrilaterals, points and normals will have 2 coordinates, while for tetrahedra or hexahedra they will have 3.

Parameters:
  • cells (npt.NDArray[np.int32]) – Array of cell ids to which the quadrature points and weights are associated to. This is a unique array. The ids are local to the MPI rank.

  • n_pts_per_entity (npt.NDArray[np.int32]) – Array indicating the number of quadrature points and weights for every cell in cells.

  • points (npt.NDArray[np.float64]) – Array of custom quadrature points for all the cells. This is a 2D array that has as many rows as points and as many columns as coordinates. It stores all the points contiguously according the ordering of cells.

  • weights (npt.NDArray[np.float64]) – Array of custom quadrature weights for all the cells. This is a 1D array with one entry per point. It stores all the weights contiguously according the ordering of cells.

  • normals (npt.NDArray[np.float64]) – Array of custom quadrature unit outer normals for all the cells. The normals are referred to the [0,1]^tdim unit domain of every cell. This is a 2D array that has as many rows as quadrature points and as many columns as coordinates. Therefore, the shape of the array must be the same as points.

Create new instance of CustomQuadUnfBoundary(cells, n_pts_per_entity, points, weights, normals)

cells: ndarray[Any, dtype[int32]]

Alias for field number 0

n_pts_per_entity: ndarray[Any, dtype[int32]]

Alias for field number 1

normals: ndarray[Any, dtype[float64]]

Alias for field number 4

points: ndarray[Any, dtype[float64]]

Alias for field number 2

weights: ndarray[Any, dtype[float64]]

Alias for field number 3

class qugar.quad.QuadGenerator(*args, **kwargs)[source]

Bases: Protocol

Protocol class for generating quadratures for unfitted domains.

create_quad_custom_cells(degree: int, cells: ndarray[Any, dtype[int32]], tag: int | None = None) CustomQuad[source]

Returns the custom quadratures for the given cells.

Among the given cells, it only creates quadratures for a subset of those (the custom cells), while no points or weights are generated for the others. Thus, the cells without custom quadratures will be listed in the returned CustomQuadInterface, but will have 0 points associated to them.

Note

This call may require 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.

  • tag (Optional[int]) – Mesh tag of the subdomain associated to the given cells. Defaults to None.

Returns:

Generated custom quadrature.

Return type:

CustomQuadInterface

create_quad_custom_facets(degree: int, facets: MeshFacets, integral_type: str, tag: int | None = None) CustomQuadFacet[source]

Returns the custom quadratures for the given facets.

Among the facets associated to tag, it only creates quadratures for a subset of those (the custom facets), while no points or weights are generated for the others. Thus, the facets without custom quadratures will be listed in the returned CustomQuadFacetInterface, but will have 0 points associated to

Among the given facets, it only creates quadratures for a subset of those (the custom facets), while no points or weights are generated for the others. Thus, the facets without custom quadratures will be listed in the returned CustomQuadFacetInterface, but will have 0 points associated to them.

Note

This call may require 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) – Mesh facets object that stores the information of the facets. The facets are local to the MPI rank.

  • integral_type (str) – Type of integral to be computed. It can be either ‘interior_facet’ or ‘exterior_facet’.

  • tag (int) – Mesh tag of the subdomain associated to the given cells.

Returns:

Generated custom facet quadratures.

Return type:

CustomQuadFacetInterface

create_quad_unf_boundaries(degree: int, cells: ndarray[Any, dtype[int32]], tag: int | None = None) CustomQuadUnfBoundary[source]

Returns the custom quadrature for unfitted boundaries for the given cells.

Warning

All the given cells associated should contain unfitted boundaries. I.e., they must be cut cells. If not, the custom coefficients generator will raise an exception.

Note

This call may require 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. It must only contain cells with unfitted boundaries.

  • tag (int) – Mesh tag of the subdomain associated to the given cells.

Returns:

Generated custom quadrature for unfitted boundaries.

Return type:

CustomQuadUnfBoundaryInterface