api.feec#

Classes#

Inheritance diagram of psydac.api.feec

DiscreteDeRham(domain_h, *spaces)

A discrete de Rham sequence built over a single-patch geometry.

MultipatchDiscreteDeRham(*, domain_h, spaces)

Represents the discrete de Rham sequence for multipatch domains.

Details#

class DiscreteDeRham(domain_h, *spaces)[source]#

Bases: BasicDiscrete

A discrete de Rham sequence built over a single-patch geometry.

Parameters:
domain_hGeometry

The discretized domain, which is a single-patch geometry.

*spaceslist of FemSpace

The discrete spaces of the de Rham sequence.

Notes

  • This constructor should not be called directly, but rather from the discretize_derham function in psydac.api.discretization.

property dim#

Dimension of the physical and logical domains, which are assumed to be the same.

property domain_h#

Discretized domain.

property spaces#

Spaces of the proper de Rham sequence (excluding Hvec).

property V0#

First space of the de Rham sequence : H1 space

property V1#

Second space of the de Rham sequence : - 1d : L2 space - 2d : either Hdiv or Hcurl space - 3d : Hcurl space

property V2#

Third space of the de Rham sequence : - 2d : L2 space - 3d : Hdiv space

property V3#

Fourth space of the de Rham sequence : L2 space in 3d

property sequence#
property H1vec#

Vector-valued H1 space built as the Cartesian product of N copies of V0, where N is the dimension of the (logical) domain.

property mapping#

The mapping from the logical space to the physical space.

property callable_mapping#

The mapping as a callable.

projectors(*, kind='global', nquads=None)[source]#

Projectors mapping callable functions of the physical coordinates to a corresponding FemField object in the de Rham sequence.

Parameters:
kindstr

