49 #ifndef __INTREPID2_HGRAD_PYR_I2_SERENDIPITY_FEM_HPP__
50 #define __INTREPID2_HGRAD_PYR_I2_SERENDIPITY_FEM_HPP__
110 typedef struct Pyramid<5> cell_topology_type;
114 template<EOperator opType>
116 template<
typename OutputViewType,
117 typename inputViewType>
118 KOKKOS_INLINE_FUNCTION
120 getValues( OutputViewType output,
121 const inputViewType input );
125 template<
typename DeviceType,
126 typename outputValueValueType,
class ...outputValueProperties,
127 typename inputPointValueType,
class ...inputPointProperties>
129 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
130 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
131 const EOperator operatorType);
136 template<
typename outputValueViewType,
137 typename inputPointViewType,
140 outputValueViewType _outputValues;
141 const inputPointViewType _inputPoints;
143 KOKKOS_INLINE_FUNCTION
144 Functor( outputValueViewType outputValues_,
145 inputPointViewType inputPoints_ )
146 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
148 KOKKOS_INLINE_FUNCTION
149 void operator()(
const ordinal_type pt)
const {
151 case OPERATOR_VALUE : {
152 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
153 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
160 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
161 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
166 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
167 opType != OPERATOR_GRAD &&
168 opType != OPERATOR_CURL &&
169 opType != OPERATOR_D2,
170 ">>> ERROR: (Intrepid2::Basis_HGRAD_PYR_I2_FEM::Serial::getValues) operator is not supported");
178 template<
typename DeviceType = void,
179 typename outputValueType = double,
180 typename pointValueType =
double>
199 getValues( OutputViewType outputValues,
200 const PointViewType inputPoints,
201 const EOperator operatorType = OPERATOR_VALUE )
const override {
202 #ifdef HAVE_INTREPID2_DEBUG
204 Intrepid2::getValues_HGRAD_Args(outputValues,
210 Impl::Basis_HGRAD_PYR_I2_FEM::
211 getValues<DeviceType>( outputValues,
218 getDofCoords( ScalarViewType dofCoords )
const override {
219 #ifdef HAVE_INTREPID2_DEBUG
221 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
222 ">>> ERROR: (Intrepid2::Basis_HGRAD_PYR_I2_FEM::getDofCoords) rank = 2 required for dofCoords array");
224 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
225 ">>> ERROR: (Intrepid2::Basis_HGRAD_PYR_I2_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
227 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
228 ">>> ERROR: (Intrepid2::Basis_HGRAD_PYR_I2_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
230 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
235 getDofCoeffs( ScalarViewType dofCoeffs )
const override {
236 #ifdef HAVE_INTREPID2_DEBUG
238 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
239 ">>> ERROR: (Intrepid2::Basis_HGRAD_PYR_I2_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
241 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
242 ">>> ERROR: (Intrepid2::Basis_HGRAD_PYR_I2_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
244 Kokkos::deep_copy(dofCoeffs, 1.0);
250 return "Intrepid2_HGRAD_PYR_I2_FEM";
261 BasisPtr<DeviceType,outputValueType,pointValueType>
263 if(subCellDim == 1) {
264 return Teuchos::rcp(
new
266 }
else if(subCellDim == 2) {
268 return Teuchos::rcp(
new
271 return Teuchos::rcp(
new
274 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
277 BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
See Intrepid2::Basis_HGRAD_PYR_I2_FEM.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Implementation of an H(grad)-compatible FEM basis of degree 2 on a Pyramid cell.
See Intrepid2::Basis_HGRAD_PYR_I2_FEM.
Header file for the Intrepid2::Basis_HGRAD_QUAD_C2_FEM class.
Definition file for FEM basis functions of degree 1 for H(grad) functions on PYR cells.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Basis_HGRAD_PYR_I2_FEM()
Constructor.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Triangle cell...
Header file for the Intrepid2::Basis_HGRAD_TRI_C2_FEM class.
See Intrepid2::Basis_HGRAD_PYR_I2_FEM.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Line 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...
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Quadrilateral cell...
virtual const char * getName() const override
Returns basis name.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Header file for the abstract base class Intrepid2::Basis.