Intrepid
test_18.cpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Intrepid Package
5 // Copyright (2007) 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 Pavel Bochev (pbboche@sandia.gov)
38 // Denis Ridzal (dridzal@sandia.gov), or
39 // Kara Peterson (kjpeter@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
44 
52 //#include "Intrepid_CubatureLineSorted.hpp"
53 #include "Intrepid_Utils.hpp"
54 #include "Teuchos_oblackholestream.hpp"
55 #include "Teuchos_RCP.hpp"
56 #include "Teuchos_RefCountPtr.hpp"
57 #include "Teuchos_GlobalMPISession.hpp"
58 
59 using namespace Intrepid;
60 std::vector<long double> alpha(1,0);
61 std::vector<long double> beta(1,0);
62 
63 template<class Scalar>
64 class StdVector {
65 private:
66  Teuchos::RefCountPtr<std::vector<Scalar> > std_vec_;
67 
68 public:
69 
70  StdVector( const Teuchos::RefCountPtr<std::vector<Scalar> > & std_vec )
71  : std_vec_(std_vec) {}
72 
73  Teuchos::RefCountPtr<StdVector<Scalar> > Create() const {
74  return Teuchos::rcp( new StdVector<Scalar>(
75  Teuchos::rcp(new std::vector<Scalar>(std_vec_->size(),0))));
76  }
77 
78  void Update( StdVector<Scalar> & s ) {
79  int dimension = (int)(std_vec_->size());
80  for (int i=0; i<dimension; i++)
81  (*std_vec_)[i] += s[i];
82  }
83 
84  void Update( Scalar alpha, StdVector<Scalar> & s ) {
85  int dimension = (int)(std_vec_->size());
86  for (int i=0; i<dimension; i++)
87  (*std_vec_)[i] += alpha*s[i];
88  }
89 
90  Scalar operator[](int i) {
91  return (*std_vec_)[i];
92  }
93 
94  void clear() {
95  std_vec_->clear();
96  }
97 
98  void resize(int n, Scalar p) {
99  std_vec_->resize(n,p);
100  }
101 
102  int size() {
103  return (int)std_vec_->size();
104  }
105 
106  void Set( Scalar alpha ) {
107  int dimension = (int)(std_vec_->size());
108  for (int i=0; i<dimension; i++)
109  (*std_vec_)[i] = alpha;
110  }
111 };
112 
113 template<class Scalar,class UserVector>
114 class ASGdata :
115  public Intrepid::AdaptiveSparseGridInterface<Scalar,UserVector> {
116 public:
117  ~ASGdata() {}
118 
119  ASGdata(int dimension,std::vector<EIntrepidBurkardt> rule1D,
120  std::vector<EIntrepidGrowth> growth1D, int maxLevel,
121  bool isNormalized) : AdaptiveSparseGridInterface<Scalar,UserVector>(
122  dimension,rule1D,growth1D,maxLevel,isNormalized) {}
123 
124  void eval_integrand(UserVector & output, std::vector<Scalar> & input) {
125  int dimension = (int)alpha.size();
126  Scalar total = 1.0;
127  Scalar point = 0;
128  for (int i=0; i<dimension; i++) {
129  point = 0.5*input[i]+0.5;
130  total *= ( 1.0/powl(alpha[i],(long double)2.0)
131  + powl(point-beta[i],(long double)2.0) );
132  }
133  output.clear(); output.resize(1,1.0/total);
134  }
135 
136  Scalar error_indicator(UserVector & input) {
137  int dimension = (int)input.size();
138  Scalar norm2 = 0.0;
139  for (int i=0; i<dimension; i++)
140  norm2 += input[i]*input[i];
141 
144  norm2 = std::sqrt(norm2)/ID;
145  return norm2;
146  }
147 };
148 
149 long double compExactInt(void) {
150  double val = 1.0;
151  int dimension = alpha.size();
152  for (int i=0; i<dimension; i++) {
153  val *= alpha[i]*( std::atan((1.0-beta[i])*alpha[i])
154  +std::atan(beta[i]*alpha[i]) );
155  }
156  return val;
157 }
158 
159 long double adaptSG(StdVector<long double> & iv,
161  problem_data,long double TOL) {
162 
163  // Construct a Container for the adapted rule
164  int dimension = problem_data.getDimension();
165  std::vector<int> index(dimension,1);
166 
167  // Initialize global error indicator
168  long double eta = 1.0;
169 
170  // Initialize the Active index set
171  std::multimap<long double,std::vector<int> > activeIndex;
172  activeIndex.insert(std::pair<long double,std::vector<int> >(eta,index));
173 
174  // Initialize the old index set
175  std::set<std::vector<int> > oldIndex;
176  /*
177  std::vector<long double> output(1,0);
178  std::vector<long double> input(dimension,0.5);
179  problem_data.eval_integrand(output,input);
180  */
181  // Perform Adaptation
182  while (eta > TOL) {
184  activeIndex,oldIndex,
185  iv,eta,
186  problem_data);
187  }
188  return eta;
189 }
190 
191 int main(int argc, char *argv[]) {
192 
193  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
194 
195  // This little trick lets us print to std::cout only if
196  // a (dummy) command-line argument is provided.
197  int iprint = argc - 1;
198  Teuchos::RCP<std::ostream> outStream;
199  Teuchos::oblackholestream bhs; // outputs nothing
200  if (iprint > 0)
201  outStream = Teuchos::rcp(&std::cout, false);
202  else
203  outStream = Teuchos::rcp(&bhs, false);
204 
205  // Save the format state of the original std::cout.
206  Teuchos::oblackholestream oldFormatState;
207  oldFormatState.copyfmt(std::cout);
208 
209  *outStream \
210  << "===============================================================================\n" \
211  << "| |\n" \
212  << "| Unit Test (AdaptiveSparseGrid) |\n" \
213  << "| |\n" \
214  << "| 1) Integrate product peaks in 5 dimensions (Genz integration test). |\n" \
215  << "| |\n" \
216  << "| Questions? Contact Drew Kouri (dpkouri@sandia.gov) or |\n" \
217  << "| Denis Ridzal (dridzal@sandia.gov). |\n" \
218  << "| |\n" \
219  << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
220  << "| Trilinos website: http://trilinos.sandia.gov |\n" \
221  << "| |\n" \
222  << "===============================================================================\n"\
223  << "| TEST 18: Integrate a product of peaks functions in 5D |\n"\
224  << "===============================================================================\n";
225 
226 
227  // internal variables:
228  int errorFlag = 0;
229  long double TOL = INTREPID_TOL;
230  int dimension = 5;
231  int maxLevel = 7;
232  bool isNormalized = true;
233 
234  std::vector<EIntrepidBurkardt> rule1D(dimension,BURK_PATTERSON);
235  std::vector<EIntrepidGrowth> growth1D(dimension,GROWTH_FULLEXP);
236 
237  alpha.resize(dimension,0); beta.resize(dimension,0);
238  for (int i=0; i<dimension; i++) {
239  alpha[i] = (long double)std::rand()/(long double)RAND_MAX;
240  beta[i] = (long double)std::rand()/(long double)RAND_MAX;
241  }
242 
244  dimension,rule1D,growth1D,
245  maxLevel,isNormalized);
246  Teuchos::RCP<std::vector<long double> > integralValue =
247  Teuchos::rcp(new std::vector<long double>(1,0.0));
248  StdVector<long double> sol(integralValue); sol.Set(0.0);
249  problem_data.init(sol);
250 
251  long double eta = adaptSG(sol,problem_data,TOL);
252 
253  long double analyticInt = compExactInt();
254  long double abstol = std::sqrt(INTREPID_TOL);
255  long double absdiff = std::abs(analyticInt-(*integralValue)[0]);
256  try {
257  *outStream << "Adaptive Sparse Grid exited with global error "
258  << std::scientific << std::setprecision(16) << eta << "\n"
259  << "Approx = " << std::scientific
260  << std::setprecision(16) << (*integralValue)[0]
261  << ", Exact = " << std::scientific
262  << std::setprecision(16) << analyticInt << "\n"
263  << "Error = " << std::scientific << std::setprecision(16)
264  << absdiff << " "
265  << "<?" << " " << abstol << "\n";
266  if (absdiff > abstol) {
267  errorFlag++;
268  *outStream << std::right << std::setw(104) << "^^^^---FAILURE!\n";
269  }
270  }
271  catch (const std::logic_error & err) {
272  *outStream << err.what() << "\n";
273  errorFlag = -1;
274  };
275 
276  if (errorFlag != 0)
277  std::cout << "End Result: TEST FAILED\n";
278  else
279  std::cout << "End Result: TEST PASSED\n";
280 
281  // reset format state of std::cout
282  std::cout.copyfmt(oldFormatState);
283 
284  return errorFlag;
285 }
Builds general adaptive sparse grid rules (Gerstner and Griebel) using the 1D cubature rules in the I...
Scalar getInitialDiff()
Return initial error indicator.
Header file for the Intrepid::AdaptiveSparseGrid class.
Intrepid utilities.
Scalar error_indicator(UserVector &input)
User defined error indicator function.
Definition: test_18.cpp:136
void eval_integrand(UserVector &output, std::vector< Scalar > &input)
Evaluate the integrand function.
Definition: test_18.cpp:124