49 #ifndef Intrepid2_SerendipityBasisFamily_h
50 #define Intrepid2_SerendipityBasisFamily_h
56 #include "Intrepid2_DerivedBasis_HDIV_QUAD.hpp"
68 template<
class FullBasis,
int numPolyOrderArgs>
74 using BasisBase =
typename FullBasis::BasisBase;
75 using BasisPtr = Teuchos::RCP<BasisBase>;
76 using DeviceType =
typename BasisBase::DeviceType;
77 using ExecutionSpace =
typename BasisBase::ExecutionSpace;
78 using OutputValueType =
typename BasisBase::OutputValueType;
79 using PointValueType =
typename BasisBase::PointValueType;
81 using OrdinalTypeArray1D =
typename BasisBase::OrdinalTypeArray1D;
82 using OrdinalTypeArray1DHost =
typename BasisBase::OrdinalTypeArray1DHost;
83 using OrdinalTypeArray2DHost =
typename BasisBase::OrdinalTypeArray2DHost;
84 using OutputViewType =
typename BasisBase::OutputViewType;
85 using PointViewType =
typename BasisBase::PointViewType;
90 SerendipityBasis<typename FullBasis::BasisBase>(Teuchos::rcp(new FullBasis(polyOrder)) )
94 template<
bool M=(numPolyOrderArgs==2)>
95 SerendipityBasisWrapper(
int polyOrder_x,
int polyOrder_y,
const EPointType pointType=POINTTYPE_DEFAULT,
typename std::enable_if<M>::type* = 0)
97 SerendipityBasis<typename FullBasis::BasisBase>(Teuchos::rcp(new FullBasis(polyOrder_x, polyOrder_y)) )
101 template<
bool M=(numPolyOrderArgs==3)>
102 SerendipityBasisWrapper(
int polyOrder_x,
int polyOrder_y,
int polyOrder_z,
const EPointType pointType=POINTTYPE_DEFAULT,
typename std::enable_if<M>::type* = 0)
104 SerendipityBasis<typename FullBasis::BasisBase>(Teuchos::rcp(new FullBasis(polyOrder_x, polyOrder_y, polyOrder_z)) )
111 template<
class LineBasisHGRAD,
class LineBasisHVOL,
class TriangleBasisFamily,
class TetrahedronBasisFamily >
115 using ExecutionSpace =
typename LineBasisHGRAD::ExecutionSpace;
116 using OutputValueType =
typename LineBasisHGRAD::OutputValueType;
117 using PointValueType =
typename LineBasisHGRAD::PointValueType;
119 using Basis =
typename LineBasisHGRAD::BasisBase;
120 using BasisPtr = Teuchos::RCP<Basis>;
124 using HGRAD_LINE = LineBasisHGRAD;
125 using HVOL_LINE = LineBasisHVOL;
140 using HGRAD_TRI =
typename TriangleBasisFamily::HGRAD;
141 using HCURL_TRI =
typename TriangleBasisFamily::HCURL;
142 using HDIV_TRI =
typename TriangleBasisFamily::HDIV;
143 using HVOL_TRI =
typename TriangleBasisFamily::HVOL;
146 using HGRAD_TET =
typename TetrahedronBasisFamily::HGRAD;
147 using HCURL_TET =
typename TetrahedronBasisFamily::HCURL;
148 using HDIV_TET =
typename TetrahedronBasisFamily::HDIV;
149 using HVOL_TET =
typename TetrahedronBasisFamily::HVOL;
155 template<
typename DeviceType,
156 typename OutputScalar = double,
157 typename PointScalar =
double>
170 template<
typename DeviceType,
171 typename OutputScalar = double,
172 typename PointScalar =
double>
SerendipityBasisWrapper(int polyOrder_x, int polyOrder_y, int polyOrder_z, const EPointType pointType=POINTTYPE_DEFAULT, typename std::enable_if< M >::type *=0)
three-argument constructor; enabled if numPolyOrderArgs is 3.
Serendipity Basis, defined as the sub-basis of a provided basis, consisting of basis elements for whi...
Serendipity basis family constructed using the DG hierarchical basis family.
SerendipityBasisWrapper(int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
single-argument constructor, for isotropic bases.
Serendipity basis family constructed in terms of arbitrary bases on the line, triangle, and tetrahedron. (These must be hierarchical bases.)
Implementation of H(div) basis on the hexahedron that is templated on H(vol) and H(grad) on the line...
Implementation of H(curl) basis on the quadrilateral that is templated on H(vol) and H(grad) on the l...
Device DeviceType
(Kokkos) Device type on which Basis is templated. Does not necessarily return true for Kokkos::is_dev...
Helper class that allows SerendipityBasis construction with poly order arguments that are passed to t...
Implementation of H(curl) basis on the hexahedron that is templated on H(vol) and H(grad) on the line...
Implementation of H(grad) basis on the quadrilateral that is templated on H(grad) on the line...
Implementation of H(vol) basis on the hexahedron that is templated on H(vol) on the line...
Basis defining Legendre basis on the line, a polynomial subspace of L^2 (a.k.a. H(vol)) on the line...
SerendipityBasisWrapper(int polyOrder_x, int polyOrder_y, const EPointType pointType=POINTTYPE_DEFAULT, typename std::enable_if< M >::type *=0)
two-argument constructor; enabled if numPolyOrderArgs is 2.
Implementation of H(grad) basis on the hexahedron that is templated on H(grad) on the line...
Implementation of H(vol) basis on the quadrilateral that is templated on H(vol) on the line...
Header file for the abstract base class Intrepid2::Basis.