49 #ifndef Intrepid2_IntegratedLegendreBasis_HGRAD_TET_h 
   50 #define Intrepid2_IntegratedLegendreBasis_HGRAD_TET_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 OutputScratchView2D = Kokkos::View<OutputScalar**,ScratchSpace,Kokkos::MemoryTraits<Kokkos::Unmanaged>>;
 
   74     using PointScratchView    = Kokkos::View<PointScalar*, ScratchSpace,Kokkos::MemoryTraits<Kokkos::Unmanaged>>;
 
   76     using TeamPolicy = Kokkos::TeamPolicy<ExecutionSpace>;
 
   77     using TeamMember = 
typename TeamPolicy::member_type;
 
   81     OutputFieldType  output_;      
 
   82     InputPointsType  inputPoints_; 
 
   85     bool defineVertexFunctions_;
 
   86     int numFields_, numPoints_;
 
   88     size_t fad_size_output_;
 
   90     static const int numVertices = 4;
 
   91     static const int numEdges    = 6;
 
   93     const int edge_start_[numEdges] = {0,1,0,0,1,2}; 
 
   94     const int edge_end_[numEdges]   = {1,2,2,3,3,3}; 
 
   96     static const int numFaces    = 4;
 
   97     const int face_vertex_0[numFaces] = {0,0,1,0}; 
 
   98     const int face_vertex_1[numFaces] = {1,1,2,2};
 
   99     const int face_vertex_2[numFaces] = {2,3,3,3};
 
  103     const int face_ordinal_of_first_edge[numFaces] = {0,0,1,2};
 
  106                                     int polyOrder, 
