49 #ifndef __INTREPID2_HGRAD_WEDGE_C1_FEM_HPP__
50 #define __INTREPID2_HGRAD_WEDGE_C1_FEM_HPP__
95 typedef struct Wedge<6> cell_topology_type;
99 template<EOperator opType>
101 template<
typename OutputViewType,
102 typename inputViewType>
103 KOKKOS_INLINE_FUNCTION
105 getValues( OutputViewType output,
106 const inputViewType input );
110 template<
typename DeviceType,
111 typename outputValueValueType,
class ...outputValueProperties,
112 typename inputPointValueType,
class ...inputPointProperties>
114 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
115 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
116 const EOperator operatorType);
121 template<
typename outputValueViewType,
122 typename inputPointViewType,
125 outputValueViewType _outputValues;
126 const inputPointViewType _inputPoints;
128 KOKKOS_INLINE_FUNCTION
129 Functor( outputValueViewType outputValues_,
130 inputPointViewType inputPoints_ )
131 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
133 KOKKOS_INLINE_FUNCTION
134 void operator()(
const ordinal_type pt)
const {
136 case OPERATOR_VALUE : {
137 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
138 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
144 case OPERATOR_MAX : {
145 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
146 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
151 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
152 opType != OPERATOR_GRAD &&
153 opType != OPERATOR_D2 &&
154 opType != OPERATOR_MAX,
155 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::Serial::getValues) operator is not supported");
164 template<
typename DeviceType = void,
165 typename outputValueType = double,
166 typename pointValueType =
double>
185 getValues( OutputViewType outputValues,
186 const PointViewType inputPoints,
187 const EOperator operatorType = OPERATOR_VALUE )
const override {
188 #ifdef HAVE_INTREPID2_DEBUG
190 Intrepid2::getValues_HGRAD_Args(outputValues,
196 Impl::Basis_HGRAD_WEDGE_C1_FEM::
197 getValues<DeviceType>( outputValues,
204 getDofCoords( ScalarViewType dofCoords )
const override {
205 #ifdef HAVE_INTREPID2_DEBUG
207 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
208 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getDofCoords) rank = 2 required for dofCoords array");
210 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
211 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
213 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
214 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
216 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
221 getDofCoeffs( ScalarViewType dofCoeffs )
const override {
222 #ifdef HAVE_INTREPID2_DEBUG
224 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
225 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
227 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
228 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
230 Kokkos::deep_copy(dofCoeffs, 1.0);
236 return "Intrepid2_HGRAD_WEDGE_C1_FEM";
247 BasisPtr<DeviceType,outputValueType,pointValueType>
249 if(subCellDim == 1) {
250 return Teuchos::rcp(
new
252 }
else if(subCellDim == 2) {
258 INTREPID2_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
"Input parameters out of bounds");
261 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.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Wedge cell.
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...
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...
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...
See Intrepid2::Basis_HGRAD_WEDGE_C1_FEM.
virtual const char * getName() const override
Returns basis name.
Definition file for FEM basis functions of degree 1 for H(grad) functions on WEDGE cells...
See Intrepid2::Basis_HGRAD_WEDGE_C1_FEM.
See Intrepid2::Basis_HGRAD_WEDGE_C1_FEM.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
Basis_HGRAD_WEDGE_C1_FEM()
Constructor.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Line cell.
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.