feec.multipatch.non_matching_operators#

This module provides utilities for constructing the conforming projections for a H1-Hcurl-L2 broken FEEC de Rham sequence.

Functions#

calculate_mass_matrix(space_1d)

Calculate the mass-matrix of a 1d spline-space.

calculate_mixed_mass_matrix(domain_space, ...)

Calculate the mixed mass-matrix of two 1d spline-spaces on the same domain.

calculate_poly_basis_integral(space_1d[, ...])

Calculate the "mixed mass-matrix" of a 1d spline-space with polynomials.

construct_extension_operator_1D(domain, codomain)

Compute the matrix of the extension operator on the interface.

construct_h1_conforming_projection(Vh[, ...])

Construct the conforming projection for a scalar space for a given regularity (0 continuous, -1 discontinuous).

construct_hcurl_conforming_projection(Vh[, ...])

Construct the conforming projection for a vector Hcurl space for a given regularity (0 continuous, -1 discontinuous).

construct_restriction_operator_1D(...[, ...])

Compute the matrix of the (moment preserving) restriction operator on the interface.

get_1d_moment_correction(space_1d[, p_moments])

Calculate the coefficients for the one-dimensional moment correction.

get_corners(domain, boundary_only)

Given the domain, extract the vertices on their respective domains with local coordinates.

get_extension_restriction(coarse_space_1d, ...)

Calculate the extension and restriction matrices for refining along an interface.

get_patch_index_from_face(domain, face)

Return the patch index of subdomain/boundary

knots_to_insert(coarse_grid, fine_grid[, tol])

knot insertion for refinement of a 1d spline space.

Classes#

Inheritance diagram of psydac.feec.multipatch.non_matching_operators

Local2GlobalIndexMap(ndim, n_patches, ...)

Details#

This module provides utilities for constructing the conforming projections for a H1-Hcurl-L2 broken FEEC de Rham sequence.

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

class Local2GlobalIndexMap(ndim, n_patches, n_components)[source]#

Bases: object

set_patch_shapes(patch_index, *shapes)[source]#
get_index(k, d, cartesian_index)[source]#

Return a global scalar index.

Parameters:
kint

The patch index.

dint

The component of a scalar field in the system of equations.

cartesian_index: tuple[int]

Multi index [i1, i2, i3 …]

Returns:
Iint

The global scalar index.

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

knot insertion for refinement of a 1d spline space.

get_corners(domain, boundary_only)[source]#

Given the domain, extract the vertices on their respective domains with local coordinates.

Parameters:
domain: <Geometry>

The discrete domain of the projector

boundary_only<bool>

Only return vertices that lie on a boundary

construct_extension_operator_1D(domain, codomain)[source]#

Compute the matrix of the extension operator on the interface.

Parameters:
domain1d spline space on the interface (coarse grid)
codomain1d spline space on the interface (fine grid)
construct_restriction_operator_1D(coarse_space_1d, fine_space_1d, E, p_moments=-1)[source]#

Compute the matrix of the (moment preserving) restriction operator on the interface.

Parameters:
coarse_space_1d1d spline space on the interface (coarse grid)
fine_space_1d1d spline space on the interface (fine grid)
EExtension matrix
p_momentsAmount of moments to be preserved
get_extension_restriction(coarse_space_1d, fine_space_1d, p_moments=-1)[source]#

Calculate the extension and restriction matrices for refining along an interface.

Parameters:
coarse_space_1dSplineSpace

Spline space of the coarse space.

fine_space_1dSplineSpace

Spline space of the fine space.

p_moments{int}

Amount of moments to be preserved.

Returns:
E_1Dnumpy array

Extension matrix.

R_1Dnumpy array

Restriction matrix.

ER_1Dnumpy array

Extension-restriction matrix.

calculate_mass_matrix(space_1d)[source]#

Calculate the mass-matrix of a 1d spline-space.

Parameters:
space_1dSplineSpace

Spline space of the fine space.

Returns:
Mass_matnumpy array

Mass matrix.

calculate_mixed_mass_matrix(domain_space, codomain_space)[source]#

Calculate the mixed mass-matrix of two 1d spline-spaces on the same domain.

Parameters:
domain_spaceSplineSpace

Spline space of the domain space.

codomain_spaceSplineSpace

Spline space of the codomain space.

Returns:
Mass_matnumpy array

Mass matrix.

calculate_poly_basis_integral(space_1d, p_moments=-1)[source]#

Calculate the “mixed mass-matrix” of a 1d spline-space with polynomials.

Parameters:
space_1dSplineSpace

Spline space of the fine space.

p_momentsInt

Amount of moments to be preserved.

Returns:
Mass_matnumpy array

Mass matrix.

get_1d_moment_correction(space_1d, p_moments=-1)[source]#

Calculate the coefficients for the one-dimensional moment correction.

Parameters:
patch_spaceSplineSpace

1d spline space.

p_momentsint

Number of moments to be preserved.

Returns:
gammaarray

Moment correction coefficients without the conformity factor.

construct_h1_conforming_projection(Vh, reg_orders=0, p_moments=-1, hom_bc=False)[source]#

Construct the conforming projection for a scalar space for a given regularity (0 continuous, -1 discontinuous).

Parameters:
VhTensorFemSpace

Finite Element Space coming from the discrete de Rham sequence.

reg_orders(int)

Regularity in each space direction -1 or 0.

p_moments(int)

Number of moments to be preserved.

hom_bc(bool)

Homogeneous boundary conditions.

Returns:
cPscipy.sparse.csr_array

Conforming projection as a sparse matrix.

construct_hcurl_conforming_projection(Vh, reg_orders=0, p_moments=-1, hom_bc=False)[source]#

Construct the conforming projection for a vector Hcurl space for a given regularity (0 continuous, -1 discontinuous).

Parameters:
VhTensorFemSpace

Finite Element Space coming from the discrete de Rham sequence.

reg_orders(int)

Regularity in each space direction -1 or 0.

p_moments(int)

Number of polynomial moments to be preserved.

hom_bc(bool)

Tangential homogeneous boundary conditions.

Returns:
cPscipy.sparse.csr_array

Conforming projection as a sparse matrix.