44#include "Teuchos_Assert.hpp"
45#include "Teuchos_BLAS.hpp"
46#include "Teuchos_TimeMonitor.hpp"
48template <
typename ordinal_type,
typename value_type>
54 bool use_pce_quad_points_,
56 bool project_integrals_,
76template <
typename ordinal_type,
typename value_type>
82template <
typename ordinal_type,
typename value_type>
86 Teuchos::Array<value_type>& quad_points,
87 Teuchos::Array<value_type>& quad_weights,
88 Teuchos::Array< Teuchos::Array<value_type> >& quad_values)
const
90#ifdef STOKHOS_TEUCHOS_TIME_MONITOR
91 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::StieltjesPCEBasis -- compute Gauss points");
104 static_cast<ordinal_type>(std::ceil((quad_order+1)/2.0));
108 if (quad_order > 2*this->
p)
109 quad_order = 2*this->
p;
116 if (quad_weights.size() < num_points) {
118 quad_weights.
resize(num_points);
119 quad_points.resize(num_points);
120 quad_values.resize(num_points);
123 quad_points[i] = quad_points[0];
124 quad_values[i].resize(this->
p+1);
130template <
typename ordinal_type,
typename value_type>
134 Teuchos::Array<value_type>&
alpha,
135 Teuchos::Array<value_type>&
beta,
136 Teuchos::Array<value_type>&
delta,
137 Teuchos::Array<value_type>&
gamma)
const
140 Teuchos::Array<value_type> nrm(n);
141 Teuchos::Array< Teuchos::Array<value_type> > vals(nqp);
157template <
typename ordinal_type,
typename value_type>
163 const Teuchos::Array< Teuchos::Array<value_type> >& quad_points =
164 quad->getQuadPoints();
190template <
typename ordinal_type,
typename value_type>
195 const Teuchos::Array<value_type>& weights,
196 const Teuchos::Array<value_type>& points,
197 Teuchos::Array<value_type>& a,
198 Teuchos::Array<value_type>& b,
199 Teuchos::Array<value_type>& nrm,
200 Teuchos::Array< Teuchos::Array<value_type> >&
phi_vals)
const
202#ifdef STOKHOS_TEUCHOS_TIME_MONITOR
203 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::StieltjesPCEBasis -- Discretized Stieltjes Procedure");
225 TEUCHOS_TEST_FOR_EXCEPTION(val1 < 0.0, std::logic_error,
226 "Stokhos::StieltjesPCEBasis::stieltjes(): "
227 <<
" Polynomial " << i <<
" out of " << nfinish
228 <<
" has norm " << val1
229 <<
"! Try increasing number of quadrature points");
232 b[i] = nrm[i]/nrm[i-1];
239template <
typename ordinal_type,
typename value_type>
243 const Teuchos::Array<value_type>& a,
244 const Teuchos::Array<value_type>& b,
245 const Teuchos::Array<value_type>& weights,
246 const Teuchos::Array<value_type>& points,
247 Teuchos::Array< Teuchos::Array<value_type> >&
phi_vals,
260template <
typename ordinal_type,
typename value_type>
264 const Teuchos::Array<value_type>& a,
265 const Teuchos::Array<value_type>& b,
266 const Teuchos::Array<value_type>& points,
267 Teuchos::Array< Teuchos::Array<value_type> >& values)
const
275 values[i][k] = points[i] - a[k-1];
279 (points[i] - a[k-1])*values[i][k-1] - b[k-1]*values[i][k-2];
282template <
typename ordinal_type,
typename value_type>
287 const Teuchos::Array<value_type>& a,
288 const Teuchos::Array<value_type>& b,
289 const Teuchos::Array<value_type>& weights,
290 const Teuchos::Array<value_type>& points,
291 Teuchos::Array< Teuchos::Array<value_type> >&
phi_vals,
296 const Teuchos::Array<value_type>&
norms =
basis->norm_squared();
316 k_it !=
Cijk->k_end(); ++k_it) {
319 j_it !=
Cijk->j_end(k_it); ++j_it) {
322 i_it !=
Cijk->i_end(j_it); ++i_it) {
331template <
typename ordinal_type,
typename value_type>
336 Teuchos::BLAS<ordinal_type, value_type> blas;
342template <
typename ordinal_type,
typename value_type>
343Teuchos::RCP<Stokhos::OneDOrthogPolyBasis<ordinal_type,value_type> >
350template <
typename ordinal_type,
typename value_type>
Class to store coefficients of a projection onto an orthogonal polynomial basis.
ordinal_type size() const
Return size.
void resize(const Teuchos::RCP< const Epetra_BlockMap > &map)
Resize to map map.
Abstract base class for quadrature methods.
Teuchos::Array< value_type > norms
Norms.
bool normalize
Normalize basis.
virtual ordinal_type size() const
Return total size of basis (given by order() + 1).
virtual void evaluateBases(const value_type &point, Teuchos::Array< value_type > &basis_pts) const
Evaluate each basis polynomial at given point point.
Teuchos::Array< value_type > alpha
Recurrence coefficients.
Teuchos::Array< value_type > beta
Recurrence coefficients.
ordinal_type p
Order of basis.
Teuchos::Array< value_type > gamma
Recurrence coefficients.
Teuchos::Array< value_type > delta
Recurrence coefficients.
virtual void setup()
Setup basis after computing recurrence coefficients.
virtual void getQuadPoints(ordinal_type quad_order, Teuchos::Array< value_type > &points, Teuchos::Array< value_type > &weights, Teuchos::Array< Teuchos::Array< value_type > > &values) const
Compute quadrature points, weights, and values of basis polynomials at given set of points points.
RecurrenceBasis(const std::string &name, ordinal_type p, bool normalize, GrowthPolicy growth=SLOW_GROWTH)
Constructor to be called by derived classes.
Data structure storing a sparse 3-tensor C(i,j,k) in a a compressed format.
kji_sparse_array::const_iterator k_iterator
Iterator for looping over k entries.
SparseArrayIterator< index_iterator, value_iterator >::value_reference value(const SparseArrayIterator< index_iterator, value_iterator > &it)
SparseArrayIterator< index_iterator, value_iterator >::value_type index(const SparseArrayIterator< index_iterator, value_iterator > &it)
j_sparse_array::const_iterator kji_iterator
Iterator for looping over i entries given k and j.
ji_sparse_array::const_iterator kj_iterator
Iterator for looping over j entries given k.
Generates three-term recurrence using the Discretized Stieltjes procedure applied to a polynomial cha...
Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > quad
Quadrature object.
Teuchos::SerialDenseMatrix< ordinal_type, value_type > fromStieltjesMat
Matrix mapping coefficients in Stieltjes basis back to original basis.
virtual void getQuadPoints(ordinal_type quad_order, Teuchos::Array< value_type > &points, Teuchos::Array< value_type > &weights, Teuchos::Array< Teuchos::Array< value_type > > &values) const
Get Gauss quadrature points, weights, and values of basis at points.
Teuchos::Array< Teuchos::Array< value_type > > phi_vals
Values of generated polynomials at PCE quadrature points.
void evaluateRecurrence(ordinal_type k, const Teuchos::Array< value_type > &a, const Teuchos::Array< value_type > &b, const Teuchos::Array< value_type > &points, Teuchos::Array< Teuchos::Array< value_type > > &values) const
Evaluate polynomials via 3-term recurrence.
void integrateBasisSquared(ordinal_type k, const Teuchos::Array< value_type > &a, const Teuchos::Array< value_type > &b, const Teuchos::Array< value_type > &weights, const Teuchos::Array< value_type > &points, Teuchos::Array< Teuchos::Array< value_type > > &phi_vals, value_type &val1, value_type &val2) const
Compute and .
Teuchos::Array< value_type > pce_vals
Values of PCE at quadrature points.
bool use_pce_quad_points
Use underlying pce's quadrature data.
Teuchos::RCP< const Cijk_type > Cijk
Triple product tensor (needed for integral projection method)
Teuchos::RCP< const Stokhos::OrthogPolyApprox< ordinal_type, value_type > > pce
PC expansion.
virtual Teuchos::RCP< OneDOrthogPolyBasis< ordinal_type, value_type > > cloneWithOrder(ordinal_type p) const
Clone this object with the option of building a higher order basis.
const Teuchos::Array< Teuchos::Array< value_type > > & basis_values
Values of PCE basis functions at quadrature points.
Teuchos::RCP< const Stokhos::OrthogPolyBasis< ordinal_type, value_type > > basis
PCE basis (needed for integral projection method)
bool project_integrals
Project Stieltjes integrals to original PCE basis.
void integrateBasisSquaredProj(ordinal_type k, const Teuchos::Array< value_type > &a, const Teuchos::Array< value_type > &b, const Teuchos::Array< value_type > &weights, const Teuchos::Array< value_type > &points, Teuchos::Array< Teuchos::Array< value_type > > &phi_vals, value_type &val1, value_type &val2) const
Compute and by projecting onto original PCE basis.
virtual void setup()
Setup basis after computing recurrence coefficients.
const Teuchos::Array< value_type > & pce_weights
PCE quadrature weights.
void stieltjes(ordinal_type nstart, ordinal_type nfinish, const Teuchos::Array< value_type > &weights, const Teuchos::Array< value_type > &points, Teuchos::Array< value_type > &a, Teuchos::Array< value_type > &b, Teuchos::Array< value_type > &nrm, Teuchos::Array< Teuchos::Array< value_type > > &phi_vals) const
Compute 3-term recurrence using Stieljtes procedure.
StieltjesPCEBasis(ordinal_type p, const Teuchos::RCP< const Stokhos::OrthogPolyApprox< ordinal_type, value_type > > &pce, const Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > &quad, bool use_pce_quad_points, bool normalize=false, bool project_integrals=false, const Teuchos::RCP< const Stokhos::Sparse3Tensor< ordinal_type, value_type > > &Cijk=Teuchos::null)
Constructor.
~StieltjesPCEBasis()
Destructor.
void transformCoeffsFromStieltjes(const value_type *in, value_type *out) const
Map expansion coefficients from this basis to original.
Teuchos::Array< value_type > phi_pce_coeffs
Array store PC expansion of generated orthogonal polynomials.
virtual bool computeRecurrenceCoefficients(ordinal_type n, Teuchos::Array< value_type > &alpha, Teuchos::Array< value_type > &beta, Teuchos::Array< value_type > &delta, Teuchos::Array< value_type > &gamma) const
Compute recurrence coefficients.