Intrepid2
Classes | Public Types | Public Member Functions | Static Public Member Functions | List of all members
Intrepid2::ProjectionTools< DeviceType > Class Template Reference

A class providing static members to perform projection-based interpolations: More...

#include <Intrepid2_ProjectionTools.hpp>

Classes

struct  ElemSystem
 Class to solve a square system A x = b on each cell A is expected to be saddle a point (KKT) matrix of the form [C B; B^T 0], where C has size nxn and B nxm, with n>0, m>=0. B^T is copied from B, so one does not have to define the B^T portion of A. b will contain the solution x. The first n-entries of x are copied into the provided basis coefficients using the provided indexing. The system is solved either with a QR factorization implemented in KokkosKernels or with Lapack GELS function. More...
 

Public Types

using ExecSpaceType = typename DeviceType::execution_space
 
using MemSpaceType = typename DeviceType::memory_space
 
using EvalPointsType = typename ProjectionStruct< DeviceType, double >::EvalPointsType
 

Public Member Functions

template<typename basisCoeffsValueType , class... basisCoeffsProperties, typename funValsValueType , class... funValsProperties, typename BasisType , typename ortValueType , class... ortProperties>
void getHVolBasisCoeffs (Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties...> basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetAtTargetEPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> orts, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct)
 

Static Public Member Functions

template<typename basisCoeffsValueType , class... basisCoeffsProperties, typename funValsValueType , class... funValsProperties, typename BasisType , typename ortValueType , class... ortProperties>
static void getL2BasisCoeffs (Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties...> basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetAtEvalPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> cellOrientations, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct)
 Computes the basis coefficients of the L2 projection of the target function. More...
 
template<typename basisCoeffsValueType , class... basisCoeffsProperties, typename funValsValueType , class... funValsProperties, typename BasisType , typename ortValueType , class... ortProperties>
static void getL2DGBasisCoeffs (Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties...> basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetAtEvalPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> cellOrientations, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct)
 Computes evaluation points for local L2 projection for broken HGRAD HCURL HDIV and HVOL spaces. More...
 
template<typename basisViewType , typename targetViewType , typename BasisType >
static void getL2DGBasisCoeffs (basisViewType basisCoeffs, const targetViewType targetAtTargetEPoints, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct)
 Computes evaluation points for local L2 projection for broken HGRAD HCURL HDIV and HVOL spaces. More...
 
template<class BasisCoeffsViewType , class TargetValueViewType , class TargetGradViewType , class BasisType , class OrientationViewType >
static void getHGradBasisCoeffs (BasisCoeffsViewType basisCoeffs, const TargetValueViewType targetAtEvalPoints, const TargetGradViewType targetGradAtGradEvalPoints, const OrientationViewType cellOrientations, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct)
 Computes the basis coefficients of the HGrad projection of the target function. More...
 
template<typename basisCoeffsValueType , class... basisCoeffsProperties, typename funValsValueType , class... funValsProperties, typename BasisType , typename ortValueType , class... ortProperties>
static void getHCurlBasisCoeffs (Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties...> basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetAtEvalPoints, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetCurlAtCurlEvalPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> cellOrientations, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct)
 Computes the basis coefficients of the HCurl projection of the target function. More...
 
template<typename basisCoeffsValueType , class... basisCoeffsProperties, typename funValsValueType , class... funValsProperties, typename BasisType , typename ortValueType , class... ortProperties>
static void getHDivBasisCoeffs (Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties...> basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetAtEvalPoints, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetDivAtDivEvalPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> cellOrientations, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct)
 Computes the basis coefficients of the HDiv projection of the target function. More...
 
template<typename basisCoeffsValueType , class... basisCoeffsProperties, typename funValsValueType , class... funValsProperties, typename BasisType , typename ortValueType , class... ortProperties>
static void getHVolBasisCoeffs (Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties...> basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetAtEvalPoints, [[maybe_unused]] const Kokkos::DynRankView< ortValueType, ortProperties...> cellOrientations, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct)
 Computes the basis coefficients of the HVol projection of the target function. More...
 
