48#ifndef Intrepid2_DerivedBasisFamily_h
49#define Intrepid2_DerivedBasisFamily_h
55#include "Intrepid2_DerivedBasis_HDIV_QUAD.hpp"
68#include "Intrepid2_SerendipityBasis.hpp"
87 template<
class LineBasisHGRAD,
class LineBasisHVOL,
class TriangleBasisFamily = EmptyBasisFamily,
class TetrahedronBasisFamily = EmptyBasisFamily,
class Pyram
idBasisFamily = EmptyBasisFamily>
91 using ExecutionSpace =
typename LineBasisHGRAD::ExecutionSpace;
92 using OutputValueType =
typename LineBasisHGRAD::OutputValueType;
93 using PointValueType =
typename LineBasisHGRAD::PointValueType;
95 using Basis =
typename LineBasisHGRAD::BasisBase;
96 using BasisPtr = Teuchos::RCP<Basis>;
100 using HGRAD_LINE = LineBasisHGRAD;
101 using HVOL_LINE = LineBasisHVOL;
106 using HDIV_QUAD = Basis_Derived_HDIV_QUAD <HGRAD_LINE, HVOL_LINE>;
107 using HVOL_QUAD = Basis_Derived_HVOL_QUAD <HVOL_LINE>;
112 using HDIV_HEX = Basis_Derived_HDIV_HEX <HGRAD_LINE, HVOL_LINE>;
113 using HVOL_HEX = Basis_Derived_HVOL_HEX <HVOL_LINE>;
116 using HGRAD_TRI =
typename TriangleBasisFamily::HGRAD;
117 using HCURL_TRI =
typename TriangleBasisFamily::HCURL;
118 using HDIV_TRI =
typename TriangleBasisFamily::HDIV;
119 using HVOL_TRI =
typename TriangleBasisFamily::HVOL;
122 using HGRAD_TET =
typename TetrahedronBasisFamily::HGRAD;
123 using HCURL_TET =
typename TetrahedronBasisFamily::HCURL;
124 using HDIV_TET =
typename TetrahedronBasisFamily::HDIV;
125 using HVOL_TET =
typename TetrahedronBasisFamily::HVOL;
130 using HDIV_WEDGE = Basis_Derived_HDIV_WEDGE < HDIV_TRI, HVOL_TRI, HGRAD_LINE, HVOL_LINE>;
131 using HVOL_WEDGE = Basis_Derived_HVOL_WEDGE < HVOL_TRI, HVOL_LINE>;
135 using HGRAD_PYR =
typename PyramidBasisFamily::HGRAD;
136 using HCURL_PYR =
typename PyramidBasisFamily::HCURL;
137 using HDIV_PYR =
typename PyramidBasisFamily::HDIV;
138 using HVOL_PYR =
typename PyramidBasisFamily::HVOL;
146 template<
class BasisFamily>
147 static typename BasisFamily::BasisPtr
getLineBasis(Intrepid2::EFunctionSpace fs,
int polyOrder,
const EPointType pointType=POINTTYPE_DEFAULT)
152 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_LINE (polyOrder,pointType));
153 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_LINE(polyOrder,pointType));
155 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
164 template<
class BasisFamily>
165 static typename BasisFamily::BasisPtr
getQuadrilateralBasis(Intrepid2::EFunctionSpace fs,
int polyOrder,
const EPointType pointType=POINTTYPE_DEFAULT)
170 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_QUAD (polyOrder,pointType));
171 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_QUAD(polyOrder,pointType));
172 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_QUAD (polyOrder,pointType));
173 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_QUAD(polyOrder,pointType));
175 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
185 template<
class BasisFamily>
186 static typename BasisFamily::BasisPtr
getQuadrilateralBasis(Intrepid2::EFunctionSpace fs,
int polyOrder_x,
int polyOrder_y,
const EPointType pointType=POINTTYPE_DEFAULT)
191 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_QUAD (polyOrder_x,polyOrder_y,pointType));
192 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_QUAD(polyOrder_x,polyOrder_y,pointType));
193 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_QUAD (polyOrder_x,polyOrder_y,pointType));
194 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_QUAD(polyOrder_x,polyOrder_y,pointType));
196 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
205 template<
class BasisFamily>
206 static typename BasisFamily::BasisPtr
getHexahedronBasis(Intrepid2::EFunctionSpace fs,
int polyOrder,
const EPointType pointType=POINTTYPE_DEFAULT)
211 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_HEX (polyOrder,pointType));
212 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_HEX(polyOrder,pointType));
213 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_HEX (polyOrder,pointType));
214 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_HEX(polyOrder,pointType));
216 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
225 template<
class BasisFamily>
226 static typename BasisFamily::BasisPtr
getHypercubeBasis_HGRAD(
int polyOrder,
int spaceDim,
const EPointType pointType=POINTTYPE_DEFAULT)
230 using BasisBase =
typename BasisFamily::HGRAD_LINE::BasisBase;
231 using BasisPtr =
typename BasisFamily::BasisPtr;
236 for (
int d=1; d<spaceDim; d++)
249 template<
class BasisFamily>
250 static typename BasisFamily::BasisPtr
getHypercubeBasis_HVOL(
int polyOrder,
int spaceDim,
const EPointType pointType=POINTTYPE_DEFAULT)
254 using BasisBase =
typename BasisFamily::HGRAD_LINE::BasisBase;
255 using BasisPtr =
typename BasisFamily::BasisPtr;
260 for (
int d=1; d<spaceDim; d++)
275 template<
class BasisFamily>
276 static typename BasisFamily::BasisPtr
getHexahedronBasis(Intrepid2::EFunctionSpace fs,
int polyOrder_x,
int polyOrder_y,
int polyOrder_z,
const EPointType pointType=POINTTYPE_DEFAULT)
281 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_HEX (polyOrder_x,polyOrder_y,polyOrder_z,pointType));
282 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_HEX(polyOrder_x,polyOrder_y,polyOrder_z,pointType));
283 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_HEX (polyOrder_x,polyOrder_y,polyOrder_z,pointType));
284 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_HEX(polyOrder_x,polyOrder_y,polyOrder_z,pointType));
286 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
294 template<
class BasisFamily>
299 using BasisBase =
typename BasisFamily::HGRAD_LINE::BasisBase;
302 return serendipityBasis;
309 template<
class BasisFamily>
314 using BasisBase =
typename BasisFamily::HGRAD_LINE::BasisBase;
317 return serendipityBasis;
325 template<
class BasisFamily>
326 static typename BasisFamily::BasisPtr
getTetrahedronBasis(Intrepid2::EFunctionSpace fs,
int polyOrder,
const EPointType pointType=POINTTYPE_DEFAULT)
331 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_TET (polyOrder,pointType));
332 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_TET(polyOrder,pointType));
333 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_TET (polyOrder,pointType));
334 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_TET(polyOrder,pointType));
336 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
345 template<
class BasisFamily>
346 static typename BasisFamily::BasisPtr
getTriangleBasis(Intrepid2::EFunctionSpace fs,
int polyOrder,
const EPointType pointType=POINTTYPE_DEFAULT)
351 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_TRI (polyOrder,pointType));
352 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_TRI(polyOrder,pointType));
353 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_TRI (polyOrder,pointType));
354 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_TRI(polyOrder,pointType));
356 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
365 template<
class BasisFamily>
366 static typename BasisFamily::BasisPtr
getWedgeBasis(Intrepid2::EFunctionSpace fs,
int polyOrder,
const EPointType pointType=POINTTYPE_DEFAULT)
371 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_WEDGE (polyOrder, pointType));
372 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_WEDGE(polyOrder, pointType));
373 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_WEDGE (polyOrder, pointType));
374 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_WEDGE(polyOrder, pointType));
376 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
385 template<
class BasisFamily>
386 static typename BasisFamily::BasisPtr
getWedgeBasis(Intrepid2::EFunctionSpace fs, ordinal_type polyOrder_xy, ordinal_type polyOrder_z,
const EPointType pointType=POINTTYPE_DEFAULT)
391 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_WEDGE (polyOrder_xy, polyOrder_z, pointType));
392 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_WEDGE(polyOrder_xy, polyOrder_z, pointType));
393 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_WEDGE (polyOrder_xy, polyOrder_z, pointType));
394 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_WEDGE(polyOrder_xy, polyOrder_z, pointType));
396 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
405 template<
class BasisFamily>
406 static typename BasisFamily::BasisPtr
getPyramidBasis(Intrepid2::EFunctionSpace fs, ordinal_type polyOrder,
const EPointType pointType=POINTTYPE_DEFAULT)
414 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_PYR(polyOrder, pointType));
416 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
429 template<
class BasisFamily>
430 static typename BasisFamily::BasisPtr
getBasis(
const shards::CellTopology &cellTopo, Intrepid2::EFunctionSpace fs,
int polyOrder,
const EPointType pointType = POINTTYPE_DEFAULT)
433 switch (cellTopo.getBaseKey())
441 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported cell topology");
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
static BasisFamily::BasisPtr getSerendipityBasis_HGRAD(int polyOrder, int spaceDim)
Factory method for isotropic HGRAD Serendipity bases on a hypercube for the given family....
static BasisFamily::BasisPtr getWedgeBasis(Intrepid2::EFunctionSpace fs, int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Factory method for isotropic wedge bases in the given family.
static BasisFamily::BasisPtr getHypercubeBasis_HVOL(int polyOrder, int spaceDim, const EPointType pointType=POINTTYPE_DEFAULT)
Factory method for isotropic HVOL bases on a hypercube for the given family. Note that this will retu...
static BasisFamily::BasisPtr getHexahedronBasis(Intrepid2::EFunctionSpace fs, int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Factory method for isotropic bases on the hexahedron in the given family.
static BasisFamily::BasisPtr getQuadrilateralBasis(Intrepid2::EFunctionSpace fs, int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Factory method for isotropic quadrilateral bases in the given family.
static BasisFamily::BasisPtr getTetrahedronBasis(Intrepid2::EFunctionSpace fs, int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Factory method for isotropic tetrahedron bases in the given family.
static BasisFamily::BasisPtr getTriangleBasis(Intrepid2::EFunctionSpace fs, int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Factory method for isotropic triangle bases in the given family.
static BasisFamily::BasisPtr getPyramidBasis(Intrepid2::EFunctionSpace fs, ordinal_type polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Factory method for pyramid bases in the given family.
static BasisFamily::BasisPtr getHypercubeBasis_HGRAD(int polyOrder, int spaceDim, const EPointType pointType=POINTTYPE_DEFAULT)
Factory method for isotropic HGRAD bases on a hypercube for the given family. Note that this will ret...
static BasisFamily::BasisPtr getLineBasis(Intrepid2::EFunctionSpace fs, int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Factory method for line bases in the given family.
static BasisFamily::BasisPtr getSerendipityBasis_HVOL(int polyOrder, int spaceDim)
Factory method for isotropic HGRAD Serendipity bases on a hypercube for the given family....
static BasisFamily::BasisPtr getBasis(const shards::CellTopology &cellTopo, Intrepid2::EFunctionSpace fs, int polyOrder, const EPointType pointType=POINTTYPE_DEFAULT)
Factory method for isotropic bases in the given family on the specified cell topology.
Implementation of H(curl) 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...
Implementation of H(curl) basis on the wedge that is templated on H(grad,tri), H(curl,...
Implementation of H(div) basis on the hexahedron that is templated on H(vol) and H(grad) on the line.
Implementation of H(div) basis on the wedge that is templated on H(div,tri), H(vol,...
Implementation of H(grad) basis on the hexahedron that is templated on 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(grad) basis on the wedge 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(vol) basis on the quadrilateral that is templated on H(vol) on the line.
Implementation of H(vol) basis on the wedge that is templated on H(grad) on the line,...
Basis defined as the tensor product of two component bases.
Device DeviceType
(Kokkos) Device type on which Basis is templated. Does not necessarily return true for Kokkos::is_dev...
A family of basis functions, constructed from H(vol) and H(grad) bases on the line.
EmptyBasisFamily allows us to set a default void family for a given topology.
Serendipity Basis, defined as the sub-basis of a provided basis, consisting of basis elements for whi...