15 #ifndef Intrepid2_DerivedBasisFamily_h
16 #define Intrepid2_DerivedBasisFamily_h
22 #include "Intrepid2_DerivedBasis_HDIV_QUAD.hpp"
35 #include "Intrepid2_SerendipityBasis.hpp"
54 template<
class LineBasisHGRAD,
class LineBasisHVOL,
class TriangleBasisFamily = EmptyBasisFamily,
class TetrahedronBasisFamily = EmptyBasisFamily,
class Pyram
idBasisFamily = EmptyBasisFamily>
58 using ExecutionSpace =
typename LineBasisHGRAD::ExecutionSpace;
59 using OutputValueType =
typename LineBasisHGRAD::OutputValueType;
60 using PointValueType =
typename LineBasisHGRAD::PointValueType;
62 using Basis =
typename LineBasisHGRAD::BasisBase;
63 using BasisPtr = Teuchos::RCP<Basis>;
67 using HGRAD_LINE = LineBasisHGRAD;
68 using HVOL_LINE = LineBasisHVOL;
83 using HGRAD_TRI =
typename TriangleBasisFamily::HGRAD;
84 using HCURL_TRI =
typename TriangleBasisFamily::HCURL;
85 using HDIV_TRI =
typename TriangleBasisFamily::HDIV;
86 using HVOL_TRI =
typename TriangleBasisFamily::HVOL;
89 using HGRAD_TET =
typename TetrahedronBasisFamily::HGRAD;
90 using HCURL_TET =
typename TetrahedronBasisFamily::HCURL;
91 using HDIV_TET =
typename TetrahedronBasisFamily::HDIV;
92 using HVOL_TET =
typename TetrahedronBasisFamily::HVOL;
101 using HGRAD_PYR =
typename PyramidBasisFamily::HGRAD;
102 using HCURL_PYR =
typename PyramidBasisFamily::HCURL;
103 using HDIV_PYR =
typename PyramidBasisFamily::HDIV;
104 using HVOL_PYR =
typename PyramidBasisFamily::HVOL;
112 template<
class BasisFamily>
113 static typename BasisFamily::BasisPtr getLineBasis(Intrepid2::EFunctionSpace fs,
int polyOrder,
const EPointType pointType=POINTTYPE_DEFAULT)
118 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_LINE (polyOrder,pointType));
119 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_LINE(polyOrder,pointType));
121 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
130 template<
class BasisFamily>
131 static typename BasisFamily::BasisPtr getQuadrilateralBasis(Intrepid2::EFunctionSpace fs,
int polyOrder,
const EPointType pointType=POINTTYPE_DEFAULT)
136 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_QUAD (polyOrder,pointType));
137 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_QUAD(polyOrder,pointType));
138 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_QUAD (polyOrder,pointType));
139 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_QUAD(polyOrder,pointType));
141 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
151 template<
class BasisFamily>
152 static typename BasisFamily::BasisPtr getQuadrilateralBasis(Intrepid2::EFunctionSpace fs,
int polyOrder_x,
int polyOrder_y,
const EPointType pointType=POINTTYPE_DEFAULT)
157 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_QUAD (polyOrder_x,polyOrder_y,pointType));
158 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_QUAD(polyOrder_x,polyOrder_y,pointType));
159 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_QUAD (polyOrder_x,polyOrder_y,pointType));
160 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_QUAD(polyOrder_x,polyOrder_y,pointType));
162 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
171 template<
class BasisFamily>
172 static typename BasisFamily::BasisPtr getHexahedronBasis(Intrepid2::EFunctionSpace fs,
int polyOrder,
const EPointType pointType=POINTTYPE_DEFAULT)
177 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_HEX (polyOrder,pointType));
178 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_HEX(polyOrder,pointType));
179 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_HEX (polyOrder,pointType));
180 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_HEX(polyOrder,pointType));
182 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
191 template<
class BasisFamily>
192 static typename BasisFamily::BasisPtr getHypercubeBasis_HGRAD(
int polyOrder,
int spaceDim,
const EPointType pointType=POINTTYPE_DEFAULT)
196 using BasisBase =
typename BasisFamily::HGRAD_LINE::BasisBase;
197 using BasisPtr =
typename BasisFamily::BasisPtr;
199 BasisPtr lineBasis = getLineBasis<BasisFamily>(FUNCTION_SPACE_HGRAD, polyOrder);
200 BasisPtr tensorBasis = lineBasis;
202 for (
int d=1; d<spaceDim; d++)
204 tensorBasis = Teuchos::rcp(
new Basis_TensorBasis<BasisBase>(tensorBasis, lineBasis, FUNCTION_SPACE_HGRAD));
215 template<
class BasisFamily>
216 static typename BasisFamily::BasisPtr getHypercubeBasis_HVOL(
int polyOrder,
int spaceDim,
const EPointType pointType=POINTTYPE_DEFAULT)
220 using BasisBase =
typename BasisFamily::HGRAD_LINE::BasisBase;
221 using BasisPtr =
typename BasisFamily::BasisPtr;
223 BasisPtr lineBasis = getLineBasis<BasisFamily>(FUNCTION_SPACE_HVOL, polyOrder);
224 BasisPtr tensorBasis = lineBasis;
226 for (
int d=1; d<spaceDim; d++)
228 tensorBasis = Teuchos::rcp(
new Basis_TensorBasis<BasisBase>(tensorBasis, lineBasis, FUNCTION_SPACE_HVOL));
241 template<
class BasisFamily>
242 static typename BasisFamily::BasisPtr getHexahedronBasis(Intrepid2::EFunctionSpace fs,
int polyOrder_x,
int polyOrder_y,
int polyOrder_z,
const EPointType pointType=POINTTYPE_DEFAULT)
247 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_HEX (polyOrder_x,polyOrder_y,polyOrder_z,pointType));
248 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_HEX(polyOrder_x,polyOrder_y,polyOrder_z,pointType));
249 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_HEX (polyOrder_x,polyOrder_y,polyOrder_z,pointType));
250 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_HEX(polyOrder_x,polyOrder_y,polyOrder_z,pointType));
252 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
260 template<
class BasisFamily>
261 static typename BasisFamily::BasisPtr getSerendipityBasis_HGRAD(
int polyOrder,
int spaceDim)
263 auto fullBasis = getHypercubeBasis_HGRAD<BasisFamily>(polyOrder, spaceDim);
265 using BasisBase =
typename BasisFamily::HGRAD_LINE::BasisBase;
267 auto serendipityBasis = Teuchos::rcp(
new SerendipityBasis<BasisBase>(fullBasis) );
268 return serendipityBasis;
275 template<
class BasisFamily>
276 static typename BasisFamily::BasisPtr getSerendipityBasis_HVOL(
int polyOrder,
int spaceDim)
278 auto fullBasis = getHypercubeBasis_HVOL<BasisFamily>(polyOrder, spaceDim);
280 using BasisBase =
typename BasisFamily::HGRAD_LINE::BasisBase;
282 auto serendipityBasis = Teuchos::rcp(
new SerendipityBasis<BasisBase>(fullBasis) );
283 return serendipityBasis;
291 template<
class BasisFamily>
292 static typename BasisFamily::BasisPtr getTetrahedronBasis(Intrepid2::EFunctionSpace fs,
int polyOrder,
const EPointType pointType=POINTTYPE_DEFAULT)
297 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_TET (polyOrder,pointType));
298 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_TET(polyOrder,pointType));
299 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_TET (polyOrder,pointType));
300 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_TET(polyOrder,pointType));
302 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
311 template<
class BasisFamily>
312 static typename BasisFamily::BasisPtr getTriangleBasis(Intrepid2::EFunctionSpace fs,
int polyOrder,
const EPointType pointType=POINTTYPE_DEFAULT)
317 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_TRI (polyOrder,pointType));
318 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_TRI(polyOrder,pointType));
319 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_TRI (polyOrder,pointType));
320 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_TRI(polyOrder,pointType));
322 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
331 template<
class BasisFamily>
332 static typename BasisFamily::BasisPtr getWedgeBasis(Intrepid2::EFunctionSpace fs,
int polyOrder,
const EPointType pointType=POINTTYPE_DEFAULT)
337 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_WEDGE (polyOrder, pointType));
338 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_WEDGE(polyOrder, pointType));
339 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_WEDGE (polyOrder, pointType));
340 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_WEDGE(polyOrder, pointType));
342 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
351 template<
class BasisFamily>
352 static typename BasisFamily::BasisPtr getWedgeBasis(Intrepid2::EFunctionSpace fs, ordinal_type polyOrder_xy, ordinal_type polyOrder_z,
const EPointType pointType=POINTTYPE_DEFAULT)
357 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_WEDGE (polyOrder_xy, polyOrder_z, pointType));
358 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_WEDGE(polyOrder_xy, polyOrder_z, pointType));
359 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_WEDGE (polyOrder_xy, polyOrder_z, pointType));
360 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_WEDGE(polyOrder_xy, polyOrder_z, pointType));
362 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
371 template<
class BasisFamily>
372 static typename BasisFamily::BasisPtr getPyramidBasis(Intrepid2::EFunctionSpace fs, ordinal_type polyOrder,
const EPointType pointType=POINTTYPE_DEFAULT)
377 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_PYR (polyOrder, pointType));
379 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_PYR (polyOrder, pointType));
380 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_PYR(polyOrder, pointType));
382 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
395 template<
class BasisFamily>
396 static typename BasisFamily::BasisPtr getBasis(
const shards::CellTopology &cellTopo, Intrepid2::EFunctionSpace fs,
int polyOrder,
const EPointType pointType = POINTTYPE_DEFAULT)
399 switch (cellTopo.getBaseKey())
401 case shards::Line<>::key:
return getLineBasis<BasisFamily>(fs,polyOrder, pointType);
402 case shards::Quadrilateral<>::key:
return getQuadrilateralBasis<BasisFamily>(fs,polyOrder,pointType);
403 case shards::Triangle<>::key:
return getTriangleBasis<BasisFamily>(fs,polyOrder,pointType);
404 case shards::Hexahedron<>::key:
return getHexahedronBasis<BasisFamily>(fs,polyOrder,pointType);
405 case shards::Tetrahedron<>::key:
return getTetrahedronBasis<BasisFamily>(fs,polyOrder,pointType);
407 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported cell topology");
EmptyBasisFamily allows us to set a default void family for a given topology.
Implementation of H(div) basis on the hexahedron that is templated on H(vol) and H(grad) on the line...
Implementation of H(vol) basis on the quadrilateral that is templated on H(vol) 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...
Implementation of H(curl) basis on the hexahedron that is templated on H(vol) and H(grad) on the line...
Implementation of H(vol) basis on the wedge that is templated on H(grad) on the line, and H(grad) on the triangle.
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...
Implementation of H(div) basis on the wedge that is templated on H(div,tri), H(vol,tri), H(grad,line), and H(vol,line).
Implementation of H(vol) basis on the quadrilateral that is templated on H(vol) on the line...
A family of basis functions, constructed from H(vol) and H(grad) bases on the line.
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...
Implementation of H(curl) basis on the wedge that is templated on H(grad,tri), H(curl,tri), H(grad,line), and H(vol,line).
Implementation of H(grad) basis on the wedge that is templated on H(grad) on the line, and H(grad) on the triangle.
Header file for the abstract base class Intrepid2::Basis.