10 #ifndef STOKHOS_CG_DIVISION_EXPANSION_STRATEGY_HPP
11 #define STOKHOS_CG_DIVISION_EXPANSION_STRATEGY_HPP
37 template <
typename ordinal_type,
typename value_type,
typename node_type>
121 template <
typename ordinal_type,
typename value_type,
typename node_type>
135 prec_iter(prec_iter_),
157 template <
typename ordinal_type,
typename value_type,
typename node_type>
166 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
195 if (pb < Cijk->num_k())
196 k_end = Cijk->find_k(pb);
202 j_it != Cijk->j_end(k_it); ++j_it) {
205 i_it != Cijk->i_end(j_it); ++i_it) {
208 (*A)(i,
j) += cijk*cb[k];
216 (*
B)(i,0) = ca[i]*basis->norm_squared(i);
231 (*A)(i,
j)=(*
A)(i,
j)/(D(i,0)*D(j,0));
237 (*B)(i,0)=(*
B)(i,0)/D(i,0);
245 pb = basis->dimension()+1;
248 if (pb < Cijk->num_k())
249 k_end = Cijk->find_k(pb);
253 j_it != Cijk->j_end(k_it); ++j_it) {
256 i_it != Cijk->i_end(j_it); ++i_it) {
259 (*M)(i,
j) += cijk*cb[k];
269 (*M)(i,
j)=(*M)(i,
j)/(D(i,0)*D(j,0));
273 CG(*
A,*X,*
B, max_it, tol, prec_iter, basis->order(), basis->dimension(), PrecNum, *M,
diag);
278 CG(*
A,*X,*
B, max_it, tol, prec_iter, basis->order(), basis->dimension(), PrecNum, *
A,
diag);
284 (*X)(i,0)=(*X)(i,0)/D(i,0);
290 cc[i] = alpha*(*X)(i,0) + beta*cc[i];
294 cc[i] = alpha*ca[i]/cb[0] + beta*cc[i];
299 template <
typename ordinal_type,
typename value_type,
typename node_type>
330 while (resid > tolerance && k < max_iter){
338 else if (PrecNum == 2){
342 else if (PrecNum == 3){
346 else if (PrecNum == 4){
359 b=rho(0,0)/oldrho(0,0);
380 #endif // STOKHOS_DIVISION_EXPANSION_STRATEGY_HPP
ScalarTraits< ScalarType >::magnitudeType normOne() const
KOKKOS_INLINE_FUNCTION PCE< Storage > sqrt(const PCE< Storage > &a)
#define TEUCHOS_FUNC_TIME_MONITOR(FUNCNAME)
void resize(ordinal_type sz)
Resize coefficient array (coefficients are preserved)
Strategy interface for computing PCE of a/b using only b[0].
virtual ordinal_type ApplyInverse(const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &Input, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &Result, ordinal_type m) const
Returns the result of a Operator inverse applied to a Teuchos::SerialDenseMatrix Input in Result...
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)
Teuchos::RCP< Teuchos::SerialDenseMatrix< ordinal_type, value_type > > M
Teuchos::RCP< Teuchos::SerialDenseMatrix< ordinal_type, value_type > > X
int scale(const ScalarType alpha)
pointer coeff()
Return coefficient array.
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.
ScalarTraits< ScalarType >::magnitudeType normFrobenius() const
virtual ordinal_type ApplyInverse(const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &Input, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &Result, ordinal_type m) const
Returns the result of a Operator inverse applied to a Teuchos::SerialDenseMatrix Input in Result...
Bi-directional iterator for traversing a sparse array.
ordinal_type prec_iter
Tolerance for CG.
Teuchos::RCP< Teuchos::SerialDenseMatrix< ordinal_type, value_type > > B
Abstract base class for multivariate orthogonal polynomials.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
virtual ~CGDivisionExpansionStrategy()
Destructor.
virtual ordinal_type ApplyInverse(const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &Input, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &Result, ordinal_type prec_iters) const
Returns the result of a Operator inverse applied to a Teuchos::SerialDenseMatrix Input in Result...
int putScalar(const ScalarType value=Teuchos::ScalarTraits< ScalarType >::zero())
Strategy interface for computing PCE of a/b.
virtual ordinal_type ApplyInverse(const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &Input, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &Result, ordinal_type m) const
Returns the result of a Operator inverse applied to a Teuchos::SerialDenseMatrix Input in Result...
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)
Class to store coefficients of a projection onto an orthogonal polynomial basis.
CGDivisionExpansionStrategy & operator=(const CGDivisionExpansionStrategy &b)
Teuchos::RCP< const Cijk_type > Cijk
Triple product.
ordinal_type size() const
Return size.
CGDivisionExpansionStrategy(const Teuchos::RCP< const Stokhos::OrthogPolyBasis< ordinal_type, value_type > > &basis_, const Teuchos::RCP< const Stokhos::Sparse3Tensor< ordinal_type, value_type > > &Cijk_, const ordinal_type prec_iter_, const value_type tol_, const ordinal_type PrecNum_, const ordinal_type max_it_, const ordinal_type linear_, const ordinal_type diag_, const ordinal_type equil_)
Constructor.
Stokhos::Sparse3Tensor< ordinal_type, value_type > Cijk_type
Short-hand for Cijk.
ordinal_type CG(const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &A, Teuchos::SerialDenseMatrix< ordinal_type, value_type > &X, const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &B, ordinal_type max_iter, value_type tolerance, ordinal_type prec_iter, ordinal_type order, ordinal_type dim, ordinal_type PrecNum, const Teuchos::SerialDenseMatrix< ordinal_type, value_type > &M, ordinal_type diag)
SerialDenseMatrix< OrdinalType, ScalarType > & assign(const SerialDenseMatrix< OrdinalType, ScalarType > &Source)
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)
OrdinalType numRows() const