bool defineVertexFunctions)
 
  107     : opType_(opType), output_(output), inputPoints_(inputPoints),
 
  108       polyOrder_(polyOrder), defineVertexFunctions_(defineVertexFunctions),
 
  109       fad_size_output_(getScalarDimensionForView(output))
 
  111       numFields_ = output.extent_int(0);
 
  112       numPoints_ = output.extent_int(1);
 
  113       INTREPID2_TEST_FOR_EXCEPTION(numPoints_ != inputPoints.extent_int(0), std::invalid_argument, 
"point counts need to match!");
 
  114       INTREPID2_TEST_FOR_EXCEPTION(numFields_ != (polyOrder_+1)*(polyOrder_+2)*(polyOrder_+3)/6, std::invalid_argument, 
"output field size does not match basis cardinality");
 
  117     KOKKOS_INLINE_FUNCTION
 
  118     void operator()( 
const TeamMember & teamMember )
 const 
  120       const int numFaceBasisFunctionsPerFace = (polyOrder_-2) * (polyOrder_-1) / 2;
 
  121       const int numInteriorBasisFunctions = (polyOrder_-3) * (polyOrder_-2) * (polyOrder_-1) / 6;
 
  123       auto pointOrdinal = teamMember.league_rank();
 
  124       OutputScratchView legendre_values1_at_point, legendre_values2_at_point;
 
  125       OutputScratchView2D jacobi_values1_at_point, jacobi_values2_at_point, jacobi_values3_at_point;
 
  126       const int numAlphaValues = (polyOrder_-1 > 1) ? (polyOrder_-1) : 1; 
 
  127       if (fad_size_output_ > 0) {
 
  128         legendre_values1_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
 
  129         legendre_values2_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1, fad_size_output_);
 
  130         jacobi_values1_at_point   = OutputScratchView2D(teamMember.team_shmem(), numAlphaValues, polyOrder_ + 1, fad_size_output_);
 
  131         jacobi_values2_at_point   = OutputScratchView2D(teamMember.team_shmem(), numAlphaValues, polyOrder_ + 1, fad_size_output_);
 
  132         jacobi_values3_at_point   = OutputScratchView2D(teamMember.team_shmem(), numAlphaValues, polyOrder_ + 1, fad_size_output_);
 
  135         legendre_values1_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
 
  136         legendre_values2_at_point = OutputScratchView(teamMember.team_shmem(), polyOrder_ + 1);
 
  137         jacobi_values1_at_point   = OutputScratchView2D(teamMember.team_shmem(), numAlphaValues, polyOrder_ + 1);
 
  138         jacobi_values2_at_point   = OutputScratchView2D(teamMember.team_shmem(), numAlphaValues, polyOrder_ + 1);
 
  139         jacobi_values3_at_point   = OutputScratchView2D(teamMember.team_shmem(), numAlphaValues, polyOrder_ + 1);
 
  142       const auto & x = inputPoints_(pointOrdinal,0);
 
  143       const auto & y = inputPoints_(pointOrdinal,1);
 
  144       const auto & z = inputPoints_(pointOrdinal,2);
 
  147       const PointScalar lambda[numVertices] = {1. - x - y - z, x, y, z};
 
  148       const PointScalar lambda_dx[numVertices] = {-1., 1., 0., 0.};
 
  149       const PointScalar lambda_dy[numVertices] = {-1., 0., 1., 0.};
 
  150       const PointScalar lambda_dz[numVertices] = {-1., 0., 0., 1.};
 
  152       const int num1DEdgeFunctions = polyOrder_ - 1;
 
  159           for (
int vertexOrdinal=0; vertexOrdinal<numVertices; vertexOrdinal++)
 
  161             output_(vertexOrdinal,pointOrdinal) = lambda[vertexOrdinal];
 
  163           if (!defineVertexFunctions_)
 
  167             output_(0,pointOrdinal) = 1.0;
 
  171           int fieldOrdinalOffset = numVertices;
 
  172           for (
int edgeOrdinal=0; edgeOrdinal<numEdges; edgeOrdinal++)
 
  174             const auto & s0 = lambda[edge_start_[edgeOrdinal]];
 
  175             const auto & s1 = lambda[  edge_end_[edgeOrdinal]];
 
  177             Polynomials::shiftedScaledIntegratedLegendreValues(legendre_values1_at_point, polyOrder_, PointScalar(s1), PointScalar(s0+s1));
 
  178             for (
int edgeFunctionOrdinal=0; edgeFunctionOrdinal<num1DEdgeFunctions; edgeFunctionOrdinal++)
 
  181               output_(edgeFunctionOrdinal+fieldOrdinalOffset,pointOrdinal) = legendre_values1_at_point(edgeFunctionOrdinal+2);
 
  183             fieldOrdinalOffset += num1DEdgeFunctions;
 
  189           for (
int faceOrdinal=0; faceOrdinal<numFaces; faceOrdinal++)
 
  191             const auto & s0 = lambda[face_vertex_0[faceOrdinal]];
 
  192             const auto & s1 = lambda[face_vertex_1[faceOrdinal]];
 
  193             const auto & s2 = lambda[face_vertex_2[faceOrdinal]];
 
  194             const PointScalar jacobiScaling = s0 + s1 + s2;
 
  197             for (
int n=2; n<=polyOrder_; n++)
 
  199               const double alpha = n*2;
 
  200               const int alphaOrdinal = n-2;
 
  201               using Kokkos::subview;
 
  203               auto jacobi_alpha = subview(jacobi_values1_at_point, alphaOrdinal, ALL);
 
  204               Polynomials::integratedJacobiValues(jacobi_alpha, alpha, polyOrder_-2, s2, jacobiScaling);
 
  207             const int edgeOrdinal = face_ordinal_of_first_edge[faceOrdinal];
 
  208             int localFaceBasisOrdinal = 0;
 
  209             for (
int totalPolyOrder=3; totalPolyOrder<=polyOrder_; totalPolyOrder++)
 
  211               for (
int i=2; i<totalPolyOrder; i++)
 
  213                 const int edgeBasisOrdinal = edgeOrdinal*num1DEdgeFunctions + i-2 + numVertices;
 
  214                 const auto & edgeValue = output_(edgeBasisOrdinal,pointOrdinal);
 
  215                 const int alphaOrdinal = i-2;
 
  217                 const int j = totalPolyOrder - i;
 
  218                 const auto & jacobiValue = jacobi_values1_at_point(alphaOrdinal,j);
 
  219                 const int fieldOrdinal = fieldOrdinalOffset + localFaceBasisOrdinal;
 
  220                 output_(fieldOrdinal,pointOrdinal) = edgeValue * jacobiValue;
 
  222                 localFaceBasisOrdinal++;
 
  225             fieldOrdinalOffset += numFaceBasisFunctionsPerFace;
 
  229           for (
int n=3; n<=polyOrder_; n++)
 
  231             const double alpha = n*2;
 
  232             const double jacobiScaling = 1.0;
 
  233             const int alphaOrdinal = n-3;
 
  234             using Kokkos::subview;
 
  236             auto jacobi_alpha = subview(jacobi_values1_at_point, alphaOrdinal, ALL);
 
  237             Polynomials::integratedJacobiValues(jacobi_alpha, alpha, polyOrder_-3, lambda[3], jacobiScaling);
 
  242           const int min_ij = min_i + min_j;
 
  243           const int min_ijk = min_ij + min_k;
 
  244           int localInteriorBasisOrdinal = 0;
 
  245           for (
int totalPolyOrder_ijk=min_ijk; totalPolyOrder_ijk <= polyOrder_; totalPolyOrder_ijk++)
 
  247             int localFaceBasisOrdinal = 0;
 
  248             for (
int totalPolyOrder_ij=min_ij; totalPolyOrder_ij <= totalPolyOrder_ijk-min_j; totalPolyOrder_ij++)
 
  250               for (
int i=2; i <= totalPolyOrder_ij-min_j; i++)
 
  252                 const int j = totalPolyOrder_ij - i;
 
  253                 const int k = totalPolyOrder_ijk - totalPolyOrder_ij;
 
  254                 const int faceBasisOrdinal = numEdges*num1DEdgeFunctions + numVertices + localFaceBasisOrdinal;
 
  255                 const auto & faceValue = output_(faceBasisOrdinal,pointOrdinal);
 
  256                 const int alphaOrdinal = (i+j)-3;
 
  257                 localFaceBasisOrdinal++;
 
  259                 const int fieldOrdinal = fieldOrdinalOffset + localInteriorBasisOrdinal;
 
  260                 const auto & jacobiValue = jacobi_values1_at_point(alphaOrdinal,k);
 
  261                 output_(fieldOrdinal,pointOrdinal) = faceValue * jacobiValue;
 
  262                 localInteriorBasisOrdinal++;
 
  266           fieldOrdinalOffset += numInteriorBasisFunctions;
 
  273           if (defineVertexFunctions_)
 
  277             output_(0,pointOrdinal,0) = -1.0;
 
  278             output_(0,pointOrdinal,1) = -1.0;
 
  279             output_(0,pointOrdinal,2) = -1.0;
 
  285             output_(0,pointOrdinal,0) = 0.0;
 
  286             output_(0,pointOrdinal,1) = 0.0;
 
  287             output_(0,pointOrdinal,2) = 0.0;
 
  290           output_(1,pointOrdinal,0) = 1.0;
 
  291           output_(1,pointOrdinal,1) = 0.0;
 
  292           output_(1,pointOrdinal,2) = 0.0;
 
  294           output_(2,pointOrdinal,0) = 0.0;
 
  295           output_(2,pointOrdinal,1) = 1.0;
 
  296           output_(2,pointOrdinal,2) = 0.0;
 
  298           output_(3,pointOrdinal,0) = 0.0;
 
  299           output_(3,pointOrdinal,1) = 0.0;
 
  300           output_(3,pointOrdinal,2) = 1.0;
 
  303           int fieldOrdinalOffset = numVertices;
 
  315           auto & P_i_minus_1 = legendre_values1_at_point;
 
  316           auto & L_i_dt      = legendre_values2_at_point;
 
  317           for (
int edgeOrdinal=0; edgeOrdinal<numEdges; edgeOrdinal++)
 
  319             const auto & s0 = lambda[edge_start_[edgeOrdinal]];
 
  320             const auto & s1 = lambda[  edge_end_[edgeOrdinal]];
 
  322             const auto & s0_dx = lambda_dx[edge_start_[edgeOrdinal]];
 
  323             const auto & s0_dy = lambda_dy[edge_start_[edgeOrdinal]];
 
  324             const auto & s0_dz = lambda_dz[edge_start_[edgeOrdinal]];
 
  325             const auto & s1_dx = lambda_dx[  edge_end_[edgeOrdinal]];
 
  326             const auto & s1_dy = lambda_dy[  edge_end_[edgeOrdinal]];
 
  327             const auto & s1_dz = lambda_dz[  edge_end_[edgeOrdinal]];
 
  329             Polynomials::shiftedScaledLegendreValues             (P_i_minus_1, polyOrder_-1, PointScalar(s1), PointScalar(s0+s1));
 
  330             Polynomials::shiftedScaledIntegratedLegendreValues_dt(L_i_dt,      polyOrder_,   PointScalar(s1), PointScalar(s0+s1));
 
  331             for (
int edgeFunctionOrdinal=0; edgeFunctionOrdinal<num1DEdgeFunctions; edgeFunctionOrdinal++)
 
  334               const int i = edgeFunctionOrdinal+2;
 
  335               output_(edgeFunctionOrdinal+fieldOrdinalOffset,pointOrdinal,0) = P_i_minus_1(i-1) * s1_dx + L_i_dt(i) * (s1_dx + s0_dx);
 
  336               output_(edgeFunctionOrdinal+fieldOrdinalOffset,pointOrdinal,1) = P_i_minus_1(i-1) * s1_dy + L_i_dt(i) * (s1_dy + s0_dy);
 
  337               output_(edgeFunctionOrdinal+fieldOrdinalOffset,pointOrdinal,2) = P_i_minus_1(i-1) * s1_dz + L_i_dt(i) * (s1_dz + s0_dz);
 
  339             fieldOrdinalOffset += num1DEdgeFunctions;
 
  360           auto & L_i            = legendre_values2_at_point;
 
  361           auto & L_2i_j_dt      = jacobi_values1_at_point;
 
  362           auto & L_2i_j         = jacobi_values2_at_point;
 
  363           auto & P_2i_j_minus_1 = jacobi_values3_at_point;
 
  365           for (
int faceOrdinal=0; faceOrdinal<numFaces; faceOrdinal++)
 
  367             const auto & s0 = lambda[face_vertex_0[faceOrdinal]];
 
  368             const auto & s1 = lambda[face_vertex_1[faceOrdinal]];
 
  369             const auto & s2 = lambda[face_vertex_2[faceOrdinal]];
 
  370             Polynomials::shiftedScaledIntegratedLegendreValues(L_i, polyOrder_, s1, s0+s1);
 
  372             const PointScalar jacobiScaling = s0 + s1 + s2;
 
  375             for (
int n=2; n<=polyOrder_; n++)
 
  377               const double alpha = n*2;
 
  378               const int alphaOrdinal = n-2;
 
  379               using Kokkos::subview;
 
  381               auto L_2i_j_dt_alpha      = subview(L_2i_j_dt,      alphaOrdinal, ALL);
 
  382               auto L_2i_j_alpha         = subview(L_2i_j,         alphaOrdinal, ALL);
 
  383               auto P_2i_j_minus_1_alpha = subview(P_2i_j_minus_1, alphaOrdinal, ALL);
 
  384               Polynomials::integratedJacobiValues_dt(L_2i_j_dt_alpha,      alpha, polyOrder_-2, s2, jacobiScaling);
 
  385               Polynomials::integratedJacobiValues   (L_2i_j_alpha,         alpha, polyOrder_-2, s2, jacobiScaling);
 
  386               Polynomials::shiftedScaledJacobiValues(P_2i_j_minus_1_alpha, alpha, polyOrder_-1, s2, jacobiScaling);
 
  389             const int edgeOrdinal = face_ordinal_of_first_edge[faceOrdinal];
 
  390             int localFaceOrdinal = 0;
 
  391             for (
int totalPolyOrder=3; totalPolyOrder<=polyOrder_; totalPolyOrder++)
 
  393               for (
int i=2; i<totalPolyOrder; i++)
 
  395                 const int edgeBasisOrdinal = edgeOrdinal*num1DEdgeFunctions + i-2 + numVertices;
 
  396                 const auto & grad_L_i_dx = output_(edgeBasisOrdinal,pointOrdinal,0);
 
  397                 const auto & grad_L_i_dy = output_(edgeBasisOrdinal,pointOrdinal,1);
 
  398                 const auto & grad_L_i_dz = output_(edgeBasisOrdinal,pointOrdinal,2);
 
  400                 const int alphaOrdinal = i-2;
 
  402                 const auto & s0_dx = lambda_dx[face_vertex_0[faceOrdinal]];
 
  403                 const auto & s0_dy = lambda_dy[face_vertex_0[faceOrdinal]];
 
  404                 const auto & s0_dz = lambda_dz[face_vertex_0[faceOrdinal]];
 
  405                 const auto & s1_dx = lambda_dx[face_vertex_1[faceOrdinal]];
 
  406                 const auto & s1_dy = lambda_dy[face_vertex_1[faceOrdinal]];
 
  407                 const auto & s1_dz = lambda_dz[face_vertex_1[faceOrdinal]];
 
  408                 const auto & s2_dx = lambda_dx[face_vertex_2[faceOrdinal]];
 
  409                 const auto & s2_dy = lambda_dy[face_vertex_2[faceOrdinal]];
 
  410                 const auto & s2_dz = lambda_dz[face_vertex_2[faceOrdinal]];
 
  412                 int j = totalPolyOrder - i;
 
  415                 auto & l_2i_j         = L_2i_j(alphaOrdinal,j);
 
  417                 auto & l_2i_j_dt      = L_2i_j_dt(alphaOrdinal,j);
 
  418                 auto & p_2i_j_minus_1 = P_2i_j_minus_1(alphaOrdinal,j-1);
 
  420                 const OutputScalar basisValue_dx = l_2i_j * grad_L_i_dx + l_i * (p_2i_j_minus_1 * s2_dx + l_2i_j_dt * (s0_dx + s1_dx + s2_dx));
 
  421                 const OutputScalar basisValue_dy = l_2i_j * grad_L_i_dy + l_i * (p_2i_j_minus_1 * s2_dy + l_2i_j_dt * (s0_dy + s1_dy + s2_dy));
 
  422                 const OutputScalar basisValue_dz = l_2i_j * grad_L_i_dz + l_i * (p_2i_j_minus_1 * s2_dz + l_2i_j_dt * (s0_dz + s1_dz + s2_dz));
 
  424                 const int fieldOrdinal = fieldOrdinalOffset + localFaceOrdinal;
 
  426                 output_(fieldOrdinal,pointOrdinal,0) = basisValue_dx;
 
  427                 output_(fieldOrdinal,pointOrdinal,1) = basisValue_dy;
 
  428                 output_(fieldOrdinal,pointOrdinal,2) = basisValue_dz;
 
  433             fieldOrdinalOffset += numFaceBasisFunctionsPerFace;
 
  452           auto & L_alpha = jacobi_values1_at_point;
 
  453           auto & P_alpha = jacobi_values2_at_point;
 
  457             const auto & s0 = lambda[0];
 
  458             const auto & s1 = lambda[1];
 
  459             const auto & s2 = lambda[2];
 
  461             Polynomials::shiftedScaledIntegratedLegendreValues(legendre_values1_at_point, polyOrder_, PointScalar(s1), PointScalar(s0+s1));
 
  464             const PointScalar jacobiScaling = s0 + s1 + s2;
 
  465             for (
int n=2; n<=polyOrder_; n++)
 
  467               const double alpha = n*2;
 
  468               const int alphaOrdinal = n-2;
 
  469               using Kokkos::subview;
 
  471               auto jacobi_alpha = subview(jacobi_values3_at_point, alphaOrdinal, ALL);
 
  472               Polynomials::integratedJacobiValues(jacobi_alpha, alpha, polyOrder_-2, s2, jacobiScaling);
 
  477           for (
int n=3; n<=polyOrder_; n++)
 
  479             const double alpha = n*2;
 
  480             const double jacobiScaling = 1.0;
 
  481             const int alphaOrdinal = n-3;
 
  482             using Kokkos::subview;
 
  486             auto L = subview(L_alpha, alphaOrdinal, ALL);
 
  487             auto P = subview(P_alpha, alphaOrdinal, ALL);
 
  488             Polynomials::integratedJacobiValues   (L, alpha, polyOrder_-3, lambda[3], jacobiScaling);
 
  489             Polynomials::shiftedScaledJacobiValues(P, alpha, polyOrder_-3, lambda[3], jacobiScaling);
 
  495           const int min_ij = min_i + min_j;
 
  496           const int min_ijk = min_ij + min_k;
 
  497           int localInteriorBasisOrdinal = 0;
 
  498           for (
int totalPolyOrder_ijk=min_ijk; totalPolyOrder_ijk <= polyOrder_; totalPolyOrder_ijk++)
 
  500             int localFaceBasisOrdinal = 0;
 
  501             for (
int totalPolyOrder_ij=min_ij; totalPolyOrder_ij <= totalPolyOrder_ijk-min_j; totalPolyOrder_ij++)
 
  503               for (
int i=2; i <= totalPolyOrder_ij-min_j; i++)
 
  505                 const int j = totalPolyOrder_ij - i;
 
  506                 const int k = totalPolyOrder_ijk - totalPolyOrder_ij;
 
  508                 const int faceBasisOrdinal = numEdges*num1DEdgeFunctions + numVertices + localFaceBasisOrdinal;
 
  510                 const auto & faceValue_dx = output_(faceBasisOrdinal,pointOrdinal,0);
 
  511                 const auto & faceValue_dy = output_(faceBasisOrdinal,pointOrdinal,1);
 
  512                 const auto & faceValue_dz = output_(faceBasisOrdinal,pointOrdinal,2);
 
  515                 OutputScalar faceValue;
 
  517                   const auto & edgeValue = legendre_values1_at_point(i);
 
  518                   const int alphaOrdinal = i-2;
 
  519                   const auto & jacobiValue = jacobi_values3_at_point(alphaOrdinal,j);
 
  520                   faceValue = edgeValue * jacobiValue;
 
  522                 localFaceBasisOrdinal++;
 
  524                 const int alphaOrdinal = (i+j)-3;
 
  526                 const int fieldOrdinal = fieldOrdinalOffset + localInteriorBasisOrdinal;
 
  527                 const auto & integratedJacobiValue = L_alpha(alphaOrdinal,k);
 
  528                 const auto & jacobiValue = P_alpha(alphaOrdinal,k-1);
 
  529                 output_(fieldOrdinal,pointOrdinal,0) = integratedJacobiValue * faceValue_dx + faceValue * jacobiValue * lambda_dx[3];
 
  530                 output_(fieldOrdinal,pointOrdinal,1) = integratedJacobiValue * faceValue_dy + faceValue * jacobiValue * lambda_dy[3];
 
  531                 output_(fieldOrdinal,pointOrdinal,2) = integratedJacobiValue * faceValue_dz + faceValue * jacobiValue * lambda_dz[3];
 
  533                 localInteriorBasisOrdinal++;
 
  537           fieldOrdinalOffset += numInteriorBasisFunctions;
 
  549           INTREPID2_TEST_FOR_ABORT( 
true,
 
  550                                    ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH::OrthPolynomialTri) Computing of second and higher-order derivatives is not currently supported");
 
  553           device_assert(
false);
 
  560     size_t team_shmem_size (
int team_size)
 const 
  567       const int numAlphaValues = std::max(polyOrder_-1, 1); 
 
  568       size_t shmem_size = 0;
 
  569       if (fad_size_output_ > 0)
 
  572         shmem_size += 2 * OutputScratchView::shmem_size(polyOrder_ + 1, fad_size_output_);
 
  574         shmem_size += 3 * OutputScratchView2D::shmem_size(numAlphaValues, polyOrder_ + 1, fad_size_output_);
 
  579         shmem_size += 2 * OutputScratchView::shmem_size(polyOrder_ + 1);
 
  581         shmem_size += 3 * OutputScratchView2D::shmem_size(numAlphaValues, polyOrder_ + 1);
 
  605   template<
