Reference documentation for deal.II version 8.5.1
solver_cg.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2017 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 
98 template <typename VectorType = Vector<double> >
99 class SolverCG : public Solver<VectorType>
100 {
101 public:
106 
111  {
117 
124 
131 
138 
144  explicit
145  AdditionalData (const bool log_coefficients,
146  const bool compute_condition_number = false,
147  const bool compute_all_condition_numbers = false,
148  const bool compute_eigenvalues = false) DEAL_II_DEPRECATED;
149 
153  AdditionalData();
154  };
155 
159  SolverCG (SolverControl &cn,
160  VectorMemory<VectorType> &mem,
161  const AdditionalData &data = AdditionalData());
162 
167  SolverCG (SolverControl &cn,
168  const AdditionalData &data=AdditionalData());
169 
173  virtual ~SolverCG ();
174 
178  template <typename MatrixType, typename PreconditionerType>
179  void
180  solve (const MatrixType &A,
181  VectorType &x,
182  const VectorType &b,
183  const PreconditionerType &precondition);
184 
191  boost::signals2::connection
193  const std_cxx11::function<void (double,double)> &slot);
194 
201  boost::signals2::connection
202  connect_condition_number_slot(const std_cxx11::function<void (double)> &slot,
203  const bool every_iteration=false);
204 
211  boost::signals2::connection
213  const std_cxx11::function<void (const std::vector<double> &)> &slot,
214  const bool every_iteration=false);
215 
216 protected:
222  virtual void print_vectors(const unsigned int step,
223  const VectorType &x,
224  const VectorType &r,
225  const VectorType &d) const;
226 
234  static void
236  const std::vector<double> &diagonal,
237  const std::vector<double> &offdiagonal,
238  const boost::signals2::signal<void (const std::vector<double> &)> &eigenvalues_signal,
239  const boost::signals2::signal<void (double)> &cond_signal,
240  const bool log_eigenvalues,
241  const bool log_cond);
242 
247  VectorType *Vr;
248  VectorType *Vp;
249  VectorType *Vz;
250 
255 
259  boost::signals2::signal<void (double,double)> coefficients_signal;
260 
265  boost::signals2::signal<void (double)> condition_number_signal;
266 
271  boost::signals2::signal<void (double)> all_condition_numbers_signal;
272 
277  boost::signals2::signal<void (const std::vector<double> &)> eigenvalues_signal;
278 
283  boost::signals2::signal<void (const std::vector<double> &)> all_eigenvalues_signal;
284 
285 private:
286  void cleanup();
287 };
288 
291 /*------------------------- Implementation ----------------------------*/
292 
293 #ifndef DOXYGEN
294 
295 template <typename VectorType>
296 inline
298 AdditionalData (const bool log_coefficients,
299  const bool compute_condition_number,
301  const bool compute_eigenvalues)
302  :
307 {}
308 
309 
310 
311 template <typename VectorType>
312 inline
315  :
316  log_coefficients (false),
317  compute_condition_number(false),
318  compute_all_condition_numbers(false),
319  compute_eigenvalues(false)
320 {}
321 
322 
323 
324 template <typename VectorType>
327  const AdditionalData &data)
328  :
329  Solver<VectorType>(cn,mem),
330  Vr(NULL),
331  Vp(NULL),
332  Vz(NULL),
333  additional_data(data)
334 {}
335 
336 
337 
338 template <typename VectorType>
340  const AdditionalData &data)
341  :
342  Solver<VectorType>(cn),
343  Vr(NULL),
344  Vp(NULL),
345  Vz(NULL),
346  additional_data(data)
347 {}
348 
349 
350 
351 template <typename VectorType>
353 {}
354 
355 
356 
357 template <typename VectorType>
358 void
360 {
361  this->memory.free(Vr);
362  this->memory.free(Vp);
363  this->memory.free(Vz);
364  deallog.pop();
365 }
366 
367 
368 
369 template <typename VectorType>
370 void
371 SolverCG<VectorType>::print_vectors(const unsigned int,
372  const VectorType &,
373  const VectorType &,
374  const VectorType &) const
375 {}
376 
377 
378 
379 template <typename VectorType>
380 inline void
382 (const std::vector<double> &diagonal,
383  const std::vector<double> &offdiagonal,
384  const boost::signals2::signal<void (const std::vector<double> &)> &eigenvalues_signal,
385  const boost::signals2::signal<void (double)> &cond_signal,
386  const bool log_eigenvalues,
387  const bool log_cond)
388 {
389  //Avoid computing eigenvalues unless they are needed.
390  if (!cond_signal.empty()|| !eigenvalues_signal.empty() || log_cond ||
391  log_eigenvalues)
392  {
393  TridiagonalMatrix<double> T(diagonal.size(), true);
394  for (size_type i=0; i<diagonal.size(); ++i)
395  {
396  T(i,i) = diagonal[i];
397  if (i< diagonal.size()-1)
398  T(i,i+1) = offdiagonal[i];
399  }
400  T.compute_eigenvalues();
401  //Need two eigenvalues to estimate the condition number.
402  if (diagonal.size()>1)
403  {
404  double condition_number=T.eigenvalue(T.n()-1)/T.eigenvalue(0);
405  cond_signal(condition_number);
406  //Send to deallog
407  if (log_cond)
408  {
409  deallog << "Condition number estimate: " <<
410  condition_number << std::endl;
411  }
412  }
413  //Avoid copying the eigenvalues of T to a vector unless a signal is
414  //connected.
415  if (!eigenvalues_signal.empty())
416  {
417  std::vector<double> eigenvalues(T.n());
418  for (unsigned int j = 0; j < T.n(); ++j)
419  {
420  eigenvalues.at(j)=T.eigenvalue(j);
421  }
422  eigenvalues_signal(eigenvalues);
423  }
424  if (log_eigenvalues)
425  {
426  for (size_type i=0; i<T.n(); ++i)
427  deallog << ' ' << T.eigenvalue(i);
428  deallog << std::endl;
429  }
430  }
431 
432 }
433 
434 
435 
436 template <typename VectorType>
437 template <typename MatrixType, typename PreconditionerType>
438 void
439 SolverCG<VectorType>::solve (const MatrixType &A,
440  VectorType &x,
441  const VectorType &b,
442  const PreconditionerType &precondition)
443 {
445 
446  deallog.push("cg");
447 
448  // Memory allocation
449  Vr = this->memory.alloc();
450  Vz = this->memory.alloc();
451  Vp = this->memory.alloc();
452  // Should we build the matrix for
453  // eigenvalue computations?
454  const bool do_eigenvalues = !condition_number_signal.empty()
455  ||!all_condition_numbers_signal.empty()
456  ||!eigenvalues_signal.empty()
457  ||!all_eigenvalues_signal.empty()
458  || additional_data.compute_condition_number
459  || additional_data.compute_all_condition_numbers
460  || additional_data.compute_eigenvalues;
461 
462  // vectors used for eigenvalue
463  // computations
464  std::vector<double> diagonal;
465  std::vector<double> offdiagonal;
466 
467  int it=0;
468  double res = -std::numeric_limits<double>::max();
469 
470  try
471  {
472  double eigen_beta_alpha = 0;
473 
474  // define some aliases for simpler access
475  VectorType &g = *Vr;
476  VectorType &d = *Vz;
477  VectorType &h = *Vp;
478  // resize the vectors, but do not set
479  // the values since they'd be overwritten
480  // soon anyway.
481  g.reinit(x, true);
482  d.reinit(x, true);
483  h.reinit(x, true);
484 
485  double gh,beta;
486 
487  // compute residual. if vector is
488  // zero, then short-circuit the
489  // full computation
490  if (!x.all_zero())
491  {
492  A.vmult(g,x);
493  g.add(-1.,b);
494  }
495  else
496  g.equ(-1.,b);
497  res = g.l2_norm();
498 
499  conv = this->iteration_status(0, res, x);
500  if (conv != SolverControl::iterate)
501  {
502  cleanup();
503  return;
504  }
505 
507  {
508  precondition.vmult(h,g);
509 
510  d.equ(-1.,h);
511 
512  gh = g*h;
513  }
514  else
515  {
516  d.equ(-1.,g);
517  gh = res*res;
518  }
519 
520  while (conv == SolverControl::iterate)
521  {
522  it++;
523  A.vmult(h,d);
524 
525  double alpha = d*h;
526  Assert(alpha != 0., ExcDivideByZero());
527  alpha = gh/alpha;
528 
529  x.add(alpha,d);
530  res = std::sqrt(g.add_and_dot(alpha, h, g));
531 
532  print_vectors(it, x, g, d);
533 
534  conv = this->iteration_status(it, res, x);
535  if (conv != SolverControl::iterate)
536  break;
537 
539  == false)
540  {
541  precondition.vmult(h,g);
542 
543  beta = gh;
544  Assert(beta != 0., ExcDivideByZero());
545  gh = g*h;
546  beta = gh/beta;
547  d.sadd(beta,-1.,h);
548  }
549  else
550  {
551  beta = gh;
552  gh = res*res;
553  beta = gh/beta;
554  d.sadd(beta,-1.,g);
555  }
556 
557  this->coefficients_signal(alpha,beta);
558  if (additional_data.log_coefficients)
559  deallog << "alpha-beta:" << alpha << '\t' << beta << std::endl;
560  // set up the vectors
561  // containing the diagonal
562  // and the off diagonal of
563  // the projected matrix.
564  if (do_eigenvalues)
565  {
566  diagonal.push_back(1./alpha + eigen_beta_alpha);
567  eigen_beta_alpha = beta/alpha;
568  offdiagonal.push_back(std::sqrt(beta)/alpha);
569  }
570  compute_eigs_and_cond(diagonal,offdiagonal,all_eigenvalues_signal,
571  all_condition_numbers_signal,false,
572  additional_data.compute_all_condition_numbers);
573  }
574  }
575  catch (...)
576  {
577  cleanup();
578  throw;
579  }
580  compute_eigs_and_cond(diagonal,offdiagonal,eigenvalues_signal,
581  condition_number_signal,
582  additional_data.compute_eigenvalues,
583  (additional_data.compute_condition_number &&
584  !additional_data.compute_all_condition_numbers));
585 
586  // Deallocate Memory
587  cleanup();
588  // in case of failure: throw exception
589  if (conv != SolverControl::success)
590  AssertThrow(false, SolverControl::NoConvergence (it, res));
591  // otherwise exit as normal
592 }
593 
594 
595 
596 template<typename VectorType>
597 boost::signals2::connection
599 (const std_cxx11::function<void(double,double)> &slot)
600 {
601  return coefficients_signal.connect(slot);
602 }
603 
604 
605 
606 template<typename VectorType>
607 boost::signals2::connection
609 (const std_cxx11::function<void(double)> &slot,
610  const bool every_iteration)
611 {
612  if (every_iteration)
613  {
614  return all_condition_numbers_signal.connect(slot);
615  }
616  else
617  {
618  return condition_number_signal.connect(slot);
619  }
620 }
621 
622 
623 
624 template<typename VectorType>
625 boost::signals2::connection
627 (const std_cxx11::function<void (const std::vector<double> &)> &slot,
628  const bool every_iteration)
629 {
630  if (every_iteration)
631  {
632  return all_eigenvalues_signal.connect(slot);
633  }
634  else
635  {
636  return eigenvalues_signal.connect(slot);
637  }
638 }
639 
640 
641 
642 #endif // DOXYGEN
643 
644 DEAL_II_NAMESPACE_CLOSE
645 
646 #endif
boost::signals2::signal< void(double, double)> coefficients_signal
Definition: solver_cg.h:259
void pop()
Definition: logstream.cc:306
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &precondition)
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, const bool log_eigenvalues, const bool log_cond)
virtual ~SolverCG()
boost::signals2::signal< void(double)> condition_number_signal
Definition: solver_cg.h:265
SolverCG(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
Continue iteration.
boost::signals2::connection connect_eigenvalues_slot(const std_cxx11::function< void(const std::vector< double > &)> &slot, const bool every_iteration=false)
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:277
STL namespace.
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
static ::ExceptionBase & ExcDivideByZero()
unsigned int global_dof_index
Definition: types.h:88
Stop iteration, goal reached.
#define Assert(cond, exc)
Definition: exceptions.h:313
AdditionalData additional_data
Definition: solver_cg.h:254
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
Definition: solver.h:325
boost::signals2::signal< void(double)> all_condition_numbers_signal
Definition: solver_cg.h:271
void push(const std::string &text)
Definition: logstream.cc:294
VectorType * Vr
Definition: solver_cg.h:247
boost::signals2::connection connect_condition_number_slot(const std_cxx11::function< void(double)> &slot, const bool every_iteration=false)
boost::signals2::signal< void(const std::vector< double > &)> all_eigenvalues_signal
Definition: solver_cg.h:283
boost::signals2::connection connect_coefficients_slot(const std_cxx11::function< void(double, double)> &slot)
types::global_dof_index size_type
Definition: solver_cg.h:105