feec.multipatch.api#

Functions#

discretize(expr, *args, **kwargs)

discretize_derham_multipatch(derham, ...)

Classes#

Inheritance diagram of psydac.feec.multipatch.api

DiscreteDerhamMultipatch(*, mapping, ...[, ...])

Represents the discrete De Rham sequence for multipatch domains.

Details#

class DiscreteDerhamMultipatch(*, mapping, domain_h, spaces, sequence=None)[source]#

Bases: DiscreteDerham

Represents the discrete De Rham sequence for multipatch domains.

It only works when the number of patches>1

Parameters:
mapping: <Mapping>

The mapping of the multipatch domain, the multipatch mapping contains the mapping of each patch

domain_h: <Geometry>

The discrete domain

spaces: <list,tuple>

The discrete spaces that are contained in the De Rham sequence

sequence: <list,tuple>

The space kind of each space in the De Rham sequence

property sequence#
property broken_derivatives_as_operators#
property broken_derivatives_as_matrices#
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: <Multipatch_Projector_H1>

Patch wise H1 projector

P1: <Multipatch_Projector_Hcurl>

Patch wise Hcurl projector

P2: <Multipatch_Projector_L2>

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

conforming_projection(space, hom_bc=False, backend_language='python', load_dir=None)[source]#

return the conforming projectors of the broken multi-patch space

Parameters:
space<str>

The space of the projector

hom_bc: <bool>

Apply homogenous boundary conditions if True

backend_language: <str>

The backend used to accelerate the code

load_dir: <str|None>

Filename for storage in sparse matrix format

Returns:
Cp: <FemLinearOperator>

The conforming projector

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

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

Parameters:
space<str>

The space of the dual dofs

f<sympy.Expr>

The function used for evaluation

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

discretize(expr, *args, **kwargs)[source]#
discretize_derham_multipatch(derham, domain_h, *args, **kwargs)[source]#