16#ifndef __INTREPID2_HGRAD_WEDGE_C1_FEM_HPP__
17#define __INTREPID2_HGRAD_WEDGE_C1_FEM_HPP__
66 template<EOperator opType>
68 template<
typename OutputViewType,
69 typename inputViewType>
70 KOKKOS_INLINE_FUNCTION
72 getValues( OutputViewType output,
73 const inputViewType input );
77 template<
typename DeviceType,
78 typename outputValueValueType,
class ...outputValueProperties,
79 typename inputPointValueType,
class ...inputPointProperties>
81 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
82 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
83 const EOperator operatorType);
88 template<
typename outputValueViewType,
89 typename inputPointViewType,
92 outputValueViewType _outputValues;
93 const inputPointViewType _inputPoints;
95 KOKKOS_INLINE_FUNCTION
96 Functor( outputValueViewType outputValues_,
97 inputPointViewType inputPoints_ )
98 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
100 KOKKOS_INLINE_FUNCTION
101 void operator()(
const ordinal_type pt)
const {
103 case OPERATOR_VALUE : {
104 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
105 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
111 case OPERATOR_MAX : {
112 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
113 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
118 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
119 opType != OPERATOR_GRAD &&
120 opType != OPERATOR_D2 &&
121 opType != OPERATOR_MAX,
122 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::Serial::getValues) operator is not supported");
131 template<
typename DeviceType = void,
132 typename outputValueType = double,
133 typename pointValueType =
double>
152 getValues( OutputViewType outputValues,
153 const PointViewType inputPoints,
154 const EOperator operatorType = OPERATOR_VALUE )
const override {
155#ifdef HAVE_INTREPID2_DEBUG
163 Impl::Basis_HGRAD_WEDGE_C1_FEM::
164 getValues<DeviceType>( outputValues,
170 getScratchSpaceSize( ordinal_type& perTeamSpaceSize,
171 ordinal_type& perThreadSpaceSize,
172 const PointViewType inputPointsconst,
173 const EOperator operatorType = OPERATOR_VALUE)
const override;
175 KOKKOS_INLINE_FUNCTION
178 OutputViewType outputValues,
179 const PointViewType inputPoints,
180 const EOperator operatorType,
181 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
182 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
183 const ordinal_type subcellDim = -1,
184 const ordinal_type subcellOrdinal = -1)
const override;
188 getDofCoords( ScalarViewType dofCoords )
const override {
189#ifdef HAVE_INTREPID2_DEBUG
191 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.rank() != 2, std::invalid_argument,
192 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getDofCoords) rank = 2 required for dofCoords array");
194 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoords.extent(0)) != this->getCardinality(), std::invalid_argument,
195 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
197 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->getBaseCellTopology().getDimension(), std::invalid_argument,
198 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
200 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
205 getDofCoeffs( ScalarViewType dofCoeffs )
const override {
206#ifdef HAVE_INTREPID2_DEBUG
208 INTREPID2_TEST_FOR_EXCEPTION( dofCoeffs.rank() != 1, std::invalid_argument,
209 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
211 INTREPID2_TEST_FOR_EXCEPTION(
static_cast<ordinal_type
>(dofCoeffs.extent(0)) != this->getCardinality(), std::invalid_argument,
212 ">>> ERROR: (Intrepid2::Basis_HGRAD_WEDGE_C1_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
214 Kokkos::deep_copy(dofCoeffs, 1.0);
220 return "Intrepid2_HGRAD_WEDGE_C1_FEM";
233 if(subCellDim == 1) {
234 return Teuchos::rcp(
new
236 }
else if(subCellDim == 2) {
242 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_QUAD_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TRI_C1_FEM class.
Definition file for FEM basis functions of degree 1 for H(grad) functions on WEDGE 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 Quadrilateral cell.
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 Wedge cell.
Basis_HGRAD_WEDGE_C1_FEM()
Constructor.
BasisPtr< DeviceType, outputValueType, pointValueType > getSubCellRefBasis(const ordinal_type subCellDim, const ordinal_type subCellOrd) const override
returns the basis associated to a subCell.
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.
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.
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.
Device DeviceType
(Kokkos) Device type on which Basis is templated. Does not necessarily return true for Kokkos::is_dev...
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.
See Intrepid2::Basis_HGRAD_WEDGE_C1_FEM.
See Intrepid2::Basis_HGRAD_WEDGE_C1_FEM.
See Intrepid2::Basis_HGRAD_WEDGE_C1_FEM.