|
Stokhos Package Browser (Single Doxygen Collection)
Version of the Day
|
Go to the documentation of this file.
46 template <
typename ordinal_type,
typename value_type>
53 bool limit_integration_order_) :
57 limit_integration_order(limit_integration_order_),
64 lanczos_vecs(nqp, p+1),
70 quad->getQuadPoints();
72 quad->getBasisAtQuadPoints();
74 pce_vals[i] = pce->evaluate(quad_points[i], basis_values[i]);
82 template <
typename ordinal_type,
typename value_type>
88 template <
typename ordinal_type,
typename value_type>
96 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
102 static_cast<ordinal_type>(std::ceil((quad_order+1)/2.0));
106 if (limit_integration_order && quad_order > 2*this->p)
107 quad_order = 2*this->p;
114 if (quad_weights.
size() < num_points) {
116 quad_weights.
resize(num_points);
117 quad_points.
resize(num_points);
118 quad_values.resize(num_points);
121 quad_points[i] = quad_points[0];
122 quad_values[i].
resize(this->p+1);
123 this->evaluateBases(quad_points[i], quad_values[i]);
128 template <
typename ordinal_type,
typename value_type>
137 template <
typename ordinal_type,
typename value_type>
145 template <
typename ordinal_type,
typename value_type>
153 value_type(1.0), fromStieltjesMat.values(), sz,
157 template <
typename ordinal_type,
typename value_type>
179 lanczos_type::computeNormalized(
n, vs,
A, u0, *lv, alpha, beta, nrm);
181 lanczos_type::compute(
n, vs,
A, u0, *lv, alpha, beta, nrm);
192 return this->normalize;
195 template <
typename ordinal_type,
typename value_type>
204 fromStieltjesMat.shape(sz, this->p+1);
205 fromStieltjesMat.putScalar(0.0);
207 quad->getBasisAtQuadPoints();
211 fromStieltjesMat(i,
j) +=
212 pce_weights[k]*lanczos_vecs(k,
j)*basis_values[k][i];
213 fromStieltjesMat(i,
j) /= pce->basis()->norm_squared(i);
218 new_pce.
resize(this->p+1);
221 u[i] = (*pce)[i]*pce->basis()->norm_squared(i);
225 new_pce[i] /= this->norms[i];
228 template <
typename ordinal_type,
typename value_type>
234 limit_integration_order(basis.limit_integration_order),
236 pce_weights(basis.pce_weights),
237 pce_vals(basis.pce_vals),
239 lanczos_vecs(nqp, p+1),
void resize(size_type new_size, const value_type &x=value_type())
RawPointerConversionTraits< Container >::Ptr_t getRawPtr(const Container &c)
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_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
virtual void setup()
Setup basis after computing recurrence coefficients.
Abstract base class for quadrature methods.
LanczosPCEBasis(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 normalize, bool limit_integration_order)
Constructor.
Generates three-term recurrence using the Lanczos procedure applied to a polynomial chaos expansion i...
~LanczosPCEBasis()
Destructor.
Implementation of OneDOrthogPolyBasis based on the general three-term recurrence relationship:
void GEMV(ETransp trans, const OrdinalType &m, const OrdinalType &n, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const x_type *x, const OrdinalType &incx, const beta_type beta, ScalarType *y, const OrdinalType &incy) const
void transformCoeffsFromLanczos(const value_type *in, value_type *out) const
Map expansion coefficients from this basis to original.
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.
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.
#define TEUCHOS_FUNC_TIME_MONITOR(FUNCNAME)
UnitTestSetup< int, double > setup
int multiply(ETransp transa, ETransp transb, ScalarType alpha, const SerialDenseMatrix< OrdinalType, ScalarType > &A, const SerialDenseMatrix< OrdinalType, ScalarType > &B, ScalarType beta)
value_type getNewCoeffs(ordinal_type i) const
Get new coefficients in this new basis.
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.