49 #ifndef Intrepid2_IntegratedLegendreBasis_HGRAD_TRI_h 
   50 #define Intrepid2_IntegratedLegendreBasis_HGRAD_TRI_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_;
 
   89     static const int numVertices = 3;
 
   90     static const int numEdges    = 3;
 
   91     const int edge_start_[numEdges] = {0,1,0}; 
 
   92     const int edge_end_[numEdges]   = {1,2,2}; 
 
   95                                     int polyOrder, 
bool defineVertexFunctions)
 
   96     : opType_(opType), output_(output), inputPoints_(inputPoints),
 
   97       polyOrder_(polyOrder), defineVertexFunctions_(defineVertexFunctions),
 
   98       fad_size_output_(getScalarDimensionForView(output))
 
  100       numFields_ = output.extent_int(0);
 
  101       numPoints_ = output.extent_int(1);
 
  102       INTREPID2_TEST_FOR_EXCEPTION(numPoints_ != inputPoints.extent_int(0), std::invalid_argument, 
"point counts need to match!");
 
  103       INTREPID2_TEST_FOR_EXCEPTION(numFields_ != (polyOrder_+1)*(polyOrder_+2)/2, std::invalid_argument, 
"output field size does not match basis cardinality");
 
  106     KOKKOS_INLINE_FUNCTION
 
  107     void operator()( 
const TeamMember & teamMember )
 const 
  109       auto pointOrdinal = teamMember.league_rank();
 
  110       OutputScratchView edge_field_values_at_point, jacobi_values_at_point, other_values_at_point, other_values2_at_point;
 
  111       if (fad_size_output_ > 0) {
 
  112         edge_field_values_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
 
  113         jacobi_values_at_point     = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
 
  114         other_values_at_point      = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
 
  115         other_values2_at_point     = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
 
  118         edge_field_values_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
 
  119         jacobi_values_at_point     = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
 
  120         other_values_at_point      = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
 
  121         other_values2_at_point     = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
 
  124       const auto & x = inputPoints_(pointOrdinal,0);
 
  125       const auto & y = inputPoints_(pointOrdinal,1);
 
  128       const PointScalar lambda[3]    = {1. - x - y, x, y};
 
  129       const PointScalar lambda_dx[3] = {-1., 1., 0.};
 
  130       const PointScalar lambda_dy[3] = {-1., 0., 1.};
 
  132       const int num1DEdgeFunctions = polyOrder_ - 1;
 
  139           for (
int vertexOrdinal=0; vertexOrdinal<numVertices; vertexOrdinal++)
 
  141             output_(vertexOrdinal,pointOrdinal) = lambda[vertexOrdinal];
 
  143           if (!defineVertexFunctions_)
 
  147             output_(0,pointOrdinal) = 1.0;
 
  151           int fieldOrdinalOffset = 3;
 
  152           for (
int edgeOrdinal=0; edgeOrdinal<numEdges; edgeOrdinal++)
 
  154             const auto & s0 = lambda[edge_start_[edgeOrdinal]];
 
  155             const auto & s1 = lambda[  edge_end_[edgeOrdinal]];
 
  157             Polynomials::shiftedScaledIntegratedLegendreValues(edge_field_values_at_point, polyOrder_, PointScalar(s1), PointScalar(s0+s1));
 
  158             for (
int edgeFunctionOrdinal=0; edgeFunctionOrdinal<num1DEdgeFunctions; edgeFunctionOrdinal++)
 
  161               output_(edgeFunctionOrdinal+fieldOrdinalOffset,pointOrdinal) = edge_field_values_at_point(edgeFunctionOrdinal+2);
 
  163             fieldOrdinalOffset += num1DEdgeFunctions;
 
  169             const double jacobiScaling = 1.0; 
 
  171             for (
int i=2; i<polyOrder_; i++)
 
  173               const int edgeBasisOrdinal = i+numVertices-2; 
 
  174               const auto & edgeValue = output_(edgeBasisOrdinal,pointOrdinal);
 
  175               const double alpha = i*2.0;
 
  177               Polynomials::integratedJacobiValues(jacobi_values_at_point, alpha, polyOrder_-2, lambda[2], jacobiScaling);
 
  178               for (
int j=1; i+j <= polyOrder_; j++)
 
  180                 const auto & jacobiValue = jacobi_values_at_point(j);
 
  181                 output_(fieldOrdinalOffset,pointOrdinal) = edgeValue * jacobiValue;
 
  182                 fieldOrdinalOffset++;
 
  192           if (defineVertexFunctions_)
 
  196             output_(0,pointOrdinal,0) = -1.0;
 
  197             output_(0,pointOrdinal,1) = -1.0;
 
  203             output_(0,pointOrdinal,0) = 0.0;
 
  204             output_(0,pointOrdinal,1) = 0.0;
 
  207           output_(1,pointOrdinal,0) = 1.0;
 
  208           output_(1,pointOrdinal,1) = 0.0;
 
  210           output_(2,pointOrdinal,0) = 0.0;
 
  211           output_(2,pointOrdinal,1) = 1.0;
 
  214           int fieldOrdinalOffset = 3;
 
  226           auto & P_i_minus_1 = edge_field_values_at_point;
 
  227           auto & L_i_dt      = jacobi_values_at_point;
 
  228           for (
int edgeOrdinal=0; edgeOrdinal<numEdges; edgeOrdinal++)
 
  230             const auto & s0 = lambda[edge_start_[edgeOrdinal]];
 
  231             const auto & s1 = lambda[  edge_end_[edgeOrdinal]];
 
  233             const auto & s0_dx = lambda_dx[edge_start_[edgeOrdinal]];
 
  234             const auto & s0_dy = lambda_dy[edge_start_[edgeOrdinal]];
 
  235             const auto & s1_dx = lambda_dx[  edge_end_[edgeOrdinal]];
 
  236             const auto & s1_dy = lambda_dy[  edge_end_[edgeOrdinal]];
 
  238             Polynomials::shiftedScaledLegendreValues             (P_i_minus_1, polyOrder_-1, PointScalar(s1), PointScalar(s0+s1));
 
  239             Polynomials::shiftedScaledIntegratedLegendreValues_dt(L_i_dt,      polyOrder_,   PointScalar(s1), PointScalar(s0+s1));
 
  240             for (
int edgeFunctionOrdinal=0; edgeFunctionOrdinal<num1DEdgeFunctions; edgeFunctionOrdinal++)
 
  243               const int i = edgeFunctionOrdinal+2;
 
  244               output_(edgeFunctionOrdinal+fieldOrdinalOffset,pointOrdinal,0) = P_i_minus_1(i-1) * s1_dx + L_i_dt(i) * (s1_dx + s0_dx);
 
  245               output_(edgeFunctionOrdinal+fieldOrdinalOffset,pointOrdinal,1) = P_i_minus_1(i-1) * s1_dy + L_i_dt(i) * (s1_dy + s0_dy);
 
  247             fieldOrdinalOffset += num1DEdgeFunctions;
 
  268           auto & P_2i_j_minus_1 = edge_field_values_at_point;
 
  269           auto & L_2i_j_dt      = jacobi_values_at_point;
 
  270           auto & L_i            = other_values_at_point;
 
  271           auto & L_2i_j         = other_values2_at_point;
 
  274             const double jacobiScaling = 1.0; 
 
  276             for (
int i=2; i<polyOrder_; i++)
 
  279               const int edgeBasisOrdinal = i+numVertices-2; 
 
  280               const auto & grad_L_i_dx = output_(edgeBasisOrdinal,pointOrdinal,0);
 
  281               const auto & grad_L_i_dy = output_(edgeBasisOrdinal,pointOrdinal,1);
 
  283               const double alpha = i*2.0;
 
  285               Polynomials::shiftedScaledIntegratedLegendreValues(L_i, polyOrder_, lambda[1], lambda[0]+lambda[1]);
 
  286               Polynomials::integratedJacobiValues_dt(     L_2i_j_dt, alpha, polyOrder_,   lambda[2], jacobiScaling);
 
  287               Polynomials::integratedJacobiValues   (        L_2i_j, alpha, polyOrder_,   lambda[2], jacobiScaling);
 
  288               Polynomials::shiftedScaledJacobiValues(P_2i_j_minus_1, alpha, polyOrder_-1, lambda[2], jacobiScaling);
 
  290               const auto & s0_dx = lambda_dx[0];
 
  291               const auto & s0_dy = lambda_dy[0];
 
  292               const auto & s1_dx = lambda_dx[1];
 
  293               const auto & s1_dy = lambda_dy[1];
 
  294               const auto & s2_dx = lambda_dx[2];
 
  295               const auto & s2_dy = lambda_dy[2];
 
  297               for (
int j=1; i+j <= polyOrder_; j++)
 
  299                 const OutputScalar basisValue_dx = L_2i_j(j) * grad_L_i_dx + L_i(i) * (P_2i_j_minus_1(j-1) * s2_dx + L_2i_j_dt(j) * (s0_dx + s1_dx + s2_dx));
 
  300                 const OutputScalar basisValue_dy = L_2i_j(j) * grad_L_i_dy + L_i(i) * (P_2i_j_minus_1(j-1) * s2_dy + L_2i_j_dt(j) * (s0_dy + s1_dy + s2_dy));
 
  302                 output_(fieldOrdinalOffset,pointOrdinal,0) = basisValue_dx;
 
  303                 output_(fieldOrdinalOffset,pointOrdinal,1) = basisValue_dy;
 
  304                 fieldOrdinalOffset++;
 
  319           INTREPID2_TEST_FOR_ABORT( 
true,
 
  320                                    ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH::OrthPolynomialTri) Computing of second and higher-order derivatives is not currently supported");
 
  323           device_assert(
false);
 
  330     size_t team_shmem_size (
int team_size)
 const 
  333       size_t shmem_size = 0;
 
  334       if (fad_size_output_ > 0)
 
  335         shmem_size += 4 * OutputScratchView::shmem_size(polyOrder_ + 1, fad_size_output_);
 
  337         shmem_size += 4 * OutputScratchView::shmem_size(polyOrder_ + 1);
 
  360   template<
typename ExecutionSpace=Kokkos::DefaultExecutionSpace,
 
  361            typename OutputScalar = double,
 
  362            typename PointScalar  = double,
 
  363            bool defineVertexFunctions = 
true>            
 
  365   : 
public Basis<ExecutionSpace,OutputScalar,PointScalar>
 
  376     bool defineVertexFunctions_; 
 
  390     polyOrder_(polyOrder)
 
  394       this->
basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Triangle<> >() );
 
  399       const int degreeLength = 1;
 
  402       int fieldOrdinalOffset = 0;
 
  405       const int numFunctionsPerVertex = 1;
 
  406       const int numVertexFunctions = numVertices * numFunctionsPerVertex;
 
  407       for (
int i=0; i<numVertexFunctions; i++)
 
  413       if (!defineVertexFunctions)
 
  417       fieldOrdinalOffset += numVertexFunctions;
 
  420       const int numFunctionsPerEdge = polyOrder - 1; 
 
  422       for (
int edgeOrdinal=0; edgeOrdinal<numEdges; edgeOrdinal++)
 
  424         for (
int i=0; i<numFunctionsPerEdge; i++)
 
  428         fieldOrdinalOffset += numFunctionsPerEdge;
 
  432       const int max_ij_sum = polyOrder;
 
  433       for (
int i=2; i<max_ij_sum; i++)
 
  435         for (
int j=1; i+j<=max_ij_sum; j++)
 
  438           fieldOrdinalOffset++;
 
  441       const int numFaces = 1;
 
  442       const int numFunctionsPerFace = ((polyOrder-1)*(polyOrder-2))/2;
 
  443       INTREPID2_TEST_FOR_EXCEPTION(fieldOrdinalOffset != this->
basisCardinality_, std::invalid_argument, 
"Internal error: basis enumeration is incorrect");
 
  450         const ordinal_type tagSize  = 4;        
 
  451         const ordinal_type posScDim = 0;        
 
  452         const ordinal_type posScOrd = 1;        
 
  453         const ordinal_type posDfOrd = 2;        
 
  455         OrdinalTypeArray1DHost tagView(
"tag view", cardinality*tagSize);
 
  456         const int vertexDim = 0, edgeDim = 1, faceDim = 2;
 
  458         if (defineVertexFunctions) {
 
  461             for (
int vertexOrdinal=0; vertexOrdinal<numVertices; vertexOrdinal++)
 
  463               for (
int functionOrdinal=0; functionOrdinal<numFunctionsPerVertex; functionOrdinal++)
 
  465                 tagView(tagNumber*tagSize+0) = vertexDim;             
 
  466                 tagView(tagNumber*tagSize+1) = vertexOrdinal;         
 
  467                 tagView(tagNumber*tagSize+2) = functionOrdinal;       
 
  468                 tagView(tagNumber*tagSize+3) = numFunctionsPerVertex; 
 
  472             for (
int edgeOrdinal=0; edgeOrdinal<numEdges; edgeOrdinal++)
 
  474               for (
int functionOrdinal=0; functionOrdinal<numFunctionsPerEdge; functionOrdinal++)
 
  476                 tagView(tagNumber*tagSize+0) = edgeDim;               
 
  477                 tagView(tagNumber*tagSize+1) = edgeOrdinal;           
 
  478                 tagView(tagNumber*tagSize+2) = functionOrdinal;       
 
  479                 tagView(tagNumber*tagSize+3) = numFunctionsPerEdge;   
 
  483             for (
int faceOrdinal=0; faceOrdinal<numFaces; faceOrdinal++)
 
  485               for (
int functionOrdinal=0; functionOrdinal<numFunctionsPerFace; functionOrdinal++)
 
  487                 tagView(tagNumber*tagSize+0) = faceDim;               
 
  488                 tagView(tagNumber*tagSize+1) = faceOrdinal;           
 
  489                 tagView(tagNumber*tagSize+2) = functionOrdinal;       
 
  490                 tagView(tagNumber*tagSize+3) = numFunctionsPerFace;   
 
  496           for (ordinal_type i=0;i<cardinality;++i) {
 
  497             tagView(i*tagSize+0) = faceDim;     
 
  498             tagView(i*tagSize+1) = 0;           
 
  499             tagView(i*tagSize+2) = i;           
 
  500             tagView(i*tagSize+3) = cardinality; 
 
  522       return "Intrepid2_IntegratedLegendreBasis_HGRAD_TRI";
 
  554     virtual void getValues( OutputViewType outputValues, 
const PointViewType inputPoints,
 
  555                            const EOperator operatorType = OPERATOR_VALUE )
 const override 
  557       auto numPoints = inputPoints.extent_int(0);
 
  561       FunctorType functor(operatorType, outputValues, inputPoints, polyOrder_, defineVertexFunctions);
 
  563       const int outputVectorSize = getVectorSizeForHierarchicalParallelism<OutputScalar>();
 
  564       const int pointVectorSize  = getVectorSizeForHierarchicalParallelism<PointScalar>();
 
  565       const int vectorSize = std::max(outputVectorSize,pointVectorSize);
 
  566       const int teamSize = 1; 
 
  568       auto policy = Kokkos::TeamPolicy<ExecutionSpace>(numPoints,teamSize,vectorSize);
 
  569       Kokkos::parallel_for( policy , functor, 
"Hierarchical_HGRAD_TRI_Functor");
 
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. ...
 
ECoordinates basisCoordinates_
The coordinate system for which the basis is defined. 
 
EFunctionSpace functionSpace_
The function space in which the basis is defined. 
 
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
 
IntegratedLegendreBasis_HGRAD_TRI(int polyOrder)
Constructor. 
 
Free functions, callable from device code, that implement various polynomials useful in basis definit...
 
ordinal_type basisDegree_
Degree of the largest complete polynomial space that can be represented by the basis. 
 
Header function for Intrepid2::Util class and other utility functions. 
 
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell. 
 
EBasis basisType_
Type of the basis. 
 
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. 
 
Functor for computing values for the IntegratedLegendreBasis_HGRAD_TRI class. 
 
OrdinalTypeArray2DHost fieldOrdinalPolynomialDegree_
Polynomial degree for each degree of freedom. Only defined for hierarchical bases right now...
 
ordinal_type getDegree() const
Returns the degree of 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. 
 
const char * getName() const override
Returns basis name. 
 
OrdinalTypeArray2DHost ordinalToTag_
"true" if tagToOrdinal_ and ordinalToTag_ have been initialized 
 
shards::CellTopology basisCellTopology_
Base topology of the cells for which the basis is defined. See the Shards package for definition of b...