Intrepid2
Intrepid2_DefaultCubatureFactoryDef.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Intrepid2 Package
4//
5// Copyright 2007 NTESS and the Intrepid2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
16#ifndef __INTREPID2_DEFAULT_CUBATURE_FACTORY_DEF_HPP__
17#define __INTREPID2_DEFAULT_CUBATURE_FACTORY_DEF_HPP__
18
19namespace Intrepid2 {
20
21 // first create method
22 template<typename DT, typename PT, typename WT>
23 Teuchos::RCP<Cubature<DT,PT,WT> >
25 create( unsigned topologyKey,
26 const std::vector<ordinal_type> &degree,
27 const EPolyType polytype,
28 const bool symmetric ) {
29
30 // Create generic cubature.
31 Teuchos::RCP<Cubature<DT,PT,WT> > r_val;
32
33 switch (topologyKey) {
34 case shards::Line<>::key:
35 case shards::Beam<>::key:
36 case shards::ShellLine<>::key: {
37 INTREPID2_TEST_FOR_EXCEPTION( degree.size() < 1, std::invalid_argument,
38 ">>> ERROR (DefaultCubatureFactory): Provided degree array is of insufficient length.");
39 if (isValidPolyType(polytype))
40 r_val = Teuchos::rcp(new CubaturePolylib<DT,PT,WT>(degree[0], polytype));
41 else
42 r_val = Teuchos::rcp(new CubatureDirectLineGauss<DT,PT,WT>(degree[0]));
43 break;
44 }
45 case shards::Triangle<>::key:
46 case shards::ShellTriangle<>::key: {
47 INTREPID2_TEST_FOR_EXCEPTION( degree.size() < 1, std::invalid_argument,
48 ">>> ERROR (DefaultCubatureFactory): Provided degree array is of insufficient length.");
49 if(symmetric || (degree[0] > 20)) {
50 r_val = Teuchos::rcp(new CubatureDirectTriSymmetric<DT,PT,WT>(degree[0])); }
51 else
52 r_val = Teuchos::rcp(new CubatureDirectTriDefault<DT,PT,WT>(degree[0]));
53 break;
54 }
55 case shards::Quadrilateral<>::key:
56 case shards::ShellQuadrilateral<>::key: {
57 INTREPID2_TEST_FOR_EXCEPTION( degree.size() < 2, std::invalid_argument,
58 ">>> ERROR (DefaultCubatureFactory): Provided degree array is of insufficient length.");
59 if (isValidPolyType(polytype)) {
60 const auto x_line = CubaturePolylib<DT,PT,WT>(degree[0], polytype);
61 const auto y_line = ( degree[1] == degree[0] ? x_line : CubaturePolylib<DT,PT,WT>(degree[1], polytype) );
62 r_val = Teuchos::rcp(new CubatureTensor<DT,PT,WT>( x_line, y_line ));
63 } else {
64 const auto x_line = CubatureDirectLineGauss<DT,PT,WT>(degree[0]);
65 const auto y_line = ( degree[1] == degree[0] ? x_line : CubatureDirectLineGauss<DT,PT,WT>(degree[1]) );
66 r_val = Teuchos::rcp(new CubatureTensor<DT,PT,WT>( x_line, y_line ));
67 }
68 break;
69 }
70 case shards::Tetrahedron<>::key: {
71 INTREPID2_TEST_FOR_EXCEPTION( degree.size() < 1, std::invalid_argument,
72 ">>> ERROR (DefaultCubatureFactory): Provided degree array is of insufficient length.");
73 if(symmetric) {
74 r_val = Teuchos::rcp(new CubatureDirectTetSymmetric<DT,PT,WT>(degree[0]));
75 } else {
76 r_val = Teuchos::rcp(new CubatureDirectTetDefault<DT,PT,WT>(degree[0]));
77 }
78 break;
79 }
80 case shards::Hexahedron<>::key: {
81 INTREPID2_TEST_FOR_EXCEPTION( degree.size() < 3, std::invalid_argument,
82 ">>> ERROR (DefaultCubatureFactory): Provided degree array is of insufficient length.");
83 if (isValidPolyType(polytype)) {
84 const auto x_line = CubaturePolylib<DT,PT,WT>(degree[0], polytype);
85 const auto y_line = ( degree[1] == degree[0] ? x_line : CubaturePolylib<DT,PT,WT>(degree[1], polytype) );
86 const auto z_line = ( degree[2] == degree[0] ? x_line :
87 degree[2] == degree[1] ? y_line : CubaturePolylib<DT,PT,WT>(degree[2], polytype) );
88
89 r_val = Teuchos::rcp(new CubatureTensor<DT,PT,WT>( x_line, y_line, z_line ));
90 } else {
91 const auto x_line = CubatureDirectLineGauss<DT,PT,WT>(degree[0]);
92 const auto y_line = ( degree[1] == degree[0] ? x_line : CubatureDirectLineGauss<DT,PT,WT>(degree[1]) );
93 const auto z_line = ( degree[2] == degree[0] ? x_line :
94 degree[2] == degree[1] ? y_line : CubatureDirectLineGauss<DT,PT,WT>(degree[2]) );
95
96 r_val = Teuchos::rcp(new CubatureTensor<DT,PT,WT>( x_line, y_line, z_line ));
97 }
98 break;
99 }
100 case shards::Wedge<>::key: {
101 INTREPID2_TEST_FOR_EXCEPTION( degree.size() < 2, std::invalid_argument,
102 ">>> ERROR (DefaultCubatureFactory): Provided degree array is of insufficient length.")
103 {
104 if(symmetric || (degree[0] > 20)) {
105 const auto xy_tri = CubatureDirectTriSymmetric<DT,PT,WT>(degree[0]);
106 const auto z_line = CubatureDirectLineGauss<DT,PT,WT>(degree[1]);
107 r_val = Teuchos::rcp(new CubatureTensor<DT,PT,WT>( xy_tri, z_line ));
108 }
109 else {
110 const auto xy_tri = CubatureDirectTriDefault<DT,PT,WT>(degree[0]);
111 const auto z_line = CubatureDirectLineGauss<DT,PT,WT>(degree[1]);
112 r_val = Teuchos::rcp(new CubatureTensor<DT,PT,WT>( xy_tri, z_line ));
113 }
114 }
115 break;
116 }
117 case shards::Pyramid<>::key: {
118 INTREPID2_TEST_FOR_EXCEPTION( degree.size() < 1, std::invalid_argument,
119 ">>> ERROR (DefaultCubatureFactory): Provided degree array is of insufficient length.");
120 {
121 // if direct gauss is used for pyramid
122 // need over-integration to account for the additional transformation
123 const auto xy_line = CubatureDirectLineGauss<DT,PT,WT>(degree[0]);
124 const auto z_line = CubatureDirectLineGaussJacobi20<DT,PT,WT>(degree[0]);
125 r_val = Teuchos::rcp(new CubatureTensorPyr<DT,PT,WT>( xy_line, xy_line, z_line ));
126 }
127 break;
128 }
129 default: {
130 INTREPID2_TEST_FOR_EXCEPTION( ( (topologyKey != shards::Line<>::key) &&
131 (topologyKey != shards::Beam<>::key) &&
132 (topologyKey != shards::ShellLine<>::key) &&
133 (topologyKey != shards::Triangle<>::key) &&
134 (topologyKey != shards::ShellTriangle<>::key) &&
135 (topologyKey != shards::Quadrilateral<>::key) &&
136 (topologyKey != shards::ShellQuadrilateral<>::key) &&
137 (topologyKey != shards::Tetrahedron<>::key) &&
138 (topologyKey != shards::Hexahedron<>::key) &&
139 (topologyKey != shards::Pyramid<>::key) &&
140 (topologyKey != shards::Wedge<>::key) ),
141 std::invalid_argument,
142 ">>> ERROR (DefaultCubatureFactory): Invalid cell topology prevents cubature creation.");
143 }
144 }
145 return r_val;
146 }
147
148 template<typename DT, typename PT, typename WT>
149 Teuchos::RCP<Cubature<DT,PT,WT> >
151 create( const shards::CellTopology cellTopology,
152 const std::vector<ordinal_type> &degree,
153 const EPolyType polytype,
154 const bool symmetric ) {
155 return create<DT,PT,WT>(cellTopology.getBaseKey(), degree, polytype, symmetric);
156 }
157
158
159 template<typename DT, typename PT, typename WT>
160 Teuchos::RCP<Cubature<DT,PT,WT> >
162 create( unsigned topologyKey,
163 const ordinal_type degree,
164 const EPolyType polytype,
165 const bool symmetric ) {
166 // uniform order for 3 axes
167 const std::vector<ordinal_type> degreeArray(3, degree);
168 return create<DT,PT,WT>(topologyKey, degreeArray, polytype, symmetric);
169 }
170
171 template<typename DT, typename PT, typename WT>
172 Teuchos::RCP<Cubature<DT,PT,WT> >
174 create( const shards::CellTopology cellTopology,
175 const ordinal_type degree,
176 const EPolyType polytype,
177 const bool symmetric ) {
178 // uniform order for 3 axes
179 const std::vector<ordinal_type> degreeArray(3, degree);
180 return create<DT,PT,WT>(cellTopology.getBaseKey(), degreeArray, polytype, symmetric);
181 }
182
183
184 // template<typename DT>
185 // template<typename cellVertexValueType, class ...cellVertexProperties>
186 // Teuchos::RCP<Cubature<DT,PT,WT> >
187 // DefaultCubatureFactory::
188 // create<DT,PT,WT>(const shards::CellTopology& cellTopology,
189 // const Kokkos::DynRankView<cellVertexValueType,cellVertexProperties> cellVertices,
190 // ordinal_type degree){
191 // return Teuchos::rcp(new CubaturePolygon<DT,PT,WT>(cellTopology,cellVertices, degree));
192 // }
193
194} // namespace Intrepid2
195
196#endif
KOKKOS_FORCEINLINE_FUNCTION bool isValidPolyType(const EPolyType polytype)
Verifies validity of a PolyType enum.
static Teuchos::RCP< Cubature< DeviceType, pointValueType, weightValueType > > create(unsigned topologyKey, const std::vector< ordinal_type > &degree, const EPolyType polytype=POLYTYPE_MAX, const bool symmetric=false)
Factory method.