template<typename dstViewType , typename dstBasisType , typename srcViewType , typename srcBasisType , typename OrientationViewType >
static void projectField (dstViewType dstCoeffs, const dstBasisType *dstBasis, const srcViewType srcCoeffs, const srcBasisType *srcBasis, const OrientationViewType cellOrientations)
 Computes the L2 projection of a finite element field into a compatible finite element space. More...
 

Detailed Description

template<typename DeviceType>
class Intrepid2::ProjectionTools< DeviceType >

A class providing static members to perform projection-based interpolations:

This class provides tools to perform projection-based interpolations of a target function $f$ that satisfy the commutativity property of the De Rham complex and that have optimal accuracy w.r.t. the corresponding norms ( $H^1$, $H^{\text{div}}$, $H^{\text{curl}}$, $L^2$).

The projection-based interpolation $\Pi_h^X$ of the target function $f$ reads

\[ \Pi_h^X(f) = \sum_i \alpha^f_i \phi_i, \quad X \in \{\text{grad},\, \text{div},\, \text{curl},\, \text{vol}\}, \]

where $\{\phi_i\}$ is the basis of the finite element, $\alpha_i^f$ are the basisCoeffs.

It also provides tools to perform a local L2 projection into HGrad, HCurl, HDiv and L2 fields. This projection does not satisfy the properties of the projection-based interpolations, but it is simpler and does not require to evaluate the derivatives of the target functions.

Use:

  1. create a ProjectionStruct object
  2. get the evaluation points where to evaluate the target function and its derivative using the ProjectionStruct methods getAllEvalPoints and getAllDerivEvalPoints
  3. Map to the physical elements the evaluation points, evaluate the target function and its derivatives at these points and map them back (inverse of pullback operator) to the reference points. Note: if the target function is defined at reference element (e.g. in case the target function is a combination of basis functions) there is no need to map points and functions between the reference and physical spaces, but one needs to simply evaluate the target function and its derivatives at the evaluation points
  4. Evaluate the basis coefficients using one of the methods getHGradBasisCoeffs, getHCURLBasisCoeffs, getHDivBasisCoeffs, getHVolBasisCoeffs, or getL2BasisCoeffs
Remarks
The projections are performed at the oriented reference element. Therefore, the target function $f$, which is contravariant, needs to be mapped back to the reference element (using inverse of pullback operator) before calling getXXXBasisCoeffs.
The projections are implemented as in "Demkowics et al., Computing with Hp-Adaptive Finite Elements, Vol. 2, Chapman & Hall/CRC 2007". However, the projections along edges for HGrad and HCurl elements have been performed on the $H^1$ seminorm and the $L^2$ norm respectively, instead of on the $L^2$ and $H^{-1}$ and norms. This requires more regularity of the target function.
Todo:
There is room for improvement. One could separate the computation of the basis function values and derivatives from the functions getXXXBasisCoeffs, so that they can be stored and reused for projecting other target functions. Similarly one could store all the QR factorizations and reuse them for other target functions.

Definition at line 145 of file Intrepid2_ProjectionTools.hpp.

Member Function Documentation

template<typename DeviceType >
template<typename basisCoeffsValueType , class... basisCoeffsProperties, typename funValsValueType , class... funValsProperties, typename BasisType , typename ortValueType , class... ortProperties>
void Intrepid2::ProjectionTools< DeviceType >::getHCurlBasisCoeffs ( Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties...>  basisCoeffs,
const Kokkos::DynRankView< funValsValueType, funValsProperties...>  targetAtEvalPoints,
const Kokkos::DynRankView< funValsValueType, funValsProperties...>  targetCurlAtCurlEvalPoints,
const Kokkos::DynRankView< ortValueType, ortProperties...>  cellOrientations,
const BasisType *  cellBasis,
ProjectionStruct< DeviceType, typename BasisType::scalarType > *  projStruct 
)
static

Computes the basis coefficients of the HCurl projection of the target function.

