16 template <
typename ordinal_type, 
typename value_type>
 
   23    bool limit_integration_order_) :
 
   25   limit_integration_order(limit_integration_order_),
 
   26   pce_sz(pce.basis()->size()),
 
   27   pce_norms(pce.basis()->norm_squared()),
 
   30   basis_vecs(pce_sz, p+1),
 
   51     pce_vals[i] = normalized_pce.
evaluate(quad_points[i], basis_values[i]);
 
   60   for (
typename Cijk_type::k_iterator k_it = Cijk.
k_begin();
 
   61        k_it != Cijk.
k_end(); ++k_it) {
 
   63     for (
typename Cijk_type::kj_iterator j_it = Cijk.
j_begin(k_it); 
 
   64    j_it != Cijk.
j_end(k_it); ++j_it) {
 
   67       for (
typename Cijk_type::kji_iterator i_it = Cijk.
i_begin(j_it);
 
   68      i_it != Cijk.
i_end(j_it); ++i_it) {
 
   99   std::cout << K << std::endl << std::endl;
 
  107          &ws_size_query, -1, &info);
 
  109          "GEQRF returned value " << info);
 
  113          &work[0], ws_size, &info);
 
  115          "GEQRF returned value " << info);
 
  119          &ws_size_query, -1, &info);
 
  121          "ORGQR returned value " << info);
 
  123   work.resize(ws_size);
 
  125          &work[0], ws_size, &info);
 
  127          "ORGQR returned value " << info);
 
  149   a[pce_sz-1] = T(pce_sz-1,pce_sz-1);
 
  157     u[i] = normalized_pce[i];
 
  164 template <
typename ordinal_type, 
typename value_type>
 
  170 template <
typename ordinal_type, 
typename value_type>
 
  178 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR 
  188   if (limit_integration_order && quad_order > 2*this->p)
 
  189     quad_order = 2*this->p;
 
  196   if (quad_weights.
size() < num_points) {
 
  198     quad_weights.
resize(num_points);
 
  199     quad_points.
resize(num_points);
 
  200     quad_values.resize(num_points);
 
  203       quad_points[i] = quad_points[0];
 
  204       quad_values[i].
resize(this->p+1);
 
  205       evaluateBases(quad_points[i], quad_values[i]);
 
  210 template <
typename ordinal_type, 
typename value_type>
 
  219 template <
typename ordinal_type, 
typename value_type>
 
  227 template <
typename ordinal_type, 
typename value_type>
 
  240     out[i] /= pce_norms[i];
 
  243 template <
typename ordinal_type, 
typename value_type>
 
  259     std::cout << 
"i = " << i << 
" alpha = " << alpha[i] << 
" beta = " << beta[i]
 
  260         << 
" gamma = " << gamma[i] << std::endl;
 
  266 template <
typename ordinal_type, 
typename value_type> 
 
  270   limit_integration_order(basis.limit_integration_order),
 
  271   pce_sz(basis.pce_sz),
 
  272   pce_norms(basis.pce_norms),
 
  275   basis_vecs(basis.basis_vecs),
 
  276   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 
 
#define TEUCHOS_FUNC_TIME_MONITOR(FUNCNAME)
 
value_type getNewCoeffs(ordinal_type i) const 
Get new coefficients in this new basis. 
 
Generates three-term recurrence using the Lanczos procedure applied to a polynomial chaos expansion i...
 
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. 
 
value_type evaluate(const Teuchos::Array< value_type > &point) const 
Evaluate polynomial approximation at a point. 
 
Data structure storing a sparse 3-tensor C(i,j,k) in a a compressed format. 
 
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 
 
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)
 
virtual ordinal_type size() const =0
Get number of quadrature points. 
 
virtual const Teuchos::Array< Teuchos::Array< value_type > > & getBasisAtQuadPoints() const =0
Get values of basis at quadrature points. 
 
matrix_type basis_vecs
Basis vectors. 
 
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. 
 
virtual const Teuchos::Array< value_type > & getQuadWeights() const =0
Get quadrature weights. 
 
virtual const Teuchos::Array< Teuchos::Array< value_type > > & getQuadPoints() const =0
Get quadrature points. 
 
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
 
Abstract base class for quadrature methods. 
 
KOKKOS_INLINE_FUNCTION PCE< Storage > ceil(const PCE< Storage > &a)
 
~MonoProjPCEBasis()
Destructor. 
 
void transformCoeffs(const value_type *in, value_type *out) const 
Map expansion coefficients from this basis to original. 
 
void resize(size_type new_size, const value_type &x=value_type())
 
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. 
 
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. 
 
ordinal_type pce_sz
Size of PC expansion. 
 
void GEQRF(const OrdinalType &m, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *TAU, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) 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. 
 
MonoProjPCEBasis(ordinal_type p, const Stokhos::OrthogPolyApprox< ordinal_type, value_type > &pce, const Stokhos::Quadrature< ordinal_type, value_type > &quad, const Stokhos::Sparse3Tensor< ordinal_type, value_type > &Cijk, bool limit_integration_order=false)
Constructor. 
 
vector_type new_pce
Projection of pce in new basis. 
 
kji_iterator i_begin(const kj_iterator &j) const 
Iterator pointing to first i entry for given j and k. 
 
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. 
 
kji_iterator i_end(const kj_iterator &j) const 
Iterator pointing to last i entry for given j and k. 
 
OrdinalType stride() const