44 #ifndef STOKHOS_DENSE_DIRECT_DIVISION_EXPANSION_STRATEGY_HPP 
   45 #define STOKHOS_DENSE_DIRECT_DIVISION_EXPANSION_STRATEGY_HPP 
   62   template <
typename ordinal_type, 
typename value_type, 
typename node_type> 
 
  114 template <
typename ordinal_type, 
typename value_type, 
typename node_type> 
 
  132 template <
typename ordinal_type, 
typename value_type, 
typename node_type> 
 
  141 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR 
  165     if (pb < Cijk->num_k())
 
  166       k_end = Cijk->find_k(pb);
 
  172      j_it != Cijk->j_end(k_it); ++j_it) {
 
  175        i_it != Cijk->i_end(j_it); ++i_it) {
 
  178     (*A)(i,
j) += cijk*cb[k];
 
  186       (*
B)(i,0) = ca[i]*basis->norm_squared(i);
 
  190     solver.setVectors(X, 
B);
 
  191     if (solver.shouldEquilibrate()) {
 
  192       solver.factorWithEquilibration(
true);
 
  193       solver.equilibrateMatrix();
 
  205       cc[i] = alpha*(*X)(i,0) + beta*cc[i];
 
  209       cc[i] = alpha*ca[i]/cb[0] + beta*cc[i];
 
  213 #endif // STOKHOS_DIVISION_EXPANSION_STRATEGY_HPP 
#define TEUCHOS_FUNC_TIME_MONITOR(FUNCNAME)
virtual void divide(Stokhos::OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const value_type &alpha, const Stokhos::OrthogPolyApprox< ordinal_type, value_type, node_type > &a, const Stokhos::OrthogPolyApprox< ordinal_type, value_type, node_type > &b, const value_type &beta)
Teuchos::RCP< const Cijk_type > Cijk
Triple product. 
void resize(ordinal_type sz)
Resize coefficient array (coefficients are preserved) 
Teuchos::SerialDenseSolver< ordinal_type, value_type > solver
Serial dense solver. 
Data structure storing a sparse 3-tensor C(i,j,k) in a a compressed format. 
pointer coeff()
Return coefficient array. 
Teuchos::RCP< Teuchos::SerialDenseMatrix< ordinal_type, value_type > > X
Bi-directional iterator for traversing a sparse array. 
Stokhos::Sparse3Tensor< ordinal_type, value_type > Cijk_type
Short-hand for Cijk. 
Abstract base class for multivariate orthogonal polynomials. 
DenseDirectDivisionExpansionStrategy & operator=(const DenseDirectDivisionExpansionStrategy &b)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
DenseDirectDivisionExpansionStrategy(const Teuchos::RCP< const Stokhos::OrthogPolyBasis< ordinal_type, value_type > > &basis_, const Teuchos::RCP< const Stokhos::Sparse3Tensor< ordinal_type, value_type > > &Cijk_)
Constructor. 
Strategy interface for computing PCE of a/b. 
Teuchos::RCP< Teuchos::SerialDenseMatrix< ordinal_type, value_type > > B
KOKKOS_INLINE_FUNCTION constexpr std::enable_if< is_view_uq_pce< view_type >::value, typename CijkType< view_type >::type >::type cijk(const view_type &view)
Teuchos::RCP< Teuchos::SerialDenseMatrix< ordinal_type, value_type > > A
Dense matrices for linear system. 
Class to store coefficients of a projection onto an orthogonal polynomial basis. 
Teuchos::RCP< const Stokhos::OrthogPolyBasis< ordinal_type, value_type > > basis
Basis. 
virtual ~DenseDirectDivisionExpansionStrategy()
Destructor. 
ordinal_type size() const 
Return size. 
Strategy interface for computing PCE of a/b using only b[0].