Type of the projection : at the moment, only global is accepted and returns geometric commuting projectors based on interpolation/histopolation for the de Rham sequence (GlobalProjector 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.

Returns:
P0, …, Pncallables

Projectors that can be called on any callable function that maps from the physical space to R (scalar case) or R^d (vector case) and returns a FemField belonging to the i-th space of the de Rham sequence

derivatives(kind='femlinop')[source]#
dirichlet_projectors(kind='femlinop')[source]#

Returns operators that apply the correct homogeneous Dirichlet BCs.

Parameters:
kindstr

The kind of the projector, can be ‘femlinop’ or ‘linop’. - ‘femlinop’ returns a psydac FemLinearOperator (default) - ‘linop’ returns a psydac LinearOperator

Returns:
d_projectorstuple

Tuple of <psydac.fem.basic.FemLinearOperator> or <psydac.fem.projectors.DirichletProjector> The Dirichlet projectors of each space and in desired form.

Notes

See examples/vector_potential_3d.py for a use case of these operators in LinearOperator form.

LST_preconditioners(*, M0=None, M1=None, M2=None, M3=None, hom_bc=False)[source]#

LST (Loli, Sangalli, Tani) preconditioners [1] are mass matrix preconditioners of the form pc = D_inv_sqrt @ D_log_sqrt @ M_log_kron_solver @ D_log_sqrt @ D_inv_sqrt, where

D_inv_sqrt is the diagonal matrix of the square roots of the inverse diagonal entries of the mass matrix M, D_log_sqrt is the diagonal matrix of the square roots of the diagonal entries of the mass matrix on the logical domain, M_log_kron_solver is the Kronecker Solver of the mass matrix on the logical domain.

These preconditioners work very well even on complex domains as numerical experiments have shown.

Upon choosing hom_bc=True, preconditioners for the modified mass matrices M{i}_0 are being returned. The preconditioner for the last mass matrix of the sequence remains identical as there are no BCs to take care of. M{i}_0 is a mass matrix of the form M{i}_0 = DP @ M{i} @ DP + (I - DP) where DP and I are the corresponding DirichletProjector and IdentityOperator. See examples/vector_potential_3d.

Parameters:
M0, M1, M2, M3psydac.linalg.stencil.StencilMatrix | psydac.linalg.block.BlockLinearOperator | None

H1, Hcurl/Hdiv, L2 (2D) or H1, Hcurl, Hdiv, L2 mass matrices or None. Returns only preconditioners for passed mass matrices.

hom_bcbool

If True, return LST preconditioners for modified M{i}_0 = DP @ M{i} @ DP + (I - DP) mass matrices (i=0,1 (2D), i=0,1,2 (3D)). The arguments M{i} in that case remain the same (M{i}, not M{i}_0). DP and I are DirichletProjector and IdentityOperator. Default: False.

Returns:
tuple

tuple of psydac.linalg.stencil.StencilMatrix and/or psydac.linalg.block.BlockLinearOperator LST preconditioner(s) for the passed mass matrices.

References

[1] Gabriele Loli, Giancarlo Sangalli, Mattia Tani. “Easy and efficient preconditioning of the isogeometric mass

matrix”. In: Computers & Mathematics with Applications 116 (2022), pp. 245–264

conforming_projectors(kind='femlinop', mom_pres=False, p_moments=-1, hom_bc=False)[source]#

return the conforming projectors of the broken multi-patch space

Parameters:
kind<str>

The kind of the projector, can be ‘femlinop’ or ‘linop’. - ‘femlinop’ returns a psydac FemLinearOperator (default) - ‘linop’ returns a psydac LinearOperator

mom_pres: <bool>

If True, preserve polynomial moments of maximal order in the projection.

p_moments: <int>

Number of polynomial moments to be preserved in the projection. (Gets overwritten if the parameter mom_pres equals True)

hom_bc: <bool>

Apply homogenous boundary conditions if True

Returns:
cP0, cP1, cP2Tuple of <psydac.fem.basic.FemLinearOperator> or <psydac.linalg.basic.LinearOperator>

The conforming projectors of each space and in desired form.

hodge_operator(space=None, dual=False, kind='femlinop', backend_language='python')[source]#

Returns the Hodge operator for the given space and specified kind.

Parameters:
spacestr or None

The space for which to return the Hodge operator, can be ‘V0’, ‘V1’, ‘V2’ or None. If None, returns a tuple with all three Hodge operators.

dualbool

If True, returns the dual Hodge operator.

kind<str>

The kind of the projector, can be ‘femlinop’ or ‘linop’. - ‘femlinop’ returns a psydac FemLinearOperator (default) - ‘linop’ returns a psydac LinearOperator

backend_languagestr

The backend used to accelerate the code, default is ‘python’.

Returns:
The Hodge operator of the space of the specified kind.
H<psydac.fem.basic.FemLinearOperator> or <psydac.linalg.basic.LinearOperator>
hodge_operators(dual=False, kind='femlinop', backend_language='python')[source]#

Returns the Hodge operators for the specified kind.

Parameters:
dualbool

If True, returns the dual Hodge operator.

kind<str>

The kind of the projector, can be ‘femlinop’ or ‘linop’. - ‘femlinop’ returns a psydac FemLinearOperator (default) - ‘linop’ returns a psydac LinearOperator

backend_languagestr

The backend used to accelerate the code, default is ‘python’.

Returns:
The Hodge operators of all spaces and of the specified kind.
class MultipatchDiscreteDeRham(*, domain_h, spaces)[source]#

Bases: DiscreteDeRham

Represents the discrete de Rham sequence for multipatch domains.

It only works when the number of patches>1.

Parameters:
domain_h: <Geometry>

The discrete domain

spaces: <list,tuple>

The discrete spaces that are contained in the de Rham sequence

property H1vec#

Vector-valued H1 space built as the Cartesian product of N copies of V0, where N is the dimension of the (logical) domain.

projectors(*, kind='global', nquads=None)[source]#

This method returns the patch-wise commuting projectors on the broken multi-patch space

Parameters:
kind: <str>

The projectors kind, can be global or local

nquads: <list,tuple>

The number of quadrature points.

Returns:
P0: <MultipatchGeometricProjector>

Patch wise H1 projector

P1: <MultipatchGeometricProjector>

Patch wise Hcurl projector

P2: <MultipatchGeometricProjector>

Patch wise L2 projector

Notes

  • when applied to smooth functions they return conforming fields

  • default ‘global projectors’ correspond to geometric interpolation/histopolation operators on Greville grids

  • here ‘global’ is a patch-level notion, as the interpolation-type problems are solved on each patch independently

dirichlet_projectors(kind='femlinop')[source]#

Returns operators that apply the correct homogeneous Dirichlet BCs.

Parameters:
kindstr

The kind of the projector, can be ‘femlinop’ or ‘linop’. - ‘femlinop’ returns a psydac FemLinearOperator (default) - ‘linop’ returns a psydac LinearOperator

Returns:
d_projectorstuple

Tuple of <psydac.fem.basic.FemLinearOperator> or <psydac.fem.projectors.MultipatchDirichletProjector> The Dirichlet projectors of each space and in desired form.

Notes

See examples/vector_potential_3d.py for a use case of these operators in LinearOperator form.