Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_IntrepidBasisFactory.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Panzer: A partial differential equation assembly
4// engine for strongly coupled complex multiphysics systems
5//
6// Copyright 2011 NTESS and the Panzer contributors.
7// SPDX-License-Identifier: BSD-3-Clause
8// *****************************************************************************
9// @HEADER
10
11#ifndef PANZER_INTREPID_BASIS_FACTORY_H
12#define PANZER_INTREPID_BASIS_FACTORY_H
13
14#include <sstream>
15#include <string>
16#include <map>
17
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"
23
24#include "Teuchos_RCP.hpp"
25#include "Intrepid2_Basis.hpp"
26
27#include "Shards_CellTopology.hpp"
28
29#include "Intrepid2_HGRAD_QUAD_C1_FEM.hpp"
30#include "Intrepid2_HGRAD_QUAD_C2_FEM.hpp"
31#include "Intrepid2_HGRAD_QUAD_Cn_FEM.hpp"
32
33#include "Intrepid2_HGRAD_HEX_C1_FEM.hpp"
34#include "Intrepid2_HGRAD_HEX_C2_FEM.hpp"
35#include "Intrepid2_HGRAD_HEX_Cn_FEM.hpp"
36
37#include "Intrepid2_HGRAD_TET_C1_FEM.hpp"
38#include "Intrepid2_HGRAD_TET_C2_FEM.hpp"
39#include "Intrepid2_HGRAD_TET_Cn_FEM.hpp"
40
41#include "Intrepid2_HGRAD_TRI_C1_FEM.hpp"
42#include "Intrepid2_HGRAD_TRI_C2_FEM.hpp"
43#include "Intrepid2_HGRAD_TRI_Cn_FEM.hpp"
44
45#include "Intrepid2_HGRAD_LINE_C1_FEM.hpp"
46#include "Intrepid2_HGRAD_LINE_Cn_FEM.hpp"
47
48#include "Intrepid2_HCURL_TRI_I1_FEM.hpp"
49#include "Intrepid2_HCURL_TRI_In_FEM.hpp"
50
51#include "Intrepid2_HCURL_TET_I1_FEM.hpp"
52#include "Intrepid2_HCURL_TET_In_FEM.hpp"
53
54#include "Intrepid2_HCURL_QUAD_I1_FEM.hpp"
55#include "Intrepid2_HCURL_QUAD_In_FEM.hpp"
56
57#include "Intrepid2_HCURL_HEX_I1_FEM.hpp"
58#include "Intrepid2_HCURL_HEX_In_FEM.hpp"
59
60#include "Intrepid2_HDIV_TRI_I1_FEM.hpp"
61#include "Intrepid2_HDIV_TRI_In_FEM.hpp"
62
63#include "Intrepid2_HDIV_QUAD_I1_FEM.hpp"
64#include "Intrepid2_HDIV_QUAD_In_FEM.hpp"
65
66#include "Intrepid2_HDIV_TET_I1_FEM.hpp"
67#include "Intrepid2_HDIV_TET_In_FEM.hpp"
68
69#include "Intrepid2_HDIV_HEX_I1_FEM.hpp"
70#include "Intrepid2_HDIV_HEX_In_FEM.hpp"
71
72
73namespace panzer {
74
75
89 template <typename ExecutionSpace, typename OutputValueType, typename PointValueType>
90 Teuchos::RCP<Intrepid2::Basis<ExecutionSpace,OutputValueType,PointValueType> >
91 createIntrepid2Basis(const std::string basis_type, int basis_order,
92 const shards::CellTopology & cell_topology)
93 {
94 // Shards supports extended topologies so the names have a "size"
95 // associated with the number of nodes. We prune the size to
96 // avoid combinatorial explosion of checks.
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);
101
102 Teuchos::RCP<Intrepid2::Basis<ExecutionSpace,OutputValueType,PointValueType> > basis;
103
104 // high order point distribution type;
105 // for now equispaced only; to get a permutation map with different orientation,
106 // orientation coeff matrix should have the same point distribution.
107 const Intrepid2::EPointType point_type = Intrepid2::POINTTYPE_EQUISPACED;
108
109 if ( (basis_type == "Const") && (basis_order == 0) )
110 basis = Teuchos::rcp( new Intrepid2::Basis_HVOL_C0_FEM<ExecutionSpace,OutputValueType,PointValueType>(cell_topology) );
111
112 else if ( (basis_type == "HVol") && (basis_order == 0) )
113 basis = Teuchos::rcp( new Intrepid2::Basis_HVOL_C0_FEM<ExecutionSpace,OutputValueType,PointValueType>(cell_topology) );
114
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) );
117
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) );
120
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) );
123
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) );
126
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> );
129
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> );
132
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) );
135
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> );
138
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) );
141
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> );
144
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) );
147
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> );
150
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> );
153
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) );
156
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> );
159
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) );
162
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> );
165
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) );
168
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> );
171
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> );
174
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) );
177
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> );
180
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) );
183
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> );
186
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) );
189
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> );
192
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> );
195
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) );
198
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> );
201
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) );
204
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> );
207
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) );
210
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> );
213
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) );
216
217 TEUCHOS_TEST_FOR_EXCEPTION(Teuchos::is_null(basis), std::runtime_error,
218 "Failed to create the requestedbasis with basis_type=\"" << basis_type <<
219 "\", basis_order=\"" << basis_order << "\", and cell_type=\"" << cell_type << "\"!\n");
220
221 // we compare that the base topologies are the same
222 // we do so using the NAME. This avoids the ugly task of getting the
223 // cell topology data and constructing a new cell topology object since you cant
224 // just get the baseCellTopology directly from a shards cell topology
225 TEUCHOS_TEST_FOR_EXCEPTION(cell_topology.getBaseName()!=basis->getBaseCellTopology().getName(),
226 std::runtime_error,
227 "Failed to create basis. Intrepid2 basis base topology does not match mesh cell base topology!");
228
229 return basis;
230 }
231
246 template <typename ExecutionSpace, typename OutputValueType, typename PointValueType>
247 Teuchos::RCP<Intrepid2::Basis<ExecutionSpace,OutputValueType,PointValueType> >
248 createIntrepid2Basis(const std::string basis_type, int basis_order,
249 const Teuchos::RCP<const shards::CellTopology> & cell_topology)
250 {
251 return createIntrepid2Basis<ExecutionSpace,OutputValueType,PointValueType>(basis_type,basis_order,*cell_topology);
252 }
253
254}
255
256
257#endif
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.