Intrepid2
|
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... | |
A class providing static members to perform projection-based interpolations:
This class provides tools to perform projection-based interpolations of a target function that satisfy the commutativity property of the De Rham complex and that have optimal accuracy w.r.t. the corresponding norms ( , , , ).
The projection-based interpolation of the target function reads
where is the basis of the finite element, 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:
Definition at line 112 of file Intrepid2_ProjectionTools.hpp.
|
static |
Computes the basis coefficients of the HCurl projection of the target function.
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 |
Definition at line 400 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().
|
static |
Computes the basis coefficients of the HDiv projection of the target function.
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 235 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().
|
static |
Computes the basis coefficients of the HGrad projection of the target function.
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 288 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().
|
static |
Computes the basis coefficients of the HVol projection of the target function.
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 |
|
static |
Computes the basis coefficients of the L2 projection of the target function.
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 |
Definition at line 358 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().
|
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
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 |
Definition at line 623 of file Intrepid2_ProjectionToolsDefL2.hpp.
References Intrepid2::OrientationTools< DeviceType >::modifyBasisByOrientationInverse().
|
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.
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 |
Definition at line 639 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().
|
inlinestatic |
Computes the L2 projection of a finite element field into a compatible finite element space.
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 366 of file Intrepid2_ProjectionTools.hpp.
References Intrepid2::ProjectionStruct< DeviceType, ValueType >::createL2ProjectionStruct(), Intrepid2::ProjectionTools< DeviceType >::getL2BasisCoeffs(), and Intrepid2::OrientationTools< DeviceType >::modifyBasisByOrientation().