Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_DenseDirectDivisionExpansionStrategy.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Stokhos Package
4 //
5 // Copyright 2009 NTESS and the Stokhos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef STOKHOS_DENSE_DIRECT_DIVISION_EXPANSION_STRATEGY_HPP
11 #define STOKHOS_DENSE_DIRECT_DIVISION_EXPANSION_STRATEGY_HPP
12 
16 
17 #include "Teuchos_TimeMonitor.hpp"
18 #include "Teuchos_RCP.hpp"
21 
22 namespace Stokhos {
23 
25 
28  template <typename ordinal_type, typename value_type, typename node_type>
30  public DivisionExpansionStrategy<ordinal_type,value_type,node_type> {
31  public:
32 
37 
40 
41  // Division operation: c = \alpha*(a/b) + beta*c
42  virtual void divide(
44  const value_type& alpha,
47  const value_type& beta);
48 
49  private:
50 
51  // Prohibit copying
54 
55  // Prohibit Assignment
58 
59  protected:
60 
63 
66 
69 
72 
75 
76  }; // class DenseDirectDivisionExpansionStrategy
77 
78 } // namespace Stokhos
79 
80 template <typename ordinal_type, typename value_type, typename node_type>
85  basis(basis_),
86  Cijk(Cijk_),
87  solver()
88 {
89  ordinal_type sz = basis->size();
91  sz, sz));
93  sz, 1));
95  sz, 1));
96 }
97 
98 template <typename ordinal_type, typename value_type, typename node_type>
99 void
102  const value_type& alpha,
105  const value_type& beta)
106 {
107 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
108  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::DenseDirectDivisionStrategy::divide()");
109 #endif
110 
111  ordinal_type sz = basis->size();
112  ordinal_type pa = a.size();
113  ordinal_type pb = b.size();
114  ordinal_type pc;
115  if (pb > 1)
116  pc = sz;
117  else
118  pc = pa;
119  if (c.size() != pc)
120  c.resize(pc);
121 
122  const value_type* ca = a.coeff();
123  const value_type* cb = b.coeff();
124  value_type* cc = c.coeff();
125 
126  if (pb > 1) {
127  // Compute A
128  A->putScalar(0.0);
129  typename Cijk_type::k_iterator k_begin = Cijk->k_begin();
130  typename Cijk_type::k_iterator k_end = Cijk->k_end();
131  if (pb < Cijk->num_k())
132  k_end = Cijk->find_k(pb);
134  ordinal_type i,j,k;
135  for (typename Cijk_type::k_iterator k_it=k_begin; k_it!=k_end; ++k_it) {
136  k = index(k_it);
137  for (typename Cijk_type::kj_iterator j_it = Cijk->j_begin(k_it);
138  j_it != Cijk->j_end(k_it); ++j_it) {
139  j = index(j_it);
140  for (typename Cijk_type::kji_iterator i_it = Cijk->i_begin(j_it);
141  i_it != Cijk->i_end(j_it); ++i_it) {
142  i = index(i_it);
143  cijk = value(i_it);
144  (*A)(i,j) += cijk*cb[k];
145  }
146  }
147  }
148 
149  // Compute B
150  B->putScalar(0.0);
151  for (ordinal_type i=0; i<pa; i++)
152  (*B)(i,0) = ca[i]*basis->norm_squared(i);
153 
154  // Setup solver
155  solver.setMatrix(A);
156  solver.setVectors(X, B);
157  if (solver.shouldEquilibrate()) {
158  solver.factorWithEquilibration(true);
159  solver.equilibrateMatrix();
160  }
161 
162  // Compute X = A^{-1}*B
163  solver.solve();
164 
165  // value_type kappa;
166  // solver.reciprocalConditionEstimate(kappa);
167  // std::cout << "kappa = " << 1.0/kappa << std::endl;
168 
169  // Compute c
170  for (ordinal_type i=0; i<pc; i++)
171  cc[i] = alpha*(*X)(i,0) + beta*cc[i];
172  }
173  else {
174  for (ordinal_type i=0; i<pc; i++)
175  cc[i] = alpha*ca[i]/cb[0] + beta*cc[i];
176  }
177 }
178 
179 #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)
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.
ordinal_type size() const
Return size.
Strategy interface for computing PCE of a/b using only b[0].