49 #ifndef __INTREPID2_HGRAD_WEDGE_C1_FEM_HPP__
50 #define __INTREPID2_HGRAD_WEDGE_C1_FEM_HPP__
93 typedef struct Wedge<6> cell_topology_type;
97 template<EOperator opType>
99 template<
typename OutputViewType,
100 typename inputViewType>
101 KOKKOS_INLINE_FUNCTION
103 getValues( OutputViewType output,
104 const inputViewType input );
108 template<
typename ExecSpaceType,
109 typename outputValueValueType,
class ...outputValueProperties,
110 typename inputPointValueType,
class ...inputPointProperties>
112 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
113 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
114 const EOperator operatorType);
119 template<
typename outputValueViewType,
120 typename inputPointViewType,
123 outputValueViewType _outputValues;
124 const inputPointViewType _inputPoints;
126 KOKKOS_INLINE_FUNCTION
127 Functor( outputValueViewType outputValues_,
128 inputPointViewType inputPoints_ )
129 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
131 KOKKOS_INLINE_FUNCTION
132 void operator()(
const ordinal_type pt)
const {
134 case OPERATOR_VALUE : {
135 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
136 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
142 case OPERATOR_MAX : {
143 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
144 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
149 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
150 opType != OPERATOR_GRAD &&
151 opType != OPERATOR_D2 &&
152 opType != OPERATOR_MAX,
153 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::Serial::getValues) operator is not supported");
162 template<
typename ExecSpaceType = void,
163 typename outputValueType = double,
164 typename pointValueType =
double>
185 const EOperator operatorType = OPERATOR_VALUE )
const {
186 #ifdef HAVE_INTREPID2_DEBUG
188 Intrepid2::getValues_HGRAD_Args(outputValues,
194 Impl::Basis_HGRAD_WEDGE_C1_FEM::
195 getValues<ExecSpaceType>( outputValues,
203 #ifdef HAVE_INTREPID2_DEBUG
205 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
206 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getDofCoords) rank = 2 required for dofCoords array");
208 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
getCardinality(), std::invalid_argument,
209 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
211 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
212 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
214 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
220 #ifdef HAVE_INTREPID2_DEBUG
222 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
223 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
225 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
226 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
228 Kokkos::deep_copy(dofCoeffs, 1.0);
234 return "Intrepid2_HGRAD_WEDGE_C1_FEM";
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Wedge cell.
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Kokkos::View< ordinal_type **, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
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...
Kokkos::DynRankView< scalarType, ExecSpaceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
See Intrepid2::Basis_HGRAD_WEDGE_C1_FEM.
ordinal_type getCardinality() const
Returns cardinality of the basis.
virtual void getDofCoeffs(ScalarViewType dofCoeffs) const
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
virtual void getDofCoords(ScalarViewType dofCoords) const
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
Definition file for FEM basis functions of degree 1 for H(grad) functions on WEDGE cells...
See Intrepid2::Basis_HGRAD_WEDGE_C1_FEM.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, ExecSpaceType > OutputViewType
View type for basis value output.
virtual const char * getName() const
Returns basis name.
Kokkos::View< ordinal_type ***, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
See Intrepid2::Basis_HGRAD_WEDGE_C1_FEM.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > PointViewType
View type for input points.
Basis_HGRAD_WEDGE_C1_FEM()
Constructor.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > ScalarViewType
View type for scalars.
Header file for the abstract base class Intrepid2::Basis.