|
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>
52 bool limit_integration_order_) :
54 limit_integration_order(limit_integration_order_),
55 pce_sz(pce.basis()->size()),
58 basis_vecs(pce_sz, p+1),
61 pce_norms = pce.
basis()->norm_squared();
67 A(0,i+1) = 1.0/pce_sz;
68 A(i+1,0) = 1.0/pce_sz;
71 for (
typename Cijk_type::k_iterator k_it = Cijk.
k_begin();
72 k_it != Cijk.
k_end(); ++k_it) {
74 for (
typename Cijk_type::kj_iterator j_it = Cijk.
j_begin(k_it);
75 j_it != Cijk.
j_end(k_it); ++j_it) {
78 for (
typename Cijk_type::kji_iterator i_it = Cijk.
i_begin(j_it);
79 i_it != Cijk.
i_end(j_it); ++i_it) {
93 lapack.SYTRD(
'L', pce_sz+1,
A.values(),
A.stride(), &a[0], &b[0], &tau[0],
94 &ws_size_query, -1, &info);
96 "SYTRD returned value " << info);
97 ws_size = static_cast<ordinal_type>(ws_size_query);
99 lapack.SYTRD(
'L', pce_sz+1,
A.values(),
A.stride(), &a[0], &b[0], &tau[0],
100 &work[0], ws_size, &info);
102 "SYTRD returned value " << info);
109 lapack.ORGQR(pce_sz+1, pce_sz+1, pce_sz-1,
A.values(),
A.stride(), &tau[0],
110 &ws_size_query, -1, &info);
112 "ORGQR returned value " << info);
113 ws_size = static_cast<ordinal_type>(ws_size_query);
114 work.resize(ws_size);
115 lapack.ORGQR(pce_sz+1, pce_sz+1, pce_sz-1,
A.values(),
A.stride(), &tau[0],
116 &work[0], ws_size, &info);
118 "ORGQR returned value " << info);
123 basis_vecs(i,
j) =
A(i+1,
j+1);
135 new_pce[i] /= this->norms[i];
137 std::cout << new_pce << std::endl;
145 std::cout <<
"orthogonalization error = " << prod.
normInf() << std::endl;
149 template <
typename ordinal_type,
typename value_type>
155 template <
typename ordinal_type,
typename value_type>
163 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
169 static_cast<ordinal_type>(std::ceil((quad_order+1)/2.0));
173 if (limit_integration_order && quad_order > 2*this->p)
174 quad_order = 2*this->p;
181 if (quad_weights.
size() < num_points) {
183 quad_weights.
resize(num_points);
184 quad_points.
resize(num_points);
185 quad_values.resize(num_points);
188 quad_points[i] = quad_points[0];
189 quad_values[i].
resize(this->p+1);
190 evaluateBases(quad_points[i], quad_values[i]);
195 template <
typename ordinal_type,
typename value_type>
205 template <
typename ordinal_type,
typename value_type>
213 template <
typename ordinal_type,
typename value_type>
222 out[i] /= pce_norms[i];
225 template <
typename ordinal_type,
typename value_type>
241 std::cout <<
"i = " << i <<
" alpha = " << alpha[i] <<
" beta = " << beta[i]
242 <<
" gamma = " << gamma[i] << std::endl;
248 template <
typename ordinal_type,
typename value_type>
252 limit_integration_order(basis.limit_integration_order),
253 pce_sz(basis.pce_sz),
256 basis_vecs(basis.basis_vecs),
257 new_pce(basis.new_pce)
kji_iterator i_begin(const kj_iterator &j) const
Iterator pointing to first i entry for given j and k.
void resize(size_type new_size, const value_type &x=value_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.
kj_iterator j_end(const k_iterator &k) const
Iterator pointing to last j entry for given k.
Generates three-term recurrence using the Lanczos procedure applied to a polynomial chaos expansion i...
value_type getNewCoeffs(ordinal_type i) const
Get new coefficients in this new basis.
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.
void transformCoeffsFromHouse(const value_type *in, value_type *out) const
Map expansion coefficients from this basis to original.
Teuchos::RCP< const Stokhos::OrthogPolyBasis< ordinal_type, value_type > > basis() const
Return basis.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
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.
kj_iterator j_begin(const k_iterator &k) const
Iterator pointing to first j entry for given k.
Data structure storing a sparse 3-tensor C(i,j,k) in a a compressed format.
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.
Stokhos::Sparse3Tensor< int, double > Cijk_type
Implementation of OneDOrthogPolyBasis based on the general three-term recurrence relationship:
KOKKOS_INLINE_FUNCTION PCE< Storage > sqrt(const PCE< Storage > &a)
ScalarTraits< ScalarType >::magnitudeType normInf() const
~HouseTriDiagPCEBasis()
Destructor.
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)
k_iterator k_begin() const
Iterator pointing to first k entry.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
kji_iterator i_end(const kj_iterator &j) const
Iterator pointing to last i entry for given j and k.
k_iterator k_end() const
Iterator pointing to last k entry.