Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_CGDivisionExpansionStrategy.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_CG_DIVISION_EXPANSION_STRATEGY_HPP
11 #define STOKHOS_CG_DIVISION_EXPANSION_STRATEGY_HPP
12 
22 
23 #include "Teuchos_TimeMonitor.hpp"
24 #include "Teuchos_RCP.hpp"
26 #include "Teuchos_BLAS.hpp"
27 #include "Teuchos_LAPACK.hpp"
28 
29 #include <iostream>
30 
31 namespace Stokhos {
32 
34 
37  template <typename ordinal_type, typename value_type, typename node_type>
39  public DivisionExpansionStrategy<ordinal_type,value_type,node_type> {
40  public:
41 
46  const ordinal_type prec_iter_,
47  const value_type tol_,
48  const ordinal_type PrecNum_,
49  const ordinal_type max_it_,
50  const ordinal_type linear_,
51  const ordinal_type diag_,
52  const ordinal_type equil_);
53 
56 
57  // Division operation: c = \alpha*(a/b) + beta*c
58  virtual void divide(
60  const value_type& alpha,
63  const value_type& beta);
64 
65  private:
66 
67  // Prohibit copying
70 
71  // Prohibit Assignment
74 
79  ordinal_type max_iter,
80  value_type tolerance,
82  ordinal_type order,
83  ordinal_type dim,
87 
88  protected:
89 
92 
95 
98 
101 
104 
106 
108 
110 
112 
114 
116 
117  }; // class CGDivisionExpansionStrategy
118 
119 } // namespace Stokhos
120 
121 template <typename ordinal_type, typename value_type, typename node_type>
126  const ordinal_type prec_iter_,
127  const value_type tol_,
128  const ordinal_type PrecNum_,
129  const ordinal_type max_it_,
130  const ordinal_type linear_,
131  const ordinal_type diag_,
132  const ordinal_type equil_):
133  basis(basis_),
134  Cijk(Cijk_),
135  prec_iter(prec_iter_),
136  tol(tol_),
137  PrecNum(PrecNum_),
138  max_it(max_it_),
139  linear(linear_),
140  diag(diag_),
141  equil(equil_)
142 
143 {
144  ordinal_type sz = basis->size();
146  sz, sz));
148  sz, 1));
150  sz, 1));
152  sz, sz));
153 
154 }
155 
156 
157 template <typename ordinal_type, typename value_type, typename node_type>
158 void
161  const value_type& alpha,
164  const value_type& beta)
165 {
166 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
167  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::CGDivisionStrategy::divide()");
168 #endif
169 
170 
171  ordinal_type sz = basis->size();
172  ordinal_type pa = a.size();
173  ordinal_type pb = b.size();
174 
175  ordinal_type pc;
176  if (pb > 1)
177  pc = sz;
178  else
179  pc = pa;
180  if (c.size() != pc)
181  c.resize(pc);
182 
183  const value_type* ca = a.coeff();
184  const value_type* cb = b.coeff();
185 
186 
187  value_type* cc = c.coeff();
188 
189  if (pb > 1) {
190  // Compute A
191  A->putScalar(0.0);
192  typename Cijk_type::k_iterator k_begin = Cijk->k_begin();
193  typename Cijk_type::k_iterator k_end = Cijk->k_end();
194 
195  if (pb < Cijk->num_k())
196  k_end = Cijk->find_k(pb);
198  ordinal_type i,j,k;
199  for (typename Cijk_type::k_iterator k_it=k_begin; k_it!=k_end; ++k_it) {
200  k = index(k_it);
201  for (typename Cijk_type::kj_iterator j_it = Cijk->j_begin(k_it);
202  j_it != Cijk->j_end(k_it); ++j_it) {
203  j = index(j_it);
204  for (typename Cijk_type::kji_iterator i_it = Cijk->i_begin(j_it);
205  i_it != Cijk->i_end(j_it); ++i_it) {
206  i = index(i_it);
207  cijk = value(i_it);
208  (*A)(i,j) += cijk*cb[k];
209  }
210  }
211  }
212 
213  // Compute B
214  B->putScalar(0.0);
215  for (ordinal_type i=0; i<pa; i++)
216  (*B)(i,0) = ca[i]*basis->norm_squared(i);
217 
219  //Equilibrate the linear system
220  if (equil == 1){
221  //Create diag mtx of max row entries
222  for (ordinal_type i=0; i<sz; i++){
224  D(i,0)=sqrt(r.normOne());
225  }
226 
227 
228  //Compute inv(D)*A*inv(D)
229  for (ordinal_type i=0; i<sz; i++){
230  for (ordinal_type j=0; j<sz; j++){
231  (*A)(i,j)=(*A)(i,j)/(D(i,0)*D(j,0));
232  }
233  }
234 
235  //Scale b by inv(D)
236  for (ordinal_type i=0; i<sz; i++){
237  (*B)(i,0)=(*B)(i,0)/D(i,0);
238  }
239 
240  }
241 
242  if (linear == 1){
243  //Compute M, the linear matrix to be used in the preconditioner
244 
245  pb = basis->dimension()+1;
246 
247  M->putScalar(0.0);
248  if (pb < Cijk->num_k())
249  k_end = Cijk->find_k(pb);
250  for (typename Cijk_type::k_iterator k_it=k_begin; k_it!=k_end; ++k_it) {
251  k = index(k_it);
252  for ( typename Cijk_type::kj_iterator j_it = Cijk->j_begin(k_it);
253  j_it != Cijk->j_end(k_it); ++j_it) {
254  j = index(j_it);
255  for ( typename Cijk_type::kji_iterator i_it = Cijk->i_begin(j_it);
256  i_it != Cijk->i_end(j_it); ++i_it) {
257  i = index(i_it);
258  cijk = value(i_it);
259  (*M)(i,j) += cijk*cb[k];
260  }
261  }
262  }
263 
264  //Scale M
265  if (equil == 1){
266  //Compute inv(D)*M*inv(D)
267  for (ordinal_type i=0; i<sz; i++){
268  for (ordinal_type j=0; j<sz; j++){
269  (*M)(i,j)=(*M)(i,j)/(D(i,0)*D(j,0));
270  }
271  }
272  }
273  CG(*A,*X,*B, max_it, tol, prec_iter, basis->order(), basis->dimension(), PrecNum, *M, diag);
274  }
275 
276  else{
277 
278  CG(*A,*X,*B, max_it, tol, prec_iter, basis->order(), basis->dimension(), PrecNum, *A, diag);
279  }
280 
281  if (equil == 1 ) {
282  //Rescale X
283  for (ordinal_type i=0; i<sz; i++){
284  (*X)(i,0)=(*X)(i,0)/D(i,0);
285  }
286  }
287 
288  // Compute c
289  for (ordinal_type i=0; i<pc; i++)
290  cc[i] = alpha*(*X)(i,0) + beta*cc[i];
291  }
292  else {
293  for (ordinal_type i=0; i<pc; i++)
294  cc[i] = alpha*ca[i]/cb[0] + beta*cc[i];
295  }
296 }
297 
298 
299 template <typename ordinal_type, typename value_type, typename node_type>
305  ordinal_type max_iter,
306  value_type tolerance,
307  ordinal_type prec_iter,
308  ordinal_type order ,
309  ordinal_type m,
310  ordinal_type PrecNum,
313 
314 {
315  ordinal_type n = A.numRows();
316  ordinal_type k=0;
317  value_type resid;
319  Ax.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,1.0, A, X, 0.0);
321  r-=Ax;
322  resid=r.normFrobenius();
328  value_type b;
329  value_type a;
330  while (resid > tolerance && k < max_iter){
332  //Solve Mz=r
333  if (PrecNum != 0){
334  if (PrecNum == 1){
336  precond.ApplyInverse(r,z,prec_iter);
337  }
338  else if (PrecNum == 2){
340  precond.ApplyInverse(r,z,2);
341  }
342  else if (PrecNum == 3){
344  precond.ApplyInverse(r,z,1);
345  }
346  else if (PrecNum == 4){
348  precond.ApplyInverse(r,z,prec_iter);
349  }
350  }
351  rho.multiply(Teuchos::TRANS,Teuchos::NO_TRANS,1.0, r, z, 0.0);
352 
353 
354  if (k==0){
355  p.assign(z);
356  rho.multiply(Teuchos::TRANS, Teuchos::NO_TRANS, 1.0, r, z, 0.0);
357  }
358  else {
359  b=rho(0,0)/oldrho(0,0);
360  p.scale(b);
361  p+=z;
362  }
363  Ap.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,1.0, A, p, 0.0);
364  pAp.multiply(Teuchos::TRANS,Teuchos::NO_TRANS,1.0, p, Ap, 0.0);
365  a=rho(0,0)/pAp(0,0);
367  scalep.scale(a);
368  X+=scalep;
369  Ap.scale(a);
370  r-=Ap;
371  oldrho.assign(rho);
372  resid=r.normFrobenius();
373  k++;
374  }
375 
376  //std::cout << "iteration count " << k << std::endl;
377  return 0;
378 }
379 
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.
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 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)
int n
OrdinalType numRows() const