29 template <
class ScalarType>
32 return 1.0/(z*z + 1.0);
42 const int EnsembleSize = 8;
63 Array< RCP<const OneDOrthogPolyBasis<int,double> > > bases(d);
64 for (
int i=0; i<d; i++) {
65 bases[i] =
rcp(
new HermiteBasis<int,double>(p));
67 RCP<const CompletePolynomialBasis<int,double> > basis =
68 rcp(
new CompletePolynomialBasis<int,double>(bases));
71 RCP<const Quadrature<int,double> > quad =
72 rcp(
new TensorProductQuadrature<int,double>(basis));
75 RCP<Sparse3Tensor<int,double> > Cijk =
76 basis->computeTripleProductTensor();
79 RCP<QuadOrthogPolyExpansion<int,double> > expn =
80 rcp(
new QuadOrthogPolyExpansion<int,double>(basis, Cijk, quad));
94 const int pce_size = basis->size();
95 const int num_quad_points = quad->size();
96 const Array<double>& quad_weights = quad->getQuadWeights();
97 const Array< Array<double> >& quad_points = quad->getQuadPoints();
98 const Array< Array<double> >& quad_values = quad->getBasisAtQuadPoints();
102 ensemble_type u_ensemble;
103 for (
int qp_block=0; qp_block<num_quad_points; qp_block+=EnsembleSize) {
104 const int qp_sz = qp_block+EnsembleSize <= num_quad_points ?
105 EnsembleSize : num_quad_points-qp_block;
108 for (
int qp=0; qp<qp_sz; ++qp)
109 u_ensemble.fastAccessCoeff(qp) =
110 u.evaluate(quad_points[qp_block+qp], quad_values[qp_block+qp]);
111 for (
int qp=qp_sz; qp<EnsembleSize; ++qp)
112 u_ensemble.fastAccessCoeff(qp) = u_ensemble.fastAccessCoeff(qp_sz-1);
118 for (
int pc=0; pc<pce_size; ++pc) {
119 const double inv_nrm_sq = 1.0 / basis->norm_squared(pc);
120 for (
int qp=0; qp<qp_sz; ++qp) {
121 const double w = quad_weights[qp_block+qp];
122 const double psi = quad_values[qp_block+qp][pc];
123 v.fastAccessCoeff(pc) +=
124 inv_nrm_sq * w * v_ensemble.fastAccessCoeff(qp) * psi;
130 std::cout <<
"\tu = ";
132 std::cout <<
"\tv = ";
136 double mean = v.mean();
137 double std_dev = v.standard_deviation();
141 for (
int i=0; i<d; i++)
143 double up = u.evaluate(pt);
145 double vp2 = v.evaluate(pt);
148 std::cout <<
"\tv mean = " << mean << std::endl;
149 std::cout <<
"\tv std. dev. = " << std_dev << std::endl;
150 std::cout <<
"\tv(0.25) (true) = " << vp << std::endl;
151 std::cout <<
"\tv(0.25) (pce) = " << vp2 << std::endl;
153 catch (std::exception& e) {
154 std::cout << e.what() << std::endl;
ScalarType simple_function(const ScalarType &u)
Hermite polynomial basis.
Sacado::ETPCE::OrthogPoly< double, Stokhos::StandardStorage< int, double > > pce_type
Data structure storing a sparse 3-tensor C(i,j,k) in a a compressed format.
Statically allocated storage class.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Abstract base class for quadrature methods.
Multivariate orthogonal polynomial basis generated from a total-order complete-polynomial tensor prod...
int main(int argc, char **argv)
Abstract base class for 1-D orthogonal polynomials.
Orthogonal polynomial expansions based on numerical quadrature.
KOKKOS_INLINE_FUNCTION PCE< Storage > log(const PCE< Storage > &a)
Defines quadrature for a tensor product basis by tensor products of 1-D quadrature rules...