Intrepid2
Static Public Member Functions | List of all members
Intrepid2::Experimental::LagrangianInterpolation< ExecSpaceType > Class Template Reference

A class providing static members to perform Lagrangian interpolation on a finite element. More...

#include <Intrepid2_LagrangianInterpolation.hpp>

Static Public Member Functions

template<typename BasisType , class... coordsProperties, class... coeffsProperties, typename ortValueType , class... ortProperties>
static void getDofCoordsAndCoeffs (Kokkos::DynRankView< typename BasisType::scalarType, coordsProperties...> dofCoords, Kokkos::DynRankView< typename BasisType::scalarType, coeffsProperties...> dofCoeffs, const BasisType *cellBasis, EPointType basisPointType, const Kokkos::DynRankView< ortValueType, ortProperties...> cellOrientations)
 Computes the points and coefficients associated with the basis DOFs for the reference oriented element. More...
 
template<typename basisCoeffsViewType , typename funcViewType , typename dofCoeffViewType >
static void getBasisCoeffs (basisCoeffsViewType basisCoeffs, const funcViewType functionAtDofCoords, const dofCoeffViewType dofCoeffs)
 Computes the basis weights of the function interpolation. More...
 

Detailed Description

template<typename ExecSpaceType>
class Intrepid2::Experimental::LagrangianInterpolation< ExecSpaceType >

A class providing static members to perform Lagrangian interpolation on a finite element.

The Lagrangian interpolant $\mathcal I_h$ is defined as

\[ \mathcal I_h(f) = \sum_i \alpha^f_i \phi_i, \quad \alpha_i^f := L_i(f), \]

where $\{\phi_i\}$ is the basis of the finite element, $\alpha_i^f$ are the basisCoeffs and $L_i$ are the DOFs of the basis. This class works only for Lagrangian elements that provide a set of orthonormal DOFs defined as

\[ L_i(f) := f(\mathbf x_i) \cdot \beta_i, \quad L_i(\phi_j) = \delta_{ij}, \]

where $\beta_i$ are referred to as dofCoeffs, and $\mathbf x_i$ are the coordinates of the basis nodes.

In order to perform the interpolation, call the member getDofCoordsAndCoeffs that returns a set of points $\{\mathbf x_i\}$ and DOF coefficients $\{\beta_i\}$, evaluate the $f$ at $\{\mathbf x_i\}$ and then obtain the basis coefficients $\alpha_i^f$ by calling the function getBasisCoeffs.

Remarks
The interpolation is performed at the oriented reference element. Therefore, the function $f$, which is contravariant, needs to mapped to the reference space, using the inverse operation of a pullback, before calling getBasisCoeffs.
Todo:
The implementation is mostly serial and needs to be improved for performance portability

Definition at line 136 of file Intrepid2_LagrangianInterpolation.hpp.

Member Function Documentation

template<typename SpT >
template<typename basisCoeffsViewType , typename funcViewType , typename dofCoeffViewType >
void Intrepid2::Experimental::LagrangianInterpolation< SpT >::getBasisCoeffs ( basisCoeffsViewType  basisCoeffs,
const funcViewType  functionAtDofCoords,
const dofCoeffViewType  dofCoeffs 
)
static

Computes the basis weights of the function interpolation.

C - num. cells
F - num. fields
D - spatial dimension
Parameters
basisCoeffs[out] - rank-2 view (C,F) that will contain the basis coefficients of the interpolation.
functionAtDofCoords[in] - variable rank view that contains the function evaluated at DOF coordinates.
dofCoeffs[in] - variable rank view that contains coefficients associated with the basis DOFs.
Remarks
The output views need to be pre-allocated. dofCoeffs and functionAtDofCoords have rank 2, (C,F) for scalar basis and 3, (C,F,D) for vector basis. functValsAtDofCoords contains the function evaluated at the dofCoords and contravariantly transformed to the reference eleemnt. dofCoeffs are as returned by getDofCoordsAndCoeffs.

Definition at line 507 of file Intrepid2_LagrangianInterpolationDef.hpp.

References Intrepid2::ArrayTools< ExecSpaceType >::dotMultiplyDataData().

template<typename SpT >
template<typename BasisType , class... coordsProperties, class... coeffsProperties, typename ortValueType , class... ortProperties>
void Intrepid2::Experimental::LagrangianInterpolation< SpT >::getDofCoordsAndCoeffs ( Kokkos::DynRankView< typename BasisType::scalarType, coordsProperties...>  dofCoords,
Kokkos::DynRankView< typename BasisType::scalarType, coeffsProperties...>  dofCoeffs,
const BasisType *  cellBasis,
EPointType  basisPointType,
const Kokkos::DynRankView< ortValueType, ortProperties...>  cellOrientations 
)
static

Computes the points and coefficients associated with the basis DOFs for the reference oriented element.

C - num. cells
F - num. fields
D - spatial dimension
Parameters
dofCoords[out] - rank-3 view (C,F,D), that will contain coordinates associated with the basis DOFs.
dofCoeffs[out] - variable rank view that will contain coefficients associated with the basis DOFs.
cellBasis[in] - pointer to the basis for the interpolation
basisPointType[in] - enum of the point type functions or for the target function
cellOrientations[in] - rank-1 view (C) containing the Orientation objects at each cell
Remarks
the output views need to be pre-allocated. dofCoeffs has rank 2, (C,F) for scalar basis and 3, (C,F,D) for vector basis.

Definition at line 275 of file Intrepid2_LagrangianInterpolationDef.hpp.

References Intrepid2::RealSpaceTools< ExecSpaceType >::clone(), Intrepid2::CellTools< ExecSpaceType >::getReferenceEdgeTangent(), Intrepid2::CellTools< ExecSpaceType >::getReferenceFaceNormal(), Intrepid2::CellTools< ExecSpaceType >::getReferenceFaceTangents(), Intrepid2::CellTools< ExecSpaceType >::getReferenceSideNormal(), and Intrepid2::CellTools< ExecSpaceType >::getSubcellParametrization().


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