51 #ifndef Intrepid2_HierarchicalBasis_HDIV_TRI_h 52 #define Intrepid2_HierarchicalBasis_HDIV_TRI_h 54 #include <Kokkos_DynRankView.hpp> 56 #include <Intrepid2_config.h> 74 template<
typename DeviceType,
75 typename OutputScalar = double,
76 typename PointScalar = double,
77 bool useCGBasis =
true>
112 return "Intrepid2_HierarchicalBasis_HDIV_TRI";
145 const EOperator operatorType = OPERATOR_VALUE )
const override 148 if (operatorType == OPERATOR_DIV)
152 else if (operatorType == OPERATOR_VALUE)
156 const ordinal_type numFields = outputValues.extent_int(0);
157 const ordinal_type numPoints = outputValues.extent_int(1);
159 auto policy = Kokkos::MDRangePolicy<ExecutionSpace,Kokkos::Rank<2>>(Kokkos::Array<int,2>{0,0},Kokkos::Array<int,2>{numFields,numPoints});
162 Kokkos::parallel_for(
"apply rotation to H(curl) values", policy,
163 KOKKOS_LAMBDA (
const int &fieldOrdinal,
const int &pointOrdinal) {
164 const auto tempValue = outputValues(fieldOrdinal,pointOrdinal,0);
165 outputValues(fieldOrdinal,pointOrdinal,0) = outputValues(fieldOrdinal,pointOrdinal,1);
166 outputValues(fieldOrdinal,pointOrdinal,1) = -tempValue;
177 using HostDeviceType =
typename Kokkos::HostSpace::device_type;
179 return Teuchos::rcp(
new HostBasisType(this->CurlBasis::polyOrder_) );
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
const char * getName() const override
Returns basis name.
EFunctionSpace functionSpace_
The function space in which the basis is defined.
virtual bool requireOrientation() const override
True if orientation is required.
Free functions, callable from device code, that implement various polynomials useful in basis definit...
H(curl) basis on the triangle using a construction involving Legendre and integrated Jacobi polynomia...
Header function for Intrepid2::Util class and other utility functions.
For mathematical details of the construction, see:
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
HierarchicalBasis_HDIV_TRI(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
EOperator
Enumeration of primitive operators available in Intrepid. Primitive operators act on reconstructed fu...
virtual void getValues(OutputViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
EPointType
Enumeration of types of point distributions in Intrepid.
ordinal_type getDegree() const
Returns the degree of the basis.
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.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
Header file for the abstract base class Intrepid2::Basis.
For mathematical details of the construction, see:
virtual BasisPtr< typename Kokkos::HostSpace::device_type, OutputScalar, PointScalar > getHostBasis() const override
Creates and returns a Basis object whose DeviceType template argument is Kokkos::HostSpace::device_ty...