fem.projectors#

Functions#

get_dual_dofs(Vh, f, domain_h[, ...])

return the dual dofs tilde_sigma_i(f) = < Lambda_i, f >_{L2} i = 1, .

knot_insertion_projection_operator(domain, ...)

Compute the projection operator based on the knot insertion technique.

knots_to_insert(coarse_grid, fine_grid[, tol])

Compute the point difference between the fine grid and coarse grid.

Classes#

Inheritance diagram of psydac.fem.projectors

DirichletProjector(fem_space, *[, bcs, ...])

A LinearOperator that applies homogeneous (unless manually given different bcs) Dirichlet boundary conditions.

MultipatchDirichletProjector(fem_space, *[, ...])

A LinearOperator (for multipatch domains) that applies homogeneous (unless manually given different bcs) Dirichlet boundary conditions.

Details#

knots_to_insert(coarse_grid, fine_grid, tol=1e-14)[source]#

Compute the point difference between the fine grid and coarse grid.

knot_insertion_projection_operator(domain, codomain)[source]#

Compute the projection operator based on the knot insertion technique.

Return a linear operator which projects an element of the domain to an element of the codomain. Domain and codomain are scalar spline spaces over a cuboid, built as the tensor product of 1D spline spaces. In particular, domain and codomain have the same multi-degree (p1, p2, …).

This function returns a LinearOperator K working at the level of the spline coefficients, which are represented by StencilVector objects.

Thanks to the tensor-product structure of the spline spaces, the projection operator is the Kronecker product of 1D projection operators K[i] operating between 1D spaces. Each 1D operators is represented by a dense matrix:

K = K[0] x K[1] x …

For each dimension i the 1D grids defined by the breakpoints of the two spaces are assumed to be identical, or one nested into the other. Let nd[i] and nc[i] be the number of cells along dimension i for domain and codomain, respectively. We then have three different cases:

  1. nd[i] == nc[i]: The two 1D grids are assumed identical, and K[i] is the identity matrix.

  2. nd[i] < nc[i]: The 1D grid of the domain is assumed nested into the 1D grid of the codomain, hence the 1D spline space of the domain is a subspace of the 1D spline space of the codomain. In this case we build K[i] using the knot insertion algorithm.

  3. nd[i] > nc[i]: The 1D grid of the codomain is assumed nested into the 1D grid of the domain, hence the 1D spline space of the codomain is a subspace of the 1D spline space of the domain. In this case we build K[i] as the transpose of the matrix obtained using the knot insertion algorithm from the codomain to the domain.

Parameters:
domainTensorFemSpace

Domain of the projection operator.

codomainTensorFemSpace

Codomain of the projection operator.

Returns:
KroneckerDenseMatrix

Matrix representation of the projection operator. This is a LinearOperator acting on the spline coefficients.

get_dual_dofs(Vh, f, domain_h, backend_language='python', return_format='stencil_array')[source]#

return the dual dofs tilde_sigma_i(f) = < Lambda_i, f >_{L2} i = 1, .. dim(Vh)) of a given function f, as a stencil array or numpy array

Parameters:
VhFemSpace

The discrete space for the dual dofs

f<sympy.Expr>

The function used for evaluation

domain_h

The discrete domain corresponding to Vh

backend_language: <str>

The backend used to accelerate the code

return_format: <str>

The format of the dofs, can be ‘stencil_array’ or ‘numpy_array’

Returns:
tilde_f: <Vector|ndarray>

The dual dofs

class DirichletProjector(fem_space, *, bcs=None, space_kind=None)[source]#

Bases: LinearOperator

A LinearOperator that applies homogeneous (unless manually given different bcs) Dirichlet boundary conditions.

Parameters:
fem_spacepsydac.fem.basic.FemSpace

fem_space.coeff_space is domain and codomain of this LO. fem_space.kind.name determines the BCs that get applied.

bcsIterable | None

Iterable of sympde.topology.Boundary objects. Allows the user to apply different kinds of BCs. Must not be passed if a space_kind argument is passed.

space_kindstr | SpaceType | None

Necessary only if fem_space.kind.name is undefined and no bcs are passed. Must not be passed if a bcs argument is passed.

Notes

See examples/vector_potential_3d.py for a use case of such an operator.

property domain#

The domain of the linear operator - an element of Vectorspace

property codomain#

The codomain of the linear operator - an element of Vectorspace

property dtype#

The data type of the coefficients of the linear operator, upon convertion to matrix.

tosparse()[source]#

Convert to a sparse matrix in any of the formats supported by scipy.sparse.

toarray()[source]#

Convert to Numpy 2D array.

dot(v, out=None)[source]#

Apply the LinearOperator self to the Vector v.

The result is written to the Vector out, if provided.

Parameters:
vVector

The vector to which the linear operator (self) is applied. It must belong to the domain of self.

outVector

The vector in which the result of the operation is stored. It must belong to the codomain of self. If out is None, a new vector is created and returned.

Returns:
Vector

The result of the operation. If out is None, a new vector is returned. Otherwise, the result is stored in out and out is returned.

transpose(conjugate=False)[source]#

Transpose the LinearOperator .

If conjugate is True, return the Hermitian transpose.

property bcs#
class MultipatchDirichletProjector(fem_space, *, bcs=None, space_kind=None)[source]#

Bases: DirichletProjector

A LinearOperator (for multipatch domains) that applies homogeneous (unless manually given different bcs) Dirichlet boundary conditions.

Parameters:
fem_spacepsydac.fem.basic.FemSpace

fem_space.coeff_space is domain and codomain of this LO. fem_space.kind.name determines the BCs that get applied.

bcsIterable | None

Iterable of sympde.topology.Boundary objects. Allows the user to apply different kinds of BCs. Must not be passed if a space_kind argument is passed.

space_kindstr | SpaceType | None

Necessary only if fem_space.kind.name is undefined and no bcs are passed. Must not be passed if a bcs argument is passed.

Notes

See examples/vector_potential_3d.py for a use case of such an operator.

dot(v, out=None)[source]#

Apply the LinearOperator self to the Vector v.

The result is written to the Vector out, if provided.

Parameters:
vVector

The vector to which the linear operator (self) is applied. It must belong to the domain of self.

outVector

The vector in which the result of the operation is stored. It must belong to the codomain of self. If out is None, a new vector is created and returned.

Returns:
Vector

The result of the operation. If out is None, a new vector is returned. Otherwise, the result is stored in out and out is returned.