pygeon.discretizations.fem.h1 module
Module for the discretizations of the H1 space.
- class pygeon.discretizations.fem.h1.Lagrange1(keyword='unitary_data')[source]
Bases:
DiscretizationClass representing the Lagrange1 finite element discretization.
- poly_order = 1
Polynomial degree of the basis functions
- tensor_order = 0
Scalar-valued discretization
- ndof(sd)[source]
Returns the number of degrees of freedom associated to the method. In this case, the number of nodes.
- Parameters:
sd – Grid, or a subclass.
- Returns:
The number of degrees of freedom.
- Return type:
ndof
- assemble_mass_matrix(sd, data=None)[source]
Returns the mass matrix for the lowest order Lagrange element
- Parameters:
sd (pg.Grid) – The grid.
data (dict | None) – Optional data for the assembly process.
- Returns:
The mass matrix obtained from the discretization.
- Return type:
sps.csc_array
- assemble_local_mass(dim)[source]
Compute the local mass matrix on an element with measure 1.
- Parameters:
dim (int) – Dimension of the matrix.
- Returns:
Local mass matrix of shape (num_nodes_of_cell, num_nodes_of_cell).
- Return type:
np.ndarray
- assemble_stiff_matrix(sd, data=None)[source]
Assembles the stiffness matrix for the finite element method.
- Parameters:
sd (pg.Grid) – The grid object representing the discretization.
data (dict) – A dictionary containing the necessary data for assembling the matrix.
- Returns:
The assembled stiffness matrix.
- Return type:
sps.csc_array
- assemble_diff_matrix(sd)[source]
Assembles the differential matrix based on the dimension of the grid.
- Parameters:
sd (pg.Grid) – The grid object.
- Returns:
The differential matrix.
- Return type:
sps.csc_array
- local_stiff(K, c_volume, coord, dim)[source]
Compute the local stiffness matrix for P1.
- Parameters:
K (np.ndarray) – Permeability of the cell of (dim, dim) shape.
c_volume (np.ndarray) – Scalar cell volume.
coord (np.ndarray) – Coordinates of the cell vertices of (dim+1, dim) shape.
dim (int) – Dimension of the problem.
- Returns:
Local stiffness matrix of (dim+1, dim+1) shape.
- Return type:
np.ndarray
- static local_grads(coord, dim)[source]
Calculate the local gradients of the finite element basis functions.
- Parameters:
coord (np.ndarray) – The coordinates of the nodes in the element.
dim (int) – The dimension of the element.
- Returns:
The local gradients of the finite element basis functions.
- Return type:
np.ndarray
- assemble_lumped_matrix(sd, data=None)[source]
Assembles the lumped mass matrix for the finite element method.
- Parameters:
sd (pg.Grid) – The grid object representing the discretization.
data (dict | None) – Optional data dictionary.
- Returns:
The assembled lumped mass matrix.
- Return type:
sps.csc_array
- proj_to_PwPolynomials(sd)[source]
Construct the matrix for projecting a Lagrangian function to a piecewise linear function.
- Parameters:
sd (pg.Grid) – The grid on which to construct the matrix.
- Returns:
The matrix representing the projection.
- Return type:
sps.csc_array
- eval_at_cell_centers(sd)[source]
Construct the matrix for evaluating a Lagrangian function at the cell centers of the given grid.
- Parameters:
sd (pg.Grid) – The grid on which to construct the matrix.
- Returns:
The matrix representing the projection at the cell centers.
- Return type:
sps.csc_array
- interpolate(sd, func)[source]
Interpolates a given function over the nodes of a 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:
An array containing the interpolated values at each node of the grid.
- Return type:
np.ndarray
- assemble_nat_bc(sd, func, b_faces)[source]
Assembles the ‘natural’ boundary condition (u, func)_Gamma with u a test function in Lagrange1
- Parameters:
sd (pg.Grid) – The grid object representing the computational domain.
func (Callable[[np.ndarray], np.ndarray]) – The function used to evaluate the ‘natural’ boundary condition.
b_faces (np.ndarray) – The array of boundary faces.
- Returns:
The assembled ‘natural’ boundary condition values.
- Return type:
np.ndarray
- get_range_discr_class(dim)[source]
Returns the appropriate range discretization class based on the dimension.
- Parameters:
dim (int) – The dimension of the problem.
- Returns:
The range discretization class.
- Return type:
- Raises:
NotImplementedError – If there’s no zero discretization in PyGeoN.
- class pygeon.discretizations.fem.h1.Lagrange2(keyword='unitary_data')[source]
Bases:
DiscretizationClass representing the Lagrange2 finite element discretization.
- poly_order = 2
Polynomial degree of the basis functions
- tensor_order = 0
Scalar-valued discretization
- ndof(sd)[source]
Returns the number of degrees of freedom associated to the method. In this case, the number of nodes plus the number of edges, where edges are one-dimensional mesh entities.
- Parameters:
sd – Grid, or a subclass.
- Returns:
The number of degrees of freedom.
- Return type:
ndof
- assemble_mass_matrix(sd, data=None)[source]
Returns the mass matrix for the second order Lagrange element
- Parameters:
sd (pg.Grid) – The grid.
data (dict | None) – Optional data for the assembly process.
- Returns:
The mass matrix.
- Return type:
sps.csc_array
- assemble_local_mass(dim)[source]
Computes the local mass matrix of the basis functions on a d-simplex with measure 1.
- Parameters:
dim (int) – The dimension of the simplex.
- Returns:
The local mass matrix.
- Return type:
np.ndarray
- assemble_barycentric_mass(expnts)[source]
Compute the inner products of all monomials up to degree 2
- Parameters:
expnts (np.ndarray) – Each column is an array of exponents alpha_i of the monomial expressed as prod_i lambda_i ^ alpha_i.
- Returns:
The inner products of the monomials on a simplex with measure 1.
- Return type:
np.ndarray
- integrate_monomial(alphas)[source]
Exact integration of products of monomials based on Vermolen and Segal (2018).
- Parameters:
alphas (np.ndarray) – Array of exponents alpha_i of the monomial expressed as prod_i lambda_i ^ alpha_i.
- Returns:
The integral of the monomial on a simplex with measure 1.
- Return type:
- num_edges_per_cell(dim)[source]
Compute the number of edges of a simplex of a given dimension.
- get_local_edge_nodes(dim)[source]
Lists the local edge-node connectivity in the cell
- Parameters:
dim (int) – Dimension.
- Returns:
Row i contains the local indices of the nodes connected to the edge with local index i.
- Return type:
np.ndarray
- eval_grads_at_nodes(dphi, e_nodes)[source]
Evaluates the gradients of the basis functions at the nodes
- Parameters:
dphi (np.ndarray) – Gradients of the P1 basis functions.
e_nodes (np.ndarray) – The local edge-node connectivity.
- Returns:
The gradient of basis function i at node j is in elements [i, 3 * (j:J + 1)].
- Return type:
np.ndarray
- get_edge_dof_indices(sd, cell, faces)[source]
Finds the indices for the edge degrees of freedom that correspond to the local numbering of the edges.
- Parameters:
sd (pg.Grid) – The grid.
cell (int) – The cell index.
faces (np.ndarray) – Face indices of the cell.
- Returns:
Indices of the edge degrees of freedom.
- Return type:
np.ndarray
- assemble_stiff_matrix(sd, data=None)[source]
Assembles the stiffness matrix for the P2 finite element method.
- Parameters:
sd (pg.Grid) – The grid object representing the discretization.
data (dict) – A dictionary containing the necessary data for assembling the matrix.
- Returns:
The stiffness matrix.
- Return type:
sps.csc_array
- assemble_diff_matrix(sd)[source]
Assembles the differential matrix based on the dimension of the grid.
- Parameters:
sd (pg.Grid) – The grid object.
- Returns:
The differential matrix.
- Return type:
sps.csc_array
- eval_at_cell_centers(sd)[source]
Construct the matrix for evaluating a P2 function at the cell centers of the given grid.
- Parameters:
sd (pg.Grid) – The grid on which to construct the matrix.
- Returns:
The matrix representing the projection at the cell centers.
- Return type:
sps.csc_array
- assemble_lumped_matrix(sd, data=None)[source]
Assembles the lumped mass matrix for the quadratic Lagrange space. This is based on the integration rule by Eggers and Radu, and is not block-diagonal for this space.
- Parameters:
sd (pg.Grid) – The grid object representing the discretization.
data (dict | None) – A dictionary containing the necessary data for assembling the matrix.
- Returns:
The lumped mass matrix.
- Return type:
sps.csc_array
- proj_to_PwPolynomials(sd)[source]
Construct the matrix for projecting a quadratic Lagrangian function to a piecewise quadratic function.
- Parameters:
sd (pg.Grid) – The grid on which to construct the matrix.
- Returns:
The matrix representing the projection.
- Return type:
sps.csc_array
- interpolate(sd, func)[source]
Interpolates a given function over the nodes of a 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:
An array containing the interpolated values at each node of the grid.
- Return type:
np.ndarray
- assemble_nat_bc(sd, func, b_faces)[source]
Assembles the ‘natural’ boundary condition (func, u)_Gamma with u a test function in Lagrange2
- Parameters:
sd (pg.Grid) – The grid object representing the computational domain.
func (Callable[[np.ndarray], np.ndarray]) – The function used to evaluate the ‘natural’ boundary condition.
b_faces (np.ndarray) – The array of boundary faces.
- Returns:
The assembled ‘natural’ boundary condition values.
- Return type:
np.ndarray
- get_range_discr_class(dim)[source]
Returns the appropriate range discretization class based on the dimension.
- Parameters:
dim (int) – The dimension of the problem.
- Returns:
The range discretization class.
- Return type:
- Raises:
NotImplementedError – There is no zero-dimensional discretization in PyGeoN.