11 #ifndef PANZER_INTREPID_BASIS_FACTORY_H
12 #define PANZER_INTREPID_BASIS_FACTORY_H
18 #include "Intrepid2_HVOL_C0_FEM.hpp"
19 #include "Intrepid2_HVOL_TRI_Cn_FEM.hpp"
20 #include "Intrepid2_HVOL_QUAD_Cn_FEM.hpp"
21 #include "Intrepid2_HVOL_HEX_Cn_FEM.hpp"
22 #include "Intrepid2_HVOL_TET_Cn_FEM.hpp"
25 #include "Intrepid2_Basis.hpp"
27 #include "Shards_CellTopology.hpp"
29 #include "Intrepid2_HGRAD_QUAD_C1_FEM.hpp"
30 #include "Intrepid2_HGRAD_QUAD_C2_FEM.hpp"
31 #include "Intrepid2_HGRAD_QUAD_Cn_FEM.hpp"
33 #include "Intrepid2_HGRAD_HEX_C1_FEM.hpp"
34 #include "Intrepid2_HGRAD_HEX_C2_FEM.hpp"
35 #include "Intrepid2_HGRAD_HEX_Cn_FEM.hpp"
37 #include "Intrepid2_HGRAD_TET_C1_FEM.hpp"
38 #include "Intrepid2_HGRAD_TET_C2_FEM.hpp"
39 #include "Intrepid2_HGRAD_TET_Cn_FEM.hpp"
41 #include "Intrepid2_HGRAD_TRI_C1_FEM.hpp"
42 #include "Intrepid2_HGRAD_TRI_C2_FEM.hpp"
43 #include "Intrepid2_HGRAD_TRI_Cn_FEM.hpp"
45 #include "Intrepid2_HGRAD_LINE_C1_FEM.hpp"
46 #include "Intrepid2_HGRAD_LINE_Cn_FEM.hpp"
48 #include "Intrepid2_HCURL_TRI_I1_FEM.hpp"
49 #include "Intrepid2_HCURL_TRI_In_FEM.hpp"
51 #include "Intrepid2_HCURL_TET_I1_FEM.hpp"
52 #include "Intrepid2_HCURL_TET_In_FEM.hpp"
54 #include "Intrepid2_HCURL_QUAD_I1_FEM.hpp"
55 #include "Intrepid2_HCURL_QUAD_In_FEM.hpp"
57 #include "Intrepid2_HCURL_HEX_I1_FEM.hpp"
58 #include "Intrepid2_HCURL_HEX_In_FEM.hpp"
60 #include "Intrepid2_HDIV_TRI_I1_FEM.hpp"
61 #include "Intrepid2_HDIV_TRI_In_FEM.hpp"
63 #include "Intrepid2_HDIV_QUAD_I1_FEM.hpp"
64 #include "Intrepid2_HDIV_QUAD_In_FEM.hpp"
66 #include "Intrepid2_HDIV_TET_I1_FEM.hpp"
67 #include "Intrepid2_HDIV_TET_In_FEM.hpp"
69 #include "Intrepid2_HDIV_HEX_I1_FEM.hpp"
70 #include "Intrepid2_HDIV_HEX_In_FEM.hpp"
89 template <
typename ExecutionSpace,
typename OutputValueType,
typename Po
intValueType>
92 const shards::CellTopology & cell_topology)
97 std::string cell_topology_type = cell_topology.getName();
98 std::size_t end_position = 0;
99 end_position = cell_topology_type.find(
"_");
100 std::string cell_type = cell_topology_type.substr(0,end_position);
107 const Intrepid2::EPointType point_type = Intrepid2::POINTTYPE_EQUISPACED;
109 if ( (basis_type ==
"Const") && (basis_order == 0) )
110 basis =
Teuchos::rcp(
new Intrepid2::Basis_HVOL_C0_FEM<ExecutionSpace,OutputValueType,PointValueType>(cell_topology) );
112 else if ( (basis_type ==
"HVol") && (basis_order == 0) )
113 basis =
Teuchos::rcp(
new Intrepid2::Basis_HVOL_C0_FEM<ExecutionSpace,OutputValueType,PointValueType>(cell_topology) );
115 else if ( (basis_type ==
"HVol") && (cell_type ==
"Quadrilateral") && (basis_order > 0) )
116 basis =
Teuchos::rcp(
new Intrepid2::Basis_HVOL_QUAD_Cn_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
118 else if ( (basis_type ==
"HVol") && (cell_type ==
"Triangle") && (basis_order > 0) )
119 basis =
Teuchos::rcp(
new Intrepid2::Basis_HVOL_TRI_Cn_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
121 else if ( (basis_type ==
"HVol") && (cell_type ==
"Hexahedron") )
122 basis =
Teuchos::rcp(
new Intrepid2::Basis_HVOL_HEX_Cn_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
124 else if ( (basis_type ==
"HVol") && (cell_type ==
"Tetrahedron") )
125 basis =
Teuchos::rcp(
new Intrepid2::Basis_HVOL_TET_Cn_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
127 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Hexahedron") && (basis_order == 1) )
128 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_HEX_C1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
130 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Hexahedron") && (basis_order == 2) )
131 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_HEX_C2_FEM<ExecutionSpace,OutputValueType,PointValueType> );
133 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Hexahedron") && (basis_order > 2) )
134 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_HEX_Cn_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
136 else if ( (basis_type ==
"HCurl") && (cell_type ==
"Hexahedron") && (basis_order == 1) )
137 basis =
Teuchos::rcp(
new Intrepid2::Basis_HCURL_HEX_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
139 else if ( (basis_type ==
"HCurl") && (cell_type ==
"Hexahedron") && (basis_order > 1) )
140 basis =
Teuchos::rcp(
new Intrepid2::Basis_HCURL_HEX_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
142 else if ( (basis_type ==
"HDiv") && (cell_type ==
"Hexahedron") && (basis_order == 1) )
143 basis =
Teuchos::rcp(
new Intrepid2::Basis_HDIV_HEX_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
145 else if ( (basis_type ==
"HDiv") && (cell_type ==
"Hexahedron") && (basis_order > 1) )
146 basis =
Teuchos::rcp(
new Intrepid2::Basis_HDIV_HEX_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
148 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Tetrahedron") && (basis_order == 1) )
149 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_TET_C1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
151 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Tetrahedron") && (basis_order == 2) )
152 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_TET_C2_FEM<ExecutionSpace,OutputValueType,PointValueType> );
154 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Tetrahedron") && (basis_order > 2) )
155 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_TET_Cn_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
157 else if ( (basis_type ==
"HCurl") && (cell_type ==
"Tetrahedron") && (basis_order == 1) )
158 basis =
Teuchos::rcp(
new Intrepid2::Basis_HCURL_TET_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
160 else if ( (basis_type ==
"HCurl") && (cell_type ==
"Tetrahedron") && (basis_order > 1) )
161 basis =
Teuchos::rcp(
new Intrepid2::Basis_HCURL_TET_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
163 else if ( (basis_type ==
"HDiv") && (cell_type ==
"Tetrahedron") && (basis_order == 1) )
164 basis =
Teuchos::rcp(
new Intrepid2::Basis_HDIV_TET_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
166 else if ( (basis_type ==
"HDiv") && (cell_type ==
"Tetrahedron") && (basis_order > 1) )
167 basis =
Teuchos::rcp(
new Intrepid2::Basis_HDIV_TET_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
169 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Quadrilateral") && (basis_order == 1) )
170 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_QUAD_C1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
172 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Quadrilateral") && (basis_order == 2) )
173 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_QUAD_C2_FEM<ExecutionSpace,OutputValueType,PointValueType> );
175 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Quadrilateral") && (basis_order > 2) )
176 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_QUAD_Cn_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
178 else if ( (basis_type ==
"HCurl") && (cell_type ==
"Quadrilateral") && (basis_order == 1) )
179 basis =
Teuchos::rcp(
new Intrepid2::Basis_HCURL_QUAD_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
181 else if ( (basis_type ==
"HCurl") && (cell_type ==
"Quadrilateral") && (basis_order > 1) )
182 basis =
Teuchos::rcp(
new Intrepid2::Basis_HCURL_QUAD_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
184 else if ( (basis_type ==
"HDiv") && (cell_type ==
"Quadrilateral") && (basis_order == 1) )
185 basis =
Teuchos::rcp(
new Intrepid2::Basis_HDIV_QUAD_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
187 else if ( (basis_type ==
"HDiv") && (cell_type ==
"Quadrilateral") && (basis_order > 1) )
188 basis =
Teuchos::rcp(
new Intrepid2::Basis_HDIV_QUAD_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
190 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Triangle") && (basis_order == 1) )
191 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_TRI_C1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
193 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Triangle") && (basis_order == 2) )
194 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_TRI_C2_FEM<ExecutionSpace,OutputValueType,PointValueType> );
196 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Triangle") && (basis_order > 2) )
197 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_TRI_Cn_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
199 else if ( (basis_type ==
"HCurl") && (cell_type ==
"Triangle") && (basis_order == 1) )
200 basis =
Teuchos::rcp(
new Intrepid2::Basis_HCURL_TRI_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
202 else if ( (basis_type ==
"HCurl") && (cell_type ==
"Triangle") && (basis_order > 1) )
203 basis =
Teuchos::rcp(
new Intrepid2::Basis_HCURL_TRI_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
205 else if ( (basis_type ==
"HDiv") && (cell_type ==
"Triangle") && (basis_order == 1) )
206 basis =
Teuchos::rcp(
new Intrepid2::Basis_HDIV_TRI_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
208 else if ( (basis_type ==
"HDiv") && (cell_type ==
"Triangle") && (basis_order > 1) )
209 basis =
Teuchos::rcp(
new Intrepid2::Basis_HDIV_TRI_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
211 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Line") && (basis_order == 1) )
212 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_LINE_C1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
214 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Line") && (basis_order >= 2) )
215 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_LINE_Cn_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
218 "Failed to create the requestedbasis with basis_type=\"" << basis_type <<
219 "\", basis_order=\"" << basis_order <<
"\", and cell_type=\"" << cell_type <<
"\"!\n");
227 "Failed to create basis. Intrepid2 basis base topology does not match mesh cell base topology!");
246 template <
typename ExecutionSpace,
typename OutputValueType,
typename Po
intValueType>
251 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.