Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_SPDDenseDirectDivisionExpansionStrategy.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_SPD_DENSE_DIRECT_DIVISION_EXPANSION_STRATEGY_HPP
11 #define STOKHOS_SPD_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 
74 
77 
78  }; // class SPDDenseDirectDivisionExpansionStrategy
79 
80 } // namespace Stokhos
81 
82 template <typename ordinal_type, typename value_type, typename node_type>
87  basis(basis_),
88  Cijk(Cijk_),
89  solver()
90 {
91  ordinal_type sz = basis->size();
93  sz, sz));
95  sz, 1));
97  sz, 1));
98 }
99 
100 template <typename ordinal_type, typename value_type, typename node_type>
101 void
104  const value_type& alpha,
107  const value_type& beta)
108 {
109 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
110  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::SPDDenseDirectDivisionStrategy::divide()");
111 #endif
112 
113  ordinal_type sz = basis->size();
114  ordinal_type pa = a.size();
115  ordinal_type pb = b.size();
116  ordinal_type pc;
117  if (pb > 1)
118  pc = sz;
119  else
120  pc = pa;
121  if (c.size() != pc)
122  c.resize(pc);
123 
124  const value_type* ca = a.coeff();
125  const value_type* cb = b.coeff();
126  value_type* cc = c.coeff();
127 
128  if (pb > 1) {
129  // Compute A
130  A->putScalar(0.0);
131  typename Cijk_type::k_iterator k_begin = Cijk->k_begin();
132  typename Cijk_type::k_iterator k_end = Cijk->k_end();
133  if (pb < Cijk->num_k())
134  k_end = Cijk->find_k(pb);
136  ordinal_type i,j,k;
137  for (typename Cijk_type::k_iterator k_it=k_begin; k_it!=k_end; ++k_it) {
138  k = index(k_it);
139  for (typename Cijk_type::kj_iterator j_it = Cijk->j_begin(k_it);
140  j_it != Cijk->j_end(k_it); ++j_it) {
141  j = index(j_it);
142  for (typename Cijk_type::kji_iterator i_it = Cijk->i_begin(j_it);
143  i_it != Cijk->i_end(j_it); ++i_it) {
144  i = index(i_it);
145  cijk = value(i_it);
146  (*A)(i,j) += cijk*cb[k];
147  }
148  }
149  }
150  A->setUpper();
151 
152  // Compute B
153  B->putScalar(0.0);
154  for (ordinal_type i=0; i<pa; i++)
155  (*B)(i,0) = ca[i]*basis->norm_squared(i);
156 
157  // Setup solver
158  solver.setMatrix(A);
159  solver.setVectors(X, B);
160  if (solver.shouldEquilibrate()) {
161  solver.factorWithEquilibration(true);
162  solver.equilibrateMatrix();
163  }
164 
165  // Compute X = A^{-1}*B
166  solver.solve();
167 
168  // value_type kappa;
169  // solver.reciprocalConditionEstimate(kappa);
170  // std::cout << "kappa = " << 1.0/kappa << std::endl;
171 
172  // Compute c
173  for (ordinal_type i=0; i<pc; i++)
174  cc[i] = alpha*(*X)(i,0) + beta*cc[i];
175  }
176  else {
177  for (ordinal_type i=0; i<pc; i++)
178  cc[i] = alpha*ca[i]/cb[0] + beta*cc[i];
179  }
180 }
181 
182 #endif // STOKHOS_DIVISION_EXPANSION_STRATEGY_HPP
SPDDenseDirectDivisionExpansionStrategy & operator=(const SPDDenseDirectDivisionExpansionStrategy &b)
#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)
Data structure storing a sparse 3-tensor C(i,j,k) in a a compressed format.
Teuchos::RCP< Teuchos::SerialSymDenseMatrix< ordinal_type, value_type > > A
Dense matrices for linear system.
pointer coeff()
Return coefficient array.
Stokhos::Sparse3Tensor< ordinal_type, value_type > Cijk_type
Short-hand for Cijk.
Bi-directional iterator for traversing a sparse array.
Abstract base class for multivariate orthogonal polynomials.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
SPDDenseDirectDivisionExpansionStrategy(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::SerialSpdDenseSolver< ordinal_type, value_type > solver
Serial dense solver.
Teuchos::RCP< Teuchos::SerialDenseMatrix< ordinal_type, value_type > > X
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.
Teuchos::RCP< Teuchos::SerialDenseMatrix< ordinal_type, value_type > > B
ordinal_type size() const
Return size.
Strategy interface for computing PCE of a/b using only b[0].
Teuchos::RCP< const Stokhos::OrthogPolyBasis< ordinal_type, value_type > > basis
Basis.