14 template <
typename ordinal_type,
typename value_type>
20 bool limit_integration_order_) :
22 limit_integration_order(limit_integration_order_),
23 pce_sz(pce.basis()->size()),
26 basis_vecs(pce_sz, p+1),
39 for (
typename Cijk_type::k_iterator k_it = Cijk.
k_begin();
40 k_it != Cijk.
k_end(); ++k_it) {
42 for (
typename Cijk_type::kj_iterator j_it = Cijk.
j_begin(k_it);
43 j_it != Cijk.
j_end(k_it); ++j_it) {
46 for (
typename Cijk_type::kji_iterator i_it = Cijk.
i_begin(j_it);
47 i_it != Cijk.
i_end(j_it); ++i_it) {
62 &ws_size_query, -1, &info);
64 "SYTRD returned value " << info);
68 &work[0], ws_size, &info);
70 "SYTRD returned value " << info);
78 &ws_size_query, -1, &info);
80 "ORGQR returned value " << info);
84 &work[0], ws_size, &info);
86 "ORGQR returned value " << info);
105 std::cout <<
new_pce << std::endl;
113 std::cout <<
"orthogonalization error = " << prod.
normInf() << std::endl;
117 template <
typename ordinal_type,
typename value_type>
123 template <
typename ordinal_type,
typename value_type>
131 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
141 if (limit_integration_order && quad_order > 2*this->p)
142 quad_order = 2*this->p;
149 if (quad_weights.
size() < num_points) {
151 quad_weights.
resize(num_points);
152 quad_points.
resize(num_points);
153 quad_values.resize(num_points);
156 quad_points[i] = quad_points[0];
157 quad_values[i].
resize(this->p+1);
158 evaluateBases(quad_points[i], quad_values[i]);
163 template <
typename ordinal_type,
typename value_type>
173 template <
typename ordinal_type,
typename value_type>
181 template <
typename ordinal_type,
typename value_type>
190 out[i] /= pce_norms[i];
193 template <
typename ordinal_type,
typename value_type>
209 std::cout <<
"i = " << i <<
" alpha = " << alpha[i] <<
" beta = " << beta[i]
210 <<
" gamma = " << gamma[i] << std::endl;
216 template <
typename ordinal_type,
typename value_type>
220 limit_integration_order(basis.limit_integration_order),
221 pce_sz(basis.pce_sz),
224 basis_vecs(basis.basis_vecs),
225 new_pce(basis.new_pce)
ScalarType * values() const
KOKKOS_INLINE_FUNCTION PCE< Storage > sqrt(const PCE< Storage > &a)
void ORGQR(const OrdinalType &m, const OrdinalType &n, const OrdinalType &k, ScalarType *A, const OrdinalType &lda, const ScalarType *TAU, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const
ordinal_type pce_sz
Size of PC expansion.
#define TEUCHOS_FUNC_TIME_MONITOR(FUNCNAME)
void SYTRD(const char &UPLO, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *D, ScalarType *E, ScalarType *TAU, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const
Teuchos::Array< value_type > norms
Norms.
Implementation of OneDOrthogPolyBasis based on the general three-term recurrence relationship: for ...
k_iterator k_begin() const
Iterator pointing to first k entry.
Data structure storing a sparse 3-tensor C(i,j,k) in a a compressed format.
int multiply(ETransp transa, ETransp transb, ScalarType alpha, const SerialDenseMatrix< OrdinalType, ScalarType > &A, const SerialDenseMatrix< OrdinalType, ScalarType > &B, ScalarType beta)
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...
kj_iterator j_begin(const k_iterator &k) const
Iterator pointing to first j entry for given k.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
value_type getNewCoeffs(ordinal_type i) const
Get new coefficients in this new basis.
Teuchos::Array< value_type > pce_norms
Basis norms.
HouseTriDiagPCEBasis(ordinal_type p, const Stokhos::OrthogPolyApprox< ordinal_type, value_type > &pce, const Stokhos::Sparse3Tensor< ordinal_type, value_type > &Cijk, bool limit_integration_order=false)
Constructor.
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.
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.
kj_iterator j_end(const k_iterator &k) const
Iterator pointing to last j entry for given k.
vector_type new_pce
Projection of pce in new basis.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
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.
KOKKOS_INLINE_FUNCTION PCE< Storage > ceil(const PCE< Storage > &a)
Generates three-term recurrence using the Lanczos procedure applied to a polynomial chaos expansion i...
Teuchos::LAPACK< ordinal_type, value_type > lapack
LAPACK routines.
Teuchos::RCP< const Stokhos::OrthogPolyBasis< ordinal_type, value_type > > basis() const
Return basis.
~HouseTriDiagPCEBasis()
Destructor.
void resize(size_type new_size, const value_type &x=value_type())
matrix_type basis_vecs
Basis vectors.
k_iterator k_end() const
Iterator pointing to last k entry.
Stokhos::Sparse3Tensor< int, double > Cijk_type
virtual void setup()
Setup basis after computing recurrence coefficients.
ScalarTraits< ScalarType >::magnitudeType normInf() const
Teuchos::Array< value_type > b
Stores full set of beta coefficients.
Teuchos::Array< value_type > a
Stores full set of alpha coefficients.
ordinal_type p
Order of basis.
void transformCoeffsFromHouse(const value_type *in, value_type *out) const
Map expansion coefficients from this basis to original.
kji_iterator i_begin(const kj_iterator &j) const
Iterator pointing to first i entry for given j and k.
kji_iterator i_end(const kj_iterator &j) const
Iterator pointing to last i entry for given j and k.
OrdinalType stride() const