43 #ifndef PANZER_INTREPID_BASIS_FACTORY_H
44 #define PANZER_INTREPID_BASIS_FACTORY_H
50 #include "Intrepid2_HVOL_C0_FEM.hpp"
53 #include "Intrepid2_Basis.hpp"
55 #include "Shards_CellTopology.hpp"
57 #include "Intrepid2_HGRAD_QUAD_C1_FEM.hpp"
58 #include "Intrepid2_HGRAD_QUAD_C2_FEM.hpp"
59 #include "Intrepid2_HGRAD_QUAD_Cn_FEM.hpp"
61 #include "Intrepid2_HGRAD_HEX_C1_FEM.hpp"
62 #include "Intrepid2_HGRAD_HEX_C2_FEM.hpp"
63 #include "Intrepid2_HGRAD_HEX_Cn_FEM.hpp"
65 #include "Intrepid2_HGRAD_TET_C1_FEM.hpp"
66 #include "Intrepid2_HGRAD_TET_C2_FEM.hpp"
67 #include "Intrepid2_HGRAD_TET_Cn_FEM.hpp"
69 #include "Intrepid2_HGRAD_TRI_C1_FEM.hpp"
70 #include "Intrepid2_HGRAD_TRI_C2_FEM.hpp"
71 #include "Intrepid2_HGRAD_TRI_Cn_FEM.hpp"
73 #include "Intrepid2_HGRAD_LINE_C1_FEM.hpp"
74 #include "Intrepid2_HGRAD_LINE_Cn_FEM.hpp"
76 #include "Intrepid2_HCURL_TRI_I1_FEM.hpp"
77 #include "Intrepid2_HCURL_TRI_In_FEM.hpp"
79 #include "Intrepid2_HCURL_TET_I1_FEM.hpp"
80 #include "Intrepid2_HCURL_TET_In_FEM.hpp"
82 #include "Intrepid2_HCURL_QUAD_I1_FEM.hpp"
83 #include "Intrepid2_HCURL_QUAD_In_FEM.hpp"
85 #include "Intrepid2_HCURL_HEX_I1_FEM.hpp"
86 #include "Intrepid2_HCURL_HEX_In_FEM.hpp"
88 #include "Intrepid2_HDIV_TRI_I1_FEM.hpp"
89 #include "Intrepid2_HDIV_TRI_In_FEM.hpp"
91 #include "Intrepid2_HDIV_QUAD_I1_FEM.hpp"
92 #include "Intrepid2_HDIV_QUAD_In_FEM.hpp"
94 #include "Intrepid2_HDIV_TET_I1_FEM.hpp"
95 #include "Intrepid2_HDIV_TET_In_FEM.hpp"
97 #include "Intrepid2_HDIV_HEX_I1_FEM.hpp"
98 #include "Intrepid2_HDIV_HEX_In_FEM.hpp"
117 template <
typename ExecutionSpace,
typename OutputValueType,
typename Po
intValueType>
120 const shards::CellTopology & cell_topology)
125 std::string cell_topology_type = cell_topology.getName();
126 std::size_t end_position = 0;
127 end_position = cell_topology_type.find(
"_");
128 std::string cell_type = cell_topology_type.substr(0,end_position);
135 const Intrepid2::EPointType point_type = Intrepid2::POINTTYPE_EQUISPACED;
137 if ( (basis_type ==
"Const") && (basis_order == 0) )
138 basis =
Teuchos::rcp(
new Intrepid2::Basis_HVOL_C0_FEM<ExecutionSpace,OutputValueType,PointValueType>(cell_topology) );
140 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Hexahedron") && (basis_order == 1) )
141 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_HEX_C1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
143 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Hexahedron") && (basis_order == 2) )
144 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_HEX_C2_FEM<ExecutionSpace,OutputValueType,PointValueType> );
146 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Hexahedron") && (basis_order > 2) )
147 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_HEX_Cn_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
149 else if ( (basis_type ==
"HCurl") && (cell_type ==
"Hexahedron") && (basis_order == 1) )
150 basis =
Teuchos::rcp(
new Intrepid2::Basis_HCURL_HEX_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
152 else if ( (basis_type ==
"HCurl") && (cell_type ==
"Hexahedron") && (basis_order > 1) )
153 basis =
Teuchos::rcp(
new Intrepid2::Basis_HCURL_HEX_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
155 else if ( (basis_type ==
"HDiv") && (cell_type ==
"Hexahedron") && (basis_order == 1) )
156 basis =
Teuchos::rcp(
new Intrepid2::Basis_HDIV_HEX_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
158 else if ( (basis_type ==
"HDiv") && (cell_type ==
"Hexahedron") && (basis_order > 1) )
159 basis =
Teuchos::rcp(
new Intrepid2::Basis_HDIV_HEX_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
161 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Tetrahedron") && (basis_order == 1) )
162 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_TET_C1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
164 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Tetrahedron") && (basis_order == 2) )
165 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_TET_C2_FEM<ExecutionSpace,OutputValueType,PointValueType> );
167 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Tetrahedron") && (basis_order > 2) )
168 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_TET_Cn_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
170 else if ( (basis_type ==
"HCurl") && (cell_type ==
"Tetrahedron") && (basis_order == 1) )
171 basis =
Teuchos::rcp(
new Intrepid2::Basis_HCURL_TET_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
173 else if ( (basis_type ==
"HCurl") && (cell_type ==
"Tetrahedron") && (basis_order > 1) )
174 basis =
Teuchos::rcp(
new Intrepid2::Basis_HCURL_TET_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
176 else if ( (basis_type ==
"HDiv") && (cell_type ==
"Tetrahedron") && (basis_order == 1) )
177 basis =
Teuchos::rcp(
new Intrepid2::Basis_HDIV_TET_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
179 else if ( (basis_type ==
"HDiv") && (cell_type ==
"Tetrahedron") && (basis_order > 1) )
180 basis =
Teuchos::rcp(
new Intrepid2::Basis_HDIV_TET_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
182 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Quadrilateral") && (basis_order == 1) )
183 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_QUAD_C1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
185 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Quadrilateral") && (basis_order == 2) )
186 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_QUAD_C2_FEM<ExecutionSpace,OutputValueType,PointValueType> );
188 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Quadrilateral") && (basis_order > 2) )
189 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_QUAD_Cn_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
191 else if ( (basis_type ==
"HCurl") && (cell_type ==
"Quadrilateral") && (basis_order == 1) )
192 basis =
Teuchos::rcp(
new Intrepid2::Basis_HCURL_QUAD_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
194 else if ( (basis_type ==
"HCurl") && (cell_type ==
"Quadrilateral") && (basis_order > 1) )
195 basis =
Teuchos::rcp(
new Intrepid2::Basis_HCURL_QUAD_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
197 else if ( (basis_type ==
"HDiv") && (cell_type ==
"Quadrilateral") && (basis_order == 1) )
198 basis =
Teuchos::rcp(
new Intrepid2::Basis_HDIV_QUAD_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
200 else if ( (basis_type ==
"HDiv") && (cell_type ==
"Quadrilateral") && (basis_order > 1) )
201 basis =
Teuchos::rcp(
new Intrepid2::Basis_HDIV_QUAD_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
203 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Triangle") && (basis_order == 1) )
204 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_TRI_C1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
206 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Triangle") && (basis_order == 2) )
207 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_TRI_C2_FEM<ExecutionSpace,OutputValueType,PointValueType> );
209 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Triangle") && (basis_order > 2) )
210 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_TRI_Cn_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
212 else if ( (basis_type ==
"HCurl") && (cell_type ==
"Triangle") && (basis_order == 1) )
213 basis =
Teuchos::rcp(
new Intrepid2::Basis_HCURL_TRI_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
215 else if ( (basis_type ==
"HCurl") && (cell_type ==
"Triangle") && (basis_order > 1) )
216 basis =
Teuchos::rcp(
new Intrepid2::Basis_HCURL_TRI_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
218 else if ( (basis_type ==
"HDiv") && (cell_type ==
"Triangle") && (basis_order == 1) )
219 basis =
Teuchos::rcp(
new Intrepid2::Basis_HDIV_TRI_I1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
221 else if ( (basis_type ==
"HDiv") && (cell_type ==
"Triangle") && (basis_order > 1) )
222 basis =
Teuchos::rcp(
new Intrepid2::Basis_HDIV_TRI_In_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
224 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Line") && (basis_order == 1) )
225 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_LINE_C1_FEM<ExecutionSpace,OutputValueType,PointValueType> );
227 else if ( (basis_type ==
"HGrad") && (cell_type ==
"Line") && (basis_order >= 2) )
228 basis =
Teuchos::rcp(
new Intrepid2::Basis_HGRAD_LINE_Cn_FEM<ExecutionSpace,OutputValueType,PointValueType>(basis_order, point_type) );
231 "Failed to create the requestedbasis with basis_type=\"" << basis_type <<
232 "\", basis_order=\"" << basis_order <<
"\", and cell_type=\"" << cell_type <<
"\"!\n");
236 "Failed to create basis. Intrepid2 basis topology does not match mesh cell topology!");
254 template <
typename ExecutionSpace,
typename OutputValueType,
typename Po
intValueType>
259 return createIntrepid2Basis<ExecutionSpace,OutputValueType,PointValueType>(basis_type,basis_order,*cell_topology);