44 #include <Teuchos_RCP.hpp> 46 #include <Shards_BasicTopologies.hpp> 48 #include "Intrepid_HGRAD_HEX_C1_FEM.hpp" 49 #include "Intrepid_HGRAD_HEX_C2_FEM.hpp" 50 #include "Intrepid_HGRAD_LINE_C1_FEM.hpp" 51 #include "Intrepid_HGRAD_PYR_C1_FEM.hpp" 52 #include "Intrepid_HGRAD_QUAD_C1_FEM.hpp" 53 #include "Intrepid_HGRAD_QUAD_C2_FEM.hpp" 54 #include "Intrepid_HGRAD_TET_C1_FEM.hpp" 55 #include "Intrepid_HGRAD_TET_C2_FEM.hpp" 56 #include "Intrepid_HGRAD_TRI_C1_FEM.hpp" 57 #include "Intrepid_HGRAD_TRI_C2_FEM.hpp" 58 #include "Intrepid_HGRAD_WEDGE_C1_FEM.hpp" 59 #include "Intrepid_HGRAD_WEDGE_C2_FEM.hpp" 60 #include <Intrepid_Basis.hpp> 70 Teuchos::RCP<Intrepid::Basis<double, Intrepid::FieldContainer<double>>>
71 IntrepidBasisFactory::create(
const shards::CellTopology &cell_topo )
74 Teuchos::RCP<Intrepid::Basis<double, Intrepid::FieldContainer<double>>>
77 switch ( cell_topo.getKey() )
80 case shards::Line<2>::key:
81 basis = Teuchos::rcp(
new Intrepid::Basis_HGRAD_LINE_C1_FEM<
82 double, Intrepid::FieldContainer<double>>() );
85 case shards::Triangle<3>::key:
86 basis = Teuchos::rcp(
new Intrepid::Basis_HGRAD_TRI_C1_FEM<
87 double, Intrepid::FieldContainer<double>>() );
90 case shards::Triangle<6>::key:
91 basis = Teuchos::rcp(
new Intrepid::Basis_HGRAD_TRI_C2_FEM<
92 double, Intrepid::FieldContainer<double>>() );
95 case shards::Quadrilateral<4>::key:
96 basis = Teuchos::rcp(
new Intrepid::Basis_HGRAD_QUAD_C1_FEM<
97 double, Intrepid::FieldContainer<double>>() );
100 case shards::Quadrilateral<9>::key:
101 basis = Teuchos::rcp(
new Intrepid::Basis_HGRAD_QUAD_C2_FEM<
102 double, Intrepid::FieldContainer<double>>() );
105 case shards::Tetrahedron<4>::key:
106 basis = Teuchos::rcp(
new Intrepid::Basis_HGRAD_TET_C1_FEM<
107 double, Intrepid::FieldContainer<double>>() );
110 case shards::Tetrahedron<10>::key:
111 basis = Teuchos::rcp(
new Intrepid::Basis_HGRAD_TET_C2_FEM<
112 double, Intrepid::FieldContainer<double>>() );
115 case shards::Hexahedron<8>::key:
116 basis = Teuchos::rcp(
new Intrepid::Basis_HGRAD_HEX_C1_FEM<
117 double, Intrepid::FieldContainer<double>>() );
120 case shards::Hexahedron<27>::key:
121 basis = Teuchos::rcp(
new Intrepid::Basis_HGRAD_HEX_C2_FEM<
122 double, Intrepid::FieldContainer<double>>() );
125 case shards::Wedge<6>::key:
126 basis = Teuchos::rcp(
new Intrepid::Basis_HGRAD_WEDGE_C1_FEM<
127 double, Intrepid::FieldContainer<double>>() );
130 case shards::Wedge<18>::key:
131 basis = Teuchos::rcp(
new Intrepid::Basis_HGRAD_WEDGE_C2_FEM<
132 double, Intrepid::FieldContainer<double>>() );
135 case shards::Pyramid<5>::key:
136 case shards::Pyramid<13>::key:
137 case shards::Pyramid<14>::key:
138 basis = Teuchos::rcp(
new Intrepid::Basis_HGRAD_PYR_C1_FEM<
139 double, Intrepid::FieldContainer<double>>() );
143 bool topology_supported =
false;
144 DTK_INSIST( topology_supported );
Assertions and Design-by-Contract for error handling.
Factory for Intrepid basis functions.