Projections into Derham#

class struphy.feec.projectors.CommutingProjector(projector_tensor: GlobalProjector, dofs_extraction_op=None, base_extraction_op=None, boundary_op=None)[source]#

Bases: object

A commuting projector of the 3d Derham diagram (can be polar).

The general structure of the inter-/histopolation problem reads: given a function \(f \in V^\alpha\) in one of the (continuous) de Rham spaces \(\alpha \in \{0,1,2,3,v\}\), find its projection \(f_h \in V_h^\alpha \subset V^\alpha\) determined by the spline coefficients \(\mathbf f \in \mathbb R^{N_\alpha}\) such that

\[\mathbb B \mathbb P\, \mathcal I\, \mathbb E^T \mathbb B^T \mathbf f = \mathbb B \mathbb P\, \mathbf d(f)\,,\]

where \(\mathbf d(f) \in \mathbb R^{N_\alpha}\) are the degrees of freedom corresponding to \(f\), and with the following linear operators:

\(\mathbb P\) and \(\mathbb E\) (and \(\mathbb B\) in case of no boundary conditions) can be identity operators, which gives the pure tensor-product Psydac GlobalProjector.

Parameters:
  • projector_tensor (GlobalProjector) – The pure tensor product projector.

  • dofs_extraction_op (PolarExtractionOperator, optional) – The degrees of freedom extraction operator mapping tensor product DOFs to polar DOFs. If not given, is set to identity.

  • base_extraction_op (PolarExtractionOperator, optional) – The basis extraction operator mapping tensor product basis functions to polar basis functions. If not given, is set to identity.

  • boundary_op (BoundaryOperator) – The boundary operator applying essential boundary conditions to a vector. If not given, is set to identity.

property projector_tensor#

Tensor product projector.

property space#

Tensor product FEM space corresponding to projector.

property dofs_extraction_op#

Degrees of freedom extraction operator (tensor product DOFs –> polar DOFs).

property base_extraction_op#

Basis functions extraction operator (tensor product basis functions –> polar basis functions).

property boundary_op#

Boundary operator setting essential boundary conditions to Stencil-/BlockVector.

property is_polar#

Whether the projector maps to polar splines (True) or pure tensor product splines.

property I#

Inter-/histopolation matrix ID * P * I * E^T * ID^T as ComposedLinearOperator (ID = IdentityOperator).

property I0#

Inter-/histopolation matrix B * P * I * E^T * B^T as ComposedLinearOperator.

property IT#

Transposed inter-/histopolation matrix ID * E * I^T * P^T * ID^T as ComposedLinearOperator (ID = IdentityOperator).

property I0T#

Transposed inter-/histopolation matrix B * E * I^T * P^T * B^T as ComposedLinearOperator.

property pc#

Preconditioner P * I^(-1) * E^T for iterative polar projections.

property pc0#

Preconditioner B * P * I^(-1) * E^T * B^T for iterative polar projections.

property pcT#

Transposed preconditioner P * I^(-T) * E^T for iterative polar projections.

property pc0T#

Transposed preconditioner B * P * I^(-T) * E^T * B^T for iterative polar projections.

solve(rhs, transposed=False, apply_bc=False, out=None, x0=None)[source]#

Solves the linear system I * x = rhs, resp. I^T * x = rhs for x, where I is the composite inter-/histopolation matrix.

Parameters:
  • rhs (psydac.linalg.basic.vector) – The right-hand side of the linear system.

  • transposed (bool, optional) – Whether to invert the transposed inter-/histopolation matrix.

  • apply_bc (bool, optional) – Whether to apply essential boundary conditions to degrees of freedom and coefficients.

  • out (psydac.linalg.basic.vector, optional) – If given, the result will be written into this vector in-place.

Returns:

x – Output vector (result of linear system).

Return type:

psydac.linalg.basic.vector

get_dofs(fun, dofs=None, apply_bc=False)[source]#

Computes the geometric degrees of freedom associated to given callable(s).

Parameters:
  • fun (callable | list) – The function for which the geometric degrees of freedom shall be computed. List of callables for vector-valued functions.

  • dofs (psydac.linalg.basic.vector, optional) – If given, the dofs will be written into this vector in-place.

  • apply_bc (bool, optional) – Whether to apply essential boundary conditions to degrees of freedom.

Returns:

dofs – The geometric degrees of freedom associated to given callable(s) “fun”.

Return type:

psydac.linalg.basic.vector