C - num. cells
F - num. fields
P1 - num. evaluation points
P2 - num. evaluation points for derivatives
D - spatial dimension
Parameters
basisCoeffs[out] - rank-2 view (C,F) containing the basis coefficients
targetAtEvalPoints[in] - rank-3 view (C,P1,D) containing the values of the target function evaluated at the evaluation points given by projStruct->getAllEvalPoints()
targetcurlAtCurlEvalPoints[in] - variable rank view containing the values of the curl of the target function evaluated at the evaluation points given by projStruct->getAllDerivEvalPoints()
cellOrientations [in] - 1-rank view (C) containing the Orientation objects at each cell
cellBasis [in] - pointer to the HCURL basis for the projection
projStruct [in] - pointer to ProjectionStruct object
Remarks
targetAtCurlEvalPoints has rank 2 (C,P2) in 2D, and rank 3 (C,P2,D) in 3D

Definition at line 433 of file Intrepid2_ProjectionToolsDefHCURL.hpp.

References Intrepid2::RealSpaceTools< DeviceType >::clone(), Intrepid2::RefSubcellParametrization< DeviceType >::get(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getAllDerivEvalPoints(), Intrepid2::Basis< DeviceType, OutputType, PointType >::getAllDofOrdinal(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getAllEvalPoints(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getBasisDerivEvalWeights(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getBasisDerivPointsRange(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getBasisEvalWeights(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getBasisPointsRange(), Intrepid2::Basis< DeviceType, OutputType, PointType >::getCardinality(), Intrepid2::Basis< DeviceType, OutputType, PointType >::getDofCount(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getNumBasisDerivEvalPoints(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getNumBasisEvalPoints(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getNumTargetDerivEvalPoints(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getNumTargetEvalPoints(), Intrepid2::CellTools< DeviceType >::getReferenceEdgeTangent(), Intrepid2::CellTools< DeviceType >::getReferenceFaceNormal(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getTargetDerivEvalWeights(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getTargetDerivPointsRange(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getTargetEvalWeights(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getTargetPointsRange(), Intrepid2::Basis< DeviceType, OutputType, PointType >::getValues(), Intrepid2::FunctionSpaceTools< DeviceType >::integrate(), Intrepid2::OrientationTools< DeviceType >::modifyBasisByOrientationInverse(), and Intrepid2::ProjectionTools< DeviceType >::ElemSystem::solve().

template<typename DeviceType >
template<typename basisCoeffsValueType , class... basisCoeffsProperties, typename funValsValueType , class... funValsProperties, typename BasisType , typename ortValueType , class... ortProperties>
void Intrepid2::ProjectionTools< DeviceType >::getHDivBasisCoeffs ( Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties...>  basisCoeffs,
const Kokkos::DynRankView< funValsValueType, funValsProperties...>  targetAtEvalPoints,
const Kokkos::DynRankView< funValsValueType, funValsProperties...>  targetDivAtDivEvalPoints,
const Kokkos::DynRankView< ortValueType, ortProperties...>  cellOrientations,
const BasisType *  cellBasis,
ProjectionStruct< DeviceType, typename BasisType::scalarType > *  projStruct 
)
static

Computes the basis coefficients of the HDiv projection of the target function.

C - num. cells
F - num. fields
P1 - num. evaluation points
P2 - num. evaluation points for derivatives
D - spatial dimension
Parameters
basisCoeffs[out] - rank-2 view (C,F) containing the basis coefficients
targetAtEvalPoints[in] - rank-3 view (C,P1,D) containing the values of the target function evaluated at the evaluation points given by projStruct->getAllEvalPoints()
targetDivAtDivEvalPoints[in] - rank-2 view (C,P2) view containing the values of the divergence of the target function evaluated at the evaluation points given by projStruct->getAllDerivEvalPoints()
cellOrientations [in] - 1-rank view (C) containing the Orientation objects at each cell
cellBasis [in] - pointer to the HDIV basis for the projection
projStruct [in] - pointer to ProjectionStruct object

Definition at line 268 of file Intrepid2_ProjectionToolsDefHDIV.hpp.

References Intrepid2::RealSpaceTools< DeviceType >::clone(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getAllDerivEvalPoints(), Intrepid2::Basis< DeviceType, OutputType, PointType >::getAllDofOrdinal(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getAllEvalPoints(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getBasisDerivEvalWeights(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getBasisDerivPointsRange(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getBasisEvalWeights(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getBasisPointsRange(), Intrepid2::Basis< DeviceType, OutputType, PointType >::getCardinality(), Intrepid2::Basis< DeviceType, OutputType, PointType >::getDofCount(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getNumBasisDerivEvalPoints(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getNumBasisEvalPoints(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getNumTargetDerivEvalPoints(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getNumTargetEvalPoints(), Intrepid2::CellTools< DeviceType >::getReferenceSideNormal(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getTargetDerivEvalWeights(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getTargetDerivPointsRange(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getTargetEvalWeights(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getTargetPointsRange(), Intrepid2::Basis< DeviceType, OutputType, PointType >::getValues(), Intrepid2::FunctionSpaceTools< DeviceType >::integrate(), Intrepid2::OrientationTools< DeviceType >::modifyBasisByOrientationInverse(), and Intrepid2::ProjectionTools< DeviceType >::ElemSystem::solve().

template<typename DeviceType >
template<class BasisCoeffsViewType , class TargetValueViewType , class TargetGradViewType , class BasisType , class OrientationViewType >
void Intrepid2::ProjectionTools< DeviceType >::getHGradBasisCoeffs ( BasisCoeffsViewType  basisCoeffs,
const TargetValueViewType  targetAtEvalPoints,
const TargetGradViewType  targetGradAtGradEvalPoints,
const OrientationViewType  cellOrientations,
const BasisType *  cellBasis,
ProjectionStruct< DeviceType, typename BasisType::scalarType > *  projStruct 
)
static

Computes the basis coefficients of the HGrad projection of the target function.

C - num. cells
F - num. fields
P1 - num. evaluation points
P2 - num. evaluation points for derivatives
D - spatial dimension
Parameters
basisCoeffs[out] - rank-2 view (C,F) containing the basis coefficients
targetAtEvalPoints[in] - rank-2 view (C,P1) containing the values of the target function evaluated at the evaluation points given by projStruct->getAllEvalPoints()
targetGradAtGradEvalPoints[in] - rank-3 view (C,P2,D) view containing the values of the gradient of the target function evaluated at the evaluation points given by projStruct->getAllDerivEvalPoints()
cellOrientations [in] - 1-rank view (C) containing the Orientation objects at each cell
cellBasis [in] - pointer to the HGRAD basis for the projection
projStruct [in] - pointer to ProjectionStruct object

Definition at line 321 of file Intrepid2_ProjectionToolsDefHGRAD.hpp.

References Intrepid2::RefSubcellParametrization< DeviceType >::get(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getAllDerivEvalPoints(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getAllEvalPoints(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getBasisDerivEvalWeights(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getBasisDerivPointsRange(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getNumBasisDerivEvalPoints(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getNumTargetDerivEvalPoints(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getNumTargetEvalPoints(), Intrepid2::CellTools< DeviceType >::getReferenceEdgeTangent(), Intrepid2::CellTools< DeviceType >::getReferenceFaceNormal(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getTargetDerivEvalWeights(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getTargetDerivPointsRange(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getTargetPointsRange(), Intrepid2::FunctionSpaceTools< DeviceType >::integrate(), Intrepid2::OrientationTools< DeviceType >::modifyBasisByOrientationInverse(), and Intrepid2::ProjectionTools< DeviceType >::ElemSystem::solve().

Referenced by Intrepid2::ProjectedGeometry< spaceDim, PointScalar, DeviceType >::projectOntoHGRADBasis().

template<typename DeviceType >
template<typename basisCoeffsValueType , class... basisCoeffsProperties, typename funValsValueType , class... funValsProperties, typename BasisType , typename ortValueType , class... ortProperties>
static void Intrepid2::ProjectionTools< DeviceType >::getHVolBasisCoeffs ( Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties...>  basisCoeffs,
const Kokkos::DynRankView< funValsValueType, funValsProperties...>  targetAtEvalPoints,
[[maybe_unused] ] const Kokkos::DynRankView< ortValueType, ortProperties...>  cellOrientations,
const BasisType *  cellBasis,
ProjectionStruct< DeviceType, typename BasisType::scalarType > *  projStruct 
)
static

Computes the basis coefficients of the HVol projection of the target function.

C - num. cells
F - num. fields
P - num. evaluation points
D - spatial dimension
Parameters
basisCoeffs[out] - rank-2 view (C,F) containing the basis coefficients
targetAtEvalPoints[in] - rank-2 view (C,P) containing the values of the target function evaluated at the evaluation points given by projStruct->getAllEvalPoints()
cellOrientations [in] - 1-rank view (C) containing the Orientation objects at each cell
cellBasis [in] - pointer to the HGRAD basis for the projection
projStruct [in] - pointer to ProjectionStruct object
template<typename DeviceType >
template<typename basisCoeffsValueType , class... basisCoeffsProperties, typename funValsValueType , class... funValsProperties, typename BasisType , typename ortValueType , class... ortProperties>
void Intrepid2::ProjectionTools< DeviceType >::getL2BasisCoeffs ( Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties...>  basisCoeffs,
const Kokkos::DynRankView< funValsValueType, funValsProperties...>  targetAtEvalPoints,
const Kokkos::DynRankView< ortValueType, ortProperties...>  cellOrientations,
const BasisType *  cellBasis,
ProjectionStruct< DeviceType, typename BasisType::scalarType > *  projStruct 
)
static

Computes the basis coefficients of the L2 projection of the target function.

C - num. cells
F - num. fields
P - num. evaluation points
D - spatial dimension
Parameters
basisCoeffs[out] - rank-2 view (C,F) containing the basis coefficients
targetAtEvalPoints[in] - variable rank view containing the values of the target function evaluated at the evaluation points given by projStruct->getAllEvalPoints()
cellOrientations [in] - 1-rank view (C) containing the Orientation objects at each cell
cellBasis [in] - pointer to the basis for the projection
projStruct [in] - pointer to ProjectionStruct object
Remarks
targetAtEvalPoints has rank 2 (C,P) for HGRAD and HVOL elements, and rank 3 (C,P,D) for HCURL and HDIV elements

Definition at line 391 of file Intrepid2_ProjectionToolsDefL2.hpp.

References Intrepid2::ProjectionStruct< DeviceType, ValueType >::getAllEvalPoints(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getBasisEvalWeights(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getBasisPointsRange(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getNumBasisEvalPoints(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getNumTargetEvalPoints(), Intrepid2::CellTools< DeviceType >::getReferenceEdgeTangent(), Intrepid2::CellTools< DeviceType >::getReferenceSideNormal(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getTargetEvalWeights(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getTargetPointsRange(), Intrepid2::FunctionSpaceTools< DeviceType >::integrate(), Intrepid2::OrientationTools< DeviceType >::modifyBasisByOrientationInverse(), and Intrepid2::ProjectionTools< DeviceType >::ElemSystem::solve().

Referenced by Intrepid2::ProjectionTools< DeviceType >::projectField().

template<typename DeviceType >
template<typename basisCoeffsValueType , class... basisCoeffsProperties, typename funValsValueType , class... funValsProperties, typename BasisType , typename ortValueType , class... ortProperties>
void Intrepid2::ProjectionTools< DeviceType >::getL2DGBasisCoeffs ( Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties...>  basisCoeffs,
const Kokkos::DynRankView< funValsValueType, funValsProperties...>  targetAtEvalPoints,
const Kokkos::DynRankView< ortValueType, ortProperties...>  cellOrientations,
const BasisType *  cellBasis,
ProjectionStruct< DeviceType, typename BasisType::scalarType > *  projStruct 
)
static

Computes evaluation points for local L2 projection for broken HGRAD HCURL HDIV and HVOL spaces.

This function simply perform an L2 projection in each cell with no guarantee of preserving continuity across cells

C - num. cells
F - num. fields
P - num. evaluation points
D - spatial dimension
Parameters
basisCoeffs[out] - rank-2 view (C,F) containing the basis coefficients
targetAtEvalPoints[in] - variable rank view containing the values of the target function evaluated at the evaluation points
cellOrientations[in] - 1-rank view (C) containing the Orientation objects at each cell
cellBasis[in] - pointer to the basis for the projection
projStruct[in] - pointer to ProjectionStruct object
Remarks
targetAtEvalPoints has rank 2 (C,P) for HGRAD and HVOL elements, and rank 3 (C,P,D) for HCURL and HDIV elements

Definition at line 656 of file Intrepid2_ProjectionToolsDefL2.hpp.

References Intrepid2::OrientationTools< DeviceType >::modifyBasisByOrientationInverse().

template<typename DeviceType >
template<typename basisViewType , typename targetViewType, typename BasisType >
void Intrepid2::ProjectionTools< DeviceType >::getL2DGBasisCoeffs ( basisViewType  basisCoeffs,
const targetViewType  targetAtTargetEPoints,
const BasisType *  cellBasis,
ProjectionStruct< DeviceType, typename BasisType::scalarType > *  projStruct 
)
static

Computes evaluation points for local L2 projection for broken HGRAD HCURL HDIV and HVOL spaces.

This function simply perform an L2 projection in each cell with no guarantee of preserving continuity across cells. It also does not account for orientation.

C - num. cells
F - num. fields
P - num. evaluation points
D - spatial dimension
Parameters
basisCoeffs[out] - rank-2 view (C,F) containing the basis coefficients
targetAtEvalPoints[in] - variable rank view containing the values of the target function evaluated at the evaluation points
cellBasis[in] - pointer to the basis for the projection
projStruct[in] - pointer to ProjectionStruct object
Remarks
targetAtEvalPoints has rank 2 (C,P) for HGRAD and HVOL elements, and rank 3 (C,P,D) for HCURL and HDIV elements

Definition at line 672 of file Intrepid2_ProjectionToolsDefL2.hpp.

References Intrepid2::RealSpaceTools< DeviceType >::clone(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getBasisEvalWeights(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getEvalPoints(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getNumBasisEvalPoints(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::getTargetEvalWeights(), Intrepid2::FunctionSpaceTools< DeviceType >::integrate(), and Intrepid2::ProjectionTools< DeviceType >::ElemSystem::solve().

template<typename DeviceType >
template<typename dstViewType , typename dstBasisType , typename srcViewType , typename srcBasisType , typename OrientationViewType >
static void Intrepid2::ProjectionTools< DeviceType >::projectField ( dstViewType  dstCoeffs,
const dstBasisType *  dstBasis,
const srcViewType  srcCoeffs,
const srcBasisType *  srcBasis,
const OrientationViewType  cellOrientations 
)
inlinestatic

Computes the L2 projection of a finite element field into a compatible finite element space.

C - num. cells
F - num. fields
P - num. evaluation points
D - spatial dimension
Parameters
dstCoeffs[out] - rank-2 view (C,F) containing the basis coefficients of the projected field over the space generated by dstBasis
dstBasis[in] - pointer to the basis we want to project to
srcCoeffs[in] - rank-2 view (C,F) containing the basis coefficients of the field w.r.t. the original basis
srcBasis[in] - original basis
cellOrientations[in] - rank-1 view (C) containing the Orientation objects at each cell

Definition at line 399 of file Intrepid2_ProjectionTools.hpp.

References Intrepid2::ProjectionStruct< DeviceType, ValueType >::createL2ProjectionStruct(), Intrepid2::ProjectionTools< DeviceType >::getL2BasisCoeffs(), and Intrepid2::OrientationTools< DeviceType >::modifyBasisByOrientation().


The documentation for this class was generated from the following files: