feec.multipatch.operators#

Functions#

allocate_interface_matrix(corners, ...)

Allocate the interface matrix for a vertex shared by two patches

get_K0_and_K0_inv(V0h[, uniform_patches])

Compute the change of basis matrices K0 and K0^{-1} in V0h.

get_K1_and_K1_inv(V1h[, uniform_patches])

Compute the change of basis matrices K1 and K1^{-1} in Hcurl space V1h.

get_interface_from_corners(corner1, corner2, ...)

Return the interface between two corners from two different patches that correspond to a single (physical) vertex.

get_patch_index_from_face(domain, face)

Return the patch index of subdomain/boundary

get_row_col_index(corner1, corner2, ...)

Return the row and column index of a corner in the StencilInterfaceMatrix

ortho_proj_Hcurl(EE, V1h, domain_h, M1[, ...])

return orthogonal projection of E on V1h, given M1 the mass matrix

Classes#

Inheritance diagram of psydac.feec.multipatch.operators

BrokenGradient_2D(V0h, V1h)

BrokenScalarCurl_2D(V1h, V2h)

BrokenTransposedGradient_2D(V0h, V1h)

BrokenTransposedScalarCurl_2D(V1h, V2h)

ConformingProjection_V0(V0h, domain_h[, ...])

Conforming projection from global broken V0 space to conforming global V0 space Defined by averaging of interface dofs

ConformingProjection_V1(V1h, domain_h[, ...])

Conforming projection from global broken V1 space to conforming V1 global space

HodgeOperator(Vh, domain_h[, metric, ...])

Change of basis operator: dual basis -> primal basis

Multipatch_Projector_H1(V0h)

to apply the H1 projection (2D) on every patch

Multipatch_Projector_Hcurl(V1h[, nquads])

to apply the Hcurl projection (2D) on every patch

Multipatch_Projector_L2(V2h[, nquads])

to apply the L2 projection (2D) on every patch

Details#

get_patch_index_from_face(domain, face)[source]#

Return the patch index of subdomain/boundary

Parameters:
domain<Sympde.topology.Domain>

The Symbolic domain

face<Sympde.topology.BasicDomain>

A patch or a boundary of a patch

Returns:
i<int>

The index of a subdomain/boundary in the multipatch domain

get_interface_from_corners(corner1, corner2, domain)[source]#

Return the interface between two corners from two different patches that correspond to a single (physical) vertex.

Parameters:
corner1<Sympde.topology.Corner>

The first corner of the 2D interface

corner2<Sympde.topology.Corner>

The second corner of the 2D interface

domain<Sympde.topology.Domain>

The Symbolic domain

Returns:
interface: <Sympde.topology.Interface|None>

The interface between two vertices

get_row_col_index(corner1, corner2, interface, axis, V1, V2)[source]#
Return the row and column index of a corner in the StencilInterfaceMatrix

for dofs of H1 type spaces

Parameters:
corner1<Sympde.topology.Corner>

The first corner of the 2D interface

corner2<Sympde.topology.Corner>

The second corner of the 2D interface

interface<Sympde.topology.Interface|None>

The interface between the two corners

axis<int|None>

Axis of the interface

V1<FemSpace>

Test Space

V2<FemSpace>

Trial Space

Returns:
index: <list>

The StencilInterfaceMatrix index of the corner, it has the form (i1, i2, k1, k2) in 2D, where (i1, i2) identifies the row and (k1, k2) the diagonal.

allocate_interface_matrix(corners, test_space, trial_space)[source]#

Allocate the interface matrix for a vertex shared by two patches

Parameters:
corners: <list>

The patch corners corresponding to the common shared vertex

test_space: <FemSpace>

The test space

trial_space: <FemSpace>

The trial space

Returns:
mat: <StencilInterfaceMatrix>

The interface matrix shared by two patches

class ConformingProjection_V0(V0h, domain_h, hom_bc=False, backend_language='python', storage_fn=None)[source]#

Bases: FemLinearOperator

Conforming projection from global broken V0 space to conforming global V0 space Defined by averaging of interface dofs

Parameters:
V0h: <FemSpace>

The discrete space

domain_h: <Geometry>

The discrete domain of the projector

hom_bc<bool>

Apply homogenous boundary conditions if True

backend_language: <str>

The backend used to accelerate the code

storage_fn:

filename to store/load the operator sparse matrix

set_homogenous_bc(boundary, rhs=None)[source]#
class ConformingProjection_V1(V1h, domain_h, hom_bc=False, backend_language='python', storage_fn=None)[source]#

