49 #ifndef Intrepid2_IntegratedLegendreBasis_HGRAD_LINE_h
50 #define Intrepid2_IntegratedLegendreBasis_HGRAD_LINE_h
52 #include <Kokkos_View.hpp>
53 #include <Kokkos_DynRankView.hpp>
55 #include <Intrepid2_config.h>
67 template<
class ExecutionSpace,
class OutputScalar,
class PointScalar,
68 class OutputFieldType,
class InputPointsType>
71 using ScratchSpace =
typename ExecutionSpace::scratch_memory_space;
72 using OutputScratchView = Kokkos::View<OutputScalar*,ScratchSpace,Kokkos::MemoryTraits<Kokkos::Unmanaged>>;
73 using PointScratchView = Kokkos::View<PointScalar*, ScratchSpace,Kokkos::MemoryTraits<Kokkos::Unmanaged>>;
75 using TeamPolicy = Kokkos::TeamPolicy<ExecutionSpace>;
76 using TeamMember =
typename TeamPolicy::member_type;
80 OutputFieldType output_;
81 InputPointsType inputPoints_;
84 bool defineVertexFunctions_;
85 int numFields_, numPoints_;
87 size_t fad_size_output_;
90 int polyOrder,
bool defineVertexFunctions)
91 : opType_(opType), output_(output), inputPoints_(inputPoints),
92 polyOrder_(polyOrder), defineVertexFunctions_(defineVertexFunctions),
95 numFields_ = output.extent_int(0);
96 numPoints_ = output.extent_int(1);
97 INTREPID2_TEST_FOR_EXCEPTION(numPoints_ != inputPoints.extent_int(0), std::invalid_argument,
"point counts need to match!");
98 INTREPID2_TEST_FOR_EXCEPTION(numFields_ != polyOrder_+1, std::invalid_argument,
"output field size does not match basis cardinality");
101 KOKKOS_INLINE_FUNCTION
102 void operator()(
const TeamMember & teamMember )
const
104 auto pointOrdinal = teamMember.league_rank();
105 OutputScratchView field_values_at_point;
106 if (fad_size_output_ > 0) {
107 field_values_at_point = OutputScratchView(teamMember.team_shmem(), numFields_, fad_size_output_);
110 field_values_at_point = OutputScratchView(teamMember.team_shmem(), numFields_);
113 const auto & input_x = inputPoints_(pointOrdinal,0);
114 const bool taking_derivative = (opType_ != OPERATOR_VALUE);
115 const bool callingShiftedScaledLegendre = (opType_ == OPERATOR_VALUE) || (opType_ == OPERATOR_GRAD) || (opType_ == OPERATOR_D1);
118 const PointScalar x = callingShiftedScaledLegendre ? PointScalar((input_x + 1.0)/2.0) : PointScalar(input_x);
119 const double legendreScaling = 1.0;
120 const double outputScaling = taking_derivative ? 0.5 : 1.0;
127 Polynomials::shiftedScaledIntegratedLegendreValues(field_values_at_point, polyOrder_, x, legendreScaling);
132 if (defineVertexFunctions_)
134 field_values_at_point(0) = 1. - x;
135 field_values_at_point(1) = x;
142 Polynomials::shiftedScaledIntegratedLegendreValues_dx(field_values_at_point, polyOrder_, x, legendreScaling);
147 if (defineVertexFunctions_)
149 field_values_at_point(0) = -1.0;
150 field_values_at_point(1) = 1.0;
164 Polynomials::legendreDerivativeValues(field_values_at_point, polyOrder_, x, derivativeOrder);
169 OutputScalar Pn_minus_one = field_values_at_point(1);
170 for (
int fieldOrdinal=2; fieldOrdinal<numFields_; fieldOrdinal++)
172 OutputScalar Pn = field_values_at_point(fieldOrdinal);
173 field_values_at_point(fieldOrdinal) = Pn_minus_one;
177 if (numFields_ >= 1) field_values_at_point(0) = 0.0;
178 if (numFields_ >= 2) field_values_at_point(1) = 0.0;
186 device_assert(
false);
190 for (
int fieldOrdinal=0; fieldOrdinal<numFields_; fieldOrdinal++)
193 output_.access(fieldOrdinal,pointOrdinal,0) = outputScaling * field_values_at_point(fieldOrdinal);
200 size_t team_shmem_size (
int team_size)
const
203 size_t shmem_size = 0;
204 if (fad_size_output_ > 0)
205 shmem_size += OutputScratchView::shmem_size(numFields_, fad_size_output_);
207 shmem_size += OutputScratchView::shmem_size(numFields_);
230 template<
typename ExecutionSpace=Kokkos::DefaultExecutionSpace,
231 typename OutputScalar = double,
232 typename PointScalar = double,
233 bool defineVertexFunctions =
true,
234 bool useMinusOneToOneReferenceElement =
true>
236 :
public Basis<ExecutionSpace,OutputScalar,PointScalar>
247 bool defineVertexFunctions_;
261 polyOrder_(polyOrder)
265 this->
basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Line<2> >() );
270 const int degreeLength = 1;
279 if (defineVertexFunctions)
289 const ordinal_type tagSize = 4;
290 const ordinal_type posScDim = 0;
291 const ordinal_type posScOrd = 1;
292 const ordinal_type posDfOrd = 2;
294 OrdinalTypeArray1DHost tagView(
"tag view", cardinality*tagSize);
296 if (defineVertexFunctions) {
298 const ordinal_type v0 = 0;
299 tagView(v0*tagSize+0) = 0;
300 tagView(v0*tagSize+1) = 0;
301 tagView(v0*tagSize+2) = 0;
302 tagView(v0*tagSize+3) = 1;
304 const ordinal_type v1 = 1;
305 tagView(v1*tagSize+0) = 0;
306 tagView(v1*tagSize+1) = 1;
307 tagView(v1*tagSize+2) = 0;
308 tagView(v1*tagSize+3) = 1;
310 const ordinal_type iend = cardinality - 2;
311 for (ordinal_type i=0;i<iend;++i) {
312 const auto e = i + 2;
313 tagView(e*tagSize+0) = 1;
314 tagView(e*tagSize+1) = 0;
315 tagView(e*tagSize+2) = i;
316 tagView(e*tagSize+3) = iend;
320 for (ordinal_type i=0;i<cardinality;++i) {
321 tagView(i*tagSize+0) = 1;
322 tagView(i*tagSize+1) = 0;
323 tagView(i*tagSize+2) = i;
324 tagView(i*tagSize+3) = cardinality;
333 this->basisCardinality_,
346 return "Intrepid2_IntegratedLegendreBasis_HGRAD_LINE";
378 virtual void getValues( OutputViewType outputValues,
const PointViewType inputPoints,
379 const EOperator operatorType = OPERATOR_VALUE )
const override
381 auto numPoints = inputPoints.extent_int(0);
385 FunctorType functor(operatorType, outputValues, inputPoints, polyOrder_, defineVertexFunctions);
387 const int outputVectorSize = getVectorSizeForHierarchicalParallelism<OutputScalar>();
388 const int pointVectorSize = getVectorSizeForHierarchicalParallelism<PointScalar>();
389 const int vectorSize = std::max(outputVectorSize,pointVectorSize);
390 const int teamSize = 1;
392 auto policy = Kokkos::TeamPolicy<ExecutionSpace>(numPoints,teamSize,vectorSize);
393 Kokkos::parallel_for( policy , functor,
"Hierarchical_HGRAD_LINE_Functor");
KOKKOS_INLINE_FUNCTION ordinal_type getOperatorOrder(const EOperator operatorType)
Returns order of an operator.
Free functions, callable from device code, that implement various polynomials useful in basis definit...
Header function for Intrepid2::Util class and other utility functions.
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 (...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, ExecSpaceType > PointViewType
View type for input points.
EBasis basisType_
Type of the basis.
ECoordinates basisCoordinates_
The coordinate system for which the basis is defined.
shards::CellTopology basisCellTopology_
Base topology of the cells for which the basis is defined. See the Shards package for definition of b...
EFunctionSpace functionSpace_
The function space in which the basis is defined.
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
OrdinalTypeArray3DHost tagToOrdinal_
DoF tag to ordinal lookup table.
ordinal_type basisCardinality_
Cardinality of the basis, i.e., the number of basis functions/degrees-of-freedom.
Kokkos::DynRankView< scalarType, Kokkos::LayoutStride, ExecSpaceType > ScalarViewType
View type for scalars.
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis.
Kokkos::View< ordinal_type **, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray2DHost
View type for 2d host array.
OrdinalTypeArray2DHost ordinalToTag_
"true" if tagToOrdinal_ and ordinalToTag_ have been initialized
Kokkos::DynRankView< OutputValueType, Kokkos::LayoutStride, ExecSpaceType > OutputViewType
View type for basis value output.
OrdinalTypeArray2DHost fieldOrdinalPolynomialDegree_
Polynomial degree for each degree of freedom. Only defined for hierarchical bases right now....
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.
Basis defining integrated Legendre basis on the line, a polynomial subspace of H(grad) on the line.
virtual bool requireOrientation() const override
True if orientation is required.
const char * getName() const override
Returns basis name.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
IntegratedLegendreBasis_HGRAD_LINE(int polyOrder)
Constructor.
Functor for computing values for the IntegratedLegendreBasis_HGRAD_LINE class.