42#include "Teuchos_Assert.hpp"
43#include "Teuchos_BLAS.hpp"
44#include "Teuchos_TimeMonitor.hpp"
46template <
typename ordinal_type,
typename value_type>
53 bool limit_integration_order_) :
78 for (
typename Cijk_type::k_iterator k_it = Cijk->k_begin();
79 k_it != Cijk->k_end(); ++k_it) {
81 for (
typename Cijk_type::kj_iterator j_it = Cijk->j_begin(k_it);
82 j_it != Cijk->j_end(k_it); ++j_it) {
85 for (
typename Cijk_type::kji_iterator i_it = Cijk->i_begin(j_it);
86 i_it != Cijk->i_end(j_it); ++i_it) {
101template <
typename ordinal_type,
typename value_type>
107template <
typename ordinal_type,
typename value_type>
111 Teuchos::Array<value_type>& quad_points,
112 Teuchos::Array<value_type>& quad_weights,
113 Teuchos::Array< Teuchos::Array<value_type> >& quad_values)
const
115#ifdef STOKHOS_TEUCHOS_TIME_MONITOR
116 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::LanczosPCEBasis -- compute Gauss points");
121 static_cast<ordinal_type>(std::ceil((quad_order+1)/2.0));
126 quad_order = 2*this->
p;
133 if (quad_weights.size() < num_points) {
135 quad_weights.
resize(num_points);
136 quad_points.resize(num_points);
137 quad_values.resize(num_points);
140 quad_points[i] = quad_points[0];
141 quad_values[i].resize(this->
p+1);
147template <
typename ordinal_type,
typename value_type>
148Teuchos::RCP<Stokhos::OneDOrthogPolyBasis<ordinal_type,value_type> >
157template <
typename ordinal_type,
typename value_type>
165template <
typename ordinal_type,
typename value_type>
171 Teuchos::BLAS<ordinal_type, value_type> blas;
172 blas.GEMV(Teuchos::NO_TRANS,
pce_sz, this->
p+1,
181template <
typename ordinal_type,
typename value_type>
185 Teuchos::Array<value_type>&
alpha,
186 Teuchos::Array<value_type>&
beta,
187 Teuchos::Array<value_type>&
delta,
188 Teuchos::Array<value_type>&
gamma)
const
190 Teuchos::Array<value_type> nrm(n);
196 Teuchos::RCP<matrix_type> lv;
234template <
typename ordinal_type,
typename value_type>
251template <
typename ordinal_type,
typename value_type>
Generates three-term recurrence using the Lanczos procedure applied to a polynomial chaos expansion i...
lanczos_type::matrix_type matrix_type
virtual void setup()
Setup basis after computing recurrence coefficients.
WeightedVectorSpace< ordinal_type, value_type > vectorspace_type
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.
lanczos_type::vector_type vector_type
vector_type weights
Weighting vector used in inner-products.
matrix_type Cijk_matrix
Triple-product matrix used in generating lanczos vectors.
void transformCoeffsFromLanczos(const value_type *in, value_type *out) const
Map expansion coefficients from this basis to original.
Teuchos::RCP< const Stokhos::OrthogPolyApprox< ordinal_type, value_type > > pce
PCE Lanczos procedure is based on.
LanczosProjPCEBasis(ordinal_type p, const Teuchos::RCP< const Stokhos::OrthogPolyApprox< ordinal_type, value_type > > &pce, const Teuchos::RCP< const Stokhos::Sparse3Tensor< ordinal_type, value_type > > &Cijk, bool normalize, bool limit_integration_order=false)
Constructor.
matrix_type lanczos_vecs
Lanczos vectors.
Teuchos::Array< value_type > pce_norms
Basis norms.
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.
DenseOperator< ordinal_type, value_type > operator_type
vector_type new_pce
Projection of pce in new basis.
~LanczosProjPCEBasis()
Destructor.
ordinal_type pce_sz
Size of PC expansion.
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.
bool limit_integration_order
Flag indicating whether to limit the integration order.
vector_type u0
Initial Lanczos vector.
value_type getNewCoeffs(ordinal_type i) const
Get new coefficients in this new basis.
static void computeNormalized(ordinal_type k, const vectorspace_type &vs, const operator_type &A, const vector_type &u_init, matrix_type &u, Teuchos::Array< value_type > &alpha, Teuchos::Array< value_type > &beta, Teuchos::Array< value_type > &nrm_sqrd)
Compute Lanczos basis.
static void compute(ordinal_type k, const vectorspace_type &vs, const operator_type &A, const vector_type &u_init, matrix_type &u, Teuchos::Array< value_type > &alpha, Teuchos::Array< value_type > &beta, Teuchos::Array< value_type > &nrm_sqrd)
Compute Lanczos basis.
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.
Teuchos::Array< value_type > norms
Norms.
bool normalize
Normalize basis.
virtual const Teuchos::Array< value_type > & norm_squared() const
Return array storing norm-squared of each basis polynomial.
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.
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)