Intrepid
test_20.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 
61 template<class Scalar>
62 class StdVector {
63 private:
64  Teuchos::RefCountPtr<std::vector<Scalar> > std_vec_;
65 
66 public:
67 
68  StdVector( const Teuchos::RefCountPtr<std::vector<Scalar> > & std_vec )
69  : std_vec_(std_vec) {}
70 
71  Teuchos::RefCountPtr<StdVector<Scalar> > Create() const {
72  return Teuchos::rcp( new StdVector<Scalar>(
73  Teuchos::rcp(new std::vector<Scalar>(std_vec_->size(),0))));
74  }
75 
76  void Update( StdVector<Scalar> & s ) {
77  int dimension = (int)(std_vec_->size());
78  for (int i=0; i<dimension; i++)
79  (*std_vec_)[i] += s[i];
80  }
81 
82  void Update( Scalar alpha, StdVector<Scalar> & s ) {
83  int dimension = (int)(std_vec_->size());
84  for (int i=0; i<dimension; i++)
85  (*std_vec_)[i] += alpha*s[i];
86  }
87 
88  Scalar operator[](int i) {
89  return (*std_vec_)[i];
90  }
91 
92  void clear() {
93  std_vec_->clear();
94  }
95 
96  void resize(int n, Scalar p) {
97  std_vec_->resize(n,p);
98  }
99 
100  int size() {
101  return (int)std_vec_->size();
102  }
103 
104  void Set( Scalar alpha ) {
105  int dimension = (int)(std_vec_->size());
106  for (int i=0; i<dimension; i++)
107  (*std_vec_)[i] = alpha;
108  }
109 };
110 
111 template<class Scalar, class UserVector>
112 class ASGdata :
113  public Intrepid::AdaptiveSparseGridInterface<Scalar,UserVector> {
114 public:
115  ~ASGdata() {}
116 
117  ASGdata(int dimension,std::vector<EIntrepidBurkardt> rule1D,
118  std::vector<EIntrepidGrowth> growth1D, int maxLevel,
119  bool isNormalized) : AdaptiveSparseGridInterface<Scalar,UserVector>(
120  dimension,rule1D,growth1D,maxLevel,isNormalized) {}
121 
122  void eval_integrand(UserVector & output, std::vector<Scalar> & input) {
123  output.clear(); output.resize(1,std::exp(-input[0]*input[0])
124  +10.0*std::exp(-input[1]*input[1]));
125  }
126  Scalar error_indicator(UserVector & input) {
127  int dimension = (int)input.size();
128  Scalar norm2 = 0.0;
129  for (int i=0; i<dimension; i++)
130  norm2 += input[i]*input[i];
131 
134  norm2 = std::sqrt(norm2)/ID;
135  return norm2;
136  }
137 };
138 
139 long double adaptSG(StdVector<long double> & iv,
141  problem_data,long double TOL) {
142 
143  // Construct a Container for the adapted rule
144  int dimension = problem_data.getDimension();
145  std::vector<int> index(dimension,1);
146 
147  // Initialize global error indicator
148  long double eta = 1.0;
149 
150  // Initialize the Active index set
151  std::multimap<long double,std::vector<int> > activeIndex;
152  activeIndex.insert(std::pair<long double,std::vector<int> >(eta,index));
153 
154  // Initialize the old index set
155  std::set<std::vector<int> > oldIndex;
156 
157  // Perform Adaptation
158  while (eta > TOL) {
160  activeIndex,oldIndex,iv,eta,problem_data);
161  }
162  return eta;
163 }
164 
165 int main(int argc, char *argv[]) {
166 
167  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
168 
169  // This little trick lets us print to std::cout only if
170  // a (dummy) command-line argument is provided.
171  int iprint = argc - 1;
172  Teuchos::RCP<std::ostream> outStream;
173  Teuchos::oblackholestream bhs; // outputs nothing
174  if (iprint > 0)
175  outStream = Teuchos::rcp(&std::cout, false);
176  else
177  outStream = Teuchos::rcp(&bhs, false);
178 
179  // Save the format state of the original std::cout.
180  Teuchos::oblackholestream oldFormatState;
181  oldFormatState.copyfmt(std::cout);
182 
183  *outStream \
184  << "===============================================================================\n" \
185  << "| |\n" \
186  << "| Unit Test (AdaptiveSparseGrid) |\n" \
187  << "| |\n" \
188  << "| 1) Integrate a sum of Gaussians in 2D (Gerstner and Griebel). |\n" \
189  << "| |\n" \
190  << "| Questions? Contact Drew Kouri (dpkouri@sandia.gov) or |\n" \
191  << "| Denis Ridzal (dridzal@sandia.gov). |\n" \
192  << "| |\n" \
193  << "| Intrepid's website: http://trilinos.sandia.gov/packages/intrepid |\n" \
194  << "| Trilinos website: http://trilinos.sandia.gov |\n" \
195  << "| |\n" \
196  << "===============================================================================\n"\
197  << "| TEST 20: Integrate an anisotropic sum of Gaussians in 2D |\n"\
198  << "===============================================================================\n";
199 
200 
201  // internal variables:
202  int errorFlag = 0;
203  long double TOL = INTREPID_TOL;
204  int dimension = 2;
205  int maxLevel = 25;
206  bool isNormalized = true;
207 
208  std::vector<EIntrepidBurkardt> rule1D(dimension,BURK_CLENSHAWCURTIS);
209  std::vector<EIntrepidGrowth> growth1D(dimension,GROWTH_FULLEXP);
210 
212  dimension,rule1D,growth1D,maxLevel,isNormalized);
213  Teuchos::RCP<std::vector<long double> > integralValue =
214  Teuchos::rcp(new std::vector<long double>(1,0.0));
215  StdVector<long double> sol(integralValue); sol.Set(0.0);
216  problem_data.init(sol);
217 
218  long double eta = adaptSG(sol,problem_data,TOL);
219 
220  long double analyticInt = (1.0+10.0)*std::sqrt(M_PI)/2.0*erff(1.0);
221  long double abstol = 1.0e1*std::sqrt(INTREPID_TOL);
222  long double absdiff = std::abs(analyticInt-sol[0]);
223  try {
224  *outStream << "Adaptive Sparse Grid exited with global error "
225  << std::scientific << std::setprecision(16) << eta << "\n"
226  << "Approx = " << std::scientific << std::setprecision(16) << sol[0]
227  << ", Exact = " << std::scientific << std::setprecision(16) << analyticInt << "\n"
228  << "Error = " << std::scientific << std::setprecision(16) << absdiff << " "
229  << "<?" << " " << abstol << "\n";
230  if (absdiff > abstol) {
231  errorFlag++;
232  *outStream << std::right << std::setw(104) << "^^^^---FAILURE!\n";
233  }
234  }
235  catch (std::logic_error err) {
236  *outStream << err.what() << "\n";
237  errorFlag = -1;
238  };
239 
240  if (errorFlag != 0)
241  std::cout << "End Result: TEST FAILED\n";
242  else
243  std::cout << "End Result: TEST PASSED\n";
244 
245  // reset format state of std::cout
246  std::cout.copyfmt(oldFormatState);
247 
248  return errorFlag;
249 }
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_20.cpp:126
void eval_integrand(UserVector &output, std::vector< Scalar > &input)
Evaluate the integrand function.
Definition: test_20.cpp:122