pygeon.discretizations.fem.hdiv module
Module for the discretizations of the H(div) space.
- class pygeon.discretizations.fem.hdiv.RT0(keyword='unitary_data')[source]
Bases:
DiscretizationDiscretization class for Raviart-Thomas of lowest order. Each degree of freedom is the integral over a mesh face.
The implementation of this class is inspired by the RT0 class in PorePy.
- poly_order = 1
Polynomial degree of the basis functions
- tensor_order = 1
Vector-valued discretization
- ndof(sd)[source]
Returns the number of faces.
- Parameters:
sd (pg.Grid) – Grid, or a subclass.
- Returns:
The number of degrees of freedom.
- Return type:
- assemble_mass_matrix(sd, data=None)[source]
Assembles the mass matrix
- Parameters:
sd (pg.Grid) – Grid object or a subclass.
data (dict | None) – Optional dictionary with physical parameters for scaling, in particular the pg.SECOND_ORDER_TENSOR that is the inverse of the diffusion tensor (permeability for porous media).
- Returns:
The mass matrix.
- Return type:
sps.csc_array
- static local_inner_product(sd)[source]
Compute the local inner product matrix for a given grid.
- Parameters:
sd (pg.Grid) – The grid object containing the discretization information.
- Returns:
Local inner product matrix.
- Return type:
np.ndarray
- static eval_basis(coord, sign, dim)[source]
Evaluate the basis functions.
- Parameters:
coord (np.ndarray) – The coordinates of the opposite node for each face.
sign (np.ndarray) – The sign associated to each of the face of the degree of freedom.
dim (int) – The dimension of the grid.
- Returns:
The value of the basis functions.
- Return type:
np.ndarray
- eval_at_cell_centers(sd)[source]
Evaluate the finite element solution at the cell centers of the given grid.
- Parameters:
sd (pg.Grid) – The grid on which to evaluate the solution.
- Returns:
The finite element solution evaluated at the cell centers.
- Return type:
sps.csc_array
- assemble_lumped_matrix(sd, data=None)[source]
Assembles the lumped mass matrix L such that B^T L^{-1} B is a TPFA method.
- Parameters:
sd (pg.Grid) – Grid object or a subclass.
data (dict | None) – Optional dictionary with physical parameters for scaling. In particular the pg.SECOND_ORDER_TENSOR that is the inverse of the diffusion tensor (permeability for porous media).
- Returns:
The lumped mass matrix.
- Return type:
sps.csc_array
- assemble_diff_matrix(sd)[source]
Assembles the matrix corresponding to the differential operator, the divergence in this case.
- Parameters:
sd (pg.Grid) – Grid object or a subclass.
- Returns:
The differential matrix.
- Return type:
sps.csc_array
- interpolate(sd, func)[source]
Interpolates a function onto the finite element space
- Parameters:
sd (pg.Grid) – Grid, or a subclass.
func (Callable[[np.ndarray], np.ndarray]) – A function that returns the function values at coordinates.
- Returns:
The values of the degrees of freedom.
- Return type:
np.ndarray
- assemble_nat_bc(sd, func, b_faces)[source]
Assembles the natural boundary condition term (n dot q, func)_Gamma
- Parameters:
sd (pg.Grid) – The grid object representing the computational domain.
func (Callable[[np.ndarray], np.ndarray]) – The function that defines the natural boundary condition.
b_faces (np.ndarray) – The array of boundary faces.
- Returns:
The assembled natural boundary condition term.
- Return type:
np.ndarray
- get_range_discr_class(dim)[source]
Returns the range discretization class for the given dimension.
- Parameters:
dim (int) – The dimension of the range space.
- Returns:
The range discretization class.
- Return type:
pg.Discretization
- proj_to_PwPolynomials(sd)[source]
Constructs the projection matrix from the current finite element space to the VecPwLinears space.
- Parameters:
sd (pg.Grid) – The grid object.
- Returns:
A sparse array in CSC format representing the projection from the current space to VecPwLinears.
- Return type:
sps.csc_array
- class pygeon.discretizations.fem.hdiv.BDM1(keyword='unitary_data')[source]
Bases:
DiscretizationBDM1 is a class that represents the BDM1 (Brezzi-Douglas-Marini) finite element method. It provides methods for assembling matrices, projecting to and from the RT0 space, evaluating the solution at cell centers, interpolating a given function onto the grid, assembling the natural boundary condition term, and more.
- poly_order = 1
Polynomial degree of the basis functions
- tensor_order = 1
Vector-valued discretization
- ndof(sd)[source]
Return the number of degrees of freedom associated to the method. In this case the number of faces times the dimension.
- Parameters:
sd (pp.Grid) – Grid object or a subclass.
- Returns:
The number of degrees of freedom.
- Return type:
- Raises:
ValueError – If the input grid is not an instance of pp.Grid.
- local_dofs_of_cell(sd, faces_loc)[source]
Compute the local degrees of freedom (DOFs) indices for a cell.
- Parameters:
sd (pp.Grid) – Grid object or a subclass.
faces_loc (np.ndarray) – Array of local face indices for the cell.
- Returns:
Array of local DOF indices associated with the cell.
- Return type:
np.ndarray
- assemble_mass_matrix(sd, data=None)[source]
Assembles the mass matrix for the given grid.
- Parameters:
sd (pg.Grid) – The grid for which the mass matrix is assembled.
data (dict | None) – Additional data for the assembly process.
- Returns:
The assembled mass matrix.
- Return type:
sps.csc_array
- eval_basis_at_node(sd, opposites, faces_loc, return_node_ind=False)[source]
Compute the local basis function for the BDM1 finite element space.
- Parameters:
sd (pg.Grid) – The grid object.
opposites (np.ndarray) – The local degrees of freedom.
cell_nodes_loc (np.ndarray) – The local nodes of the cell.
faces_loc (np.ndarray) – The local faces.
return_node_ind (bool) – Whether to return the local indexing of the nodes, used in assemble_lumped_matrix.
- Returns:
The local mass matrix.
- Return type:
np.ndarray
- static local_inner_product(dim)[source]
Compute the local inner product matrix for the given dimension.
- Parameters:
dim (int) – The dimension of the matrix.
- Returns:
The computed local inner product matrix.
- Return type:
np.ndarray
- proj_to_RT0(sd)[source]
Project the function space to the lowest order Raviart-Thomas (RT0) space.
- Parameters:
sd (pg.Grid) – The grid object representing the computational domain.
- Returns:
The projection matrix to the RT0 space.
- Return type:
sps.csc_array
- proj_from_RT0(sd)[source]
Project the RT0 finite element space onto the faces of the given grid.
- Parameters:
sd (pg.Grid) – The grid on which the projection is performed.
- Returns:
The projection matrix.
- Return type:
sps.csc_array
- assemble_diff_matrix(sd)[source]
Assembles the matrix corresponding to the differential operator.
- Parameters:
sd (pg.Grid) – Grid object or a subclass.
- Returns:
The differential matrix.
- Return type:
sps.csc_array
- eval_at_cell_centers(sd)[source]
Evaluate the finite element solution at the cell centers of the given grid.
- Parameters:
sd (pg.Grid) – The grid on which to evaluate the solution.
- Returns:
The finite element solution evaluated at the cell centers.
- Return type:
sps.csc_array
- interpolate(sd, func)[source]
Interpolates a given function onto the grid.
- Parameters:
sd (pg.Grid) – The grid on which to interpolate the function.
func (Callable[[np.ndarray], np.ndarray]) – The function to be interpolated.
- Returns:
The interpolated values on the grid.
- Return type:
np.ndarray
- assemble_nat_bc(sd, func, b_faces)[source]
Assembles the natural boundary condition term (n dot q, func)_Gamma
- Parameters:
sd (pg.Grid) – The grid object representing the computational domain.
func (Callable[[np.ndarray], np.ndarray]) – The function that defines the natural boundary condition.
b_faces (np.ndarray) – The array of boundary faces.
- Returns:
The assembled natural boundary condition term.
- Return type:
np.ndarray
- get_range_discr_class(dim)[source]
Returns the range discretization class for the given dimension.
- Parameters:
dim (int) – The dimension of the range space.
- Returns:
The range discretization class.
- Return type:
pg.Discretization
- assemble_lumped_matrix(sd, data=None)[source]
Assembles the lumped matrix for the given grid.
- Parameters:
sd (pg.Grid) – The grid object.
data (dict | None) – Optional data dictionary.
- Returns:
The assembled lumped matrix.
- Return type:
sps.csc_array
- proj_to_PwPolynomials(sd)[source]
Constructs the projection matrix from the current finite element space to the VecPwLinears space.
- Parameters:
sd (pg.Grid) – The grid object.
- Returns:
A sparse array in CSC format representing the projection from the current space to VecPwLinears.
- Return type:
sps.csc_array
- class pygeon.discretizations.fem.hdiv.RT1(keyword='unitary_data')[source]
Bases:
DiscretizationRT1 Discretization class for H(div) finite element method.
This class implements the Raviart-Thomas elements of order 1 (RT1) for discretizing vector fields in H(div) space. It provides methods for assembling mass matrices, differential matrices, evaluating basis functions, and interpolating functions onto the finite element space.
- poly_order = 2
Polynomial degree of the basis functions
- tensor_order = 1
Vector-valued discretization
- ndof(sd)[source]
Returns the number of degrees of freedom.
- Parameters:
sd (pg.Grid) – Grid, or a subclass.
- Returns:
The number of degrees of freedom.
- Return type:
- local_dofs_of_cell(sd, faces_loc, c)[source]
Compute the local degrees of freedom (DOFs) indices for a cell.
- Parameters:
sd (pp.Grid) – Grid object or a subclass.
faces_loc (np.ndarray) – Array of local face indices for the cell.
c (int) – Cell index.
- Returns:
Array of local DOF indices associated with the cell.
- Return type:
np.ndarray
- assemble_mass_matrix(sd, data=None)[source]
Assembles the mass matrix
- Parameters:
sd (pg.Grid) – Grid object or a subclass.
data (dict | None) – Optional dictionary with physical parameters for scaling, in particular the pg.SECOND_ORDER_TENSOR that is the inverse of the diffusion tensor (permeability for porous media).
- Returns:
The mass matrix.
- Return type:
sps.csc_array
- local_inner_product(dim)[source]
Assembles the local inner products based on the Lagrange2 element
- Parameters:
dim (int) – Dimension of the grid.
- Returns:
The local mass matrix.
- Return type:
np.ndarray
- reorder_faces(cell_faces, opposite_nodes, cell)[source]
Reorders the local nodes, faces, and corresponding cell-face orientations
- Parameters:
cell_faces (sps.csc_array) – Cell_face connectivity of the grid.
opposite_nodes (sps.csc_array) – Opposite nodes for each face.
cell (int) – Cell index.
- Returns:
The reordered local node indices np.ndarray: The reordered local face indices np.ndarray: The reordered cell-face orientation signs
- Return type:
np.ndarray
- eval_basis_functions(sd, nodes_loc, signs_loc, volume)[source]
Evaluates the basis functions at the nodes and edges of a cell.
- Parameters:
sd (pg.Grid) – The grid.
nodes_loc (np.ndarray) – Nodes of the cell.
signs_loc (np.ndarray) – Cell-face orientation signs.
volume (float) – Cell volume.
- Returns:
An array Psi in which [i, 3j : 3(j + 1)] contains the values of basis function phi_i at evaluation point j
- Return type:
np.ndarray
- eval_basis_functions_at_center(sd, nodes_loc, volume)[source]
Evaluates the basis functions at the center of a cell.
- Parameters:
sd (pg.Grid) – The grid.
nodes_loc (np.ndarray) – Nodes of the cell.
volume (float) – Cell volume.
- Returns:
A (3 x dim) array with the values of the cell-based basis functions at the cell center.
- Return type:
np.ndarray
- eval_at_cell_centers(sd)[source]
Evaluate the finite element solution at the cell centers of the given grid.
- Parameters:
sd (pg.Grid) – The grid on which to evaluate the solution.
- Returns:
The finite element solution evaluated at the cell centers.
- Return type:
sps.csc_array
- assemble_diff_matrix(sd)[source]
Assembles the matrix corresponding to the differential operator, the divergence in this case.
- Parameters:
sd (pg.Grid) – Grid object or a subclass.
- Returns:
The differential matrix.
- Return type:
sps.csc_array
- compute_local_div_matrix(dim)[source]
Assembles the local divergence matrix using local node and face ordering
- Parameters:
dim (int) – Dimension of the grid.
- Returns:
The local divergence matrix
- Return type:
np.ndarray
- interpolate(sd, func)[source]
Interpolates a function onto the finite element space
- Parameters:
sd (pg.Grid) – Grid, or a subclass.
func (Callable[[np.ndarray], np.ndarray]) – A function that returns the function values at coordinates.
- Returns:
The values of the degrees of freedom.
- Return type:
np.ndarray
- assemble_nat_bc(sd, func, b_faces)[source]
Assembles the natural boundary condition term (n dot q, func)_Gamma
- Parameters:
sd (pg.Grid) – The grid object representing the computational domain.
func (Callable[[np.ndarray], np.ndarray]) – The function that defines the natural boundary condition.
b_faces (np.ndarray) – The array of boundary faces.
- Returns:
The assembled natural boundary condition term.
- Return type:
np.ndarray
- get_range_discr_class(dim)[source]
Returns the range discretization class for the given dimension.
- Parameters:
dim (int) – The dimension of the range space.
- Returns:
The range discretization class.
- Return type:
pg.Discretization
- assemble_lumped_matrix(sd, data=None)[source]
Assembles the lumped matrix for the given grid, using the integration rule from Egger & Radu (2020)
- Parameters:
sd (pg.Grid) – The grid object.
data (dict | None) – Optional data dictionary.
- Returns:
The assembled lumped matrix.
- Return type:
sps.csc_array
- proj_to_PwPolynomials(sd)[source]
Constructs the projection matrix from the current finite element space to the VecPwQuadratics space.
- Parameters:
sd (pg.Grid) – The grid object.
- Returns:
A sparse array in CSC format representing the projection from the current space to VecPwQuadratics.
- Return type:
sps.csc_array