Intrepid2
|
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... | |
A class providing static members to perform Lagrangian interpolation on a finite element.
The Lagrangian interpolant is defined as
where is the basis of the finite element, are the basisCoeffs and are the DOFs of the basis. This class works only for Lagrangian elements that provide a set of orthonormal DOFs defined as
where are referred to as dofCoeffs, and are the coordinates of the basis nodes.
In order to perform the interpolation, call the member getDofCoordsAndCoeffs that returns a set of points and DOF coefficients , evaluate the at and then obtain the basis coefficients by calling the function getBasisCoeffs.
Definition at line 136 of file Intrepid2_LagrangianInterpolation.hpp.
|
static |
Computes the basis weights of the function interpolation.
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. |
Definition at line 507 of file Intrepid2_LagrangianInterpolationDef.hpp.
References Intrepid2::ArrayTools< ExecSpaceType >::dotMultiplyDataData().
|
static |
Computes the points and coefficients associated with the basis DOFs for the reference oriented element.
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 |
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().