49 #ifndef __INTREPID2_DEFAULT_CUBATURE_FACTORY_DEF_HPP__
50 #define __INTREPID2_DEFAULT_CUBATURE_FACTORY_DEF_HPP__
55 template<
typename SpT,
typename PT,
typename WT>
56 Teuchos::RCP<Cubature<SpT,PT,WT> >
58 create(
unsigned topologyKey,
59 const std::vector<ordinal_type> °ree,
60 const EPolyType polytype ) {
63 Teuchos::RCP<Cubature<SpT,PT,WT> > r_val;
65 switch (topologyKey) {
66 case shards::Line<>::key: {
67 INTREPID2_TEST_FOR_EXCEPTION( degree.size() < 1, std::invalid_argument,
68 ">>> ERROR (DefaultCubatureFactory): Provided degree array is of insufficient length.");
69 if (isValidPolyType(polytype))
70 r_val = Teuchos::rcp(
new CubaturePolylib<SpT,PT,WT>(degree[0], polytype));
72 r_val = Teuchos::rcp(
new CubatureDirectLineGauss<SpT,PT,WT>(degree[0]));
75 case shards::Triangle<>::key: {
76 INTREPID2_TEST_FOR_EXCEPTION( degree.size() < 1, std::invalid_argument,
77 ">>> ERROR (DefaultCubatureFactory): Provided degree array is of insufficient length.");
78 r_val = Teuchos::rcp(
new CubatureDirectTriDefault<SpT,PT,WT>(degree[0]));
81 case shards::Quadrilateral<>::key:
82 case shards::ShellQuadrilateral<>::key: {
83 INTREPID2_TEST_FOR_EXCEPTION( degree.size() < 2, std::invalid_argument,
84 ">>> ERROR (DefaultCubatureFactory): Provided degree array is of insufficient length.");
85 if (isValidPolyType(polytype)) {
86 const auto x_line = CubaturePolylib<SpT,PT,WT>(degree[0], polytype);
87 const auto y_line = ( degree[1] == degree[0] ? x_line : CubaturePolylib<SpT,PT,WT>(degree[1], polytype) );
88 r_val = Teuchos::rcp(
new CubatureTensor<SpT,PT,WT>( x_line, y_line ));
90 const auto x_line = CubatureDirectLineGauss<SpT,PT,WT>(degree[0]);
91 const auto y_line = ( degree[1] == degree[0] ? x_line : CubatureDirectLineGauss<SpT,PT,WT>(degree[1]) );
92 r_val = Teuchos::rcp(
new CubatureTensor<SpT,PT,WT>( x_line, y_line ));
96 case shards::Tetrahedron<>::key: {
97 INTREPID2_TEST_FOR_EXCEPTION( degree.size() < 1, std::invalid_argument,
98 ">>> ERROR (DefaultCubatureFactory): Provided degree array is of insufficient length.");
99 r_val = Teuchos::rcp(
new CubatureDirectTetDefault<SpT,PT,WT>(degree[0]));
102 case shards::Hexahedron<>::key: {
103 INTREPID2_TEST_FOR_EXCEPTION( degree.size() < 3, std::invalid_argument,
104 ">>> ERROR (DefaultCubatureFactory): Provided degree array is of insufficient length.");
105 if (isValidPolyType(polytype)) {
106 const auto x_line = CubaturePolylib<SpT,PT,WT>(degree[0], polytype);
107 const auto y_line = ( degree[1] == degree[0] ? x_line : CubaturePolylib<SpT,PT,WT>(degree[1], polytype) );
108 const auto z_line = ( degree[2] == degree[0] ? x_line :
109 degree[2] == degree[1] ? y_line : CubaturePolylib<SpT,PT,WT>(degree[2], polytype) );
111 r_val = Teuchos::rcp(
new CubatureTensor<SpT,PT,WT>( x_line, y_line, z_line ));
113 const auto x_line = CubatureDirectLineGauss<SpT,PT,WT>(degree[0]);
114 const auto y_line = ( degree[1] == degree[0] ? x_line : CubatureDirectLineGauss<SpT,PT,WT>(degree[1]) );
115 const auto z_line = ( degree[2] == degree[0] ? x_line :
116 degree[2] == degree[1] ? y_line : CubatureDirectLineGauss<SpT,PT,WT>(degree[2]) );
118 r_val = Teuchos::rcp(
new CubatureTensor<SpT,PT,WT>( x_line, y_line, z_line ));
122 case shards::Wedge<>::key: {
123 INTREPID2_TEST_FOR_EXCEPTION( degree.size() < 2, std::invalid_argument,
124 ">>> ERROR (DefaultCubatureFactory): Provided degree array is of insufficient length.")
126 const auto xy_tri = CubatureDirectTriDefault<SpT,PT,WT>(degree[0]);
127 const auto z_line = CubatureDirectLineGauss<SpT,PT,WT>(degree[1]);
128 r_val = Teuchos::rcp(
new CubatureTensor<SpT,PT,WT>( xy_tri, z_line ));
132 case shards::Pyramid<>::key: {
133 INTREPID2_TEST_FOR_EXCEPTION( degree.size() < 1, std::invalid_argument,
134 ">>> ERROR (DefaultCubatureFactory): Provided degree array is of insufficient length.");
138 const auto xy_line = CubatureDirectLineGauss<SpT,PT,WT>(degree[0]);
139 const auto z_line = CubatureDirectLineGaussJacobi20<SpT,PT,WT>(degree[0]);
140 r_val = Teuchos::rcp(
new CubatureTensorPyr<SpT,PT,WT>( xy_line, xy_line, z_line ));
145 INTREPID2_TEST_FOR_EXCEPTION( ( (topologyKey != shards::Line<>::key) &&
146 (topologyKey != shards::Triangle<>::key) &&
147 (topologyKey != shards::Quadrilateral<>::key) &&
148 (topologyKey != shards::ShellQuadrilateral<>::key) &&
149 (topologyKey != shards::Tetrahedron<>::key) &&
150 (topologyKey != shards::Hexahedron<>::key) &&
151 (topologyKey != shards::Pyramid<>::key) &&
152 (topologyKey != shards::Wedge<>::key) ),
153 std::invalid_argument,
154 ">>> ERROR (DefaultCubatureFactory): Invalid cell topology prevents cubature creation.");
160 template<
typename SpT,
typename PT,
typename WT>
161 Teuchos::RCP<Cubature<SpT,PT,WT> >
163 create(
const shards::CellTopology cellTopology,
164 const std::vector<ordinal_type> °ree,
165 const EPolyType polytype) {
166 return create<SpT,PT,WT>(cellTopology.getBaseKey(), degree, polytype);
170 template<
typename SpT,
typename PT,
typename WT>
171 Teuchos::RCP<Cubature<SpT,PT,WT> >
173 create(
unsigned topologyKey,
174 const ordinal_type degree,
175 const EPolyType polytype ) {
177 const std::vector<ordinal_type> degreeArray(3, degree);
178 return create<SpT,PT,WT>(topologyKey, degreeArray, polytype);
181 template<
typename SpT,
typename PT,
typename WT>
182 Teuchos::RCP<Cubature<SpT,PT,WT> >
184 create(
const shards::CellTopology cellTopology,
185 const ordinal_type degree,
186 const EPolyType polytype ) {
188 const std::vector<ordinal_type> degreeArray(3, degree);
189 return create<SpT,PT,WT>(cellTopology.getBaseKey(), degreeArray, polytype);
static Teuchos::RCP< Cubature< ExecSpaceType, pointValueType, weightValueType > > create(unsigned topologyKey, const std::vector< ordinal_type > °ree, const EPolyType polytype=POLYTYPE_MAX)
Factory method.