Reference documentation for deal.II version 9.0.0
solver_cg.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2018 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_solver_cg_h
17 #define dealii_solver_cg_h
18 
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/lac/tridiagonal_matrix.h>
22 #include <deal.II/lac/solver.h>
23 #include <deal.II/lac/solver_control.h>
24 #include <deal.II/base/exceptions.h>
25 #include <deal.II/base/logstream.h>
26 #include <deal.II/base/subscriptor.h>
27 #include <cmath>
28 
29 DEAL_II_NAMESPACE_OPEN
30 
31 // forward declaration
33 
34 
37 
92 template <typename VectorType = Vector<double> >
93 class SolverCG : public Solver<VectorType>
94 {
95 public:
100 
106  struct AdditionalData {};
107 
111  SolverCG (SolverControl &cn,
113  const AdditionalData &data = AdditionalData());
114 
119  SolverCG (SolverControl &cn,
120  const AdditionalData &data=AdditionalData());
121 
125  virtual ~SolverCG () = default;
126 
130  template <typename MatrixType, typename PreconditionerType>
131  void
132  solve (const MatrixType &A,
133  VectorType &x,
134  const VectorType &b,
135  const PreconditionerType &preconditioner);
136 
143  boost::signals2::connection
145  const std::function<void (double,double)> &slot);
146 
153  boost::signals2::connection
154  connect_condition_number_slot(const std::function<void (double)> &slot,
155  const bool every_iteration=false);
156 
163  boost::signals2::connection
165  const std::function<void (const std::vector<double> &)> &slot,
166  const bool every_iteration=false);
167 
168 protected:
174  virtual void print_vectors(const unsigned int step,
175  const VectorType &x,
176  const VectorType &r,
177  const VectorType &d) const;
178 
184  static void
186  const std::vector<double> &diagonal,
187  const std::vector<double> &offdiagonal,
188  const boost::signals2::signal<void (const std::vector<double> &)> &eigenvalues_signal,
189  const boost::signals2::signal<void (double)> &cond_signal);
190 
195 
199  boost::signals2::signal<void (double,double)> coefficients_signal;
200 
205  boost::signals2::signal<void (double)> condition_number_signal;
206 
211  boost::signals2::signal<void (double)> all_condition_numbers_signal;
212 
217  boost::signals2::signal<void (const std::vector<double> &)> eigenvalues_signal;
218 
223  boost::signals2::signal<void (const std::vector<double> &)> all_eigenvalues_signal;
224 };
225 
228 /*------------------------- Implementation ----------------------------*/
229 
230 #ifndef DOXYGEN
231 
232 template <typename VectorType>
235  const AdditionalData &data)
236  :
237  Solver<VectorType>(cn,mem),
238  additional_data(data)
239 {}
240 
241 
242 
243 template <typename VectorType>
245  const AdditionalData &data)
246  :
247  Solver<VectorType>(cn),
248  additional_data(data)
249 {}
250 
251 
252 
253 template <typename VectorType>
254 void
255 SolverCG<VectorType>::print_vectors(const unsigned int,
256  const VectorType &,
257  const VectorType &,
258  const VectorType &) const
259 {}
260 
261 
262 
263 template <typename VectorType>
264 inline void
266 (const std::vector<double> &diagonal,
267  const std::vector<double> &offdiagonal,
268  const boost::signals2::signal<void (const std::vector<double> &)> &eigenvalues_signal,
269  const boost::signals2::signal<void (double)> &cond_signal)
270 {
271  //Avoid computing eigenvalues unless they are needed.
272  if (!cond_signal.empty()|| !eigenvalues_signal.empty())
273  {
274  TridiagonalMatrix<double> T(diagonal.size(), true);
275  for (size_type i=0; i<diagonal.size(); ++i)
276  {
277  T(i,i) = diagonal[i];
278  if (i< diagonal.size()-1)
279  T(i,i+1) = offdiagonal[i];
280  }
281  T.compute_eigenvalues();
282  //Need two eigenvalues to estimate the condition number.
283  if (diagonal.size()>1)
284  {
285  double condition_number=T.eigenvalue(T.n()-1)/T.eigenvalue(0);
286  cond_signal(condition_number);
287  }
288  //Avoid copying the eigenvalues of T to a vector unless a signal is
289  //connected.
290  if (!eigenvalues_signal.empty())
291  {
292  std::vector<double> eigenvalues(T.n());
293  for (unsigned int j = 0; j < T.n(); ++j)
294  {
295  eigenvalues.at(j)=T.eigenvalue(j);
296  }
297  eigenvalues_signal(eigenvalues);
298  }
299  }
300 
301 }
302 
303 
304 
305 template <typename VectorType>
306 template <typename MatrixType, typename PreconditionerType>
307 void
308 SolverCG<VectorType>::solve (const MatrixType &A,
309  VectorType &x,
310  const VectorType &b,
311  const PreconditionerType &preconditioner)
312 {
314 
315  LogStream::Prefix prefix("cg");
316 
317  // Memory allocation
318  typename VectorMemory<VectorType>::Pointer g_pointer(this->memory);
319  typename VectorMemory<VectorType>::Pointer d_pointer(this->memory);
320  typename VectorMemory<VectorType>::Pointer h_pointer(this->memory);
321 
322  // define some aliases for simpler access
323  VectorType &g = *g_pointer;
324  VectorType &d = *d_pointer;
325  VectorType &h = *h_pointer;
326 
327  // Should we build the matrix for eigenvalue computations?
328  const bool do_eigenvalues = !condition_number_signal.empty()
329  ||!all_condition_numbers_signal.empty()
330  ||!eigenvalues_signal.empty()
331  ||!all_eigenvalues_signal.empty();
332 
333  // vectors used for eigenvalue
334  // computations
335  std::vector<double> diagonal;
336  std::vector<double> offdiagonal;
337 
338  int it=0;
339  double res = -std::numeric_limits<double>::max();
340 
341  double eigen_beta_alpha = 0;
342 
343  // resize the vectors, but do not set
344  // the values since they'd be overwritten
345  // soon anyway.
346  g.reinit(x, true);
347  d.reinit(x, true);
348  h.reinit(x, true);
349 
350  double gh,beta;
351 
352  // compute residual. if vector is
353  // zero, then short-circuit the
354  // full computation
355  if (!x.all_zero())
356  {
357  A.vmult(g,x);
358  g.add(-1.,b);
359  }
360  else
361  g.equ(-1.,b);
362  res = g.l2_norm();
363 
364  conv = this->iteration_status(0, res, x);
365  if (conv != SolverControl::iterate)
366  return;
367 
368  if (std::is_same<PreconditionerType,PreconditionIdentity>::value == false)
369  {
370  preconditioner.vmult(h,g);
371 
372  d.equ(-1.,h);
373 
374  gh = g*h;
375  }
376  else
377  {
378  d.equ(-1.,g);
379  gh = res*res;
380  }
381 
382  while (conv == SolverControl::iterate)
383  {
384  it++;
385  A.vmult(h,d);
386 
387  double alpha = d*h;
388  Assert(alpha != 0., ExcDivideByZero());
389  alpha = gh/alpha;
390 
391  x.add(alpha,d);
392  res = std::sqrt(g.add_and_dot(alpha, h, g));
393 
394  print_vectors(it, x, g, d);
395 
396  conv = this->iteration_status(it, res, x);
397  if (conv != SolverControl::iterate)
398  break;
399 
400  if (std::is_same<PreconditionerType,PreconditionIdentity>::value
401  == false)
402  {
403  preconditioner.vmult(h,g);
404 
405  beta = gh;
406  Assert(beta != 0., ExcDivideByZero());
407  gh = g*h;
408  beta = gh/beta;
409  d.sadd(beta,-1.,h);
410  }
411  else
412  {
413  beta = gh;
414  gh = res*res;
415  beta = gh/beta;
416  d.sadd(beta,-1.,g);
417  }
418 
419  this->coefficients_signal(alpha,beta);
420  // set up the vectors
421  // containing the diagonal
422  // and the off diagonal of
423  // the projected matrix.
424  if (do_eigenvalues)
425  {
426  diagonal.push_back(1./alpha + eigen_beta_alpha);
427  eigen_beta_alpha = beta/alpha;
428  offdiagonal.push_back(std::sqrt(beta)/alpha);
429  }
430  compute_eigs_and_cond(diagonal,offdiagonal,all_eigenvalues_signal,
431  all_condition_numbers_signal);
432  }
433 
434  compute_eigs_and_cond(diagonal,offdiagonal,eigenvalues_signal,
435  condition_number_signal);
436 
437  // in case of failure: throw exception
438  if (conv != SolverControl::success)
439  AssertThrow(false, SolverControl::NoConvergence (it, res));
440  // otherwise exit as normal
441 }
442 
443 
444 
445 template <typename VectorType>
446 boost::signals2::connection
448 (const std::function<void(double,double)> &slot)
449 {
450  return coefficients_signal.connect(slot);
451 }
452 
453 
454 
455 template <typename VectorType>
456 boost::signals2::connection
458 (const std::function<void(double)> &slot,
459  const bool every_iteration)
460 {
461  if (every_iteration)
462  {
463  return all_condition_numbers_signal.connect(slot);
464  }
465  else
466  {
467  return condition_number_signal.connect(slot);
468  }
469 }
470 
471 
472 
473 template <typename VectorType>
474 boost::signals2::connection
476 (const std::function<void (const std::vector<double> &)> &slot,
477  const bool every_iteration)
478 {
479  if (every_iteration)
480  {
481  return all_eigenvalues_signal.connect(slot);
482  }
483  else
484  {
485  return eigenvalues_signal.connect(slot);
486  }
487 }
488 
489 
490 
491 #endif // DOXYGEN
492 
493 DEAL_II_NAMESPACE_CLOSE
494 
495 #endif
boost::signals2::signal< void(double, double)> coefficients_signal
Definition: solver_cg.h:199
boost::signals2::signal< void(double)> condition_number_signal
Definition: solver_cg.h:205
SolverCG(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
Continue iteration.
virtual void print_vectors(const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const
boost::signals2::signal< void(const std::vector< double > &)> eigenvalues_signal
Definition: solver_cg.h:217
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
static ::ExceptionBase & ExcDivideByZero()
virtual ~SolverCG()=default
Matrix is diagonal.
boost::signals2::connection connect_condition_number_slot(const std::function< void(double)> &slot, const bool every_iteration=false)
unsigned int global_dof_index
Definition: types.h:88
Stop iteration, goal reached.
#define Assert(cond, exc)
Definition: exceptions.h:1142
AdditionalData additional_data
Definition: solver_cg.h:194
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
boost::signals2::connection connect_coefficients_slot(const std::function< void(double, double)> &slot)
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
Definition: solver.h:322
boost::signals2::signal< void(double)> all_condition_numbers_signal
Definition: solver_cg.h:211
boost::signals2::connection connect_eigenvalues_slot(const std::function< void(const std::vector< double > &)> &slot, const bool every_iteration=false)
boost::signals2::signal< void(const std::vector< double > &)> all_eigenvalues_signal
Definition: solver_cg.h:223
std::array< Number, 1 > eigenvalues(const SymmetricTensor< 2, 1, Number > &T)
static void compute_eigs_and_cond(const std::vector< double > &diagonal, const std::vector< double > &offdiagonal, const boost::signals2::signal< void(const std::vector< double > &)> &eigenvalues_signal, const boost::signals2::signal< void(double)> &cond_signal)
types::global_dof_index size_type
Definition: solver_cg.h:99