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 INTREPID2_TEST_FOR_EXCEPTION( degree.size() < 1, std::invalid_argument,
36 ">>> ERROR (DefaultCubatureFactory): Provided degree array is of insufficient length.");
37 if (isValidPolyType(polytype))
38 r_val = Teuchos::rcp(new CubaturePolylib<DT,PT,WT>(degree[0], polytype));
39 else
40 r_val = Teuchos::rcp(new CubatureDirectLineGauss<DT,PT,WT>(degree[0]));
41 break;
42 }
43 case shards::Triangle<>::key: {
44 INTREPID2_TEST_FOR_EXCEPTION( degree.size() < 1, std::invalid_argument,
45 ">>> ERROR (DefaultCubatureFactory): Provided degree array is of insufficient length.");
46 if(symmetric || (degree[0] > 20)) {
47 r_val = Teuchos::rcp(new CubatureDirectTriSymmetric<DT,PT,WT>(degree[0])); }
48 else
49 r_val = Teuchos::rcp(new CubatureDirectTriDefault<DT,PT,WT>(degree[0]));
50 break;
51 }
52 case shards::Quadrilateral<>::key:
53 case shards::ShellQuadrilateral<>::key: {
54 INTREPID2_TEST_FOR_EXCEPTION( degree.size() < 2, std::invalid_argument,
55 ">>> ERROR (DefaultCubatureFactory): Provided degree array is of insufficient length.");
56 if (isValidPolyType(polytype)) {
57 const auto x_line = CubaturePolylib<DT,PT,WT>(degree[0], polytype);
58 const auto y_line = ( degree[1] == degree[0] ? x_line : CubaturePolylib<DT,PT,WT>(degree[1], polytype) );
59 r_val = Teuchos::rcp(new CubatureTensor<DT,PT,WT>( x_line, y_line ));
60 } else {
61 const auto x_line = CubatureDirectLineGauss<DT,PT,WT>(degree[0]);
62 const auto y_line = ( degree[1] == degree[0] ? x_line : CubatureDirectLineGauss<DT,PT,WT>(degree[1]) );
63 r_val = Teuchos::rcp(new CubatureTensor<DT,PT,WT>( x_line, y_line ));
64 }
65 break;
66 }
67 case shards::Tetrahedron<>::key: {
68 INTREPID2_TEST_FOR_EXCEPTION( degree.size() < 1, std::invalid_argument,
69 ">>> ERROR (DefaultCubatureFactory): Provided degree array is of insufficient length.");
70 if(symmetric) {
71 r_val = Teuchos::rcp(new CubatureDirectTetSymmetric<DT,PT,WT>(degree[0]));
72 } else {
73 r_val = Teuchos::rcp(new CubatureDirectTetDefault<DT,PT,WT>(degree[0]));
74 }
75 break;
76 }
77 case shards::Hexahedron<>::key: {
78 INTREPID2_TEST_FOR_EXCEPTION( degree.size() < 3, std::invalid_argument,
79 ">>> ERROR (DefaultCubatureFactory): Provided degree array is of insufficient length.");
80 if (isValidPolyType(polytype)) {
81 const auto x_line = CubaturePolylib<DT,PT,WT>(degree[0], polytype);
82 const auto y_line = ( degree[1] == degree[0] ? x_line : CubaturePolylib<DT,PT,WT>(degree[1], polytype) );
83 const auto z_line = ( degree[2] == degree[0] ? x_line :
84 degree[2] == degree[1] ? y_line : CubaturePolylib<DT,PT,WT>(degree[2], polytype) );
85
86 r_val = Teuchos::rcp(new CubatureTensor<DT,PT,WT>( x_line, y_line, z_line ));
87 } else {
88 const auto x_line = CubatureDirectLineGauss<DT,PT,WT>(degree[0]);
89 const auto y_line = ( degree[1] == degree[0] ? x_line : CubatureDirectLineGauss<DT,PT,WT>(degree[1]) );
90 const auto z_line = ( degree[2] == degree[0] ? x_line :
91 degree[2] == degree[1] ? y_line : CubatureDirectLineGauss<DT,PT,WT>(degree[2]) );
92
93 r_val = Teuchos::rcp(new CubatureTensor<DT,PT,WT>( x_line, y_line, z_line ));
94 }
95 break;
96 }
97 case shards::Wedge<>::key: {
98 INTREPID2_TEST_FOR_EXCEPTION( degree.size() < 2, std::invalid_argument,
99 ">>> ERROR (DefaultCubatureFactory): Provided degree array is of insufficient length.")
100 {
101 if(symmetric || (degree[0] > 20)) {
102 const auto xy_tri = CubatureDirectTriSymmetric<DT,PT,WT>(degree[0]);
103 const auto z_line = CubatureDirectLineGauss<DT,PT,WT>(degree[1]);
104 r_val = Teuchos::rcp(new CubatureTensor<DT,PT,WT>( xy_tri, z_line ));
105 }
106 else {
107 const auto xy_tri = CubatureDirectTriDefault<DT,PT,WT>(degree[0]);
108 const auto z_line = CubatureDirectLineGauss<DT,PT,WT>(degree[1]);
109 r_val = Teuchos::rcp(new CubatureTensor<DT,PT,WT>( xy_tri, z_line ));
110 }
111 }
112 break;
113 }
114 case shards::Pyramid<>::key: {
115 INTREPID2_TEST_FOR_EXCEPTION( degree.size() < 1, std::invalid_argument,
116 ">>> ERROR (DefaultCubatureFactory): Provided degree array is of insufficient length.");
117 {
118 // if direct gauss is used for pyramid
119 // need over-integration to account for the additional transformation
120 const auto xy_line = CubatureDirectLineGauss<DT,PT,WT>(degree[0]);
121 const auto z_line = CubatureDirectLineGaussJacobi20<DT,PT,WT>(degree[0]);
122 r_val = Teuchos::rcp(new CubatureTensorPyr<DT,PT,WT>( xy_line, xy_line, z_line ));
123 }
124 break;
125 }
126 default: {
127 INTREPID2_TEST_FOR_EXCEPTION( ( (topologyKey != shards::Line<>::key) &&
128 (topologyKey != shards::Triangle<>::key) &&
129 (topologyKey != shards::Quadrilateral<>::key) &&
130 (topologyKey != shards::ShellQuadrilateral<>::key) &&
131 (topologyKey != shards::Tetrahedron<>::key) &&
132 (topologyKey != shards::Hexahedron<>::key) &&
133 (topologyKey != shards::Pyramid<>::key) &&
134 (topologyKey != shards::Wedge<>::key) ),
135 std::invalid_argument,
136 ">>> ERROR (DefaultCubatureFactory): Invalid cell topology prevents cubature creation.");
137 }
138 }
139 return r_val;
140 }
141
142 template<typename DT, typename PT, typename WT>
143 Teuchos::RCP<Cubature<DT,PT,WT> >
145 create( const shards::CellTopology cellTopology,
146 const std::vector<ordinal_type> &degree,
147 const EPolyType polytype,
148 const bool symmetric ) {
149 return create<DT,PT,WT>(cellTopology.getBaseKey(), degree, polytype, symmetric);
150 }
151
152
153 template<typename DT, typename PT, typename WT>
154 Teuchos::RCP<Cubature<DT,PT,WT> >
156 create( unsigned topologyKey,
157 const ordinal_type degree,
158 const EPolyType polytype,
159 const bool symmetric ) {
160 // uniform order for 3 axes
161 const std::vector<ordinal_type> degreeArray(3, degree);
162 return create<DT,PT,WT>(topologyKey, degreeArray, polytype, symmetric);
163 }
164
165 template<typename DT, typename PT, typename WT>
166 Teuchos::RCP<Cubature<DT,PT,WT> >
168 create( const shards::CellTopology cellTopology,
169 const ordinal_type degree,
170 const EPolyType polytype,
171 const bool symmetric ) {
172 // uniform order for 3 axes
173 const std::vector<ordinal_type> degreeArray(3, degree);
174 return create<DT,PT,WT>(cellTopology.getBaseKey(), degreeArray, polytype, symmetric);
175 }
176
177
178 // template<typename DT>
179 // template<typename cellVertexValueType, class ...cellVertexProperties>
180 // Teuchos::RCP<Cubature<DT,PT,WT> >
181 // DefaultCubatureFactory::
182 // create<DT,PT,WT>(const shards::CellTopology& cellTopology,
183 // const Kokkos::DynRankView<cellVertexValueType,cellVertexProperties> cellVertices,
184 // ordinal_type degree){
185 // return Teuchos::rcp(new CubaturePolygon<DT,PT,WT>(cellTopology,cellVertices, degree));
186 // }
187
188} // namespace Intrepid2
189
190#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.