58 double operator() (
const double& a,
const double& b)
const {
65 mutable Teuchos::Array<double>
vec;
71 double operator() (
const Teuchos::Array<double>& x)
const {
80 double operator() (
const Teuchos::Array<double>& x)
const {
87 return std::abs(a-b)/std::abs(b);
94 const unsigned int d = 2;
95 const unsigned int pmin = 1;
96 const unsigned int pmax = 10;
97 const unsigned int np = pmax-pmin+1;
98 bool use_pce_quad_points =
false;
99 bool normalize =
false;
100 bool sparse_grid =
true;
101#ifndef HAVE_STOKHOS_DAKOTA
104 Teuchos::Array<double> mean(np), mean_st(np), std_dev(np), std_dev_st(np);
105 Teuchos::Array<double> pt(np), pt_st(np);
107 Teuchos::Array< Teuchos::RCP<const Stokhos::OneDOrthogPolyBasis<int,double> > > bases(d);
108 Teuchos::Array<double> eval_pt(d, 0.5);
113 for (
unsigned int p=pmin; p<=pmax; p++) {
115 std::cout <<
"p = " << p << std::endl;
118 for (
unsigned int i=0; i<d; i++)
120 Teuchos::RCP<const Stokhos::CompletePolynomialBasis<int,double> > basis =
126 for (
unsigned int i=0; i<d; i++) {
130 double x_pt = x.evaluate(eval_pt);
131 pt_true = std::sin(x_pt)/(1.0 + exp(x_pt));
134 Teuchos::RCP<const Stokhos::Quadrature<int,double> > quad;
135#ifdef HAVE_STOKHOS_DAKOTA
138 Teuchos::rcp(
new Stokhos::SparseGridQuadrature<int,double>(basis, p));
145 Teuchos::RCP<Stokhos::Sparse3Tensor<int,double> > Cijk =
146 basis->computeTripleProductTensor();
158 Teuchos::Array< Teuchos::RCP<const Stokhos::OneDOrthogPolyBasis<int,double> > > st_bases(2);
159 Teuchos::RCP<sin_func> sf = Teuchos::rcp(
new sin_func(x));
160 Teuchos::RCP<exp_func> ef = Teuchos::rcp(
new exp_func(x));
163 p, sf, quad, use_pce_quad_points, normalize));
166 p, ef, quad, use_pce_quad_points, normalize));
167 Teuchos::RCP<const Stokhos::CompletePolynomialBasis<int,double> >
174 u_st.term(0, 0) = u.mean();
175 u_st.term(0, 1) = 1.0;
176 v_st.term(0, 0) = v.mean();
177 v_st.term(1, 1) = 1.0;
180 Teuchos::RCP<Stokhos::Sparse3Tensor<int,double> > st_Cijk =
181 st_basis->computeTripleProductTensor();
184 Teuchos::RCP<const Stokhos::Quadrature<int,double> > st_quad;
185 if (!use_pce_quad_points) {
186#ifdef HAVE_STOKHOS_DAKOTA
188 st_quad = Teuchos::rcp(
new Stokhos::SparseGridQuadrature<int,double>(st_basis, p));
194 Teuchos::Array<double> st_points_0;
195 Teuchos::Array<double> st_weights_0;
196 Teuchos::Array< Teuchos::Array<double> > st_values_0;
197 st_bases[0]->getQuadPoints(p+1, st_points_0, st_weights_0, st_values_0);
198 Teuchos::Array<double> st_points_1;
199 Teuchos::Array<double> st_weights_1;
200 Teuchos::Array< Teuchos::Array<double> > st_values_1;
201 st_bases[1]->getQuadPoints(p+1, st_points_1, st_weights_1, st_values_1);
202 Teuchos::RCP< Teuchos::Array< Teuchos::Array<double> > > st_points =
203 Teuchos::rcp(
new Teuchos::Array< Teuchos::Array<double> >(st_points_0.size()));
204 for (
int i=0; i<st_points_0.size(); i++) {
205 (*st_points)[i].resize(2);
206 (*st_points)[i][0] = st_points_0[i];
207 (*st_points)[i][1] = st_points_1[i];
209 Teuchos::RCP< Teuchos::Array<double> > st_weights =
210 Teuchos::rcp(
new Teuchos::Array<double>(st_weights_0));
211 Teuchos::RCP< const Stokhos::OrthogPolyBasis<int,double> > st_b =
225 st_quad_exp.
divide(w_st, u_st, v_st);
236 mean_st[n] = w2.
mean();
237 std_dev[n] = w.standard_deviation();
239 pt[n] = w.evaluate(eval_pt);
246 std::cout <<
"Statistical error:" << std::endl;
248 << std::setw(wi) <<
"mean" <<
" "
249 << std::setw(wi) <<
"mean_st" <<
" "
250 << std::setw(wi) <<
"std_dev" <<
" "
251 << std::setw(wi) <<
"std_dev_st" <<
" "
252 << std::setw(wi) <<
"point" <<
" "
253 << std::setw(wi) <<
"point_st" << std::endl;
254 for (
unsigned int p=pmin; p<pmax; p++) {
255 std::cout.precision(3);
256 std::cout.setf(std::ios::scientific);
257 std::cout << p <<
" "
258 << std::setw(wi) <<
rel_err(mean[n], mean[np-1]) <<
" "
259 << std::setw(wi) <<
rel_err(mean_st[n], mean[np-1]) <<
" "
260 << std::setw(wi) <<
rel_err(std_dev[n], std_dev[np-1]) <<
" "
261 << std::setw(wi) <<
rel_err(std_dev_st[n], std_dev[np-1])
263 << std::setw(wi) <<
rel_err(pt[n], pt_true) <<
" "
264 << std::setw(wi) <<
rel_err(pt_st[n], pt_true)
270 catch (std::exception& e) {
271 std::cout << e.what() << std::endl;
Legendre polynomial basis using Clenshaw-Curtis quadrature points.
Multivariate orthogonal polynomial basis generated from a total-order complete-polynomial tensor prod...
Class to store coefficients of a projection onto an orthogonal polynomial basis.
value_type mean() const
Compute mean of expansion.
value_type standard_deviation() const
Compute standard deviation of expansion.
value_type evaluate(const Teuchos::Array< value_type > &point) const
Evaluate polynomial approximation at a point.
Abstract base class for multivariate orthogonal polynomials.
void plusEqual(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &x)
Orthogonal polynomial expansions based on numerical quadrature.
void exp(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void divide(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a, const OrthogPolyApprox< ordinal_type, value_type, node_type > &b)
void binary_op(const FuncT &func, OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a, const OrthogPolyApprox< ordinal_type, value_type, node_type > &b)
Nonlinear binary function.
void sin(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
Generates three-term recurrence using the Discretized Stieltjes procedure applied to a functional map...
Defines quadrature for a tensor product basis by tensor products of 1-D quadrature rules.
int main(int argc, char **argv)
double rel_err(double a, double b)
Stokhos::ClenshawCurtisLegendreBasis< int, double > basis_type
const Stokhos::OrthogPolyApprox< int, double > & pce
double operator()(const Teuchos::Array< double > &x) const
exp_func(const Stokhos::OrthogPolyApprox< int, double > &pce_)
const Stokhos::OrthogPolyApprox< int, double > & pce
double operator()(const double &a, const double &b) const
const Stokhos::OrthogPolyBasis< int, double > & basis
pce_quad_func(const Stokhos::OrthogPolyApprox< int, double > &pce_, const Stokhos::OrthogPolyBasis< int, double > &basis_)
Teuchos::Array< double > vec
sin_func(const Stokhos::OrthogPolyApprox< int, double > &pce_)
double operator()(const Teuchos::Array< double > &x) const
const Stokhos::OrthogPolyApprox< int, double > & pce