49 #ifndef Intrepid2_DerivedBasisFamily_h 
   50 #define Intrepid2_DerivedBasisFamily_h 
   56 #include "Intrepid2_DerivedBasis_HDIV_QUAD.hpp" 
   81   template<
class LineBasisHGRAD, 
class LineBasisHVOL, 
class TriangleBasisFamily = EmptyBasisFamily, 
class TetrahedronBasisFamily = EmptyBasisFamily>
 
   85     using ExecutionSpace  = 
typename LineBasisHGRAD::ExecutionSpace;
 
   86     using OutputValueType = 
typename LineBasisHGRAD::OutputValueType;
 
   87     using PointValueType  = 
typename LineBasisHGRAD::PointValueType;
 
   90     using BasisPtr = Teuchos::RCP<Basis>;
 
   93     using HGRAD_LINE = LineBasisHGRAD;
 
   94     using HVOL_LINE  = LineBasisHVOL;
 
  109     using HGRAD_TRI = 
typename TriangleBasisFamily::HGRAD;
 
  112     using HGRAD_TET = 
typename TetrahedronBasisFamily::HGRAD;
 
  119   template<
class BasisFamily>
 
  120   static typename BasisFamily::BasisPtr getLineBasis(Intrepid2::EFunctionSpace fs, 
int polyOrder)
 
  125       case FUNCTION_SPACE_HVOL:  
return rcp(
new typename BasisFamily::HVOL_LINE (polyOrder));
 
  126       case FUNCTION_SPACE_HGRAD: 
return rcp(
new typename BasisFamily::HGRAD_LINE(polyOrder));
 
  128         INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument, 
"Unsupported function space");
 
  136   template<
class BasisFamily>
 
  137   static typename BasisFamily::BasisPtr getQuadrilateralBasis(Intrepid2::EFunctionSpace fs, 
int polyOrder)
 
  142       case FUNCTION_SPACE_HVOL:  
return rcp(
new typename BasisFamily::HVOL_QUAD (polyOrder));
 
  143       case FUNCTION_SPACE_HCURL: 
return rcp(
new typename BasisFamily::HCURL_QUAD(polyOrder));
 
  144       case FUNCTION_SPACE_HDIV:  
return rcp(
new typename BasisFamily::HDIV_QUAD (polyOrder));
 
  145       case FUNCTION_SPACE_HGRAD: 
return rcp(
new typename BasisFamily::HGRAD_QUAD(polyOrder));
 
  147         INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument, 
"Unsupported function space");
 
  156   template<
class BasisFamily>
 
  157   static typename BasisFamily::BasisPtr getQuadrilateralBasis(Intrepid2::EFunctionSpace fs, 
int polyOrder_x, 
int polyOrder_y)
 
  162       case FUNCTION_SPACE_HVOL:  
return rcp(
new typename BasisFamily::HVOL_QUAD (polyOrder_x,polyOrder_y));
 
  163       case FUNCTION_SPACE_HCURL: 
return rcp(
new typename BasisFamily::HCURL_QUAD(polyOrder_x,polyOrder_y));
 
  164       case FUNCTION_SPACE_HDIV:  
return rcp(
new typename BasisFamily::HDIV_QUAD (polyOrder_x,polyOrder_y));
 
  165       case FUNCTION_SPACE_HGRAD: 
return rcp(
new typename BasisFamily::HGRAD_QUAD(polyOrder_x,polyOrder_y));
 
  167         INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument, 
"Unsupported function space");
 
  175   template<
class BasisFamily>
 
  176   static typename BasisFamily::BasisPtr getHexahedronBasis(Intrepid2::EFunctionSpace fs, 
int polyOrder)
 
  181       case FUNCTION_SPACE_HVOL:  
return rcp(
new typename BasisFamily::HVOL_HEX (polyOrder));
 
  182       case FUNCTION_SPACE_HCURL: 
return rcp(
new typename BasisFamily::HCURL_HEX(polyOrder));
 
  183       case FUNCTION_SPACE_HDIV:  
return rcp(
new typename BasisFamily::HDIV_HEX (polyOrder));
 
  184       case FUNCTION_SPACE_HGRAD: 
return rcp(
new typename BasisFamily::HGRAD_HEX(polyOrder));
 
  186         INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument, 
"Unsupported function space");
 
  196   template<
class BasisFamily>
 
  197   static typename BasisFamily::BasisPtr getHexahedronBasis(Intrepid2::EFunctionSpace fs, 
int polyOrder_x, 
int polyOrder_y, 
int polyOrder_z)
 
  202       case FUNCTION_SPACE_HVOL:  
return rcp(
new typename BasisFamily::HVOL_HEX (polyOrder_x,polyOrder_y,polyOrder_z));
 
  203       case FUNCTION_SPACE_HCURL: 
return rcp(
new typename BasisFamily::HCURL_HEX(polyOrder_x,polyOrder_y,polyOrder_z));
 
  204       case FUNCTION_SPACE_HDIV:  
return rcp(
new typename BasisFamily::HDIV_HEX (polyOrder_x,polyOrder_y,polyOrder_z));
 
  205       case FUNCTION_SPACE_HGRAD: 
return rcp(
new typename BasisFamily::HGRAD_HEX(polyOrder_x,polyOrder_y,polyOrder_z));
 
  207         INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument, 
"Unsupported function space");
 
  215   template<
class BasisFamily>
 
  216   static typename BasisFamily::BasisPtr getTetrahedronBasis(Intrepid2::EFunctionSpace fs, 
int polyOrder)
 
  224       case FUNCTION_SPACE_HGRAD: 
return rcp(
new typename BasisFamily::HGRAD_TET(polyOrder));
 
  226         INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument, 
"Unsupported function space");
 
  234   template<
class BasisFamily>
 
  235   static typename BasisFamily::BasisPtr getTriangleBasis(Intrepid2::EFunctionSpace fs, 
int polyOrder)
 
  243       case FUNCTION_SPACE_HGRAD: 
return rcp(
new typename BasisFamily::HGRAD_TRI(polyOrder));
 
  245         INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument, 
"Unsupported function space");
 
  257   template<
class BasisFamily>
 
  258   static typename BasisFamily::BasisPtr getBasis(shards::CellTopology &cellTopo, Intrepid2::EFunctionSpace fs, 
int polyOrder)
 
  261     switch (cellTopo.getBaseKey())
 
  263       case shards::Line<>::key:          
return getLineBasis<BasisFamily>(fs,polyOrder);
 
  264       case shards::Quadrilateral<>::key: 
return getQuadrilateralBasis<BasisFamily>(fs,polyOrder);
 
  265       case shards::Triangle<>::key:      
return getTriangleBasis<BasisFamily>(fs,polyOrder);
 
  266       case shards::Hexahedron<>::key:    
return getHexahedronBasis<BasisFamily>(fs,polyOrder);
 
  267       case shards::Tetrahedron<>::key:   
return getTetrahedronBasis<BasisFamily>(fs,polyOrder);
 
  269         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...
 
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
 
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...
 
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...
 
Header file for the abstract base class Intrepid2::Basis.