feec.global_projectors#
Functions#
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
Classes#
|
Projects callable functions to some scalar or vector FEM space. |
|
Projector from H1 to an H1-conforming finite element space (i.e. a finite dimensional subspace of H1) constructed with tensor-product B-splines in 1, 2 or 3 dimensions. |
|
Projector from H1^3 = H1 x H1 x H1 to a conforming finite element space, i.e. a finite dimensional subspace of H1^3, constructed with tensor-product B-splines in 2 or 3 dimensions. |
|
Projector from H(curl) to an H(curl)-conforming finite element space, i.e. a finite dimensional subspace of H(curl), constructed with tensor-product B- and M-splines in 2 or 3 dimensions. |
|
Projector from H(div) to an H(div)-conforming finite element space, i.e. a finite dimensional subspace of H(div), constructed with tensor-product B- and M-splines in 2 or 3 dimensions. |
|
Projector from L2 to an L2-conforming finite element space (i.e. a finite dimensional subspace of L2) constructed with tensor-product M-splines in 1, 2 or 3 dimensions. |
Details#
- class GlobalProjector(space, nquads=None)[source]#
Bases:
object
Projects callable functions to some scalar or vector FEM space.
A global projector is constructed over a tensor-product grid in the logical domain. The vertices of this grid are obtained as the tensor product of the 1D splines’ Greville points along each direction.
This projector matches the “geometric” degrees of freedom of discrete n-forms (where n depends on the underlying space). This is done by projecting each component of the vector field independently, by combining 1D histopolation with 1D interpolation.
This class cannot be instantiated directly (use a subclass instead).
- Parameters:
- spaceVectorFemSpace | TensorFemSpace
Some finite element space, codomain of the projection operator. The exact structure where to use histopolation and where interpolation has to be given by a subclass of the GlobalProjector class. As of now, it is implicitly assumed for a VectorFemSpace, that for each direction that all spaces with interpolation are the same, and all spaces with histopolation are the same (i.e. yield the same quadrature/interpolation points etc.); so use with care on an arbitrary VectorFemSpace. It is right now only intended to be used with VectorFemSpaces or TensorFemSpaces from DeRham complex objects.
- nquadslist(int) | tuple(int)
Number of quadrature points along each direction, to be used in Gauss quadrature rule for computing the (approximated) degrees of freedom. This parameter is ignored, if the projector only uses interpolation (and no histopolation).
- property space#
The space to which this Projector projects.
- property dim#
The dimension of the underlying TensorFemSpaces.
- property blockcount#
The number of blocks. In case that self.space is a TensorFemSpace, this is 1, otherwise it denotes the number of blocks in the VectorFemSpace.
- property grid_x#
The local interpolation/histopolation grids which are used; it denotes the position of the interpolation/quadrature points. All the grids are stored inside a two-dimensional array; the outer dimension denotes the block, the inner the tensor space direction.
- property grid_w#
The local interpolation/histopolation grids which are used; it denotes the weights of the quadrature points (in the case of interpolation, this will return the weight 1 for the given positions). All the grids are stored inside a two-dimensional array; the outer dimension denotes the block, the inner the tensor space direction.
- property func#
The function which is used for projecting a given callable (or list thereof) to the DOFs in the target space.
- property solver#
The solver used for transforming the DOFs in the target space into spline coefficients.
- class Projector_H1(space, nquads=None)[source]#
Bases:
GlobalProjector
Projector from H1 to an H1-conforming finite element space (i.e. a finite dimensional subspace of H1) constructed with tensor-product B-splines in 1, 2 or 3 dimensions.
This is a global projector based on interpolation over a tensor-product grid in the logical domain. The interpolation grid is the tensor product of the 1D splines’ Greville points along each direction.
- Parameters:
- H1SplineSpace or TensorFemSpace
H1-conforming finite element space, codomain of the projection operator
- class Projector_Hcurl(space, nquads=None)[source]#
Bases:
GlobalProjector
Projector from H(curl) to an H(curl)-conforming finite element space, i.e. a finite dimensional subspace of H(curl), constructed with tensor-product B- and M-splines in 2 or 3 dimensions.
This is a global projector constructed over a tensor-product grid in the logical domain. The vertices of this grid are obtained as the tensor product of the 1D splines’ Greville points along each direction.
The H(curl) projector matches the “geometric” degrees of freedom of discrete 1-forms, which are the line integrals of a vector field along cell edges. To achieve this, each component of the vector field is projected independently, by combining 1D histopolation along the direction of the edges with 1D interpolation along the other directions.
- Parameters:
- HcurlVectorFemSpace
H(curl)-conforming finite element space, codomain of the projection operator.
- nquadslist(int) | tuple(int)
Number of quadrature points along each direction, to be used in Gauss quadrature rule for computing the (approximated) degrees of freedom.
- class Projector_Hdiv(space, nquads=None)[source]#
Bases:
GlobalProjector
Projector from H(div) to an H(div)-conforming finite element space, i.e. a finite dimensional subspace of H(div), constructed with tensor-product B- and M-splines in 2 or 3 dimensions.
This is a global projector constructed over a tensor-product grid in the logical domain. The vertices of this grid are obtained as the tensor product of the 1D splines’ Greville points along each direction.
The H(div) projector matches the “geometric” degrees of freedom of discrete (N-1)-forms in N dimensions, which are the integrated flux of a vector field through cell faces (in 3D) or cell edges (in 2D).
To achieve this, each component of the vector field is projected independently, by combining histopolation along the direction(s) tangential to the face (in 3D) or edge (in 2D), with 1D interpolation along the normal direction.
- Parameters:
- HdivVectorFemSpace
H(div)-conforming finite element space, codomain of the projection operator.
- nquadslist(int) | tuple(int)
Number of quadrature points along each direction, to be used in Gauss quadrature rule for computing the (approximated) degrees of freedom.
- class Projector_L2(space, nquads=None)[source]#
Bases:
GlobalProjector
Projector from L2 to an L2-conforming finite element space (i.e. a finite dimensional subspace of L2) constructed with tensor-product M-splines in 1, 2 or 3 dimensions.
This is a global projector constructed over a tensor-product grid in the logical domain. The vertices of this grid are obtained as the tensor product of the 1D splines’ Greville points along each direction.
The L2 projector matches the “geometric” degrees of freedom of discrete N-forms in N dimensions, which are line/surface/volume integrals of a scalar field over an edge/face/cell in 1/2/3 dimension(s). To this end histopolation is used along each direction.
- Parameters:
- L2SplineSpace
L2-conforming finite element space, codomain of the projection operator
- nquadslist(int) | tuple(int)
Number of quadrature points along each direction, to be used in Gauss quadrature rule for computing the (approximated) degrees of freedom.
- evaluate_dofs_2d_1form_hcurl(intp_x1, intp_x2, quad_x1, quad_x2, quad_w1, quad_w2, F1, F2, f1, f2)[source]#
- evaluate_dofs_2d_1form_hdiv(intp_x1, intp_x2, quad_x1, quad_x2, quad_w1, quad_w2, F1, F2, f1, f2)[source]#
- evaluate_dofs_3d_1form(intp_x1, intp_x2, intp_x3, quad_x1, quad_x2, quad_x3, quad_w1, quad_w2, quad_w3, F1, F2, F3, f1, f2, f3)[source]#