44#ifndef STOKHOS_DENSE_DIRECT_DIVISION_EXPANSION_STRATEGY_HPP
45#define STOKHOS_DENSE_DIRECT_DIVISION_EXPANSION_STRATEGY_HPP
51#include "Teuchos_TimeMonitor.hpp"
52#include "Teuchos_RCP.hpp"
53#include "Teuchos_SerialDenseMatrix.hpp"
54#include "Teuchos_SerialDenseSolver.hpp"
62 template <
typename ordinal_type,
typename value_type,
typename node_type>
78 const value_type& alpha,
81 const value_type& beta);
96 Teuchos::RCP<const Stokhos::OrthogPolyBasis<ordinal_type, value_type> >
basis;
102 Teuchos::RCP<const Cijk_type>
Cijk;
105 Teuchos::RCP< Teuchos::SerialDenseMatrix<ordinal_type,value_type> >
A,
X,
B;
108 Teuchos::SerialDenseSolver<ordinal_type,value_type>
solver;
114template <
typename ordinal_type,
typename value_type,
typename node_type>
123 ordinal_type sz =
basis->size();
124 A = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<ordinal_type,value_type>(
126 B = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<ordinal_type,value_type>(
128 X = Teuchos::rcp(
new Teuchos::SerialDenseMatrix<ordinal_type,value_type>(
132template <
typename ordinal_type,
typename value_type,
typename node_type>
136 const value_type& alpha,
139 const value_type& beta)
141#ifdef STOKHOS_TEUCHOS_TIME_MONITOR
142 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::DenseDirectDivisionStrategy::divide()");
145 ordinal_type sz = basis->size();
146 ordinal_type pa = a.
size();
147 ordinal_type pb = b.
size();
156 const value_type* ca = a.
coeff();
157 const value_type* cb = b.
coeff();
158 value_type* cc = c.
coeff();
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];
185 for (ordinal_type i=0; i<pa; i++)
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();
204 for (ordinal_type i=0; i<pc; i++)
205 cc[i] = alpha*(*X)(i,0) + beta*cc[i];
208 for (ordinal_type i=0; i<pc; i++)
209 cc[i] = alpha*ca[i]/cb[0] + beta*cc[i];
Strategy interface for computing PCE of a/b using only b[0].
Teuchos::RCP< const Stokhos::OrthogPolyBasis< ordinal_type, value_type > > basis
Basis.
Teuchos::RCP< Teuchos::SerialDenseMatrix< ordinal_type, value_type > > A
Dense matrices for linear system.
Teuchos::RCP< const Cijk_type > Cijk
Triple product.
DenseDirectDivisionExpansionStrategy(const DenseDirectDivisionExpansionStrategy &)
DenseDirectDivisionExpansionStrategy & operator=(const DenseDirectDivisionExpansionStrategy &b)
Teuchos::SerialDenseSolver< ordinal_type, value_type > solver
Serial dense solver.
Teuchos::RCP< Teuchos::SerialDenseMatrix< ordinal_type, value_type > > B
Stokhos::Sparse3Tensor< ordinal_type, value_type > Cijk_type
Short-hand for Cijk.
virtual ~DenseDirectDivisionExpansionStrategy()
Destructor.
Teuchos::RCP< Teuchos::SerialDenseMatrix< ordinal_type, value_type > > X
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)
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.
Class to store coefficients of a projection onto an orthogonal polynomial basis.
void resize(ordinal_type sz)
Resize coefficient array (coefficients are preserved)
pointer coeff()
Return coefficient array.
ordinal_type size() const
Return size.
Abstract base class for multivariate orthogonal polynomials.
Data structure storing a sparse 3-tensor C(i,j,k) in a a compressed format.
kji_sparse_array::const_iterator k_iterator
Iterator for looping over k entries.
j_sparse_array::const_iterator kji_iterator
Iterator for looping over i entries given k and j.
Top-level namespace for Stokhos classes and functions.
Bi-directional iterator for traversing a sparse array.