43 #ifndef PANZER_INTREPID_BASIS_FACTORY_H 
   44 #define PANZER_INTREPID_BASIS_FACTORY_H 
   50 #include "Intrepid2_HVOL_C0_FEM.hpp" 
   51 #include "Intrepid2_HVOL_TRI_Cn_FEM.hpp" 
   52 #include "Intrepid2_HVOL_QUAD_Cn_FEM.hpp" 
   55 #include "Intrepid2_Basis.hpp" 
   57 #include "Shards_CellTopology.hpp" 
   59 #include "Intrepid2_HGRAD_QUAD_C1_FEM.hpp" 
   60 #include "Intrepid2_HGRAD_QUAD_C2_FEM.hpp" 
   61 #include "Intrepid2_HGRAD_QUAD_Cn_FEM.hpp" 
   63 #include "Intrepid2_HGRAD_HEX_C1_FEM.hpp" 
   64 #include "Intrepid2_HGRAD_HEX_C2_FEM.hpp" 
   65 #include "Intrepid2_HGRAD_HEX_Cn_FEM.hpp" 
   67 #include "Intrepid2_HGRAD_TET_C1_FEM.hpp" 
   68 #include "Intrepid2_HGRAD_TET_C2_FEM.hpp" 
   69 #include "Intrepid2_HGRAD_TET_Cn_FEM.hpp" 
   71 #include "Intrepid2_HGRAD_TRI_C1_FEM.hpp" 
   72 #include "Intrepid2_HGRAD_TRI_C2_FEM.hpp" 
   73 #include "Intrepid2_HGRAD_TRI_Cn_FEM.hpp" 
   75 #include "Intrepid2_HGRAD_LINE_C1_FEM.hpp" 
   76 #include "Intrepid2_HGRAD_LINE_Cn_FEM.hpp" 
   78 #include "Intrepid2_HCURL_TRI_I1_FEM.hpp" 
   79 #include "Intrepid2_HCURL_TRI_In_FEM.hpp" 
   81 #include "Intrepid2_HCURL_TET_I1_FEM.hpp" 
   82 #include "Intrepid2_HCURL_TET_In_FEM.hpp" 
   84 #include "Intrepid2_HCURL_QUAD_I1_FEM.hpp" 
   85 #include "Intrepid2_HCURL_QUAD_In_FEM.hpp" 
   87 #include "Intrepid2_HCURL_HEX_I1_FEM.hpp" 
   88 #include "Intrepid2_HCURL_HEX_In_FEM.hpp" 
   90 #include "Intrepid2_HDIV_TRI_I1_FEM.hpp" 
   91 #include "Intrepid2_HDIV_TRI_In_FEM.hpp" 
   93 #include "Intrepid2_HDIV_QUAD_I1_FEM.hpp" 
   94 #include "Intrepid2_HDIV_QUAD_In_FEM.hpp" 
   96 #include "Intrepid2_HDIV_TET_I1_FEM.hpp" 
   97 #include "Intrepid2_HDIV_TET_In_FEM.hpp" 
   99 #include "Intrepid2_HDIV_HEX_I1_FEM.hpp" 
  100 #include "Intrepid2_HDIV_HEX_In_FEM.hpp" 
  119   template <
