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