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>
102 CurlBasis(polyOrder,pointType)
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_) );
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
H(curl) basis on the triangle using a construction involving Legendre and integrated Jacobi polynomia...
Free functions, callable from device code, that implement various polynomials useful in basis definit...
Header function for Intrepid2::Util class and other utility functions.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, DeviceType > OutputViewType
ordinal_type getDegree() const
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
virtual void getValues(OutputViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
typename DeviceType::execution_space ExecutionSpace
EFunctionSpace functionSpace_
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
HierarchicalBasis_HCURL_TRI(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
const char * getName() const override
Returns basis name.
HierarchicalBasis_HDIV_TRI(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
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...
virtual bool requireOrientation() const override
True if orientation is required.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.