Intrepid2
Classes | Public Types | Public Member Functions | Static Public Member Functions | List of all members
Intrepid2::Experimental::ProjectionTools< ExecSpaceType > 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 EvalPointsType = typename ProjectionStruct< ExecSpaceType, double >::EvalPointsType
 

Public Member Functions

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)
 

Static Public Member Functions

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...
 

Detailed Description

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 $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. allocate views for storing the points where to evaluate the target function and its derivatives
  3. evaluate the points/weights using one of the methods getHGradEvaluationPoints, getHCURLEvaluationPoints, getHDivEvaluationPoints, getHVolEvaluationPoints, or getL2EvaluationPoints
  4. 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
  5. 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 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.

Member Function Documentation

template<typename ExecSpaceType >
template<typename basisCoeffsValueType , class... basisCoeffsProperties, typename funValsValueType , class... funValsProperties, typename BasisType , typename ortValueType , class... ortProperties>
static void Intrepid2::Experimental::ProjectionTools< ExecSpaceType >::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 
)
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
targetcurlAtCurlEvalPoints[in] - variable rank view containing the values of the curl of the target function evaluated at the evaluation points
evaluationPoints[in] - rank-3 view (C,P1,D) containing the coordinates of the evaluation points for the projection at each cell
curlEvalPoints[in] - rank-3 view (C,P2,D) containing the coordinates of the points where to evaluate the function curls, at each cell
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
template<typename ExecSpaceType >
template<typename BasisType , typename ortValueType , class... ortProperties>
static void Intrepid2::Experimental::ProjectionTools< ExecSpaceType >::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 
)
static

Computes evaluation points for HCurl projection.

C - num. cells
P1 - num. evaluation points
P2 - num. evaluation points for derivatives
D - spatial dimension
Parameters
evaluationPoints[out] - rank-3 view (C,P1,D) containing the coordinates of the evaluation points for the projection at each cell
curlEvalPoints[in] - rank-3 view (C,P2,D) containing the coordinates of the points where to evaluate the function curls, at each cell
cellOrientations[in] - rank-1 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
evalPointType[in] - enum selecting whether the points should be computed for the basis functions or for the target function
template<typename ExecSpaceType >
template<typename basisCoeffsValueType , class... basisCoeffsProperties, typename funValsValueType , class... funValsProperties, typename BasisType , typename ortValueType , class... ortProperties>
static void Intrepid2::Experimental::ProjectionTools< ExecSpaceType >::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 
)
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
targetDivAtDivEvalPoints[in] - rank-2 view (C,P2) view containing the values of the divergence of the target function evaluated at the evaluation points
evaluationPoints[in] - rank-3 view (C,P1,D) containing the coordinates of the evaluation points, at each cell
divEvalPoints[in] - rank-3 view (C,P2,D) containing the coordinates of the points where to evaluate the function divergence, at each cell
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
template<typename ExecSpaceType >
template<typename BasisType , typename ortValueType , class... ortProperties>
static void Intrepid2::Experimental::ProjectionTools< ExecSpaceType >::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 
)
static

Computes evaluation points for HDiv projection.

C - num. cells
P1 - num. evaluation points
P2 - num. evaluation points for derivatives
D - spatial dimension
Parameters
evaluationPoints[out] - rank-3 view (C,P1,D) containing the coordinates of the evaluation points for the projection at each cell
divEvalPoints[in] - rank-3 view (C,P2,D) containing the coordinates of the points where to evaluate the function divergence, at each cell
cellOrientations[in] - rank-1 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
evalPointType[in] - enum selecting whether the points should be computed for the basis functions or for the target function
template<typename ExecSpaceType >
template<typename basisCoeffsValueType , class... basisCoeffsProperties, typename funValsValueType , class... funValsProperties, typename BasisType , typename ortValueType , class... ortProperties>
static void Intrepid2::Experimental::ProjectionTools< ExecSpaceType >::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 
)
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
targetGradAtGradEvalPoints[in] - rank-3 view (C,P2,D) view containing the values of the gradient of the target function evaluated at the evaluation points
evaluationPoints[in] - rank-3 view (C,P1,D) containing the coordinates of the evaluation points, at each cell
gradEvalPoints[in] - rank-3 view (C,P2,D) containing the coordinates of the points where to evaluate the function gradients, at each cell
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 ExecSpaceType >
template<typename BasisType , typename ortValueType , class... ortProperties>
static void Intrepid2::Experimental::ProjectionTools< ExecSpaceType >::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 
)
static

Computes evaluation points for HGrad projection.

C - num. cells
P1 - num. evaluation points
P2 - num. evaluation points for derivatives
D - spatial dimension
Parameters
evaluationPoints[out] - rank-3 view (C,P1,D) containing the coordinates of the evaluation points, at each cell
gradEvalPoints[in] - rank-3 view (C,P2,D) containing the coordinates of the points where to evaluate the function gradients, at each cell
cellOrientations[in] - rank-1 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
evalPointType[in] - enum selecting whether the points should be computed for the basis functions or for the target function
template<typename ExecSpaceType >
template<typename basisCoeffsValueType , class... basisCoeffsProperties, typename funValsValueType , class... funValsProperties, typename BasisType , typename ortValueType , class... ortProperties>
static void Intrepid2::Experimental::ProjectionTools< ExecSpaceType >::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 
)
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
evaluationPoints[in] - rank-3 view (C,P,D) containing the coordinates of the evaluation points, at each cell
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 ExecSpaceType >
template<typename BasisType , typename ortValueType , class... ortProperties>
static void Intrepid2::Experimental::ProjectionTools< ExecSpaceType >::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 
)
static

Computes evaluation points for HVol projection.

C - num. cells
P - num. evaluation points
D - spatial dimension
Parameters
evaluationPoints[out] - rank-3 view (C,P,D) containing the coordinates of the evaluation points, at each cell
cellOrientations[in] - rank-1 view (C) containing the Orientation objects at each cell
cellBasis[in] - pointer to the HVOL basis for the projection
projStruct[in] - pointer to ProjectionStruct object
evalPointType[in] - enum selecting whether the points should be computed for the basis functions or for the target function
template<typename ExecSpaceType >
template<typename basisCoeffsValueType , class... basisCoeffsProperties, typename funValsValueType , class... funValsProperties, typename BasisType , typename ortValueType , class... ortProperties>
static void Intrepid2::Experimental::ProjectionTools< ExecSpaceType >::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 
)
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
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
template<typename ExecSpaceType >
template<typename BasisType , typename ortValueType , class... ortProperties>
static void Intrepid2::Experimental::ProjectionTools< ExecSpaceType >::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 
)
static

Computes evaluation points for L2 projection.

C - num. cells
P - num. evaluation points
D - spatial dimension
Parameters
evaluationPoints[out] - rank-3 view (C,P,D) containing the coordinates of the evaluation points for the projection at each cell
cellOrientations[in] - rank-1 view (C) containing the Orientation objects at each cell
cellBasis[in] - pointer to the basis for the projection
projStruct[in] - pointer to ProjectionStruct object
evalPointType[in] - enum selecting whether the points should be computed for the basis functions or for the target function

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