16#ifndef __INTREPID2_HGRAD_TRI_C1_FEM_HPP__
17#define __INTREPID2_HGRAD_TRI_C1_FEM_HPP__
58 template<EOperator opType>
60 template<
typename OutputViewType,
61 typename inputViewType>
62 KOKKOS_INLINE_FUNCTION
64 getValues( OutputViewType output,
65 const inputViewType input );
69 template<
typename DeviceType,
70 typename outputValueValueType,
class ...outputValueProperties,
71 typename inputPointValueType,
class ...inputPointProperties>
73 getValues(
const typename DeviceType::execution_space& space,
74 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
75 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
76 const EOperator operatorType);
81 template<
typename outputValueViewType,
82 typename inputPointViewType,
85 outputValueViewType _outputValues;
86 const inputPointViewType _inputPoints;
88 KOKKOS_INLINE_FUNCTION
89 Functor( outputValueViewType outputValues_,
90 inputPointViewType inputPoints_ )
91 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
93 KOKKOS_INLINE_FUNCTION
94 void operator()(
const ordinal_type pt)
const {
96 case OPERATOR_VALUE : {
97 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
98 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
104 case OPERATOR_MAX : {
105 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
106 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
111 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
112 opType != OPERATOR_GRAD &&
113 opType != OPERATOR_CURL &&
114 opType != OPERATOR_MAX,
115 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C1_FEM::Serial::getValues) operator is not supported");
123 template<
typename DeviceType = void,
124 typename outputValueType = double,
125 typename pointValueType =
double>
149 const EOperator operatorType = OPERATOR_VALUE )
const override {
150#ifdef HAVE_INTREPID2_DEBUG
158 Impl::Basis_HGRAD_TRI_C1_FEM::
159 getValues<DeviceType>(space,
167 ordinal_type& perThreadSpaceSize,
169 const EOperator operatorType = OPERATOR_VALUE)
const override;
171 KOKKOS_INLINE_FUNCTION
176 const EOperator operatorType,
177 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
178 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
179 const ordinal_type subcellDim = -1,
180 const ordinal_type subcellOrdinal = -1)
const override;
185#ifdef HAVE_INTREPID2_DEBUG
187 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
188 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C1_FEM::getDofCoords) rank = 2 required for dofCoords array");
190 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
191 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
193 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
194 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
196 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
202#ifdef HAVE_INTREPID2_DEBUG
204 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
205 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C1_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
207 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
208 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C1_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
210 Kokkos::deep_copy(dofCoeffs, 1.0);
216 return "Intrepid2_HGRAD_TRI_C1_FEM";
229 if(subCellDim == 1) {
230 return Teuchos::rcp(
new
233 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
void getValues_HGRAD_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 HGRAD-conforming FEM basis....
Header file for the Intrepid2::Basis_HGRAD_LINE_C1_FEM class.
Definition file for FEM basis functions of degree 1 for H(grad) functions on TRI cells.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Line cell.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Triangle cell.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
virtual void getValues(const ExecutionSpace &space, OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
BasisPtr< typename Kokkos::HostSpace::device_type, outputValueType, pointValueType > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...
virtual const char * getName() const override
Returns basis name.
Basis_HGRAD_TRI_C1_FEM()
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 void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
virtual void getScratchSpaceSize(ordinal_type &perTeamSpaceSize, ordinal_type &perThreadSpaceSize, const PointViewType inputPointsconst, const EOperator operatorType=OPERATOR_VALUE) const override
Return the size of the scratch space, in bytes, needed for using the team-level implementation of get...
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.
virtual KOKKOS_INLINE_FUNCTION void getValues(OutputViewType, const PointViewType, const EOperator, const typename Kokkos::TeamPolicy< ExecutionSpace >::member_type &teamMember, const typename ExecutionSpace::scratch_memory_space &scratchStorage, const ordinal_type subcellDim=-1, const ordinal_type subcellOrdinal=-1) const
Team-level evaluation of basis functions on a reference cell.
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 getCardinality() const
Returns cardinality of the basis.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
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.
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.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
See Intrepid2::Basis_HGRAD_TRI_C1_FEM.
See Intrepid2::Basis_HGRAD_TRI_C1_FEM.
See Intrepid2::Basis_HGRAD_TRI_C1_FEM.
Triangle topology, 3 nodes.