typename ExecutionSpace=Kokkos::DefaultExecutionSpace,
 
  606            typename OutputScalar = double,
 
  607            typename PointScalar  = double,
 
  608            bool defineVertexFunctions = 
true>  
 
  610   : 
public Basis<ExecutionSpace,OutputScalar,PointScalar>
 
  621     bool defineVertexFunctions_;  
 
  635     polyOrder_(polyOrder)
 
  639       this->
basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Tetrahedron<> >() );
 
  644       const int degreeLength = 1;
 
  647       int fieldOrdinalOffset = 0;
 
  650       const int numFunctionsPerVertex = 1;
 
  651       const int numVertexFunctions = numVertices * numFunctionsPerVertex;
 
  652       for (
int i=0; i<numVertexFunctions; i++)
 
  658       if (!defineVertexFunctions)
 
  662       fieldOrdinalOffset += numVertexFunctions;
 
  665       const int numFunctionsPerEdge = polyOrder - 1; 
 
  667       for (
int edgeOrdinal=0; edgeOrdinal<numEdges; edgeOrdinal++)
 
  669         for (
int i=0; i<numFunctionsPerEdge; i++)
 
  673         fieldOrdinalOffset += numFunctionsPerEdge;
 
  677       const int numFunctionsPerFace   = ((polyOrder-1)*(polyOrder-2))/2;
 
  678       const int numFaces = 4;
 
  679       for (
int faceOrdinal=0; faceOrdinal<numFaces; faceOrdinal++)
 
  681         for (
int totalPolyOrder=3; totalPolyOrder<=polyOrder_; totalPolyOrder++)
 
  683           const int totalFaceDofs         = (totalPolyOrder-2)*(totalPolyOrder-1)/2;
 
  684           const int totalFaceDofsPrevious = (totalPolyOrder-3)*(totalPolyOrder-2)/2;
 
  685           const int faceDofsForPolyOrder  = totalFaceDofs - totalFaceDofsPrevious;
 
  686           for (
int i=0; i<faceDofsForPolyOrder; i++)
 
  689             fieldOrdinalOffset++;
 
  695       const int numFunctionsPerVolume = ((polyOrder-1)*(polyOrder-2)*(polyOrder-3))/6;
 
  696       const int numVolumes = 1; 
 
  697       for (
int volumeOrdinal=0; volumeOrdinal<numVolumes; volumeOrdinal++)
 
  699         for (
int totalPolyOrder=4; totalPolyOrder<=polyOrder_; totalPolyOrder++)
 
  701           const int totalInteriorDofs         = (totalPolyOrder-3)*(totalPolyOrder-2)*(totalPolyOrder-1)/6;
 
  702           const int totalInteriorDofsPrevious = (totalPolyOrder-4)*(totalPolyOrder-3)*(totalPolyOrder-2)/6;
 
  703           const int interiorDofsForPolyOrder  = totalInteriorDofs - totalInteriorDofsPrevious;
 
  705           for (
int i=0; i<interiorDofsForPolyOrder; i++)
 
  708             fieldOrdinalOffset++;
 
  713       INTREPID2_TEST_FOR_EXCEPTION(fieldOrdinalOffset != this->
basisCardinality_, std::invalid_argument, 
"Internal error: basis enumeration is incorrect");
 
  720         const int intrepid2FaceOrdinals[4] {3,0,1,2}; 
 
  725         const ordinal_type tagSize  = 4;        
 
  726         const ordinal_type posScDim = 0;        
 
  727         const ordinal_type posScOrd = 1;        
 
  728         const ordinal_type posDfOrd = 2;        
 
  730         OrdinalTypeArray1DHost tagView(
"tag view", cardinality*tagSize);
 
  731         const int vertexDim = 0, edgeDim = 1, faceDim = 2, volumeDim = 3;
 
  733         if (defineVertexFunctions) {
 
  736             for (
int vertexOrdinal=0; vertexOrdinal<numVertices; vertexOrdinal++)
 
  738               for (
int functionOrdinal=0; functionOrdinal<numFunctionsPerVertex; functionOrdinal++)
 
  740                 tagView(tagNumber*tagSize+0) = vertexDim;             
 
  741                 tagView(tagNumber*tagSize+1) = vertexOrdinal;         
 
  742                 tagView(tagNumber*tagSize+2) = functionOrdinal;       
 
  743                 tagView(tagNumber*tagSize+3) = numFunctionsPerVertex; 
 
  747             for (
int edgeOrdinal=0; edgeOrdinal<numEdges; edgeOrdinal++)
 
  749               for (
int functionOrdinal=0; functionOrdinal<numFunctionsPerEdge; functionOrdinal++)
 
  751                 tagView(tagNumber*tagSize+0) = edgeDim;               
 
  752                 tagView(tagNumber*tagSize+1) = edgeOrdinal;           
 
  753                 tagView(tagNumber*tagSize+2) = functionOrdinal;       
 
  754                 tagView(tagNumber*tagSize+3) = numFunctionsPerEdge;   
 
  758             for (
int faceOrdinalESEAS=0; faceOrdinalESEAS<numFaces; faceOrdinalESEAS++)
 
  760               int faceOrdinalIntrepid2 = intrepid2FaceOrdinals[faceOrdinalESEAS];
 
  761               for (
int functionOrdinal=0; functionOrdinal<numFunctionsPerFace; functionOrdinal++)
 
  763                 tagView(tagNumber*tagSize+0) = faceDim;               
 
  764                 tagView(tagNumber*tagSize+1) = faceOrdinalIntrepid2;  
 
  765                 tagView(tagNumber*tagSize+2) = functionOrdinal;       
 
  766                 tagView(tagNumber*tagSize+3) = numFunctionsPerFace;   
 
  770             for (
int volumeOrdinal=0; volumeOrdinal<numVolumes; volumeOrdinal++)
 
  772               for (
int functionOrdinal=0; functionOrdinal<numFunctionsPerVolume; functionOrdinal++)
 
  774                 tagView(tagNumber*tagSize+0) = volumeDim;               
 
  775                 tagView(tagNumber*tagSize+1) = volumeOrdinal;           
 
  776                 tagView(tagNumber*tagSize+2) = functionOrdinal;         
 
  777                 tagView(tagNumber*tagSize+3) = numFunctionsPerVolume;   
 
  781             INTREPID2_TEST_FOR_EXCEPTION(tagNumber != this->
basisCardinality_, std::invalid_argument, 
"Internal error: basis tag enumeration is incorrect");
 
  784           for (ordinal_type i=0;i<cardinality;++i) {
 
  785             tagView(i*tagSize+0) = volumeDim;   
 
  786             tagView(i*tagSize+1) = 0;           
 
  787             tagView(i*tagSize+2) = i;           
 
  788             tagView(i*tagSize+3) = cardinality; 
 
  810       return "Intrepid2_IntegratedLegendreBasis_HGRAD_TET";
 
  842     virtual void getValues( OutputViewType outputValues, 
const PointViewType  inputPoints,
 
  843                            const EOperator operatorType = OPERATOR_VALUE )
 const override 
  845       auto numPoints = inputPoints.extent_int(0);
 
  849       FunctorType functor(operatorType, outputValues, inputPoints, polyOrder_, defineVertexFunctions);
 
  851       const int outputVectorSize = getVectorSizeForHierarchicalParallelism<OutputScalar>();
 
  852       const int pointVectorSize  = getVectorSizeForHierarchicalParallelism<PointScalar>();
 
  853       const int vectorSize = std::max(outputVectorSize,pointVectorSize);
 
  854       const int teamSize = 1; 
 
  856       auto policy = Kokkos::TeamPolicy<ExecutionSpace>(numPoints,teamSize,vectorSize);
 
  857       Kokkos::parallel_for( policy , functor, 
"Hierarchical_HGRAD_TET_Functor");
 
const char * getName() const override
Returns basis name. 
 
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. ...
 
virtual bool requireOrientation() const override
True if orientation is required. 
 
ECoordinates basisCoordinates_
The coordinate system for which the basis is defined. 
 
Basis defining integrated Legendre basis on the line, a polynomial subspace of H(grad) on the line...
 
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 (...
 
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. 
 
EBasis basisType_
Type of the basis. 
 
IntegratedLegendreBasis_HGRAD_TET(int polyOrder)
Constructor. 
 
OrdinalTypeArray2DHost fieldOrdinalPolynomialDegree_
Polynomial degree for each degree of freedom. Only defined for hierarchical bases right now...
 
Functor for computing values for the IntegratedLegendreBasis_HGRAD_TET class. 
 
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. 
 
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell. 
 
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...