51 template <
class Scalar,
class ArrayPo
int,
class ArrayWeight>
53 int degree, EIntrepidBurkardt rule,
bool isNormalized ) {
54 TEUCHOS_TEST_FOR_EXCEPTION((degree < 0),std::out_of_range,
55 ">>> ERROR (CubatureLineSorted): No rule implemented for desired polynomial degree.");
59 if (rule==BURK_CHEBYSHEV1||rule==BURK_CHEBYSHEV2||
60 rule==BURK_LAGUERRE ||rule==BURK_LEGENDRE ||
62 numPoints_ = (degree+1)/2+1;
64 else if (rule==BURK_CLENSHAWCURTIS||rule==BURK_FEJER2) {
65 numPoints_ = degree+1;
67 else if (rule==BURK_TRAPEZOIDAL) {
70 else if (rule==BURK_PATTERSON) {
71 int l = 0, o = (degree-0.5)/1.5;
72 for (
int i=0; i<8; i++) {
73 l = (int)pow(2.0,(
double)i+1.0)-1;
80 else if (rule==BURK_GENZKEISTER) {
81 int o_ghk[8] = {1,3,9,19,35,37,41,43};
82 int o = (degree-0.5)/1.5;
83 for (
int i=0; i<8; i++) {
85 numPoints_ = o_ghk[i];
91 Teuchos::Array<Scalar> nodes(numPoints_), weights(numPoints_);
93 if (rule==BURK_CHEBYSHEV1) {
94 IntrepidBurkardtRules::chebyshev1_compute<Scalar>(numPoints_,
95 nodes.getRawPtr(),weights.getRawPtr());
97 else if (rule==BURK_CHEBYSHEV2) {
98 IntrepidBurkardtRules::chebyshev2_compute<Scalar>(numPoints_,
99 nodes.getRawPtr(),weights.getRawPtr());
101 else if (rule==BURK_CLENSHAWCURTIS) {
102 IntrepidBurkardtRules::clenshaw_curtis_compute<Scalar>(numPoints_,
103 nodes.getRawPtr(),weights.getRawPtr());
105 else if (rule==BURK_FEJER2) {
106 IntrepidBurkardtRules::fejer2_compute<Scalar>(numPoints_,
107 nodes.getRawPtr(),weights.getRawPtr());
109 else if (rule==BURK_LEGENDRE) {
110 IntrepidBurkardtRules::legendre_compute<Scalar>(numPoints_,
111 nodes.getRawPtr(),weights.getRawPtr());
113 else if (rule==BURK_PATTERSON) {
114 IntrepidBurkardtRules::patterson_lookup<Scalar>(numPoints_,
115 nodes.getRawPtr(),weights.getRawPtr());
117 else if (rule==BURK_TRAPEZOIDAL) {
118 IntrepidBurkardtRules::trapezoidal_compute<Scalar>(numPoints_,
119 nodes.getRawPtr(),weights.getRawPtr());
121 else if (rule==BURK_HERMITE) {
122 IntrepidBurkardtRules::hermite_compute<Scalar>(numPoints_,
123 nodes.getRawPtr(),weights.getRawPtr());
125 else if (rule==BURK_GENZKEISTER) {
126 IntrepidBurkardtRules::hermite_genz_keister_lookup<Scalar>(numPoints_,
127 nodes.getRawPtr(),weights.getRawPtr());
128 if (numPoints_>=37) {
129 for (
int i=0; i<numPoints_; i++) {
130 weights[i] *= sqrt(M_PI);
134 else if (rule==BURK_LAGUERRE) {
135 IntrepidBurkardtRules::laguerre_compute<Scalar>(numPoints_,
136 nodes.getRawPtr(),weights.getRawPtr());
141 for (
int i=0; i<numPoints_; i++) {
144 for (
int i=0; i<numPoints_; i++) {
149 points_.clear(); weights_.clear();
150 typename std::map<Scalar,int>::iterator it(points_.begin());
151 for (
int i=0; i<numPoints_; i++) {
152 points_.insert(it,std::pair<Scalar,int>(nodes[i],i));
153 weights_.push_back(weights[i]);
158 template <
class Scalar,
class ArrayPo
int,
class ArrayWeight>
160 EIntrepidBurkardt rule,
int numPoints,
bool isNormalized ) {
161 TEUCHOS_TEST_FOR_EXCEPTION((numPoints < 0),std::out_of_range,
162 ">>> ERROR (CubatureLineSorted): No rule implemented for desired number of points.");
163 numPoints_ = numPoints;
166 Teuchos::Array<Scalar> nodes(numPoints_), weights(numPoints_);
167 if (rule==BURK_CHEBYSHEV1) {
168 degree_ = 2*numPoints-1;
169 IntrepidBurkardtRules::chebyshev1_compute<Scalar>(numPoints_,
170 nodes.getRawPtr(),weights.getRawPtr());
172 else if (rule==BURK_CHEBYSHEV2) {
173 degree_ = 2*numPoints-1;
174 IntrepidBurkardtRules::chebyshev2_compute<Scalar>(numPoints_,
175 nodes.getRawPtr(),weights.getRawPtr());
177 else if (rule==BURK_CLENSHAWCURTIS) {
178 degree_ = numPoints-1;
179 IntrepidBurkardtRules::clenshaw_curtis_compute<Scalar>(numPoints_,
180 nodes.getRawPtr(),weights.getRawPtr());
182 else if (rule==BURK_FEJER2) {
183 degree_ = numPoints-1;
184 IntrepidBurkardtRules::fejer2_compute<Scalar>(numPoints_,
185 nodes.getRawPtr(),weights.getRawPtr());
187 else if (rule==BURK_LEGENDRE) {
188 degree_ = 2*numPoints-1;
189 IntrepidBurkardtRules::legendre_compute<Scalar>(numPoints_,
190 nodes.getRawPtr(),weights.getRawPtr());
192 else if (rule==BURK_PATTERSON) {
193 bool correctNumPoints =
false;
194 for (
int i=0; i<8; i++) {
195 int l = (int)pow(2.0,(
double)i+1.0)-1;
197 correctNumPoints =
true;
201 TEUCHOS_TEST_FOR_EXCEPTION((correctNumPoints==
false),std::out_of_range,
202 ">>> ERROR (CubatureLineSorted): Number of points must be numPoints = 1, 3, 7, 15, 31, 63, 127, 255.");
203 Scalar degree = 1.5*(double)numPoints+0.5;
204 degree_ = (int)degree;
205 IntrepidBurkardtRules::patterson_lookup<Scalar>(numPoints_,
206 nodes.getRawPtr(),weights.getRawPtr());
208 else if (rule==BURK_TRAPEZOIDAL) {
210 IntrepidBurkardtRules::trapezoidal_compute<Scalar>(numPoints_,
211 nodes.getRawPtr(),weights.getRawPtr());
213 else if (rule==BURK_HERMITE) {
214 degree_ = 2*numPoints-1;
215 IntrepidBurkardtRules::hermite_compute<Scalar>(numPoints_,
216 nodes.getRawPtr(),weights.getRawPtr());
218 else if (rule==BURK_GENZKEISTER) {
219 bool correctNumPoints =
false;
220 int o_ghk[8] = {1,3,9,19,35,37,41,43};
221 for (
int i=0; i<8; i++) {
222 if (o_ghk[i]==numPoints) {
223 correctNumPoints =
true;
227 TEUCHOS_TEST_FOR_EXCEPTION((correctNumPoints==
false),std::out_of_range,
228 ">>> ERROR (CubatureLineSorted): Number of points must be numPoints = 1, 3, 9, 35, 37, 41, 43.");
229 Scalar degree = 1.5*(double)numPoints+0.5;
230 degree_ = (int)degree;
231 IntrepidBurkardtRules::hermite_genz_keister_lookup<Scalar>(numPoints_,
232 nodes.getRawPtr(),weights.getRawPtr());
233 if (numPoints_>=37) {
234 for (
int i=0; i<numPoints_; i++) {
235 weights[i] *= sqrt(M_PI);
239 else if (rule==BURK_LAGUERRE) {
240 degree_ = 2*numPoints-1;
241 IntrepidBurkardtRules::laguerre_compute<Scalar>(numPoints_,
242 nodes.getRawPtr(),weights.getRawPtr());
247 for (
int i=0; i<numPoints_; i++) {
250 for (
int i=0; i<numPoints_; i++) {
254 points_.clear(); weights_.clear();
255 typename std::map<Scalar,int>::iterator it(points_.begin());
256 for (
int i=0; i<numPoints; i++) {
257 points_.insert(it,std::pair<Scalar,int>(nodes[i],i));
258 weights_.push_back(weights[i]);
263 template <
class Scalar,
class ArrayPo
int,
class ArrayWeight>
265 std::vector<Scalar> & points, std::vector<Scalar> & weights) {
267 int size = (int)weights.size();
268 TEUCHOS_TEST_FOR_EXCEPTION(((
int)points.size()!=size),std::out_of_range,
269 ">>> ERROR (CubatureLineSorted): Input dimension mismatch.");
270 points_.clear(); weights.clear();
271 for (
int loc=0; loc<size; loc++) {
272 points_.insert(std::pair<Scalar,int>(points[loc],loc));
273 weights_.push_back(weights[loc]);
278 template <
class Scalar,
class ArrayPo
int,
class ArrayWeight>
280 return cubature_name_;
283 template <
class Scalar,
class ArrayPo
int,
class ArrayWeight>
288 template <
class Scalar,
class ArrayPo
int,
class ArrayWeight>
290 std::vector<int> & accuracy)
const {
291 accuracy.assign(1, degree_);
294 template <
class Scalar,
class ArrayPo
int,
class ArrayWeight>
297 template <
class Scalar,
class ArrayPo
int,
class ArrayWeight>
299 ArrayPoint & cubPoints, ArrayWeight & cubWeights)
const {
300 typename std::map<Scalar,int>::const_iterator it;
302 for (it = points_.begin(); it!=points_.end(); it++) {
303 cubPoints(i) = it->first;
304 cubWeights(i) = weights_[it->second];
309 template <
class Scalar,
class ArrayPo
int,
class ArrayWeight>
311 ArrayWeight& cubWeights,
312 ArrayPoint& cellCoords)
const
314 TEUCHOS_TEST_FOR_EXCEPTION( (
true), std::logic_error,
315 ">>> ERROR (CubatureLineSorted): Cubature defined in reference space calling method for physical space cubature.");
318 template <
class Scalar,
class ArrayPo
int,
class ArrayWeight>
323 template <
class Scalar,
class ArrayPo
int,
class ArrayWeight>
325 typename std::map<Scalar,int>::iterator it) {
329 template <
class Scalar,
class ArrayPo
int,
class ArrayWeight>
331 return weights_[node];
334 template <
class Scalar,
class ArrayPo
int,
class ArrayWeight>
337 return weights_[points_[point]];
341 template <
class Scalar,
class ArrayPo
int,
class ArrayWeight>
343 return points_.begin();
346 template <
class Scalar,
class ArrayPo
int,
class ArrayWeight>
348 return points_.end();
351 template <
class Scalar,
class ArrayPo
int,
class ArrayWeight>
356 typename std::map<Scalar,int>::iterator it;
359 typename std::map<Scalar,int> newPoints;
360 std::vector<Scalar> newWeights;
366 typename std::map<Scalar,int> inter;
367 std::set_intersection(points_.begin(),points_.end(),
369 inserter(inter,inter.begin()),inter.value_comp());
370 for (it=inter.begin(); it!=inter.end(); it++) {
372 newWeights.push_back( alpha1*weights_[it->second]
374 newPoints.insert(std::pair<Scalar,int>(node,loc));
377 int isize = inter.size();
380 int size = weights_.size();
382 typename std::map<Scalar,int> diff1;
383 std::set_difference(points_.begin(),points_.end(),
385 inserter(diff1,diff1.begin()),diff1.value_comp());
386 for (it=diff1.begin(); it!=diff1.end(); it++) {
388 newWeights.push_back(alpha1*weights_[it->second]);
389 newPoints.insert(std::pair<Scalar,int>(node,loc));
397 typename std::map<Scalar,int> diff2;
398 std::set_difference(cubRule2.
begin(),cubRule2.
end(),
399 points_.begin(),points_.end(),
400 inserter(diff2,diff2.begin()),diff2.value_comp());
401 for (it=diff2.begin(); it!=diff2.end(); it++) {
403 newWeights.push_back(alpha2*cubRule2.
getWeight(it->second));
404 newPoints.insert(std::pair<Scalar,int>(node,loc));
409 points_.clear(); points_.insert(newPoints.begin(),newPoints.end());
410 weights_.clear(); weights_.assign(newWeights.begin(),newWeights.end());
411 numPoints_ = (int)points_.size();
414 int growthRule1D(
int index, EIntrepidGrowth growth, EIntrepidBurkardt rule) {
432 if (rule==BURK_CLENSHAWCURTIS) {
433 if (growth==GROWTH_SLOWLIN) {
436 else if (growth==GROWTH_SLOWLINODD) {
437 return 2*((level+1)/2)+1;
439 else if (growth==GROWTH_MODLIN) {
442 else if (growth==GROWTH_SLOWEXP) {
454 else if (growth==GROWTH_MODEXP||growth==GROWTH_DEFAULT) {
460 while (o<4*level+1) {
466 else if (growth==GROWTH_FULLEXP) {
471 return (
int)pow(2.0,(
double)level)+1;
475 else if (rule==BURK_FEJER2) {
476 if (growth==GROWTH_SLOWLIN) {
479 else if (growth==GROWTH_SLOWLINODD) {
480 return 2*((level+1)/2)+1;
482 else if (growth==GROWTH_MODLIN) {
485 else if (growth==GROWTH_SLOWEXP) {
487 while (o<2*level+1) {
492 else if (growth==GROWTH_MODEXP||growth==GROWTH_DEFAULT) {
494 while (o<4*level+1) {
499 else if (growth==GROWTH_FULLEXP) {
500 return (
int)pow(2.0,(
double)level+1.0)-1;
504 else if (rule==BURK_PATTERSON) {
505 if (growth==GROWTH_SLOWLIN||
506 growth==GROWTH_SLOWLINODD||
507 growth==GROWTH_MODLIN) {
508 std::cout <<
"Specified Growth Rule Not Allowed!\n";
511 else if (growth==GROWTH_SLOWEXP) {
518 while (p<2*level+1) {
525 else if (growth==GROWTH_MODEXP||growth==GROWTH_DEFAULT) {
532 while (p<4*level+1) {
539 else if (growth==GROWTH_FULLEXP) {
540 return (
int)pow(2.0,(
double)level+1.0)-1;
544 else if (rule==BURK_LEGENDRE) {
545 if (growth==GROWTH_SLOWLIN) {
548 else if (growth==GROWTH_SLOWLINODD) {
549 return 2*((level+1)/2)+1;
551 else if (growth==GROWTH_MODLIN||growth==GROWTH_DEFAULT) {
554 else if (growth==GROWTH_SLOWEXP) {
556 while (2*o-1<2*level+1) {
561 else if (growth==GROWTH_MODEXP) {
563 while (2*o-1<4*level+1) {
568 else if (growth==GROWTH_FULLEXP) {
569 return (
int)pow(2.0,(
double)level+1.0)-1;
573 else if (rule==BURK_HERMITE) {
574 if (growth==GROWTH_SLOWLIN) {
577 else if (growth==GROWTH_SLOWLINODD) {
578 return 2*((level+1)/2)+1;
580 else if (growth==GROWTH_MODLIN||growth==GROWTH_DEFAULT) {
583 else if (growth==GROWTH_SLOWEXP) {
585 while (2*o-1<2*level+1) {
590 else if (growth==GROWTH_MODEXP) {
592 while (2*o-1<4*level+1) {
597 else if (growth==GROWTH_FULLEXP) {
598 return (
int)pow(2.0,(
double)level+1.0)-1;
602 else if (rule==BURK_LAGUERRE) {
603 if (growth==GROWTH_SLOWLIN) {
606 else if (growth==GROWTH_SLOWLINODD) {
607 return 2*((level+1)/2)+1;
609 else if (growth==GROWTH_MODLIN||growth==GROWTH_DEFAULT) {
612 else if (growth==GROWTH_SLOWEXP) {
614 while (2*o-1<2*level+1) {
619 else if (growth==GROWTH_MODEXP) {
621 while (2*o-1<4*level+1) {
626 else if (growth==GROWTH_FULLEXP) {
627 return (
int)pow(2.0,(
double)level+1.0)-1;
631 else if (rule==BURK_CHEBYSHEV1) {
632 if (growth==GROWTH_SLOWLIN) {
635 else if (growth==GROWTH_SLOWLINODD) {
636 return 2*((level+1)/2)+1;
638 else if (growth==GROWTH_MODLIN||growth==GROWTH_DEFAULT) {
641 else if (growth==GROWTH_SLOWEXP) {
643 while (2*o-1<2*level+1) {
648 else if (growth==GROWTH_MODEXP) {
650 while (2*o-1<4*level+1) {
655 else if (growth==GROWTH_FULLEXP) {
656 return (
int)pow(2.0,(
double)level+1.0)-1;
661 else if (rule==BURK_CHEBYSHEV2) {
662 if (growth==GROWTH_SLOWLIN) {
665 else if (growth==GROWTH_SLOWLINODD) {
666 return 2*((level+1)/2)+1;
668 else if (growth==GROWTH_MODLIN||growth==GROWTH_DEFAULT) {
671 else if (growth==GROWTH_SLOWEXP) {
673 while (2*o-1<2*level+1) {
678 else if (growth==GROWTH_MODEXP) {
680 while (2*o-1<4*level+1) {
685 else if (growth==GROWTH_FULLEXP) {
686 return (
int)pow(2.0,(
double)level+1.0)-1;
690 else if (rule==BURK_GENZKEISTER) {
691 static int o_hgk[5] = { 1, 3, 9, 19, 35 };
692 static int p_hgk[5] = { 1, 5, 15, 29, 51 };
693 if (growth==GROWTH_SLOWLIN||
694 growth==GROWTH_SLOWLINODD||
695 growth==GROWTH_MODLIN) {
696 std::cout <<
"Specified Growth Rule Not Allowed!\n";
699 else if (growth==GROWTH_SLOWEXP) {
700 int l = 0, p = p_hgk[l], o = o_hgk[l];
701 while (p<2*level+1 && l<4) {
708 else if (growth==GROWTH_MODEXP||growth==GROWTH_DEFAULT) {
709 int l = 0, p = p_hgk[l], o = o_hgk[l];
710 while (p<4*level+1 && l<4) {
717 else if (growth==GROWTH_FULLEXP) {
718 int l = level; l = std::max(l,0); l = std::min(l,4);
723 else if (rule==BURK_TRAPEZOIDAL) {
724 if (growth==GROWTH_SLOWLIN) {
727 else if (growth==GROWTH_SLOWLINODD) {
728 return 2*((level+1)/2)+1;
730 else if (growth==GROWTH_MODLIN) {
733 else if (growth==GROWTH_SLOWEXP) {
745 else if (growth==GROWTH_MODEXP||growth==GROWTH_DEFAULT) {
751 while (o<4*level+1) {
757 else if (growth==GROWTH_FULLEXP) {
762 return (
int)pow(2.0,(
double)level)+1;
771 #if defined(Intrepid_SHOW_DEPRECATED_WARNINGS)
773 #warning "The Intrepid package is deprecated"
CubatureLineSorted(int degree=0, EIntrepidBurkardt rule=BURK_CLENSHAWCURTIS, bool isNormalized=false)
Constructor.
Utilizes cubature (integration) rules contained in the library sandia_rules (John Burkardt...
std::map< Scalar, int >::iterator begin(void)
Initiate iterator at the beginning of data.
void update(Scalar alpha2, CubatureLineSorted< Scalar > &cubRule2, Scalar alpha1)
Replace CubatureLineSorted values with "this = alpha1*this+alpha2*cubRule2".
void getAccuracy(std::vector< int > &accuracy) const
Returns max. degree of polynomials that are integrated exactly. The return vector has size 1...
void getCubature(ArrayPoint &cubPoints, ArrayWeight &cubWeights) const
Returns cubature points and weights (return arrays must be pre-sized/pre-allocated).
int getNumPoints() const
Returns the number of cubature points.
int getDimension() const
Returns dimension of domain of integration.
Scalar getNode(typename std::map< Scalar, int >::iterator it)
Get a specific node described by the iterator location.
const char * getName() const
Returns cubature name.
Scalar getWeight(int weight)
Get a specific weight described by the integer location.
std::map< Scalar, int >::iterator end(void)
Initiate iterator at the end of data.