16 #ifndef __INTREPID2_HGRAD_PYR_C1_FEM_HPP__
17 #define __INTREPID2_HGRAD_PYR_C1_FEM_HPP__
60 typedef struct Pyramid<5> cell_topology_type;
64 template<EOperator opType>
66 template<
typename OutputViewType,
67 typename inputViewType>
68 KOKKOS_INLINE_FUNCTION
70 getValues( OutputViewType output,
71 const inputViewType input );
75 template<
typename DeviceType,
76 typename outputValueValueType,
class ...outputValueProperties,
77 typename inputPointValueType,
class ...inputPointProperties>
79 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
80 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
81 const EOperator operatorType);
86 template<
typename outputValueViewType,
87 typename inputPointViewType,
90 outputValueViewType _outputValues;
91 const inputPointViewType _inputPoints;
93 KOKKOS_INLINE_FUNCTION
94 Functor( outputValueViewType outputValues_,
95 inputPointViewType inputPoints_ )
96 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
98 KOKKOS_INLINE_FUNCTION
99 void operator()(
const ordinal_type pt)
const {
101 case OPERATOR_VALUE : {
102 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
103 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
110 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
111 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
116 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
117 opType != OPERATOR_GRAD &&
118 opType != OPERATOR_CURL &&
119 opType != OPERATOR_D2,
120 ">>> ERROR: (Intrepid2::Basis_HGRAD_PYR_C1_FEM::Serial::getValues) operator is not supported");
128 template<
typename DeviceType = void,
129 typename outputValueType = double,
130 typename pointValueType =
double>
149 getValues( OutputViewType outputValues,
150 const PointViewType inputPoints,
151 const EOperator operatorType = OPERATOR_VALUE )
const override {
152 #ifdef HAVE_INTREPID2_DEBUG
154 Intrepid2::getValues_HGRAD_Args(outputValues,
160 Impl::Basis_HGRAD_PYR_C1_FEM::
161 getValues<DeviceType>( outputValues,
167 getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
168 ordinal_type& perThreadSpaceSize,
169 const PointViewType inputPointsconst,
170 const EOperator operatorType = OPERATOR_VALUE)
const override;
172 KOKKOS_INLINE_FUNCTION
175 OutputViewType outputValues,
176 const PointViewType inputPoints,
177 const EOperator operatorType,
178 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
179 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
180 const ordinal_type subcellDim = -1,
181 const ordinal_type subcellOrdinal = -1)
const override;
185 getDofCoords( ScalarViewType dofCoords )
const override {
186 #ifdef HAVE_INTREPID2_DEBUG
188 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
189 ">>> ERROR: (Intrepid2::Basis_HGRAD_PYR_C1_FEM::getDofCoords) rank = 2 required for dofCoords array");
191 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
192 ">>> ERROR: (Intrepid2::Basis_HGRAD_PYR_C1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
194 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
195 ">>> ERROR: (Intrepid2::Basis_HGRAD_PYR_C1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
197 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
202 getDofCoeffs( ScalarViewType dofCoeffs )
const override {
203 #ifdef HAVE_INTREPID2_DEBUG
205 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
206 ">>> ERROR: (Intrepid2::Basis_HGRAD_PYR_C1_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
208 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
209 ">>> ERROR: (Intrepid2::Basis_HGRAD_PYR_C1_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
211 Kokkos::deep_copy(dofCoeffs, 1.0);
217 return "Intrepid2_HGRAD_PYR_C1_FEM";
228 BasisPtr<DeviceType,outputValueType,pointValueType>
230 if(subCellDim == 1) {
231 return Teuchos::rcp(
new
233 }
else if(subCellDim == 2) {
235 return Teuchos::rcp(
new
238 return Teuchos::rcp(
new
241 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
244 BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Quadrilateral cell...
Header file for the Intrepid2::Basis_HGRAD_TRI_C1_FEM class.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Header file for the Intrepid2::Basis_HGRAD_QUAD_C1_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Triangle cell...
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Pyramid cell...
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
See Intrepid2::Basis_HGRAD_PYR_C1_FEM.
Basis_HGRAD_PYR_C1_FEM()
Constructor.
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...
Definition file for FEM basis functions of degree 1 for H(grad) functions on PYR cells.
See Intrepid2::Basis_HGRAD_PYR_C1_FEM.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Line cell.
See Intrepid2::Basis_HGRAD_PYR_C1_FEM.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
virtual const char * getName() const override
Returns basis name.
Header file for the abstract base class Intrepid2::Basis.