ROL
ROL_SparseGridGeneratorDef.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 #ifndef ROL_SPARSEGRIDGENERATOR_DEF_HPP
45 #define ROL_SPARSEGRIDGENERATOR_DEF_HPP
46 
47 namespace ROL {
48 
49 template <class Real>
51  SparseGridInfo &info, bool adaptive )
52  : SampleGenerator<Real>(bman), adaptive_(adaptive), info_(info) {
53  if ( adaptive ) {
54  //index_.resize(info.dim,1);
55  index_.clear();
56  direction_ = 0;
57  error_ = 0.0;
58  activeIndex_.clear();
59  oldIndex_.clear();
60  grid_ = Teuchos::rcp(new Quadrature<Real>(info.dim));
61  }
62  else {
63  int dim_num = info.dim;
64  int level_max = info.maxLevel - 1;
65  std::vector<EROLBurkardt> rule1D = info.rule1D;
66  std::vector<EROLGrowth> growth1D = info.growth1D;
67  bool isNormalized = info.normalize;
68  bool useSandia = info.useSandia;
69  grid_ = Teuchos::rcp(new Quadrature<Real>(dim_num,level_max,rule1D,growth1D,isNormalized,useSandia));
70  // Split points and weights across processors
71  this->setSamples(true);
72  }
73 }
74 
75 template <class Real>
77  const char* SGinfo,const char* SGdata,
78  bool isNormalized) : SampleGenerator<Real>(bman),
79  adaptive_(false) {
80  grid_ = Teuchos::rcp(new Quadrature<Real>(SGinfo,SGdata,isNormalized));
81  // Split points and weights across processors
82  this->setSamples(true);
83 }
84 
85 template<class Real>
86 bool SparseGridGenerator<Real>::checkMaxLevel(std::vector<int> &index) {
87  int level = 0;
88  bool useMax = true;
89  for ( unsigned l = 0; l < index.size(); l++ ) {
90  if ( useMax ) {
91  level = std::max(level,index[l]);
92  }
93  else {
94  level += index[l];
95  }
96  }
97  if ( level >= (this->info_).maxLevel ) {
98  return true;
99  }
100  return false;
101 }
102 
103 template<class Real>
104 bool SparseGridGenerator<Real>::isAdmissible(std::vector<int> &index, int direction) {
105  for ( int i = 0; i < (this->info_).dim; i++ ) {
106  if ( index[i] > 1 && i != direction ) {
107  index[i]--;
108  if ( !((this->oldIndex_).count(index)) ) {
109  return false;
110  }
111  index[i]++;
112  }
113  }
114  return true;
115 }
116 
117 template<class Real>
118 void SparseGridGenerator<Real>::buildDiffRule(Quadrature<Real> &outRule, std::vector<int> &index) {
119  int numPoints = 0;
120  for ( int i = 0; i < (this->info_).dim; i++ ) {
121  numPoints = growthRule1D(index[i],((this->info_).growth1D)[i],((this->info_).rule1D)[i]);
122  Quadrature<Real> diffRule(((this->info_).rule1D)[i],numPoints,((this->info_).normalize));
123  if ( numPoints != 1 ) {
124  numPoints = growthRule1D(index[i]-1,((this->info_).growth1D)[i],((this->info_).rule1D)[i]);
125  Quadrature<Real> rule(((this->info_).rule1D)[i],numPoints,((this->info_).normalize));
126  diffRule.update(-1.0,rule);
127  }
128  outRule = kron_prod<Real>(outRule,diffRule);
129  }
130 }
131 
132 template<class Real>
135  if ( this->adaptive_ ) {
136  //(this->index_).resize((this->info_).dim,1);
137  (this->index_).clear();
138  this->direction_ = 0;
139  this->error_ = 0.0;
140  (this->activeIndex_).clear();
141  (this->oldIndex_).clear();
142  this->grid_ = Teuchos::rcp(new Quadrature<Real>((this->info_).dim));
143  }
144 }
145 
146 template<class Real>
147 Real SparseGridGenerator<Real>::computeError(std::vector<Real> &vals){
148  if ( this->adaptive_ ) {
149  Real myerror = 0.0;
150  for ( int i = 0; i < SampleGenerator<Real>::numMySamples(); i++ ) {
151  myerror += vals[i]*SampleGenerator<Real>::getMyWeight(i);
152  }
153  Real error = 0.0;
154  SampleGenerator<Real>::sumAll(&myerror,&error,1);
155  error = std::abs(error);
156  // Update global error and index sets
157  this->error_ += error;
158  if ( (this->activeIndex_).end() != (this->activeIndex_).begin() ) {
159  (this->activeIndex_).insert((this->activeIndex_).end()--,
160  std::pair<Real,std::vector<int> >(error,this->search_index_));
161  }
162  else {
163  (this->activeIndex_).insert(std::pair<Real,std::vector<int> >(error,this->search_index_));
164  }
165  // Return
166  return this->error_;
167  }
168  else {
169  return 0.0;
170  }
171 }
172 
173 template<class Real>
174 Real SparseGridGenerator<Real>::computeError(std::vector<Teuchos::RCP<Vector<Real> > > &vals, const Vector<Real> &x ){
175  if ( this->adaptive_ ) {
176  Teuchos::RCP<Vector<Real> > mydiff = x.clone();
177  mydiff->zero();
178  for ( int i = 0; i < SampleGenerator<Real>::numMySamples(); i++ ) {
179  mydiff->axpy(SampleGenerator<Real>::getMyWeight(i),(*vals[i]));
180  }
181  Teuchos::RCP<Vector<Real> > diff = mydiff->clone();
182  SampleGenerator<Real>::sumAll(*mydiff,*diff);
183  Real error = diff->norm();
184  // Update global error and index sets
185  this->error_ += error;
186  if ( (this->activeIndex_).end() != (this->activeIndex_).begin() ) {
187  (this->activeIndex_).insert((this->activeIndex_).end()--,
188  std::pair<Real,std::vector<int> >(error,this->search_index_));
189  }
190  else {
191  (this->activeIndex_).insert(std::pair<Real,std::vector<int> >(error,this->search_index_));
192  }
193  // Return
194  return this->error_;
195  }
196  else {
197  return 0.0;
198  }
199 }
200 
201 template<class Real>
203  if ( this->adaptive_ ) {
204  // Select index to investigate
205  if ( !((this->index_).empty()) && ((this->direction_) < ((this->info_).dim-1)) ) {
206  // Select index corresponding to next direction
207  this->search_index_ = (this->index_);
208  this->direction_++;
209  (this->search_index_)[this->direction_]++;
210  }
211  else {
212  // Select index corresponding to largest error
213  if ( (this->index_).empty() ) {
214  (this->index_).resize((this->info_).dim,1);
215  this->search_index_ = this->index_;
216  this->direction_ = (this->info_).dim;
217  }
218  else {
219  typename std::multimap<Real,std::vector<int> >::iterator it = (this->activeIndex_).end();
220  it--;
221  this->error_ -= it->first;
222  this->index_ = it->second;
223  (this->activeIndex_).erase(it);
224  (this->oldIndex_).insert((this->oldIndex_).end(),(this->index_));
225  this->search_index_ = (this->index_);
226  this->direction_ = 0;
227  (this->search_index_)[this->direction_]++;
228  }
229  }
230  // Check to see if maxLevel is violated
231  Quadrature<Real> rule((this->info_).dim);
232  if ( !(this->checkMaxLevel(this->search_index_)) ) {
233  // Check admissibility of index
234  if ( this->isAdmissible((this->search_index_),this->direction_) ) {
235  // Build difference rule
236  this->buildDiffRule(rule,this->search_index_);
237  }
238  }
239  // Set values of difference rule as points and weights
240  this->updateSamples(rule);
241  }
242 }
243 
244 template<class Real>
245 void SparseGridGenerator<Real>::splitSamples(std::vector<std::vector<Real> > &mypts, std::vector<Real> &mywts) {
246  // Get global points and weights
247  std::vector<std::vector<Real> > pts;
248  std::vector<Real> wts;
249  this->grid_->getCubature(pts,wts);
250  // Split global points and weights across processors
251  int frac = pts.size()/SampleGenerator<Real>::numBatches();
252  int rem = pts.size()%SampleGenerator<Real>::numBatches();
253  int npts = frac;
254  if ( SampleGenerator<Real>::batchID() < rem ) {
255  npts++;
256  }
257  mypts.resize(npts);
258  mywts.resize(npts);
259  int index = 0;
260  for ( int i = 0; i < npts; i++ ) {
262  mywts[i] = wts[index];
263  mypts[i] = pts[index];
264  }
265 }
266 
267 template<class Real>
268 void SparseGridGenerator<Real>::setSamples(bool inConstructor) {
269  if ( !(SampleGenerator<Real>::batchID()) && !inConstructor) {
270  std::cout << "Number of Global Points: " << (this->grid_)->getNumPoints() << "\n";
271  }
272  if ( this->adaptive_ || inConstructor ) {
273  // Split samples based on PID
274  std::vector<std::vector<Real> > mypts;
275  std::vector<Real> mywts;
276  this->splitSamples(mypts,mywts);
277  // Set local points and weights
280  }
281 }
282 
283 template<class Real>
285  // Add increment to stored sparse grid
286  this->grid_->update(1.0,grid);
287  // Split global stored points and weights across processors
288  std::vector<std::vector<Real> > mypts;
289  std::vector<Real> mywts;
290  this->splitSamples(mypts,mywts);
291  // Get global incremental points and weights
292  std::vector<std::vector<Real> > pts_inc;
293  std::vector<Real> wts_inc;
294  grid.getCubature(pts_inc,wts_inc);
295  // STL set intesection of global incremental points and local stored points
296  std::vector<std::vector<Real> > pts_inter;
297  typename std::vector<std::vector<Real> >::iterator mypts_it;
298  typename std::vector<std::vector<Real> >::iterator pts_inc_it;
299  for ( pts_inc_it=pts_inc.begin(); pts_inc_it!=pts_inc.end(); pts_inc_it++ ) {
300  for ( mypts_it=mypts.begin(); mypts_it!=mypts.end(); mypts_it++ ) {
301  if ( *mypts_it == *pts_inc_it ) {
302  pts_inter.push_back(*mypts_it);
303  break;
304  }
305  }
306  }
307  // Get intersection weights
308  std::vector<Real> wts_inter;
309  for ( unsigned i = 0; i < pts_inter.size(); i++ ) {
310  wts_inter.push_back(grid.getWeight(pts_inter[i]));
311  }
312  // Set local points and weights
315 }
316 
317 } // namespace ROL
318 
319 #endif
virtual void update(const Vector< Real > &x)
void splitSamples(std::vector< std::vector< Real > > &mypts, std::vector< Real > &mywts)
virtual void getCubature(std::vector< std::vector< Real > > &cubPoints, std::vector< Real > &cubWeights) const
virtual Teuchos::RCP< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual Real getMyWeight(const int i) const
bool checkMaxLevel(std::vector< int > &index)
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:72
void buildDiffRule(Quadrature< Real > &outRule, std::vector< int > &index)
void sumAll(Real *input, Real *output, int dim) const
virtual void update(Real alpha, Quadrature< Real > &rule)
bool isAdmissible(std::vector< int > &index, int direction)
std::multimap< Real, std::vector< int > > activeIndex_
std::vector< EROLBurkardt > rule1D
Real computeError(std::vector< Real > &vals)
std::set< std::vector< int > > oldIndex_
SparseGridGenerator(Teuchos::RCP< BatchManager< Real > > &bman, SparseGridInfo &info, bool adaptive=false)
void updateSamples(Quadrature< Real > &grid)
virtual Real getWeight(int node)
Teuchos::RCP< Quadrature< Real > > grid_
void update(const Vector< Real > &x)
int growthRule1D(int index, EROLGrowth growth, EROLBurkardt rule)
void setPoints(std::vector< std::vector< Real > > &p)
void setSamples(bool inConstructor=false)
void setWeights(std::vector< Real > &w)
std::vector< EROLGrowth > growth1D