48#ifndef __INTREPID2_HVOL_TRI_CN_FEM_HPP__
49#define __INTREPID2_HVOL_TRI_CN_FEM_HPP__
54#include "Teuchos_LAPACK.hpp"
83 template<EOperator opType>
99 getWorkSizePerPoint(ordinal_type
order) {
117 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...>
outputValues,
118 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...>
inputPoints,
119 const Kokkos::DynRankView<vinvValueType, vinvProperties...>
vinv,
143 _vinv(vinv_), _work(
work_) {}
146 void operator()(
const size_type
iter)
const {
151 const auto input = Kokkos::subview( _inputPoints,
ptRange, Kokkos::ALL() );
153 typename workViewType::pointer_type
ptr = _work.data() + _work.extent(0)*
ptBegin*get_dimension_scalar(_work);
155 auto vcprop = Kokkos::common_view_alloc_prop(_work);
159 case OPERATOR_VALUE : {
160 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(),
ptRange );
166 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(),
ptRange, Kokkos::ALL() );
171 INTREPID2_TEST_FOR_ABORT(
true,
172 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::Functor) operator is not supported");
181 template<
typename DeviceType = void,
182 typename outputValueType = double,
183 typename pointValueType =
double>
185 :
public Basis<DeviceType,outputValueType,pointValueType> {
196 const EPointType pointType = POINTTYPE_EQUISPACED);
210 const PointViewType inputPoints,
211 const EOperator operatorType = OPERATOR_VALUE)
const override {
212#ifdef HAVE_INTREPID2_DEBUG
220 Impl::Basis_HVOL_TRI_Cn_FEM::
221 getValues<DeviceType,numPtsPerEval>( outputValues,
230#ifdef HAVE_INTREPID2_DEBUG
232 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
233 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getDofCoords) rank = 2 required for dofCoords array");
235 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
236 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
238 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
239 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
241 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
247#ifdef HAVE_INTREPID2_DEBUG
249 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
250 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
252 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
253 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
255 Kokkos::deep_copy(dofCoeffs, 1.0);
259 getVandermondeInverse( ScalarViewType vinv )
const {
261 Kokkos::deep_copy(vinv, this->
vinv_);
267 return "Intrepid2_HVOL_TRI_Cn_FEM";
285 Kokkos::DynRankView<scalarType,DeviceType>
vinv_;
286 EPointType pointType_;
Header file for the abstract base class Intrepid2::Basis.
BasisPtr< typename Kokkos::HostSpace::device_type, OutputType, PointType > HostBasisPtr
Pointer to a Basis whose device type is on the host (Kokkos::HostSpace::device_type),...
void getValues_HVOL_Args(const outputValueViewType outputValues, const inputPointViewType inputPoints, const EOperator operatorType, const shards::CellTopology cellTopo, const ordinal_type basisCard)
Runtime check of the arguments for the getValues method in an HVOL-conforming FEM basis....
Definition file for FEM basis functions of degree n for H(vol) functions on TRI.
Implementation of the default HVOL-compatible Lagrange basis of arbitrary degree on Triangle cell.
virtual const char * getName() const override
Returns basis name.
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
virtual HostBasisPtr< outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
Basis_HVOL_TRI_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
virtual void getDofCoeffs(ScalarViewType dofCoeffs) const override
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
virtual bool requireOrientation() const override
True if orientation is required.
Kokkos::DynRankView< scalarType, DeviceType > vinv_
inverse of Generalized Vandermonde matrix, whose columns store the expansion coefficients of the noda...
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Device DeviceType
(Kokkos) Device type on which Basis is templated. Does not necessarily return true for Kokkos::is_dev...
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
See Intrepid2::Basis_HVOL_TRI_Cn_FEM.
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
See Intrepid2::Basis_HVOL_TRI_Cn_FEM.
See Intrepid2::Basis_HVOL_TRI_Cn_FEM.
Triangle topology, 3 nodes.