49 #ifndef Intrepid2_LegendreBasis_HVOL_LINE_h
50 #define Intrepid2_LegendreBasis_HVOL_LINE_h
52 #include <Kokkos_View.hpp>
53 #include <Kokkos_DynRankView.hpp>
55 #include <Intrepid2_config.h>
58 #ifdef HAVE_INTREPID2_SACADO
59 #include <KokkosExp_View_Fad.hpp>
72 template<
class DeviceType,
class OutputScalar,
class PointScalar,
73 class OutputFieldType,
class InputPointsType>
76 using ExecutionSpace =
typename DeviceType::execution_space;
77 using ScratchSpace =
typename ExecutionSpace::scratch_memory_space;
78 using OutputScratchView = Kokkos::View<OutputScalar*,ScratchSpace,Kokkos::MemoryTraits<Kokkos::Unmanaged>>;
79 using PointScratchView = Kokkos::View<PointScalar*, ScratchSpace,Kokkos::MemoryTraits<Kokkos::Unmanaged>>;
81 using TeamPolicy = Kokkos::TeamPolicy<ExecutionSpace>;
82 using TeamMember =
typename TeamPolicy::member_type;
84 OutputFieldType output_;
85 InputPointsType inputPoints_;
88 int numFields_, numPoints_;
92 size_t fad_size_output_;
95 int polyOrder, EOperator op)
96 : output_(output), inputPoints_(inputPoints), polyOrder_(polyOrder), op_(op),
99 numFields_ = output.extent_int(0);
100 numPoints_ = output.extent_int(1);
101 INTREPID2_TEST_FOR_EXCEPTION(numPoints_ != inputPoints.extent_int(0), std::invalid_argument,
"point counts need to match!");
102 INTREPID2_TEST_FOR_EXCEPTION(numFields_ != polyOrder_+1, std::invalid_argument,
"output field size does not match basis cardinality");
105 KOKKOS_INLINE_FUNCTION
106 void operator()(
const TeamMember & teamMember )
const
108 auto pointOrdinal = teamMember.league_rank();
109 OutputScratchView field_values_at_point;
110 if (fad_size_output_ > 0) {
111 field_values_at_point = OutputScratchView(teamMember.team_shmem(), numFields_, fad_size_output_);
114 field_values_at_point = OutputScratchView(teamMember.team_shmem(), numFields_);
117 const PointScalar x = inputPoints_(pointOrdinal,0);
122 Polynomials::legendreValues(field_values_at_point, polyOrder_, x);
140 Polynomials::legendreDerivativeValues(field_values_at_point, polyOrder_, x, derivativeOrder);
145 device_assert(
false);
148 for (
int fieldOrdinal=0; fieldOrdinal<numFields_; fieldOrdinal++)
150 output_.access(fieldOrdinal,pointOrdinal,0) = field_values_at_point(fieldOrdinal);
157 size_t team_shmem_size (
int team_size)
const
160 size_t shmem_size = 0;
161 if (fad_size_output_ > 0)
162 shmem_size += OutputScratchView::shmem_size(numFields_, fad_size_output_);
164 shmem_size += OutputScratchView::shmem_size(numFields_);
181 template<
typename DeviceType,
182 typename OutputScalar = double,
183 typename PointScalar =
double>
185 :
public Basis<DeviceType,OutputScalar,PointScalar>
199 EPointType pointType_;
211 polyOrder_(polyOrder),
212 pointType_(pointType)
214 INTREPID2_TEST_FOR_EXCEPTION(pointType!=POINTTYPE_DEFAULT,std::invalid_argument,
"PointType not supported");
217 this->
basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Line<2> >() );
222 const int degreeLength = 1;
236 const ordinal_type tagSize = 4;
237 const ordinal_type posScDim = 0;
238 const ordinal_type posScOrd = 1;
239 const ordinal_type posDfOrd = 2;
241 OrdinalTypeArray1DHost tagView(
"tag view", cardinality*tagSize);
243 for (ordinal_type i=0;i<cardinality;++i) {
244 tagView(i*tagSize+0) = 1;
245 tagView(i*tagSize+1) = 0;
246 tagView(i*tagSize+2) = i;
247 tagView(i*tagSize+3) = cardinality;
255 this->basisCardinality_,
275 return "Intrepid2_LegendreBasis_HVOL_LINE";
296 virtual void getValues( OutputViewType outputValues,
const PointViewType inputPoints,
297 const EOperator operatorType = OPERATOR_VALUE )
const override
299 auto numPoints = inputPoints.extent_int(0);
303 FunctorType functor(outputValues, inputPoints, polyOrder_, operatorType);
305 const int outputVectorSize = getVectorSizeForHierarchicalParallelism<OutputScalar>();
306 const int pointVectorSize = getVectorSizeForHierarchicalParallelism<PointScalar>();
307 const int vectorSize = std::max(outputVectorSize,pointVectorSize);
308 const int teamSize = 1;
312 auto policy = Kokkos::TeamPolicy<ExecutionSpace>(numPoints,teamSize,vectorSize);
313 Kokkos::parallel_for( policy , functor,
"Hierarchical_HVOL_LINE_Functor");
322 using HostDeviceType =
typename Kokkos::HostSpace::device_type;
324 return Teuchos::rcp(
new HostBasisType(polyOrder_, pointType_) );
KOKKOS_INLINE_FUNCTION ordinal_type getOperatorOrder(const EOperator operatorType)
Returns order of an operator.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
Implementation of an assert that can safely be called from device code.
Free functions, callable from device code, that implement various polynomials useful in basis definit...
constexpr KOKKOS_INLINE_FUNCTION unsigned getScalarDimensionForView(const ViewType &view)
Returns the size of the Scalar dimension for the View. This is 0 for non-AD types....
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
ECoordinates basisCoordinates_
The coordinate system for which the basis is defined.
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.
EBasis basisType_
Type of the basis.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
void setOrdinalTagData(OrdinalTypeView3D &tagToOrdinal, OrdinalTypeView2D &ordinalToTag, const OrdinalTypeView1D tags, const ordinal_type basisCard, const ordinal_type tagSize, const ordinal_type posScDim, const ordinal_type posScOrd, const ordinal_type posDfOrd)
Fills ordinalToTag_ and tagToOrdinal_ by basis-specific tag data.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, DeviceType > ScalarViewType
View type for scalars.
OrdinalTypeArray2DHost ordinalToTag_
"true" if tagToOrdinal_ and ordinalToTag_ have been initialized
Kokkos::View< ordinal_type **, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom.
OrdinalTypeArray3DHost tagToOrdinal_
DoF tag to ordinal lookup table.
virtual void getValues(OutputViewType, const PointViewType, const EOperator=OPERATOR_VALUE) const
Evaluation of a FEM basis on a reference cell.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
shards::CellTopology basisCellTopology_
Base topology of the cells for which the basis is defined. See the Shards package for definition of b...
typename DeviceType::execution_space ExecutionSpace
(Kokkos) Execution space for basis.
OrdinalTypeArray2DHost fieldOrdinalPolynomialDegree_
Polynomial degree for each degree of freedom. Only defined for hierarchical bases right now....
EFunctionSpace functionSpace_
The function space in which the basis is defined.
Basis defining Legendre basis on the line, a polynomial subspace of L^2 (a.k.a. H(vol)) on the line.
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 void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
virtual const char * getName() const override
Returns basis name.
LegendreBasis_HVOL_LINE(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Constructor.
Functor for computing values for the LegendreBasis_HVOL_LINE class.