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;
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;
134 using HGRAD_PYR =
typename PyramidBasisFamily::HGRAD;
135 using HCURL_PYR =
typename PyramidBasisFamily::HCURL;
136 using HDIV_PYR =
typename PyramidBasisFamily::HDIV;
137 using HVOL_PYR =
typename PyramidBasisFamily::HVOL;
145 template<
class BasisFamily>
146 static typename BasisFamily::BasisPtr getLineBasis(Intrepid2::EFunctionSpace fs,
int polyOrder,
const EPointType pointType=POINTTYPE_DEFAULT)
151 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_LINE (polyOrder,pointType));
152 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_LINE(polyOrder,pointType));
154 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
163 template<
class BasisFamily>
164 static typename BasisFamily::BasisPtr getQuadrilateralBasis(Intrepid2::EFunctionSpace fs,
int polyOrder,
const EPointType pointType=POINTTYPE_DEFAULT)
169 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_QUAD (polyOrder,pointType));
170 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_QUAD(polyOrder,pointType));
171 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_QUAD (polyOrder,pointType));
172 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_QUAD(polyOrder,pointType));
174 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
184 template<
class BasisFamily>
185 static typename BasisFamily::BasisPtr getQuadrilateralBasis(Intrepid2::EFunctionSpace fs,
int polyOrder_x,
int polyOrder_y,
const EPointType pointType=POINTTYPE_DEFAULT)
190 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_QUAD (polyOrder_x,polyOrder_y,pointType));
191 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_QUAD(polyOrder_x,polyOrder_y,pointType));
192 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_QUAD (polyOrder_x,polyOrder_y,pointType));
193 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_QUAD(polyOrder_x,polyOrder_y,pointType));
195 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
204 template<
class BasisFamily>
205 static typename BasisFamily::BasisPtr getHexahedronBasis(Intrepid2::EFunctionSpace fs,
int polyOrder,
const EPointType pointType=POINTTYPE_DEFAULT)
210 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_HEX (polyOrder,pointType));
211 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_HEX(polyOrder,pointType));
212 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_HEX (polyOrder,pointType));
213 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_HEX(polyOrder,pointType));
215 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
224 template<
class BasisFamily>
225 static typename BasisFamily::BasisPtr getHypercubeBasis_HGRAD(
int polyOrder,
int spaceDim,
const EPointType pointType=POINTTYPE_DEFAULT)
229 using BasisBase =
typename BasisFamily::HGRAD_LINE::BasisBase;
230 using BasisPtr =
typename BasisFamily::BasisPtr;
232 BasisPtr lineBasis = getLineBasis<BasisFamily>(FUNCTION_SPACE_HGRAD, polyOrder);
233 BasisPtr tensorBasis = lineBasis;
235 for (
int d=1; d<spaceDim; d++)
237 tensorBasis = Teuchos::rcp(
new Basis_TensorBasis<BasisBase>(tensorBasis, lineBasis, FUNCTION_SPACE_HGRAD));
248 template<
class BasisFamily>
249 static typename BasisFamily::BasisPtr getHypercubeBasis_HVOL(
int polyOrder,
int spaceDim,
const EPointType pointType=POINTTYPE_DEFAULT)
253 using BasisBase =
typename BasisFamily::HGRAD_LINE::BasisBase;
254 using BasisPtr =
typename BasisFamily::BasisPtr;
256 BasisPtr lineBasis = getLineBasis<BasisFamily>(FUNCTION_SPACE_HVOL, polyOrder);
257 BasisPtr tensorBasis = lineBasis;
259 for (
int d=1; d<spaceDim; d++)
261 tensorBasis = Teuchos::rcp(
new Basis_TensorBasis<BasisBase>(tensorBasis, lineBasis, FUNCTION_SPACE_HVOL));
274 template<
class BasisFamily>
275 static typename BasisFamily::BasisPtr getHexahedronBasis(Intrepid2::EFunctionSpace fs,
int polyOrder_x,
int polyOrder_y,
int polyOrder_z,
const EPointType pointType=POINTTYPE_DEFAULT)
280 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_HEX (polyOrder_x,polyOrder_y,polyOrder_z,pointType));
281 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_HEX(polyOrder_x,polyOrder_y,polyOrder_z,pointType));
282 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_HEX (polyOrder_x,polyOrder_y,polyOrder_z,pointType));
283 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_HEX(polyOrder_x,polyOrder_y,polyOrder_z,pointType));
285 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
293 template<
class BasisFamily>
294 static typename BasisFamily::BasisPtr getSerendipityBasis_HGRAD(
int polyOrder,
int spaceDim)
296 auto fullBasis = getHypercubeBasis_HGRAD<BasisFamily>(polyOrder, spaceDim);
298 using BasisBase =
typename BasisFamily::HGRAD_LINE::BasisBase;
300 auto serendipityBasis = Teuchos::rcp(
new SerendipityBasis<BasisBase>(fullBasis) );
301 return serendipityBasis;
308 template<
class BasisFamily>
309 static typename BasisFamily::BasisPtr getSerendipityBasis_HVOL(
int polyOrder,
int spaceDim)
311 auto fullBasis = getHypercubeBasis_HVOL<BasisFamily>(polyOrder, spaceDim);
313 using BasisBase =
typename BasisFamily::HGRAD_LINE::BasisBase;
315 auto serendipityBasis = Teuchos::rcp(
new SerendipityBasis<BasisBase>(fullBasis) );
316 return serendipityBasis;
324 template<
class BasisFamily>
325 static typename BasisFamily::BasisPtr getTetrahedronBasis(Intrepid2::EFunctionSpace fs,
int polyOrder,
const EPointType pointType=POINTTYPE_DEFAULT)
330 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_TET (polyOrder,pointType));
331 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_TET(polyOrder,pointType));
332 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_TET (polyOrder,pointType));
333 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_TET(polyOrder,pointType));
335 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
344 template<
class BasisFamily>
345 static typename BasisFamily::BasisPtr getTriangleBasis(Intrepid2::EFunctionSpace fs,
int polyOrder,
const EPointType pointType=POINTTYPE_DEFAULT)
350 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_TRI (polyOrder,pointType));
351 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_TRI(polyOrder,pointType));
352 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_TRI (polyOrder,pointType));
353 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_TRI(polyOrder,pointType));
355 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
364 template<
class BasisFamily>
365 static typename BasisFamily::BasisPtr getWedgeBasis(Intrepid2::EFunctionSpace fs,
int polyOrder,
const EPointType pointType=POINTTYPE_DEFAULT)
370 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_WEDGE (polyOrder, pointType));
371 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_WEDGE(polyOrder, pointType));
372 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_WEDGE (polyOrder, pointType));
373 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_WEDGE(polyOrder, pointType));
375 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
384 template<
class BasisFamily>
385 static typename BasisFamily::BasisPtr getWedgeBasis(Intrepid2::EFunctionSpace fs, ordinal_type polyOrder_xy, ordinal_type polyOrder_z,
const EPointType pointType=POINTTYPE_DEFAULT)
390 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_WEDGE (polyOrder_xy, polyOrder_z, pointType));
391 case FUNCTION_SPACE_HCURL:
return rcp(
new typename BasisFamily::HCURL_WEDGE(polyOrder_xy, polyOrder_z, pointType));
392 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_WEDGE (polyOrder_xy, polyOrder_z, pointType));
393 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_WEDGE(polyOrder_xy, polyOrder_z, pointType));
395 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
404 template<
class BasisFamily>
405 static typename BasisFamily::BasisPtr getPyramidBasis(Intrepid2::EFunctionSpace fs, ordinal_type polyOrder,
const EPointType pointType=POINTTYPE_DEFAULT)
410 case FUNCTION_SPACE_HVOL:
return rcp(
new typename BasisFamily::HVOL_PYR (polyOrder, pointType));
412 case FUNCTION_SPACE_HDIV:
return rcp(
new typename BasisFamily::HDIV_PYR (polyOrder, pointType));
413 case FUNCTION_SPACE_HGRAD:
return rcp(
new typename BasisFamily::HGRAD_PYR(polyOrder, pointType));
415 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
"Unsupported function space");
428 template<
class BasisFamily>
429 static typename BasisFamily::BasisPtr getBasis(
const shards::CellTopology &cellTopo, Intrepid2::EFunctionSpace fs,
int polyOrder,
const EPointType pointType = POINTTYPE_DEFAULT)
432 switch (cellTopo.getBaseKey())
434 case shards::Line<>::key:
return getLineBasis<BasisFamily>(fs,polyOrder, pointType);
435 case shards::Quadrilateral<>::key:
return getQuadrilateralBasis<BasisFamily>(fs,polyOrder,pointType);
436 case shards::Triangle<>::key:
return getTriangleBasis<BasisFamily>(fs,polyOrder,pointType);
437 case shards::Hexahedron<>::key:
return getHexahedronBasis<BasisFamily>(fs,polyOrder,pointType);
438 case shards::Tetrahedron<>::key:
return getTetrahedronBasis<BasisFamily>(fs,polyOrder,pointType);
440 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.