class struphy.feec.projectors.CommutingProjectorLocal(space_id: str, space_key: str, fem_space: FemSpace, pts: list, wts: list, wij: list, whij: list, fem_space_B: TensorFemSpace, fem_space_D: TensorFemSpace)[source]#

Bases: object

A commuting projector of the 3d Derham diagram, based on local quasi-inter/histopolation.

We shall describe the algortihm by means of 1d inter- and histopolation, which is then combined to give the 3d projections.

For interpolation, given a knot vector \(\hat{T} = \{ \eta_i \}_{0 \leq i \leq n+2p}\), we perform the following steps to obtain the i-th coefficient \(\lambda_i(f)\) of the quasi-interpolant

\[I^p f := \sum_{i=0}^{\hat{n}_N -1} \lambda_i(f) N_i^p\,.\]
  1. For \(i\) fixed, choose \(\nu - \mu +p\) equidistant interpolation points \(\{ x^i_j \}_{0 \leq j < 2p -1}\) in the sub-interval \(Q = [\eta_\mu , \eta_\nu]\) given by:

    • Clamped:

    \[\begin{split}Q = \left\{\begin{array}{lr} [\eta_p = 0, \eta_{p+1}], & i = 0 \,,\\ {[\eta_p = 0, \eta_{p+i}]}, & 0 < i < p-1\,,\\ {[\eta_{i+1}, \eta_{i+p}]}, & p-1 \leq i \leq \hat{n}_N - p\,,\\ {[\eta_{i+1}, \eta_{\hat{n}_N} = 1]}, & \hat{n}_N - p < i < \hat{n}_N -1\,,\\ {[\eta_{\hat{n}_N -1}, \eta_{\hat{n}_N} = 1]}, & i = \hat{n}_N -1 \,. \end{array} \; \right .\end{split}\]
    • Periodic:

    \[Q = [\eta_{i + 1}, \eta_{i + p}] \:\:\:\:\: \forall \:\: i.\]
    • In the periodic case the point set \(\{ x^i_j \}_{0 \leq j < 2p -1}\) is then the union of the \(p\) knots in \(Q\) plus their \(p-1\) mid-points.

2. Determine the “local coefficients” \((f_k)_{k=\mu-p}^{\nu-1} \in \mathbb R^{\nu - \mu +p}\) by solving the local interpolation problem

\[\sum_{k = \mu - p}^{\nu -1} f_k N^p_k(x^i_j) = f(x^i_j),\qquad \forall j \in \{0, ..., \nu - \mu +p -1\} .\]
  1. Set \(\lambda_i(f) = f_i\).

Solving the local interpolation problem in step 2 means that \(\lambda_i(f)\) can be written as a linear combination of \(f(x^i_j)\)

\[\lambda_i(f) = \sum_{j=0}^{\nu - \mu +p-1}\omega^i_j f(x^i_j)\,,\]

where \(\omega^i\) is the \(i\)-th line of the inverse collocation matrix \(\omega = C^{-1} \in \mathbb R^{(\nu - \mu +p)\times (\nu - \mu +p)}\) with \(C_{jk} = N^p_k(x^i_j)\).

On the other hand, the histopolation operator is defined by

\[H^{p-1}f := \sum_{i=0}^{\hat{n}_N -1} \tilde{\lambda}_i(f) D^{p-1}_i\,.\]

For the periodic case the FEEC coefficients \(\tilde{\lambda}_i(f)\) are computed with

\[\tilde{\lambda}_i(f) = \sum_{j=0}^{2p-1}\tilde{\omega}^i_j \int_{x^i_j}^{x^i_{j+1}}f(t)dt,\]

with the weights given by

\[\begin{split}\tilde{\omega}^i_j = \left\{\begin{array}{lr} \omega^i_0, & j=0 \,, \\ \omega^i_0 + \omega^i_1, & j = 1 \,, \\ \sum_{q=0}^{j}\omega^i_q - \sum_{q=0}^{j-2}\omega^{i+1}_q, & j = 2,...,2p-2 \,, \\ \sum_{q=0}^{2p-2}\omega^i_q - \sum_{q=0}^{2p-3}\omega^{i+1}_q, & j= 2p-1 \,. \end{array} \; \right .\end{split}\]

For the clamped case the FEEC coefficients \(\tilde{\lambda}_i(f)\) are computed with

\[\tilde{\lambda}_i(f) = \sum_{j=0}^{4p-5}\tilde{\omega}^i_j \int_{x^i_j}^{x^i_{j+1}}f(t)dt,\]

if \(i=0, \hat{n}_N -2\), the weights are given by

\[\begin{split}\tilde{\omega}^i_j = \left\{\begin{array}{lr} \sum_{q=0}^j (\omega^i_q - \omega^{i+1}_q), & j=0,...,p-1 \,, \\ 0, & j= p, ..., 4p-5 \,. \end{array} \; \right .\end{split}\]

For \(0<i<p-1\), they are obtained from

\[\begin{split}\tilde{\omega}^i_j = \left\{\begin{array}{lr} -\sum_{k=j+1}^{p+i-1} \omega^i_k, & j\leq p+i-2 \,, \\ 0, & j = p+i-1 \,, \\ \sum_{k=j-p-i+1}^{p+i}\omega^{i+1}_k, & p+i \leq j \leq 2p+2i-1 \,, \\ 0, & 2p+2i-1<j \leq 4p-5 \,. \end{array} \; \right .\end{split}\]

For \(p-1 \leq i < \hat{n}_N-p\), we use

\[\begin{split}\tilde{\omega}^i_j = \left\{\begin{array}{lr} \omega^i_0, & j=0 \,, \\ \omega^i_0 + \omega^i_1, & j = 1 \,, \\ \sum_{q=0}^{j}\omega^i_q - \sum_{q=0}^{j-2}\omega^{i+1}_q, & j = 2,...,2p-2 \,, \\ \sum_{q=0}^{2p-2}\omega^i_q - \sum_{q=0}^{2p-3}\omega^{i+1}_q, & j= 2p-1 \,, \\ 0, & j= 2p, ..., 4p-5 \,. \end{array} \; \right .\end{split}\]

Finally, for \(\hat{n}_N-p\leq i < \hat{n}_N-2\), we have

\[\begin{split}\tilde{\omega}^i_j = \left\{\begin{array}{lr} \sum_{k=0}^{j}\omega^i_k, & j\leq \hat{n}_N+p-i-3 \,, \\ 0, & j = \hat{n}_N+p-i-2 \,, \\ -\sum_{k=0}^{j-\hat{n}_N-p+i+1}\omega^{i+1}_k, & \hat{n}_N+p-i-1 \leq j \leq 2\hat{n}_N+2p-2i-5 \,, \\ 0, & j= 2\hat{n}_N+2p-2i-4, ..., 4p-5 \,. \end{array} \; \right .\end{split}\]

Furthermore, in the particular case \(p=1\), the weights are given by

\[\tilde{\omega}^i_j = \omega^i_0 , \:\:\:\: j=0,1.\]

The integration points \(\{x^i_j\}\) are the same as the quasi-interpolation points. Except in the case \(p=1, n=1\) with periodic boundary conditions, in which they are \(\{0,0.5,1\}\).

Parameters:#

space_idstr

One of “H1”, “Hcurl”, “Hdiv”, “L2” or “H1vec”.

space_keystr

One of “0”, “1”, “2”, “3” or “v”.

fem_spaceFemSpace

FEEC space into which the functions shall be projected.

ptslist of xp.array

3-list (or nested 3-list[3-list] for BlockVectors) of 2D arrays with the quasi-interpolation points (or Gauss-Legendre quadrature points for histopolation). In format [spatial direction](B-spline index, point) for StencilVector spaces or [vector component][spatial direction](B-spline index, point) for BlockVector spaces.

wtslist of xp.array

3D (4D for BlockVectors) list of 2D array with the Gauss-Legendre quadrature weights (full of ones for interpolation). In format [spatial direction](B-spline index, point) for StencilVector spaces or [vector component][spatial direction](B-spline index, point) for BlockVector spaces.

wijlist of xp.array

List of 2D arrays for the coefficients \(\omega_j^i\) obtained by inverting the local collocation matrix. Use for obtaining the FE coefficients of a function via interpolation. In format [spatial direction](B-spline index, point).

whijlist of xp.array

List of 2D arrays for the coefficients \(\hat{\omega}_j^i\) obtained from the \(\omega_j^i\). Use for obtaining the FE coefficients of a function via histopolation. In format [spatial direction](D-spline index, point).

fem_space_BTensorFemSpace

FEEC space for the zero forms.

fem_space_DTensorFemSpace

FEEC space for the three forms.

property space_id#

The ID of the space (H1, Hcurl, Hdiv, L2 or H1vec).

property space_key#

The key of the space (0, 1, 2, 3 or v).

property fem_space#

The Finite Elements spline space

property coeff_space#

The vector space underlying the FEM space.

property pts#

3D (4D for BlockVectors) list of 2D array with the quasi-interpolation points (or Gauss-Legendre quadrature points for histopolation). In format (ns, nb, np) = (spatial direction, B-spline index, point) for StencilVector spaces or (nv,ns, nb, np) = (vector entry,spatial direction, B-spline index, point) for BlockVector spaces.

property wts#

3D (4D for BlockVectors) list of 2D array with the Gauss-Legendre quadrature points (full of ones for interpolation). In format (ns, nb, np) = (spatial direction, B-spline index, point) for StencilVector spaces or (nv,ns, nb, np) = (vector entry,spatial direction, B-spline index, point) for BlockVector spaces.

property wij#

List of 2D arrays for the coefficients \(\omega_j^i\) obtained by inverting the local collocation matrix. Use for obtaining the FE coefficients of a function via interpolation. In format (ns, nb, np) = (spatial direction, B-spline index, point).

property whij#

List of 2D arrays for the coefficients \(\hat{\omega}_j^i\) obtained from the \(\omega_j^i\). Use for obtaining the FE coefficients of a function via histopolation. In format (ns, nb, np) = (spatial direction, D-spline index, point).

solve(rhs, out=None)[source]#

Solves

Parameters:
  • rhs (numpy array) – The right-hand side of the linear system.

  • out (psydac.linalg.basic.vector, optional) – If given, the result will be written into this vector in-place.

Returns:

out – Output vector (result of linear system).

Return type:

psydac.linalg.basic.vector

solve_weighted(rhs, out=None)[source]#

Solves

Parameters:
  • rhs (numpy array) – The right-hand side of the linear system.

  • out (numpy array, optional) – If given, the result will be written into this array in-place.

Returns:

out – 3d numpy array containing the output vector (result of linear system), only the portion that corresponds to the current MPI rank.

Return type:

numpy array

get_dofs(fun, dofs=None)[source]#

Builds 3D numpy array with the evaluation of the right-hand-side

get_dofs_weighted(fun, dofs=None, first_go=True, pre_computed_dofs=None)[source]#

Builds 3D numpy array with the evaluation of the right-hand-side.

get_translation_b(i, h)[source]#

Selects the correct translation index value. The only real functionality of this function is to make the code easier to read by hiding from the user the intricate (but necessary) way of accessing this data.

Parameters:
  • i (int) – Determines for which of three spatial directions we want to get the translation index

  • h (int) – Only for BlockVector spaces, determines the blockvector index.

Returns:

self._translation_indices_block_B_or_D_splines[i][h][self._B_or_D[i]][self._basis_indices[i]] – index of self._Basis_functions_indices_B(or D) where the B(or D) spline with the label self._basis_indices[i] is stored. This applies for the first i-th spatial direction and the h blockvector entry.

Return type:

int

get_rowstarts(i, h=None)[source]#

Selects the correct rowstarts array. The only real functionality of this function is to make the code easier to read by hiding from the user the intricate (but necessary) way of accessing this data.

Parameters:
  • i (int) – Determines for which of three spatial directions we want to get the rowstarts

  • h (int) – Only for BlockVector spaces, determines the blockvector index.

Returns:

self._rows_B_or_D_splines_i[self._B_or_D[i]][self._translation_indices_B_or_D_splines_i[self._B_or_D[i]][self._basis_indices[i]]] – Array that tell us for which rows the basis function in the i-th direction produces non-zero entries in the BasisProjectionOperatorLocal matrix. This array contains the start indices of said regions.

Return type:

1d int array

get_rowends(i, h=None)[source]#

Selects the correct rowends array. The only real functionality of this function is to make the code easier to read by hiding from the user the intricate (but necessary) way of accessing this data.

Parameters:
  • i (int) – intiger that determines for which of three spatial directions we want to get the rowends

  • h (int) – Only for BlockVector spaces, determines the blockvector index.

Returns:

self._rowe_B_or_D_splines_i[self._B_or_D[i]][self._translation_indices_B_or_D_splines_i[self._B_or_D[i]][self._basis_indices[i]]] – Array that tell us for which rows the basis function in the i-th direction produces non-zero entries in the BasisProjectionOperatorLocal matrix. This array contains the end indices of said regions.

Return type:

1d int array

get_values(i, h=None)[source]#

Returns array with the evaluated basis function for the i-th direction. The only real functionality of this function is to make the code easier to read by hiding from the user the intricate (but necessary) way of accessing this data.

Parameters:
  • i (int) – intiger that determines for which of three spatial directions we want to get the values.

  • h (int) – Only for BlockVector spaces, determines the blockvector index.

Returns:

self._values_B_or_D_splines_i[self._B_or_D[i]][self._translation_indices_B_or_D_splines_i[self._B_or_D[i]][self._basis_indices[i]]] – Array with the evaluated basis function for the i-th direction.

Return type:

1d float array

get_are_zero(i, h=None)[source]#

Selects the correct are_zero array. The only real functionality of this function is to make the code easier to read by hiding from the user the intricate (but necessary) way of accessing this data.

Parameters:
  • i (int) – intiger that determines for which of three spatial directions we want to get the arezero

  • h (int) – Only for BlockVector spaces, determines the blockvector index.

Returns:

self._are_zero_B_or_D_splines_i[self._B_or_D[i]][self._translation_indices_B_or_D_splines_i[self._B_or_D[i]][self._basis_indices[i]]] – Array of zeros or ones. A one at index j means that for the set of quadrature points found in self._localpts[i][j] the basis function is not zero for at least one of them.

Return type:

1d int array

class struphy.feec.projectors.L2Projector(space_id, mass_ops, **params)[source]#

Bases: object

An orthogonal projection into a discrete Derham space based on the L2-scalar product.

It solves the following system for the FE-coefficients \(\mathbf f = (f_{lmn}) \in \mathbb R^{N_\alpha}\):

\[\mathbb M^\alpha_{ijk, lmn} f_{lmn} = (f^\alpha, \Lambda^\alpha_{ijk})_{L^2}\,,\]

where \(\mathbb M^\alpha\) denotes the mass matrix of space \(\alpha \in \{0,1,2,3,v\}\) and \(f^\alpha\) is a \(\alpha\)-form proxy function.

Parameters:#

space_idstr

One of “H1”, “Hcurl”, “Hdiv”, “L2” or “H1vec”.

mass_opsstruphy.mass.WeighteMassOperators

Mass operators object, see mass_ops.

paramsdict

Keyword arguments for the solver parameters.

property mass_ops#

Struphy mass operators object, see mass_ops..

property space_id#

The ID of the space (H1, Hcurl, Hdiv, L2 or H1vec).

property space_key#

The key of the space (0, 1, 2, 3 or v).

property space#

The Derham finite element space (from Derham.Vh_fem).

property params#

Parameters for the iterative solver.

property Mmat#

The mass matrix of space.

property quad_grid_pts#

List of quadrature points in each direction for integration over grid cells in format (ni, nq) = (cell, quadrature point).

property quad_grid_mesh#

Mesh grids of quad_grid_pts.

property geom_weights#

Geometric coefficients (e.g. Jacobians) evaluated at quad_grid_mesh, stored as list[list] either 1x1 or 3x3.

solve(rhs, out=None)[source]#

Solves the linear system M * x = rhs, where M is the mass matrix.

Parameters:
  • rhs (psydac.linalg.basic.vector) – The right-hand side of the linear system.

  • out (psydac.linalg.basic.vector, optional) – If given, the result will be written into this vector in-place.

Returns:

out – Output vector (result of linear system).

Return type:

psydac.linalg.basic.vector

get_dofs(fun, dofs=None, apply_bc=False, clear=True)[source]#

Assembles (in 3d) the Stencil-/BlockVector

\[V_{ijk} = \int f * w_\textrm{geom} * \Lambda^\alpha_{ijk}\,\textrm d \boldsymbol \eta = \left( f\,, \Lambda^\alpha_{ijk}\right)_{L^2}\,,\]

where \(\Lambda^\alpha_{ijk}\) are the basis functions of \(V_h^\alpha\), \(f\) is an \(\alpha\)-form proxy function and \(w_\textrm{geom}\) stand for metric coefficients.

Note that any geometric terms (e.g. Jacobians) in the L2 scalar product are automatically assembled into \(w_\textrm{geom}\), depending on the space of \(\alpha\)-forms.

The integration is performed with Gauss-Legendre quadrature over the whole logical domain.

Parameters:
  • fun (callable | list) – Weight function(s) (callables or xp.ndarrays) in a 1d list of shape corresponding to number of components.

  • dofs (StencilVector | BlockVector, optional) – The vector for the output.

  • apply_bc (bool, optional) – Whether to apply essential boundary conditions to degrees of freedom.

  • clear (bool) – Whether to first set all data to zero before assembly. If False, the new contributions are added to existing ones in vec.

class struphy.feec.projectors.ProjectorPreconditioner(projector, transposed=False, apply_bc=False)[source]#

Bases: LinearOperator

Preconditioner for approximately inverting a (polar) 3d inter-/histopolation matrix via

\[(B * P * I * E^T * B^T)^{-1} \approx B * P * I^{-1} * E^T * B^T.\]

In case that $P$ and $E$ are identity operators, the solution is exact (pure tensor product case).

Parameters:
  • projector (CommutingProjector) – The global commuting projector for which the inter-/histopolation matrix shall be inverted.

  • transposed (bool, optional) – Whether to invert the transposed inter-/histopolation matrix.

  • apply_bc (bool, optional) – Whether to include the boundary operators.

property space#

Stencil-/BlockVectorSpace or PolarDerhamSpace.

property solver#

KroneckerLinearSolver for exactly inverting tensor product inter-histopolation matrix.

property transposed#

Whether to invert the transposed inter-/histopolation matrix.

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.

transpose(conjugate=False)[source]#

Returns the transposed operator.

solve(rhs, out=None)[source]#

Computes (B * P * I^(-1) * E^T * B^T) * rhs, resp. (B * P * I^(-T) * E^T * B^T) * rhs (transposed=True) as an approximation for an inverse inter-/histopolation matrix.

Parameters:
  • rhs (psydac.linalg.basic.Vector) – The right-hand side vector.

  • out (psydac.linalg.basic.Vector, optional) – If given, the output vector will be written into this vector in-place.

Returns:

out – The result of (B * E * M^(-1) * E^T * B^T) * rhs, resp. (B * P * I^(-T) * E^T * B^T) * rhs (transposed=True).

Return type:

psydac.linalg.basic.Vector

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

Apply linear operator to Vector v. Result is written to Vector out, if provided.

class psydac.feec.global_projectors.GlobalProjector(space, nquads=None)[source]#

Bases: object

Projects callable functions to some scalar or vector FEM space.

A global projector is constructed over a tensor-product grid in the logical domain. The vertices of this grid are obtained as the tensor product of the 1D splines’ Greville points along each direction.

This projector matches the “geometric” degrees of freedom of discrete n-forms (where n depends on the underlying space). This is done by projecting each component of the vector field independently, by combining 1D histopolation with 1D interpolation.

This class cannot be instantiated directly (use a subclass instead).

Parameters:
  • space (VectorFemSpace | TensorFemSpace) – Some finite element space, codomain of the projection operator. The exact structure where to use histopolation and where interpolation has to be given by a subclass of the GlobalProjector class. As of now, it is implicitly assumed for a VectorFemSpace, that for each direction that all spaces with interpolation are the same, and all spaces with histopolation are the same (i.e. yield the same quadrature/interpolation points etc.); so use with care on an arbitrary VectorFemSpace. It is right now only intended to be used with VectorFemSpaces or TensorFemSpaces from DeRham complex objects.

  • nquads (list(int) | tuple(int)) – Number of quadrature points along each direction, to be used in Gauss quadrature rule for computing the (approximated) degrees of freedom. This parameter is ignored, if the projector only uses interpolation (and no histopolation).

property space#

The space to which this Projector projects.

property dim#

The dimension of the underlying TensorFemSpaces.

property blockcount#

The number of blocks. In case that self.space is a TensorFemSpace, this is 1, otherwise it denotes the number of blocks in the VectorFemSpace.

property grid_x#

The local interpolation/histopolation grids which are used; it denotes the position of the interpolation/quadrature points. All the grids are stored inside a two-dimensional array; the outer dimension denotes the block, the inner the tensor space direction.

property grid_w#

The local interpolation/histopolation grids which are used; it denotes the weights of the quadrature points (in the case of interpolation, this will return the weight 1 for the given positions). All the grids are stored inside a two-dimensional array; the outer dimension denotes the block, the inner the tensor space direction.

property func#

The function which is used for projecting a given callable (or list thereof) to the DOFs in the target space.

property solver#

The solver used for transforming the DOFs in the target space into spline coefficients.

property imat_kronecker#

Inter-/Histopolation matrix in distributed format.