ROL
ROL_Quadrature1dConstructors.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Rapid Optimization Library (ROL) Package
5 // Copyright (2014) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact lead developers:
38 // Drew Kouri (dpkouri@sandia.gov) and
39 // Denis Ridzal (dridzal@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
44 
50 namespace ROL {
51 
52 template <class Real>
53 Quadrature<Real>::Quadrature( int degree, EROLBurkardt rule, bool isNormalized )
54  : dimension_(1) {
55  TEUCHOS_TEST_FOR_EXCEPTION((degree < 0),std::out_of_range,
56  ">>> ERROR (Quadrature): No rule implemented for desired polynomial degree.");
57  accuracy_.clear();
58  accuracy_.push_back(degree);
59 
60  int numPoints = (degree+1)/2+1;
61  if (rule==BURK_CLENSHAWCURTIS||rule==BURK_FEJER2) {
62  numPoints = degree+1;
63  }
64  else if (rule==BURK_TRAPEZOIDAL) {
65  numPoints = 2;
66  }
67  else if (rule==BURK_PATTERSON) {
68  int l = 0, o = (degree-0.5)/1.5;
69  for (int i=0; i<8; i++) {
70  l = (int)pow(2.0,(double)i+1.0)-1;
71  if (l>=o) {
72  numPoints = l;
73  break;
74  }
75  }
76  }
77  else if (rule==BURK_GENZKEISTER) {
78  int o_ghk[8] = {1,3,9,19,35,37,41,43};
79  int o = (degree-0.5)/1.5;
80  for (int i=0; i<8; i++) {
81  if (o_ghk[i]>=o) {
82  numPoints = o_ghk[i];
83  break;
84  }
85  }
86  }
87 
88  std::vector<Real> points(numPoints), weights(numPoints);
89  switch(rule) {
90  case BURK_CHEBYSHEV1: // Gauss-Chebyshev Type 1
91  ROLBurkardtRules::chebyshev1_compute<Real>(numPoints,&points[0],&weights[0]); break;
92  case BURK_CHEBYSHEV2: // Gauss-Chebyshev Type 2
93  ROLBurkardtRules::chebyshev2_compute<Real>(numPoints,&points[0],&weights[0]); break;
94  case BURK_CLENSHAWCURTIS: // Clenshaw-Curtis
95  ROLBurkardtRules::clenshaw_curtis_compute<Real>(numPoints,&points[0],&weights[0]); break;
96  case BURK_FEJER2: // Fejer Type 2
97  ROLBurkardtRules::fejer2_compute<Real>(numPoints,&points[0],&weights[0]); break;
98  case BURK_LEGENDRE: // Gauss-Legendre
99  ROLBurkardtRules::legendre_compute<Real>(numPoints,&points[0],&weights[0]); break;
100  case BURK_PATTERSON: // Gauss-Patterson
101  ROLBurkardtRules::patterson_lookup<Real>(numPoints,&points[0],&weights[0]); break;
102  case BURK_TRAPEZOIDAL: // Trapezoidal Rule
103  ROLBurkardtRules::trapezoidal_compute<Real>(numPoints,&points[0],&weights[0]); break;
104  case BURK_HERMITE: // Gauss-Hermite
105  ROLBurkardtRules::hermite_compute<Real>(numPoints,&points[0],&weights[0]); break;
106  case BURK_GENZKEISTER: // Hermite-Genz-Keister
107  ROLBurkardtRules::hermite_genz_keister_lookup<Real>(numPoints,&points[0],&weights[0]); break;
108  case BURK_LAGUERRE: // Gauss-Laguerre
109  ROLBurkardtRules::laguerre_compute<Real>(numPoints,&points[0],&weights[0]); break;
110  default:
111  TEUCHOS_TEST_FOR_EXCEPTION(true,std::out_of_range,
112  ">>> ERROR (Quadrature): No such 1D quadrature rule."); break;
113  }
114 
115  if (isNormalized) {
116  Real sum = 0.0;
117  for (int i=0; i<numPoints; i++) {
118  sum += weights[i];
119  }
120  for (int i=0; i<numPoints; i++) {
121  weights[i] /= sum;
122  }
123  }
124  else {
125  if (rule==BURK_GENZKEISTER && numPoints >= 37) {
126  for (int i=0; i<numPoints; i++) {
127  weights[i] *= sqrt(M_PI);
128  }
129  }
130  }
131 
132  std::vector<Real> point(1,0.0);
133  for (int i=0; i<numPoints; i++) {
134  point[0] = points[i];
135  addPointAndWeight(point,weights[i],i);
136  }
137 } // end constructor
138 
139 template <class Real>
140 Quadrature<Real>::Quadrature( EROLBurkardt rule, int numPoints, bool isNormalized )
141  : dimension_(1) {
142  TEUCHOS_TEST_FOR_EXCEPTION((numPoints < 0),std::out_of_range,
143  ">>> ERROR (Quadrature): No rule implemented for desired number of points.");
144 
145  accuracy_.clear();
146  std::vector<Real> points(numPoints), weights(numPoints);
147  if (rule==BURK_CHEBYSHEV1) { // Gauss-Chebyshev Type 1
148  accuracy_.push_back(2*numPoints-1);
149  ROLBurkardtRules::chebyshev1_compute<Real>(numPoints,&points[0],&weights[0]);
150  }
151  else if (rule==BURK_CHEBYSHEV2) { // Gauss-Chebyshev Type 2
152  accuracy_.push_back(2*numPoints-1);
153  ROLBurkardtRules::chebyshev2_compute<Real>(numPoints,&points[0],&weights[0]);
154  }
155  else if (rule==BURK_CLENSHAWCURTIS) { // Clenshaw-Curtis
156  accuracy_.push_back(numPoints-1);
157  ROLBurkardtRules::clenshaw_curtis_compute<Real>(numPoints,&points[0],&weights[0]);
158  }
159  else if (rule==BURK_FEJER2) { // Fejer Type 2
160  accuracy_.push_back(numPoints-1);
161  ROLBurkardtRules::fejer2_compute<Real>(numPoints,&points[0],&weights[0]);
162  }
163  else if (rule==BURK_LEGENDRE) { // Gauss-Legendre
164  accuracy_.push_back(2*numPoints-1);
165  ROLBurkardtRules::legendre_compute<Real>(numPoints,&points[0],&weights[0]);
166  }
167  else if (rule==BURK_PATTERSON) { // Gauss-Patterson
168  bool correctNumPoints = false;
169  for (int i=0; i<8; i++) {
170  int l = (int)pow(2.0,(double)i+1.0)-1;
171  if (numPoints==l) {
172  correctNumPoints = true;
173  break;
174  }
175  }
176  TEUCHOS_TEST_FOR_EXCEPTION((correctNumPoints==false),std::out_of_range,
177  ">>> ERROR (Quadrature): Number of points must be numPoints = 1, 3, 7, 15, 31, 63, 127, 255.");
178  Real degree = 1.5*(double)numPoints+0.5;
179  accuracy_.push_back((int)degree);
180  ROLBurkardtRules::patterson_lookup<Real>(numPoints,&points[0],&weights[0]);
181  }
182  else if (rule==BURK_TRAPEZOIDAL) { // Trapezoidal Rule
183  accuracy_.push_back(2);
184  ROLBurkardtRules::trapezoidal_compute<Real>(numPoints,&points[0],&weights[0]);
185  }
186  else if (rule==BURK_HERMITE) { // Gauss-Hermite
187  accuracy_.push_back(2*numPoints-1);
188  ROLBurkardtRules::hermite_compute<Real>(numPoints,&points[0],&weights[0]);
189  }
190  else if (rule==BURK_GENZKEISTER) { // Hermite-Genz-Keister
191  bool correctNumPoints = false;
192  int o_ghk[8] = {1,3,9,19,35,37,41,43};
193  for (int i=0; i<8; i++) {
194  if (o_ghk[i]==numPoints) {
195  correctNumPoints = true;
196  break;
197  }
198  }
199  TEUCHOS_TEST_FOR_EXCEPTION((correctNumPoints==false),std::out_of_range,
200  ">>> ERROR (Quadrature): Number of points must be numPoints = 1, 3, 9, 35, 37, 41, 43.");
201  Real degree = 1.5*(double)numPoints+0.5;
202  accuracy_.push_back((int)degree);
203  ROLBurkardtRules::hermite_genz_keister_lookup<Real>(numPoints,&points[0],&weights[0]);
204  }
205  else if (rule==BURK_LAGUERRE) { // Gauss-Laguerre
206  accuracy_.push_back(2*numPoints-1);
207  ROLBurkardtRules::laguerre_compute<Real>(numPoints,&points[0],&weights[0]);
208  }
209 
210  if (isNormalized) {
211  Real sum = 0.0;
212  for (int i=0; i<numPoints; i++) {
213  sum += weights[i];
214  }
215  for (int i=0; i<numPoints; i++) {
216  weights[i] /= sum;
217  }
218  }
219  else {
220  if (rule==BURK_GENZKEISTER && numPoints >= 37) {
221  for (int i=0; i<numPoints; i++) {
222  weights[i] *= sqrt(M_PI);
223  }
224  }
225  }
226 
227  std::vector<Real> point(1,0.0);
228  for (int i=0; i<numPoints; i++) {
229  point[0] = points[i];
230  addPointAndWeight(point,weights[i],i);
231  }
232 } // end constructor
233 
234 template <class Real>
235 Quadrature<Real>::Quadrature( std::vector<Real>& points, std::vector<Real>& weights)
236  : dimension_(1) {
237  int size = (int)weights.size();
238  TEUCHOS_TEST_FOR_EXCEPTION(((int)points.size()!=size),std::out_of_range,
239  ">>> ERROR (Quadrature): Input dimension mismatch.");
240  std::vector<Real> point(1,0.0);
241  for (int i=0; i<size; i++) {
242  point[0] = points[i];
243  addPointAndWeight(point,weights[i],i);
244  }
245 }
246 
247 } // ROL namespace
const double weights[4][5]
Definition: ROL_Types.hpp:1001
EROLBurkardt
Enumeration of integration rules provided by John Burkardt.
void addPointAndWeight(std::vector< Real > point, Real weight, int loc)
Quadrature(int dimension=1)
std::vector< int > accuracy_