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_) :
54 RecurrenceBasis<ordinal_type, value_type>(
"Lanczos-proj PCE", p, normalize),
56 limit_integration_order(limit_integration_order_),
57 pce_sz(pce->basis()->size()),
58 Cijk_matrix(pce_sz,pce_sz),
60 const_cast<value_type*>(pce->basis()->norm_squared().getRawPtr()),
63 lanczos_vecs(pce_sz, p+1),
66 u0[0] = value_type(1);
69 for (ordinal_type i=0; i<
pce_sz; i++) {
78 for (
typename Cijk_type::k_iterator k_it = Cijk->k_begin();
79 k_it != Cijk->k_end(); ++k_it) {
80 ordinal_type k = index(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) {
83 ordinal_type
j = index(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) {
87 ordinal_type i = index(i_it);
88 value_type c = value(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");
120 ordinal_type num_points =
121 static_cast<ordinal_type
>(std::ceil((quad_order+1)/2.0));
125 if (limit_integration_order && quad_order > 2*this->p)
126 quad_order = 2*this->p;
133 if (quad_weights.size() < num_points) {
134 ordinal_type old_size = quad_weights.size();
135 quad_weights.resize(num_points);
136 quad_points.resize(num_points);
137 quad_values.resize(num_points);
138 for (ordinal_type i=old_size; i<num_points; i++) {
139 quad_weights[i] = value_type(0);
140 quad_points[i] = quad_points[0];
141 quad_values[i].resize(this->p+1);
142 this->evaluateBases(quad_points[i], quad_values[i]);
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,
173 value_type(1.0), lanczos_vecs.values(), pce_sz,
174 in, ordinal_type(1), value_type(0.0), out, ordinal_type(1));
177 for (ordinal_type i=0; i<pce_sz; i++)
178 out[i] /= pce_norms[i];
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;
198 lv = Teuchos::rcp(&lanczos_vecs,
false);
203 lanczos_type::computeNormalized(n, vs, A, u0, *lv, alpha, beta, nrm);
205 lanczos_type::compute(n, vs, A, u0, *lv, alpha, beta, nrm);
207 for (ordinal_type i=0; i<n; i++) {
208 delta[i] = value_type(1.0);
213 for (ordinal_type i=0; i<n; i++)
214 gamma[i] = value_type(1.0);
231 return this->normalize;
234template <
typename ordinal_type,
typename value_type>
243 for (ordinal_type i=0; i<pce_sz; i++)
244 u[i] = (*pce)[i]*pce_norms[i];
245 new_pce.multiply(Teuchos::TRANS, Teuchos::NO_TRANS, 1.0, lanczos_vecs, u,
247 for (ordinal_type i=0; i<=this->p; i++)
248 new_pce[i] /= this->norms[i];
251template <
typename ordinal_type,
typename value_type>
254 RecurrenceBasis<ordinal_type, value_type>(
"Lanczos-proj PCE", p, false),
256 limit_integration_order(basis.limit_integration_order),
257 pce_sz(basis.pce_sz),
258 pce_norms(basis.pce_norms),
259 Cijk_matrix(basis.Cijk_matrix),
260 weights(basis.weights),
262 lanczos_vecs(pce_sz, p+1),
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.
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.
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.
~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.
vector_type u0
Initial Lanczos vector.
value_type getNewCoeffs(ordinal_type i) const
Get new coefficients in this new basis.
Class to store coefficients of a projection onto an orthogonal polynomial basis.
Implementation of OneDOrthogPolyBasis based on the general three-term recurrence relationship:
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.
Data structure storing a sparse 3-tensor C(i,j,k) in a a compressed format.