typename ExecutionSpace, 
typename OutputValueType, 
typename Po
intValueType>
 
  122                        const shards::CellTopology & cell_topology)
 
  127     std::string cell_topology_type = cell_topology.getName();
 
  128     std::size_t end_position = 0;
 
  129     end_position = cell_topology_type.find(
"_");
 
  130     std::string cell_type = cell_topology_type.substr(0,end_position);
 
  137     const Intrepid2::EPointType point_type = Intrepid2::POINTTYPE_EQUISPACED;
 
  139     if ( (basis_type == 
"Const") && (basis_order == 0) )
 
  140       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HVOL_C0_FEM<ExecutionSpace,OutputValueType,PointValueType>(cell_topology) );
 
  142     else if ( (basis_type == 
"HVol") && (basis_order == 0) )
 
  143       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HVOL_C0_FEM<ExecutionSpace,OutputValueType,PointValueType>(cell_topology) );
 
  145     else if ( (basis_type == 
"HVol") && (cell_type == 
"Quadrilateral") && (basis_order > 0) )
 
  146       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HVOL_QUAD_Cn_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
 
  148     else if ( (basis_type == 
"HVol") && (cell_type == 
"Triangle") && (basis_order > 0) )
 
  149       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HVOL_TRI_Cn_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
 
  151     else if ( (basis_type == 
"HGrad") && (cell_type == 
"Hexahedron") && (basis_order == 1) )
 
  152       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HGRAD_HEX_C1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
 
  154     else if ( (basis_type == 
"HGrad") && (cell_type == 
"Hexahedron") && (basis_order == 2) )
 
  155       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HGRAD_HEX_C2_FEM<ExecutionSpace,OutputValueType,PointValueType> );
 
  157     else if ( (basis_type == 
"HGrad") && (cell_type == 
"Hexahedron") && (basis_order > 2) )
 
  158       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HGRAD_HEX_Cn_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
 
  160     else if ( (basis_type == 
"HCurl") && (cell_type == 
"Hexahedron") && (basis_order == 1) )
 
  161       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HCURL_HEX_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
 
  163     else if ( (basis_type == 
"HCurl") && (cell_type == 
"Hexahedron") && (basis_order > 1) )
 
  164       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HCURL_HEX_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
 
  166     else if ( (basis_type == 
"HDiv") && (cell_type == 
"Hexahedron") && (basis_order == 1) )
 
  167       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HDIV_HEX_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
 
  169     else if ( (basis_type == 
"HDiv") && (cell_type == 
"Hexahedron") && (basis_order > 1) )
 
  170       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HDIV_HEX_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
 
  172     else if ( (basis_type == 
"HGrad") && (cell_type == 
"Tetrahedron") && (basis_order == 1) )
 
  173       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HGRAD_TET_C1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
 
  175     else if ( (basis_type == 
"HGrad") && (cell_type == 
"Tetrahedron") && (basis_order == 2) )
 
  176       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HGRAD_TET_C2_FEM<ExecutionSpace,OutputValueType,PointValueType> );
 
  178     else if ( (basis_type == 
"HGrad") && (cell_type == 
"Tetrahedron") && (basis_order > 2) )
 
  179       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HGRAD_TET_Cn_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
 
  181     else if ( (basis_type == 
"HCurl") && (cell_type == 
"Tetrahedron") && (basis_order == 1) )
 
  182       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HCURL_TET_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
 
  184     else if ( (basis_type == 
"HCurl") && (cell_type == 
"Tetrahedron") && (basis_order > 1) )
 
  185       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HCURL_TET_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
 
  187     else if ( (basis_type == 
"HDiv") && (cell_type == 
"Tetrahedron") && (basis_order == 1) )
 
  188       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HDIV_TET_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> ); 
 
  190     else if ( (basis_type == 
"HDiv") && (cell_type == 
"Tetrahedron") && (basis_order > 1) )
 
  191       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HDIV_TET_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) ); 
 
  193     else if ( (basis_type == 
"HGrad") && (cell_type == 
"Quadrilateral") && (basis_order == 1) )
 
  194       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HGRAD_QUAD_C1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
 
  196     else if ( (basis_type == 
"HGrad") && (cell_type == 
"Quadrilateral") && (basis_order == 2) )
 
  197       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HGRAD_QUAD_C2_FEM<ExecutionSpace,OutputValueType,PointValueType> );
 
  199     else if ( (basis_type == 
"HGrad") && (cell_type == 
"Quadrilateral") && (basis_order > 2) )
 
  200       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HGRAD_QUAD_Cn_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
 
  202     else if ( (basis_type == 
"HCurl") && (cell_type == 
"Quadrilateral") && (basis_order == 1) )
 
  203       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HCURL_QUAD_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
 
  205     else if ( (basis_type == 
"HCurl") && (cell_type == 
"Quadrilateral") && (basis_order > 1) )
 
  206       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HCURL_QUAD_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
 
  208     else if ( (basis_type == 
"HDiv") && (cell_type == 
"Quadrilateral") && (basis_order == 1) )
 
  209       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HDIV_QUAD_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
 
  211     else if ( (basis_type == 
"HDiv") && (cell_type == 
"Quadrilateral") && (basis_order > 1) )
 
  212       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HDIV_QUAD_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
 
  214     else if ( (basis_type == 
"HGrad") && (cell_type == 
"Triangle") && (basis_order == 1) )
 
  215       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HGRAD_TRI_C1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
 
  217     else if ( (basis_type == 
"HGrad") && (cell_type == 
"Triangle") && (basis_order == 2) )
 
  218       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HGRAD_TRI_C2_FEM<ExecutionSpace,OutputValueType,PointValueType> );
 
  220     else if ( (basis_type == 
"HGrad") && (cell_type == 
"Triangle") && (basis_order > 2) )
 
  221       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HGRAD_TRI_Cn_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
 
  223     else if ( (basis_type == 
"HCurl") && (cell_type == 
"Triangle") && (basis_order == 1) )
 
  224       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HCURL_TRI_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
 
  226     else if ( (basis_type == 
"HCurl") && (cell_type == 
"Triangle") && (basis_order > 1) )
 
  227       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HCURL_TRI_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
 
  229     else if ( (basis_type == 
"HDiv") && (cell_type == 
"Triangle") && (basis_order == 1) )
 
  230       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HDIV_TRI_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
 
  232     else if ( (basis_type == 
"HDiv") && (cell_type == 
"Triangle") && (basis_order > 1) )
 
  233       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HDIV_TRI_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
 
  235     else if ( (basis_type == 
"HGrad") && (cell_type == 
"Line") && (basis_order == 1) )
 
  236       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HGRAD_LINE_C1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
 
  238     else if ( (basis_type == 
"HGrad") && (cell_type == 
"Line") && (basis_order >= 2) )
 
  239       basis = 
Teuchos::rcp( 
new Intrepid2::Basis_HGRAD_LINE_Cn_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
 
  242                                "Failed to create the requestedbasis with basis_type=\"" << basis_type << 
 
  243                                "\", basis_order=\"" << basis_order << 
"\", and cell_type=\"" << cell_type << 
"\"!\n");
 
  247                                "Failed to create basis.  Intrepid2 basis topology does not match mesh cell topology!");
 
  265   template <
typename ExecutionSpace, 
typename OutputValueType, 
typename Po
intValueType>
 
  270     return createIntrepid2Basis<ExecutionSpace,OutputValueType,PointValueType>(basis_type,basis_order,*cell_topology);
 
bool is_null(const std::shared_ptr< T > &p)
 
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
 
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
 
Teuchos::RCP< Intrepid2::Basis< ExecutionSpace, OutputValueType, PointValueType > > createIntrepid2Basis(const std::string basis_type, int basis_order, const shards::CellTopology &cell_topology)
Creates an Intrepid2::Basis object given the basis, order and cell topology.