16 #ifndef __INTREPID2_HGRAD_LINE_C1_FEM_HPP__
17 #define __INTREPID2_HGRAD_LINE_C1_FEM_HPP__
52 typedef struct Line<2> cell_topology_type;
56 template<EOperator opType>
58 template<
typename OutputViewType,
59 typename inputViewType>
60 KOKKOS_INLINE_FUNCTION
62 getValues( OutputViewType output,
63 const inputViewType input );
67 template<
typename DeviceType,
68 typename outputValueValueType,
class ...outputValueProperties,
69 typename inputPointValueType,
class ...inputPointProperties>
71 getValues(
const typename DeviceType::execution_space& space,
72 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
73 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
74 const EOperator operatorType);
79 template<
typename outputValueViewType,
80 typename inputPointViewType,
83 outputValueViewType _outputValues;
84 const inputPointViewType _inputPoints;
86 KOKKOS_INLINE_FUNCTION
87 Functor( outputValueViewType outputValues_,
88 inputPointViewType inputPoints_ )
89 : _outputValues(outputValues_), _inputPoints(inputPoints_) {}
91 KOKKOS_INLINE_FUNCTION
92 void operator()(
const ordinal_type pt)
const {
94 case OPERATOR_VALUE : {
95 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt );
96 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
101 case OPERATOR_MAX : {
102 auto output = Kokkos::subview( _outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
103 const auto input = Kokkos::subview( _inputPoints, pt, Kokkos::ALL() );
108 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
109 opType != OPERATOR_GRAD &&
110 opType != OPERATOR_MAX,
111 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C1_FEM::Serial::getValues) operator is not supported");
120 template<
typename DeviceType = void,
121 typename outputValueType = double,
122 typename pointValueType =
double>
124 :
public Basis<DeviceType,outputValueType,pointValueType> {
148 const EOperator operatorType = OPERATOR_VALUE )
const override {
149 #ifdef HAVE_INTREPID2_DEBUG
151 Intrepid2::getValues_HGRAD_Args(outputValues,
157 Impl::Basis_HGRAD_LINE_C1_FEM::
158 getValues<DeviceType>(space,
166 ordinal_type& perThreadSpaceSize,
168 const EOperator operatorType = OPERATOR_VALUE)
const override;
170 KOKKOS_INLINE_FUNCTION
175 const EOperator operatorType,
176 const typename Kokkos::TeamPolicy<typename DeviceType::execution_space>::member_type& team_member,
177 const typename DeviceType::execution_space::scratch_memory_space & scratchStorage,
178 const ordinal_type subcellDim = -1,
179 const ordinal_type subcellOrdinal = -1)
const override;
184 #ifdef HAVE_INTREPID2_DEBUG
186 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoords) != 2, std::invalid_argument,
187 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C1_FEM::getDofCoords) rank = 2 required for dofCoords array");
189 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoords.extent(0)) != this->
basisCardinality_, std::invalid_argument,
190 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C1_FEM::getDofCoords) mismatch in number of dof and 0th dimension of dofCoords array");
192 INTREPID2_TEST_FOR_EXCEPTION( dofCoords.extent(1) != this->
getBaseCellTopology().getDimension(), std::invalid_argument,
193 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C1_FEM::getDofCoords) incorrect reference cell (1st) dimension in dofCoords array");
195 Kokkos::deep_copy(dofCoords, this->
dofCoords_);
201 #ifdef HAVE_INTREPID2_DEBUG
203 INTREPID2_TEST_FOR_EXCEPTION( rank(dofCoeffs) != 1, std::invalid_argument,
204 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C1_FEM::getdofCoeffs) rank = 1 required for dofCoeffs array");
206 INTREPID2_TEST_FOR_EXCEPTION( static_cast<ordinal_type>(dofCoeffs.extent(0)) != this->
getCardinality(), std::invalid_argument,
207 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C1_FEM::getdofCoeffs) mismatch in number of dof and 0th dimension of dofCoeffs array");
209 Kokkos::deep_copy(dofCoeffs, 1.0);
215 return "Intrepid2_HGRAD_LINE_C1_FEM";
218 BasisPtr<typename Kokkos::HostSpace::device_type,outputValueType,pointValueType>
See Intrepid2::Basis_HGRAD_LINE_C1_FEM.
virtual void getDofCoeffs(ScalarViewType dofCoeffs) const override
Coefficients for computing degrees of freedom for Lagrangian basis If P is an element of the space sp...
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Definition file for FEM basis functions of degree 1 for H(grad) functions on a Line.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
virtual void getDofCoords(ScalarViewType dofCoords) const override
Returns spatial locations (coordinates) of degrees of freedom on the reference cell.
shards::CellTopology getBaseCellTopology() const
Returns the base cell topology for which the basis is defined. See Shards documentation https://trili...
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
virtual KOKKOS_INLINE_FUNCTION void getValues(OutputViewType, const PointViewType, const EOperator, const typename Kokkos::TeamPolicy< ExecutionSpace >::member_type &teamMember, const typename ExecutionSpace::scratch_memory_space &scratchStorage, const ordinal_type subcellDim=-1, const ordinal_type subcellOrdinal=-1) const
Team-level evaluation of basis functions on a reference cell.
virtual void getValues(const ExecutionSpace &space, OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference 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...
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
View type for basis value output.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom. ...
See Intrepid2::Basis_HGRAD_LINE_C1_FEM.
Kokkos::View< ordinal_type ***, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray3DHost
View type for 3d host array.
virtual const char * getName() const override
Returns basis name.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Line cell.
virtual void getScratchSpaceSize(ordinal_type &perTeamSpaceSize, ordinal_type &perThreadSpaceSize, const PointViewType inputPointsconst, const EOperator operatorType=OPERATOR_VALUE) const override
Return the size of the scratch space, in bytes, needed for using the team-level implementation of get...
Basis_HGRAD_LINE_C1_FEM()
Constructor.
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
Coordinates of degrees-of-freedom for basis functions defined in physical space.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Header file for the abstract base class Intrepid2::Basis.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
See Intrepid2::Basis_HGRAD_LINE_C1_FEM.