|
template<typename BasisType , typename ortValueType , class... ortProperties> |
void | getHCurlEvaluationPoints (typename BasisType::ScalarViewType targetEPoints, typename BasisType::ScalarViewType targetCurlEPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> orts, const BasisType *cellBasis, ProjectionStruct< SpT, typename BasisType::scalarType > *projStruct, const EvalPointsType evalPointType) |
|
template<typename basisCoeffsValueType , class... basisCoeffsProperties, typename funValsValueType , class... funValsProperties, typename BasisType , typename ortValueType , class... ortProperties> |
void | getHCurlBasisCoeffs (Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties...> basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetAtTargetEPoints, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetCurlAtTargetCurlEPoints, const typename BasisType::ScalarViewType targetEPoints, const typename BasisType::ScalarViewType targetCurlEPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> orts, const BasisType *cellBasis, ProjectionStruct< SpT, typename BasisType::scalarType > *projStruct) |
|
template<typename BasisType , typename ortValueType , class... ortProperties> |
void | getHDivEvaluationPoints (typename BasisType::ScalarViewType targetEPoints, typename BasisType::ScalarViewType targetDivEPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> orts, const BasisType *cellBasis, ProjectionStruct< SpT, typename BasisType::scalarType > *projStruct, const EvalPointsType evalPointType) |
|
template<typename basisCoeffsValueType , class... basisCoeffsProperties, typename funValsValueType , class... funValsProperties, typename BasisType , typename ortValueType , class... ortProperties> |
void | getHDivBasisCoeffs (Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties...> basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetAtEPoints, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetDivAtDivEPoints, const typename BasisType::ScalarViewType targetEPoints, const typename BasisType::ScalarViewType targetDivEPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> orts, const BasisType *cellBasis, ProjectionStruct< SpT, typename BasisType::scalarType > *projStruct) |
|
template<typename BasisType , typename ortValueType , class... ortProperties> |
void | getHGradEvaluationPoints (typename BasisType::ScalarViewType ePoints, typename BasisType::ScalarViewType gradEPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> orts, const BasisType *cellBasis, ProjectionStruct< SpT, typename BasisType::scalarType > *projStruct, const EvalPointsType evalPointType) |
|
template<typename basisCoeffsValueType , class... basisCoeffsProperties, typename funValsValueType , class... funValsProperties, typename BasisType , typename ortValueType , class... ortProperties> |
void | getHGradBasisCoeffs (Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties...> basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetAtTargetEPoints, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetGradAtTargetGradEPoints, const typename BasisType::ScalarViewType targetEPoints, const typename BasisType::ScalarViewType targetGradEPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> orts, const BasisType *cellBasis, ProjectionStruct< SpT, typename BasisType::scalarType > *projStruct) |
|
template<typename BasisType , typename ortValueType , class... ortProperties> |
void | getHVolEvaluationPoints (typename BasisType::ScalarViewType ePoints, const Kokkos::DynRankView< ortValueType, ortProperties...>, const BasisType *cellBasis, ProjectionStruct< SpT, typename BasisType::scalarType > *projStruct, const EvalPointsType ePointType) |
|
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 typename BasisType::ScalarViewType targetEPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> orts, const BasisType *cellBasis, ProjectionStruct< SpT, typename BasisType::scalarType > *projStruct) |
|
template<typename BasisType , typename ortValueType , class... ortProperties> |
void | getL2EvaluationPoints (typename BasisType::ScalarViewType ePoints, const Kokkos::DynRankView< ortValueType, ortProperties...> orts, const BasisType *cellBasis, ProjectionStruct< SpT, typename BasisType::scalarType > *projStruct, const EvalPointsType ePointType) |
|
template<typename basisCoeffsValueType , class... basisCoeffsProperties, typename funValsValueType , class... funValsProperties, typename BasisType , typename ortValueType , class... ortProperties> |
void | getL2BasisCoeffs (Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties...> basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetAtTargetEPoints, const typename BasisType::ScalarViewType targetEPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> orts, const BasisType *cellBasis, ProjectionStruct< SpT, typename BasisType::scalarType > *projStruct) |
|
|
template<typename BasisType , typename ortValueType , class... ortProperties> |
static void | getL2EvaluationPoints (typename BasisType::ScalarViewType evaluationPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> cellOrientations, const BasisType *cellBasis, ProjectionStruct< ExecSpaceType, typename BasisType::scalarType > *projStruct, const EvalPointsType evalPointType=EvalPointsType::TARGET) |
| Computes evaluation points for L2 projection. More...
|
|
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 typename BasisType::ScalarViewType evaluationPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> cellOrientations, const BasisType *cellBasis, ProjectionStruct< ExecSpaceType, typename BasisType::scalarType > *projStruct) |
| Computes the basis coefficients of the L2 projection of the target function. More...
|
|
template<typename BasisType , typename ortValueType , class... ortProperties> |
static void | getHGradEvaluationPoints (typename BasisType::ScalarViewType evaluationPoints, typename BasisType::ScalarViewType gradEvalPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> cellOrientations, const BasisType *cellBasis, ProjectionStruct< ExecSpaceType, typename BasisType::scalarType > *projStruct, const EvalPointsType evalPointType=EvalPointsType::TARGET) |
| Computes evaluation points for HGrad projection. More...
|
|
template<typename basisCoeffsValueType , class... basisCoeffsProperties, typename funValsValueType , class... funValsProperties, typename BasisType , typename ortValueType , class... ortProperties> |
static void | getHGradBasisCoeffs (Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties...> basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetAtEvalPoints, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetGradAtGradEvalPoints, const typename BasisType::ScalarViewType evaluationPoints, const typename BasisType::ScalarViewType gradEvalPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> cellOrientations, const BasisType *cellBasis, ProjectionStruct< ExecSpaceType, typename BasisType::scalarType > *projStruct) |
| Computes the basis coefficients of the HGrad projection of the target function. More...
|
|
template<typename BasisType , typename ortValueType , class... ortProperties> |
static void | getHCurlEvaluationPoints (typename BasisType::ScalarViewType evaluationPoints, typename BasisType::ScalarViewType curlEvalPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> cellOrientations, const BasisType *cellBasis, ProjectionStruct< ExecSpaceType, typename BasisType::scalarType > *projStruct, const EvalPointsType evalPointType=EvalPointsType::TARGET) |
| Computes evaluation points for HCurl projection. 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 typename BasisType::ScalarViewType evaluationPoints, const typename BasisType::ScalarViewType curlEvalPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> cellOrientations, const BasisType *cellBasis, ProjectionStruct< ExecSpaceType, typename BasisType::scalarType > *projStruct) |
| Computes the basis coefficients of the HCurl projection of the target function. More...
|
|
template<typename BasisType , typename ortValueType , class... ortProperties> |
static void | getHDivEvaluationPoints (typename BasisType::ScalarViewType evaluationPoints, typename BasisType::ScalarViewType divEvalPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> cellOrientations, const BasisType *cellBasis, ProjectionStruct< ExecSpaceType, typename BasisType::scalarType > *projStruct, const EvalPointsType evalPointType=EvalPointsType::TARGET) |
| Computes evaluation points for HDiv projection. 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 typename BasisType::ScalarViewType evaluationPoints, const typename BasisType::ScalarViewType divEvalPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> cellOrientations, const BasisType *cellBasis, ProjectionStruct< ExecSpaceType, typename BasisType::scalarType > *projStruct) |
| Computes the basis coefficients of the HDiv projection of the target function. More...
|
|
template<typename BasisType , typename ortValueType , class... ortProperties> |
static void | getHVolEvaluationPoints (typename BasisType::ScalarViewType evaluationPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> cellOrientations, const BasisType *cellBasis, ProjectionStruct< ExecSpaceType, typename BasisType::scalarType > *projStruct, const EvalPointsType evalPointType=EvalPointsType::TARGET) |
| Computes evaluation points for HVol projection. 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, const typename BasisType::ScalarViewType evaluationPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> cellOrientations, const BasisType *cellBasis, ProjectionStruct< ExecSpaceType, typename BasisType::scalarType > *projStruct) |
| Computes the basis coefficients of the HVol projection of the target function. More...
|
|
template<typename ExecSpaceType>
class Intrepid2::Experimental::ProjectionTools< ExecSpaceType >
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:
- create a ProjectionStruct object
- allocate views for storing the points where to evaluate the target function and its derivatives
- evaluate the points/weights using one of the methods getHGradEvaluationPoints, getHCURLEvaluationPoints, getHDivEvaluationPoints, getHVolEvaluationPoints, or getL2EvaluationPoints
- 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
- Evaluate the basis coefficients using one of the methods getHGradBasisCoeffs, getHCURLBasisCoeffs, getHDivBasisCoeffs, getHVolBasisCoeffs, or getL2BasisCoeffs
- Todo:
- There is room for significant 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. For internal evaluation points (that are not affected by orientation) one could compute the QR factorization on the reference cell and then use on all the cells.
Note: Other algorithmic improvements could be enabled by accessing the implementation of the orientation tools,
however, we preferred the projections to work with any orientation, and assuming only that internal basis functions are not affected by
the orientation.
Definition at line 183 of file Intrepid2_ProjectionTools.hpp.