ROL
ROL_QuadratureTPConstructors.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 dimension, std::vector<int> numPoints1D, std::vector<EROLBurkardt> rule1D,
54  bool isNormalized)
55  : dimension_(dimension) {
56  /*
57  This constructor builds a tensor product rule according to quadInfo myRule.
58  */
59  TEUCHOS_TEST_FOR_EXCEPTION((dimension!=(int)numPoints1D.size()||
60  dimension!=(int)rule1D.size()),std::out_of_range,
61  ">>> ERROR (Quadrature): Dimension mismatch for inputs.");
62  accuracy_.clear();
63  std::vector<int> degree(1,0);
64  Quadrature<Real> newRule(0,1);
65  for (int i=0; i<dimension; i++) {
66  // Compute 1D rules
67  Quadrature<Real> rule1(rule1D[i],numPoints1D[i],isNormalized);
68  rule1.getAccuracy(degree);
69  accuracy_.push_back(degree[0]);
70  // Build Tensor Rule
71  newRule = kron_prod<Real>(newRule,rule1);
72  }
73  typename std::map<std::vector<Real>,int>::iterator it;
74  int loc = 0;
75  std::vector<Real> node(dimension,0.0);
76  for (it=newRule.begin(); it!=newRule.end(); it++) {
77  node = it->first;
78  addPointAndWeight(node,newRule.getWeight(node),loc);
79  loc++;
80  }
81  if (isNormalized) {
82  normalize();
83  }
84 }
85 
86 template <class Real>
87 Quadrature<Real>::Quadrature(int dimension, std::vector<int> numPoints1D, std::vector<EROLBurkardt> rule1D,
88  std::vector<EROLGrowth> growth1D, bool isNormalized)
89  : dimension_(dimension) {
90  /*
91  This constructor builds a tensor product rule according to quadInfo myRule.
92  */
93  TEUCHOS_TEST_FOR_EXCEPTION((dimension!=(int)numPoints1D.size()||
94  dimension!=(int)rule1D.size()||
95  dimension!=(int)growth1D.size()),std::out_of_range,
96  ">>> ERROR (Quadrature): Dimension mismatch for inputs.");
97  accuracy_.clear();
98  accuracy_.resize(dimension);
99  std::vector<int> degree(1);
100  Quadrature<Real> newRule(0,1);
101  for (int i=0; i<dimension; i++) {
102  // Compute 1D rules
103  int numPoints = growthRule1D(numPoints1D[i],growth1D[i],rule1D[i]);
104  Quadrature<Real> rule1(rule1D[i],numPoints,isNormalized);
105  rule1.getAccuracy(degree);
106  accuracy_.push_back(degree[0]);
107  // Build Tensor Rule
108  newRule = kron_prod<Real>(newRule,rule1);
109  }
110  typename std::map<std::vector<Real>,int>::iterator it;
111  int loc = 0;
112  std::vector<Real> node;
113  for (it=newRule.begin(); it!=newRule.end(); it++) {
114  node = it->first;
115  addPointAndWeight(node,newRule.getWeight(node),loc);
116  loc++;
117  }
118  if (isNormalized) {
119  normalize();
120  }
121 }
122 
123 template <class Real>
124 Quadrature<Real>::Quadrature(int dimension, int maxNumPoints, std::vector<EROLBurkardt> rule1D,
125  std::vector<EROLGrowth> growth1D, bool isNormalized)
126  : dimension_(dimension) {
127  /*
128  This constructor builds a tensor product rule according to quadInfo myRule.
129  */
130  TEUCHOS_TEST_FOR_EXCEPTION((dimension!=(int)rule1D.size()||
131  dimension!=(int)growth1D.size()),std::out_of_range,
132  ">>> ERROR (Quadrature): Dimension mismatch for inputs.");
133  accuracy_.clear();
134  accuracy_.resize(dimension);
135  std::vector<int> degree(1);
136  Quadrature<Real> newRule(0,1);
137  for (int i=0; i<dimension; i++) {
138  // Compute 1D rules
139  int numPoints = growthRule1D(maxNumPoints,growth1D[i],rule1D[i]);
140  Quadrature<Real> rule1(rule1D[i],numPoints,isNormalized);
141  rule1.getAccuracy(degree);
142  accuracy_.push_back(degree[0]);
143  // Build Tensor Rule
144  newRule = kron_prod<Real>(newRule,rule1);
145  }
146  typename std::map<std::vector<Real>,int>::iterator it;
147  int loc = 0;
148  std::vector<Real> node;
149  for (it=newRule.begin(); it!=newRule.end(); it++) {
150  node = it->first;
151  addPointAndWeight(node,newRule.getWeight(node),loc);
152  loc++;
153  }
154  if (isNormalized) {
155  normalize();
156  }
157 }
158 
159 typedef void ( *GWPointer ) ( int order, int np, double p[], double w[] );
160 typedef int ( SandiaRules::*GWPointer2 ) ( int level, int growth );
161 typedef void ( SandiaRules2::*GWPointer1 ) ( int order, int dim, double w[] );
162 
163 template <class Real>
164 Quadrature<Real>::Quadrature(int dim, int maxLevel, std::vector<EROLBurkardt> rule1D,
165  std::vector<EROLGrowth> growth1D, bool isNormalized,
166  bool useSandia) : dimension_(dim) {
167  std::vector<Real> sparse_point;
168  std::vector<Real> sparse_weight;
169  int point_num = 0;
170 
171  if (!useSandia) {
172  SandiaRules SR; sgmga webbur;
173 
174  Real tol = std::sqrt(SR.r8_epsilon());
175 
176  GWPointer *gw_compute_points;
177  gw_compute_points = new GWPointer[dim];
178  GWPointer *gw_compute_weights;
179  gw_compute_weights = new GWPointer[dim];
180 
181  std::vector<int> growth(dim,0);
182  std::vector<int> rule(dim,0);
183  for( int i = 0; i < dim; i++ ) {
184  switch( rule1D[i] ) {
185  case BURK_CLENSHAWCURTIS: rule[i] = 1; break;
186  case BURK_FEJER2: rule[i] = 2; break;
187  case BURK_PATTERSON: rule[i] = 3; break;
188  case BURK_LEGENDRE: rule[i] = 4; break;
189  case BURK_HERMITE: rule[i] = 5; break;
190  case BURK_LAGUERRE: rule[i] = 7; break;
191  case BURK_GENZKEISTER: rule[i] = 10; break;
192  default: rule[i] = 1;
193  }
194  switch( growth1D[i] ) {
195  case GROWTH_DEFAULT: growth[i] = 0; break;
196  case GROWTH_SLOWLIN: growth[i] = 1; break;
197  case GROWTH_SLOWLINODD: growth[i] = 2; break;
198  case GROWTH_MODLIN: growth[i] = 3; break;
199  case GROWTH_SLOWEXP: growth[i] = 4; break;
200  case GROWTH_MODEXP: growth[i] = 5; break;
201  case GROWTH_FULLEXP: growth[i] = 6; break;
202  default: growth[i] = 0;
203  }
204  }
205 
206  std::vector<Real> level_weight(dim,1.0);
207  std::vector<int> np(dim,0);
208  int np_sum = SR.i4vec_sum(dim,&np[0]);
209  std::vector<Real> p(np_sum,0.0);
210 
211  // Compute necessary data.
212  int point_total_num = webbur.sgmga_size_total(dim,&level_weight[0],maxLevel,&rule[0],&growth[0]);
213 
214  point_num = webbur.sgmga_size(dim,&level_weight[0],maxLevel,&rule[0],&np[0],&p[0],
215  gw_compute_points,tol,&growth[0]);
216 
217  std::vector<int> sparse_unique_index( point_total_num );
218  webbur.sgmga_unique_index(dim,&level_weight[0],maxLevel,&rule[0],&np[0],&p[0],gw_compute_points,
219  tol,point_num,point_total_num,&growth[0],&sparse_unique_index[0]);
220 
221  std::vector<int> sparse_order(dim*point_num,0);
222  std::vector<int> sparse_index(dim*point_num,0);
223  webbur.sgmga_index(dim,&level_weight[0],maxLevel,&rule[0],point_num,point_total_num,
224  &sparse_unique_index[0],&growth[0],&sparse_order[0],&sparse_index[0]);
225 
226  // Compute points and weights.
227  sparse_point.resize(dim*point_num,0.0);
228  webbur.sgmga_point(dim,&level_weight[0],maxLevel,&rule[0],&np[0],&p[0],gw_compute_points,point_num,
229  &sparse_order[0],&sparse_index[0],&growth[0],&sparse_point[0]);
230 
231  sparse_weight.resize(point_num,0.0);
232  webbur.sgmga_weight(dim,&level_weight[0],maxLevel,&rule[0],&np[0],&p[0],gw_compute_weights,
233  point_num,point_total_num,&sparse_unique_index[0],&growth[0],&sparse_weight[0]);
234 
235  delete [] gw_compute_points;
236  delete [] gw_compute_weights;
237 
238  /***********************************************************************/
239  /******** END BUILD SPARSE GRID USING SGMGA ****************************/
240  /***********************************************************************/
241  }
242  else {
243  SandiaRules2 SR; SandiaSGMGA webbur;
244 
245  Real tol = std::sqrt(SR.r8_epsilon());
246 
247  GWPointer1 *gw_compute_points;
248  gw_compute_points = new GWPointer1[dim];
249  GWPointer1 *gw_compute_weights;
250  gw_compute_weights = new GWPointer1[dim];
251  GWPointer2 *gw_compute_order;
252  gw_compute_order = new GWPointer2[dim];
253 
254  int growth = 0;
255 
256  for (int i=0; i<dim; i++) {
257  if(rule1D[i] == BURK_CLENSHAWCURTIS) {
258  gw_compute_points[i] = &SandiaRules2::clenshaw_curtis_points;
259  gw_compute_weights[i] = &SandiaRules2::clenshaw_curtis_weights;
260  gw_compute_order[i] = &SandiaRules::level_to_order_exp_cc;
261  }
262  else if (rule1D[i] == BURK_FEJER2) {
263  gw_compute_points[i] = &SandiaRules2::fejer2_points;
264  gw_compute_weights[i] = &SandiaRules2::fejer2_weights;
265  gw_compute_order[i] = &SandiaRules::level_to_order_exp_f2;
266  }
267  else if (rule1D[i] == BURK_PATTERSON) {
268  gw_compute_points[i] = &SandiaRules2::patterson_points;
269  gw_compute_weights[i] = &SandiaRules2::patterson_weights;
270  gw_compute_order[i] = &SandiaRules::level_to_order_exp_gp;
271  }
272  else if (rule1D[i] == BURK_LEGENDRE) {
273  gw_compute_points[i] = &SandiaRules2::legendre_points;
274  gw_compute_weights[i] = &SandiaRules2::legendre_weights;
275  gw_compute_order[i] = &SandiaRules::level_to_order_exp_gauss;
276  }
277  else if (rule1D[i] == BURK_HERMITE) {
278  gw_compute_points[i] = &SandiaRules2::hermite_points;
279  gw_compute_weights[i] = &SandiaRules2::hermite_weights;
280  gw_compute_order[i] = &SandiaRules::level_to_order_exp_gauss;
281  }
282  else if (rule1D[i] == BURK_LAGUERRE) {
283  gw_compute_points[i] = &SandiaRules2::laguerre_points;
284  gw_compute_weights[i] = &SandiaRules2::laguerre_weights;
285  gw_compute_order[i] = &SandiaRules::level_to_order_exp_gauss;
286  }
287  else if (rule1D[i] == BURK_GENZKEISTER) {
288  gw_compute_points[i] = &SandiaRules2::hermite_genz_keister_points;
289  gw_compute_weights[i] = &SandiaRules2::hermite_genz_keister_weights;
290  gw_compute_order[i] = &SandiaRules::level_to_order_exp_hgk;
291  }
292 
293  if( growth1D[i] == GROWTH_DEFAULT ||
294  growth1D[i] == GROWTH_FULLEXP ) {
295  growth = 2;
296  }
297  else if ( growth1D[i] == GROWTH_SLOWLIN ||
298  growth1D[i] == GROWTH_SLOWLINODD ||
299  growth1D[i] == GROWTH_SLOWEXP ) {
300  growth = 0;
301  }
302  else if ( growth1D[i] == GROWTH_MODLIN ||
303  growth1D[i] == GROWTH_MODEXP ) {
304  growth = 1;
305  }
306  }
307 
308  std::vector<Real> level_weight(dim,1.0);
309  std::vector<int> np(dim,0);
310  int np_sum = SR.i4vec_sum(dim,&np[0]);
311  std::vector<Real> p(np_sum,0.0);
312 
313  // Compute necessary data.
314  int point_total_num = webbur.sgmga_size_total(dim,&level_weight[0],maxLevel,growth,
315  gw_compute_order);
316 
317  point_num = webbur.sgmga_size(dim,&level_weight[0],maxLevel,gw_compute_points,tol,growth,
318  gw_compute_order);
319 
320  std::vector<int> sparse_unique_index(point_total_num);
321  webbur.sgmga_unique_index(dim,&level_weight[0],maxLevel,gw_compute_points,tol,point_num,
322  point_total_num,growth,gw_compute_order,&sparse_unique_index[0]);
323 
324  std::vector<int> sparse_order(dim*point_num,0);
325  std::vector<int> sparse_index(dim*point_num,0);
326  webbur.sgmga_index(dim,&level_weight[0],maxLevel,point_num,point_total_num,&sparse_unique_index[0],
327  growth,gw_compute_order,&sparse_order[0],&sparse_index[0]);
328 
329  // Compute points and weights.
330  sparse_point.resize(dim*point_num,0.0);
331  webbur.sgmga_point(dim,&level_weight[0],maxLevel,gw_compute_points,point_num,&sparse_order[0],
332  &sparse_index[0],growth,gw_compute_order,&sparse_point[0]);
333 
334  sparse_weight.resize(point_num,0.0);
335  webbur.sgmga_weight(dim,&level_weight[0],maxLevel,gw_compute_weights,point_num,point_total_num,
336  &sparse_unique_index[0],growth,gw_compute_order,&sparse_weight[0]);
337 
338  delete [] gw_compute_points;
339  delete [] gw_compute_weights;
340  delete [] gw_compute_order;
341  /***********************************************************************/
342  /******** END BUILD SPARSE GRID USING SANDIA_SGMGA *********************/
343  /***********************************************************************/
344  }
345 
346  typename std::map<std::vector<Real>,int>::iterator it;
347  Real weight = 0.0;
348  std::vector<Real> node(dim,0.0);;
349  for (int i = 0; i < point_num; i++) {
350  weight = sparse_weight[i];
351  for (int j = 0; j < dim; j++) {
352  node[j] = sparse_point[j+i*dim];
353  }
354  addPointAndWeight(node,weight,i);
355  }
356  if (isNormalized) {
357  normalize();
358  }
359 }
360 
361 template<class Real>
362 Quadrature<Real>::Quadrature(const char* SGinfo, const char* SGdata, bool isNormalized) {
363  // Open SGdata
364  std::fstream info; info.open(SGinfo);
365  Real buf;
366 
367  // Get Dimension Number
368  info >> buf;
369  int dim_num = (int)buf;
370  dimension_ = dim_num;
371 
372  // Get Number of Cubature Points/Weights
373  info >> buf;
374  int point_num = (int)buf;
375 
376  // Close SGinfo
377  info.close();
378 
379  // Open SGdata File
380  std::fstream file; file.open(SGdata);
381 
382  // Get Cubature Points/Weights
383  Real weight = 0.0;
384  std::vector<Real> node(dim_num, 0.0);
385  int loc = 0;
386  for ( int i = 0; i < point_num; i++ ) {
387  // Get Cubature Point
388  for ( int j = 0; j < dim_num; j++ ) {
389  file >> buf;
390  node[j] = buf;
391  }
392  // Get Cubature Weight
393  file >> buf;
394  weight = buf;
395  // Insert Point/Weight into Sparse Grid Storage
396  addPointAndWeight(node,weight,loc);
397  loc++;
398  // Move to Next Line in SGdata
399  file.ignore(256,'\n');
400  }
401  file.close();
402  if (isNormalized) {
403  normalize();
404  }
405 }
406 } // end ROL namespace
void sgmga_index(int dim_num, double level_weight[], int level_max, int rule[], int point_num, int point_total_num, int sparse_unique_index[], int growth[], int sparse_order[], int sparse_index[])
void fejer2_points(int n, int dim, double x[])
void sgmga_point(int dim_num, double level_weight[], int level_max, int rule[], int np[], double p[], void(*gw_compute_points[])(int order, int np, double p[], double x[]), int point_num, int sparse_order[], int sparse_index[], int growth[], double sparse_point[])
void sgmga_weight(int dim_num, double level_weight[], int level_max, void(SandiaRules2::*gw_compute_weights[])(int order, int dim, double w[]), int point_num, int point_total_num, int sparse_unique_index[], int growth, int(SandiaRules::*gw_compute_order[])(int level, int growth), double sparse_weight[])
void patterson_points(int n, int dim, double x[])
int sgmga_size_total(int dim_num, double level_weight[], int level_max, int growth, int(SandiaRules::*gw_compute_order[])(int level, int growth))
void laguerre_weights(int n, int dim, double w[])
void hermite_points(int n, int dim, double x[])
virtual void getAccuracy(std::vector< int > &accuracy) const
int sgmga_size_total(int dim_num, double level_weight[], int level_max, int rule[], int growth[])
void laguerre_points(int n, int dim, double x[])
void sgmga_point(int dim_num, double level_weight[], int level_max, void(SandiaRules2::*gw_compute_points[])(int order, int dim, double x[]), int point_num, int sparse_order[], int sparse_index[], int growth, int(SandiaRules::*gw_compute_order[])(int level, int growth), double sparse_point[])
virtual void normalize()
void clenshaw_curtis_weights(int n, int dim, double w[])
int level_to_order_exp_hgk(int level, int growth)
void clenshaw_curtis_points(int n, int dim, double x[])
void sgmga_weight(int dim_num, double level_weight[], int level_max, int rule[], int np[], double p[], void(*gw_compute_weights[])(int order, int np, double p[], double w[]), int point_num, int point_total_num, int sparse_unique_index[], int growth[], double sparse_weight[])
void sgmga_unique_index(int dim_num, double level_weight[], int level_max, int rule[], int np[], double p[], void(*gw_compute_points[])(int order, int np, double p[], double x[]), double tol, int point_num, int point_total_num, int growth[], int sparse_unique_index[])
void hermite_weights(int n, int dim, double w[])
int sgmga_size(int dim_num, double level_weight[], int level_max, void(SandiaRules2::*gw_compute_points[])(int order, int dim, double x[]), double tol, int growth, int(SandiaRules::*gw_compute_order[])(int level, int growth))
virtual std::map< std::vector< Real >, int >::iterator end()
void addPointAndWeight(std::vector< Real > point, Real weight, int loc)
void(* GWPointer)(int order, int np, double p[], double w[])
Quadrature(int dimension=1)
int i4vec_sum(int n, int a[])
int level_to_order_exp_f2(int level, int growth)
void sgmga_unique_index(int dim_num, double level_weight[], int level_max, void(SandiaRules2::*gw_compute_points[])(int order, int dim, double x[]), double tol, int point_num, int point_total_num, int growth, int(SandiaRules::*gw_compute_order[])(int level, int growth), int sparse_unique_index[])
virtual Real getWeight(int node)
int(SandiaRules::* GWPointer2)(int level, int growth)
std::vector< int > accuracy_
void patterson_weights(int n, int dim, double w[])
void(SandiaRules2::* GWPointer1)(int order, int dim, double w[])
int level_to_order_exp_gauss(int level, int growth)
void legendre_points(int n, int dim, double x[])
void hermite_genz_keister_points(int n, int dim, double x[])
virtual std::map< std::vector< Real >, int >::iterator begin()
int level_to_order_exp_gp(int level, int growth)
void hermite_genz_keister_weights(int n, int dim, double w[])
void sgmga_index(int dim_num, double level_weight[], int level_max, int point_num, int point_total_num, int sparse_unique_index[], int growth, int(SandiaRules::*gw_compute_order[])(int level, int growth), int sparse_order[], int sparse_index[])
int sgmga_size(int dim_num, double level_weight[], int level_max, int rule[], int np[], double p[], void(*gw_compute_points[])(int order, int np, double p[], double x[]), double tol, int growth[])
int growthRule1D(int index, EROLGrowth growth, EROLBurkardt rule)
int level_to_order_exp_cc(int level, int growth)
void legendre_weights(int n, int dim, double w[])
void fejer2_weights(int n, int dim, double w[])