Reference documentation for deal.II version 9.0.0
precondition_block_base.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 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_precondition_block_base_h
17 #define dealii_precondition_block_base_h
18 
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/exceptions.h>
22 #include <deal.II/base/subscriptor.h>
23 #include <deal.II/base/smartpointer.h>
24 #include <deal.II/base/memory_consumption.h>
25 #include <deal.II/lac/householder.h>
26 #include <deal.II/lac/lapack_full_matrix.h>
27 
28 #include <vector>
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 template <typename number> class FullMatrix;
33 template <typename number> class Vector;
34 
53 template <typename number>
55 {
56 public:
61 
66  enum Inversion
67  {
81  };
82 
87  Inversion method = gauss_jordan);
88 
92  ~PreconditionBlockBase() = default;
93 
98  void clear();
99 
104  void reinit(unsigned int nblocks, size_type blocksize, bool compress,
105  Inversion method = gauss_jordan);
106 
110  void inverses_computed(bool are_they);
111 
115  bool same_diagonal () const;
116 
120  bool store_diagonals() const;
121 
125  bool inverses_ready () const;
126 
130  unsigned int size() const;
131 
135  template <typename number2>
136  void inverse_vmult(size_type i, Vector<number2> &dst, const Vector<number2> &src) const;
137 
141  template <typename number2>
142  void inverse_Tvmult(size_type i, Vector<number2> &dst, const Vector<number2> &src) const;
143 
148 
153 
158 
162  const FullMatrix<number> &inverse (size_type i) const;
163 
168 
173 
178 
182  const FullMatrix<number> &diagonal (size_type i) const;
183 
189  void log_statistics () const;
190 
195  std::size_t memory_consumption () const;
196 
202 
209 
210 protected:
215 
216 private:
220  unsigned int n_diagonal_blocks;
221 
228  std::vector<FullMatrix<number> > var_inverse_full;
229 
236  std::vector<Householder<number> > var_inverse_householder;
237 
244  std::vector<LAPACKFullMatrix<number> > var_inverse_svd;
245 
251  std::vector<FullMatrix<number> > var_diagonal;
252 
253 
258 
263 
269 };
270 
271 //----------------------------------------------------------------------//
272 
273 template <typename number>
274 inline
276  bool store, Inversion method)
277  :
278  inversion(method),
279  n_diagonal_blocks(0),
280  var_store_diagonals(store),
281  var_same_diagonal(false),
282  var_inverses_ready(false)
283 {}
284 
285 
286 template <typename number>
287 inline
288 void
290 {
291  if (var_inverse_full.size()!=0)
292  var_inverse_full.erase(var_inverse_full.begin(), var_inverse_full.end());
293  if (var_inverse_householder.size()!=0)
294  var_inverse_householder.erase(var_inverse_householder.begin(), var_inverse_householder.end());
295  if (var_inverse_svd.size()!=0)
296  var_inverse_svd.erase(var_inverse_svd.begin(), var_inverse_svd.end());
297  if (var_diagonal.size()!=0)
298  var_diagonal.erase(var_diagonal.begin(), var_diagonal.end());
299  var_same_diagonal = false;
300  var_inverses_ready = false;
301  n_diagonal_blocks = 0;
302 }
303 
304 template <typename number>
305 inline
306 void
307 PreconditionBlockBase<number>::reinit(unsigned int n, size_type b, bool compress,
308  Inversion method)
309 {
310  inversion = method;
311  var_same_diagonal = compress;
312  var_inverses_ready = false;
313  n_diagonal_blocks = n;
314 
315  if (compress)
316  {
317  switch (inversion)
318  {
319  case gauss_jordan:
320  var_inverse_full.resize(1);
321  var_inverse_full[0].reinit(b,b);
322  break;
323  case householder:
324  var_inverse_householder.resize(1);
325  break;
326  case svd:
327  var_inverse_svd.resize(1);
328  var_inverse_svd[0].reinit(b,b);
329  break;
330  default:
331  Assert(false, ExcNotImplemented());
332  }
333 
334  if (store_diagonals())
335  {
336  var_diagonal.resize(1);
337  var_diagonal[0].reinit(b,b);
338  }
339  }
340  else
341  {
342  // set the arrays to the right
343  // size. we could do it like this:
344  // var_inverse = vector<>(nblocks,FullMatrix<>())
345  // but this would involve copying many
346  // FullMatrix objects.
347  //
348  // the following is a neat trick which
349  // avoids copying
350  if (store_diagonals())
351  {
352  std::vector<FullMatrix<number> >
353  tmp(n, FullMatrix<number>(b));
354  var_diagonal.swap (tmp);
355  }
356 
357  switch (inversion)
358  {
359  case gauss_jordan:
360  {
361  std::vector<FullMatrix<number> >
362  tmp(n, FullMatrix<number>(b));
363  var_inverse_full.swap (tmp);
364  break;
365  }
366  case householder:
367  var_inverse_householder.resize(n);
368  break;
369  case svd:
370  {
371  std::vector<LAPACKFullMatrix<number> >
372  tmp(n, LAPACKFullMatrix<number>(b));
373  var_inverse_svd.swap (tmp);
374  break;
375  }
376  default:
377  Assert(false, ExcNotImplemented());
378  }
379  }
380 }
381 
382 
383 template <typename number>
384 inline
385 unsigned int
387 {
388  return n_diagonal_blocks;
389 }
390 
391 template <typename number>
392 template <typename number2>
393 inline
394 void
396  size_type i, Vector<number2> &dst, const Vector<number2> &src) const
397 {
398  const size_type ii = same_diagonal() ? 0U : i;
399 
400  switch (inversion)
401  {
402  case gauss_jordan:
403  AssertIndexRange (ii, var_inverse_full.size());
404  var_inverse_full[ii].vmult(dst, src);
405  break;
406  case householder:
407  AssertIndexRange (ii, var_inverse_householder.size());
408  var_inverse_householder[ii].vmult(dst, src);
409  break;
410  case svd:
411  AssertIndexRange (ii, var_inverse_svd.size());
412  var_inverse_svd[ii].vmult(dst, src);
413  break;
414  default:
415  Assert(false, ExcNotImplemented());
416  }
417 }
418 
419 
420 template <typename number>
421 template <typename number2>
422 inline
423 void
425  size_type i, Vector<number2> &dst, const Vector<number2> &src) const
426 {
427  const size_type ii = same_diagonal() ? 0U : i;
428 
429  switch (inversion)
430  {
431  case gauss_jordan:
432  AssertIndexRange (ii, var_inverse_full.size());
433  var_inverse_full[ii].Tvmult(dst, src);
434  break;
435  case householder:
436  AssertIndexRange (ii, var_inverse_householder.size());
437  var_inverse_householder[ii].Tvmult(dst, src);
438  break;
439  case svd:
440  AssertIndexRange (ii, var_inverse_svd.size());
441  var_inverse_svd[ii].Tvmult(dst, src);
442  break;
443  default:
444  Assert(false, ExcNotImplemented());
445  }
446 }
447 
448 
449 template <typename number>
450 inline
451 const FullMatrix<number> &
453 {
454  if (same_diagonal())
455  return var_inverse_full[0];
456 
457  Assert (i < var_inverse_full.size(), ExcIndexRange(i,0,var_inverse_full.size()));
458  return var_inverse_full[i];
459 }
460 
461 
462 template <typename number>
463 inline
464 const Householder<number> &
466 {
467  if (same_diagonal())
468  return var_inverse_householder[0];
469 
470  AssertIndexRange (i, var_inverse_householder.size());
471  return var_inverse_householder[i];
472 }
473 
474 
475 template <typename number>
476 inline
479 {
480  if (same_diagonal())
481  return var_inverse_svd[0];
482 
483  AssertIndexRange (i, var_inverse_svd.size());
484  return var_inverse_svd[i];
485 }
486 
487 
488 template <typename number>
489 inline
490 const FullMatrix<number> &
492 {
493  Assert(store_diagonals(), ExcDiagonalsNotStored());
494 
495  if (same_diagonal())
496  return var_diagonal[0];
497 
498  Assert (i < var_diagonal.size(), ExcIndexRange(i,0,var_diagonal.size()));
499  return var_diagonal[i];
500 }
501 
502 
503 template <typename number>
504 inline
507 {
508  Assert(var_inverse_full.size() != 0, ExcInverseNotAvailable());
509 
510  if (same_diagonal())
511  return var_inverse_full[0];
512 
513  Assert (i < var_inverse_full.size(), ExcIndexRange(i,0,var_inverse_full.size()));
514  return var_inverse_full[i];
515 }
516 
517 
518 template <typename number>
519 inline
522 {
523  Assert(var_inverse_householder.size() != 0, ExcInverseNotAvailable());
524 
525  if (same_diagonal())
526  return var_inverse_householder[0];
527 
528  AssertIndexRange (i, var_inverse_householder.size());
529  return var_inverse_householder[i];
530 }
531 
532 
533 template <typename number>
534 inline
537 {
538  Assert(var_inverse_svd.size() != 0, ExcInverseNotAvailable());
539 
540  if (same_diagonal())
541  return var_inverse_svd[0];
542 
543  AssertIndexRange (i, var_inverse_svd.size());
544  return var_inverse_svd[i];
545 }
546 
547 
548 template <typename number>
549 inline
552 {
553  Assert(store_diagonals(), ExcDiagonalsNotStored());
554 
555  if (same_diagonal())
556  return var_diagonal[0];
557 
558  Assert (i < var_diagonal.size(), ExcIndexRange(i,0,var_diagonal.size()));
559  return var_diagonal[i];
560 }
561 
562 
563 template <typename number>
564 inline bool
566 {
567  return var_same_diagonal;
568 }
569 
570 
571 template <typename number>
572 inline bool
574 {
575  return var_store_diagonals;
576 }
577 
578 
579 template <typename number>
580 inline void
582 {
583  var_inverses_ready = x;
584 }
585 
586 
587 template <typename number>
588 inline bool
590 {
591  return var_inverses_ready;
592 }
593 
594 
595 template <typename number>
596 inline void
598 {
599  deallog << "PreconditionBlockBase: " << size() << " blocks; ";
600 
601  if (inversion == svd)
602  {
603  unsigned int kermin = 100000000, kermax = 0;
604  double sigmin = 1.e300, sigmax= -1.e300;
605  double kappamin = 1.e300, kappamax= -1.e300;
606 
607  for (size_type b=0; b<size(); ++b)
608  {
609  const LAPACKFullMatrix<number> &matrix = inverse_svd(b);
610  size_type k=1;
611  while (k <= matrix.n_cols() && matrix.singular_value(matrix.n_cols()-k) == 0)
612  ++k;
613  const double s0 = matrix.singular_value(0);
614  const double sm = matrix.singular_value(matrix.n_cols()-k);
615  const double co = sm/s0;
616 
617  if (kermin > k) kermin = k-1;
618  if (kermax < k) kermax = k-1;
619  if (s0 < sigmin) sigmin = s0;
620  if (sm > sigmax) sigmax = sm;
621  if (co < kappamin) kappamin = co;
622  if (co > kappamax) kappamax = co;
623  }
624  deallog << "dim ker [" << kermin << ':' << kermax
625  << "] sigma [" << sigmin << ':' << sigmax
626  << "] kappa [" << kappamin << ':' << kappamax << ']' << std::endl;
627 
628  }
629  else if (inversion == householder)
630  {
631  }
632  else if (inversion == gauss_jordan)
633  {
634  }
635  else
636  {
637  Assert(false, ExcNotImplemented());
638  }
639 }
640 
641 
642 template <typename number>
643 inline
644 std::size_t
646 {
647  std::size_t mem = sizeof(*this);
648  for (size_type i=0; i<var_inverse_full.size(); ++i)
649  mem += MemoryConsumption::memory_consumption(var_inverse_full[i]);
650  for (size_type i=0; i<var_diagonal.size(); ++i)
651  mem += MemoryConsumption::memory_consumption(var_diagonal[i]);
652  return mem;
653 }
654 
655 
656 DEAL_II_NAMESPACE_CLOSE
657 
658 #endif
types::global_dof_index size_type
#define AssertIndexRange(index, range)
Definition: exceptions.h:1284
~PreconditionBlockBase()=default
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
void inverses_computed(bool are_they)
void inverse_vmult(size_type i, Vector< number2 > &dst, const Vector< number2 > &src) const
std::vector< FullMatrix< number > > var_diagonal
FullMatrix< number > & inverse(size_type i)
PreconditionBlockBase(bool store_diagonals=false, Inversion method=gauss_jordan)
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:1142
#define DeclException0(Exception0)
Definition: exceptions.h:323
std::vector< FullMatrix< number > > var_inverse_full
void reinit(unsigned int nblocks, size_type blocksize, bool compress, Inversion method=gauss_jordan)
void inverse_Tvmult(size_type i, Vector< number2 > &dst, const Vector< number2 > &src) const
static ::ExceptionBase & ExcDiagonalsNotStored()
std::vector< LAPACKFullMatrix< number > > var_inverse_svd
LAPACKFullMatrix< number > & inverse_svd(size_type i)
Householder< number > & inverse_householder(size_type i)
static ::ExceptionBase & ExcNotImplemented()
std::size_t memory_consumption() const
std::vector< Householder< number > > var_inverse_householder
FullMatrix< number > & diagonal(size_type i)
static ::ExceptionBase & ExcInverseNotAvailable()
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)