Bases: FemLinearOperator

Conforming projection from global broken V1 space to conforming V1 global space

proj.dot(v) returns the conforming projection of v, computed by solving linear system

Parameters:
V1h: <FemSpace>

The discrete space

domain_h: <Geometry>

The discrete domain of the projector

hom_bc<bool>

Apply homogenous boundary conditions if True

backend_language: <str>

The backend used to accelerate the code

storage_fn:

filename to store/load the operator sparse matrix

set_homogenous_bc(boundary)[source]#
get_K0_and_K0_inv(V0h, uniform_patches=False)[source]#

Compute the change of basis matrices K0 and K0^{-1} in V0h.

With K0_ij = sigma^0_i(B_j) = B_jx(n_ix) * B_jy(n_iy) where sigma_i is the geometric (interpolation) dof and B_j is the tensor-product B-spline

get_K1_and_K1_inv(V1h, uniform_patches=False)[source]#

Compute the change of basis matrices K1 and K1^{-1} in Hcurl space V1h.

With K1_ij = sigma^1_i(B_j) = int_{e_ix}(M_jx) * B_jy(n_iy) if i = horizontal edge [e_ix, n_iy] and j = (M_jx o B_jy) x-oriented MoB spline or = B_jx(n_ix) * int_{e_iy}(M_jy) if i = vertical edge [n_ix, e_iy] and j = (B_jx o M_jy) y-oriented BoM spline (above, ‘o’ denotes tensor-product for functions)

class HodgeOperator(Vh, domain_h, metric='identity', backend_language='python', load_dir=None, load_space_index='')[source]#

Bases: FemLinearOperator

Change of basis operator: dual basis -> primal basis

self._matrix: matrix of the primal Hodge = this is the mass matrix ! self.dual_Hodge_matrix: this is the INVERSE mass matrix

Parameters:
Vh: <FemSpace>

The discrete space

domain_h: <Geometry>

The discrete domain of the projector

metric<str>

the metric of the de Rham complex

backend_language: <str>

The backend used to accelerate the code

load_dir: <str>

storage files for the primal and dual Hodge sparse matrice

load_space_index: <str>

the space index in the derham sequence

Notes

Either we use a storage, or these matrices are only computed on demand # todo: we compute the sparse matrix when to_sparse_matrix is called – but never the stencil matrix (should be fixed…) We only support the identity metric, this implies that the dual Hodge is the inverse of the primal one. # todo: allow for non-identity metrics

to_sparse_matrix()[source]#

the Hodge matrix is the patch-wise multi-patch mass matrix it is not stored by default but assembled on demand

assemble_primal_Hodge_matrix()[source]#

the Hodge matrix is the patch-wise multi-patch mass matrix it is not stored by default but assembled on demand

get_dual_Hodge_sparse_matrix()[source]#
assemble_dual_Hodge_matrix()[source]#

the dual Hodge matrix is the patch-wise inverse of the multi-patch mass matrix it is not stored by default but computed on demand, by local (patch-wise) inversion of the mass matrix

class BrokenGradient_2D(V0h, V1h)[source]#

Bases: FemLinearOperator

transpose(conjugate=False)[source]#

Transpose the LinearOperator .

If conjugate is True, return the Hermitian transpose.

class BrokenTransposedGradient_2D(V0h, V1h)[source]#

Bases: FemLinearOperator

transpose(conjugate=False)[source]#

Transpose the LinearOperator .

If conjugate is True, return the Hermitian transpose.

class BrokenScalarCurl_2D(V1h, V2h)[source]#

Bases: FemLinearOperator

transpose(conjugate=False)[source]#

Transpose the LinearOperator .

If conjugate is True, return the Hermitian transpose.

class BrokenTransposedScalarCurl_2D(V1h, V2h)[source]#

Bases: FemLinearOperator

transpose(conjugate=False)[source]#

Transpose the LinearOperator .

If conjugate is True, return the Hermitian transpose.

ortho_proj_Hcurl(EE, V1h, domain_h, M1, backend_language='python')[source]#

return orthogonal projection of E on V1h, given M1 the mass matrix

class Multipatch_Projector_H1(V0h)[source]#

Bases: object

to apply the H1 projection (2D) on every patch

class Multipatch_Projector_Hcurl(V1h, nquads=None)[source]#

Bases: object

to apply the Hcurl projection (2D) on every patch

class Multipatch_Projector_L2(V2h, nquads=None)[source]#

Bases: object

to apply the L2 projection (2D) on every patch