16 #ifndef dealii__mg_smoother_h 17 #define dealii__mg_smoother_h 20 #include <deal.II/base/config.h> 21 #include <deal.II/base/smartpointer.h> 22 #include <deal.II/lac/pointer_matrix.h> 23 #include <deal.II/lac/vector_memory.h> 24 #include <deal.II/multigrid/mg_base.h> 25 #include <deal.II/base/mg_level_object.h> 28 DEAL_II_NAMESPACE_OPEN
45 template <
typename VectorType>
82 void set_debug (
const unsigned int level);
132 template <
typename VectorType>
141 virtual void smooth (
const unsigned int level,
143 const VectorType &rhs)
const;
144 virtual void clear ();
179 template<
class RelaxationType,
typename VectorType>
199 template <
typename MatrixType2>
201 const typename RelaxationType::AdditionalData &additional_data
202 =
typename RelaxationType::AdditionalData());
211 template <
typename MatrixType2,
class DATA>
223 virtual void smooth (
const unsigned int level,
225 const VectorType &rhs)
const;
265 template<
typename MatrixType,
class RelaxationType,
typename VectorType>
285 template <
typename MatrixType2>
287 const typename RelaxationType::AdditionalData &additional_data
288 =
typename RelaxationType::AdditionalData());
298 template <
typename MatrixType2,
class DATA>
311 template <
typename MatrixType2,
class DATA>
313 const DATA &additional_data,
314 const unsigned int block_row,
315 const unsigned int block_col);
326 template <
typename MatrixType2,
class DATA>
329 const unsigned int block_row,
330 const unsigned int block_col);
340 virtual void smooth (
const unsigned int level,
342 const VectorType &rhs)
const;
395 template<
typename MatrixType,
typename PreconditionerType,
typename VectorType>
416 template <
typename MatrixType2>
418 const typename PreconditionerType::AdditionalData &additional_data =
typename PreconditionerType::AdditionalData());
429 template <
typename MatrixType2,
class DATA>
443 template <
typename MatrixType2,
class DATA>
445 const DATA &additional_data,
446 const unsigned int block_row,
447 const unsigned int block_col);
459 template <
typename MatrixType2,
class DATA>
462 const unsigned int block_row,
463 const unsigned int block_col);
473 virtual void smooth (
const unsigned int level,
475 const VectorType &rhs)
const;
502 template <
typename VectorType>
506 const VectorType &)
const 509 template <
typename VectorType>
516 template <
typename VectorType>
520 const bool symmetric,
525 symmetric(symmetric),
531 template <
typename VectorType>
539 template <
typename VectorType>
547 template <
typename VectorType>
555 template <
typename VectorType>
563 template <
typename VectorType>
574 template <
class RelaxationType,
typename VectorType>
577 (
const unsigned int steps,
579 const bool symmetric,
580 const bool transpose)
581 :
MGSmoother<VectorType>(steps, variable, symmetric, transpose)
585 template <
class RelaxationType,
typename VectorType>
587 SmootherRelaxation<RelaxationType, VectorType>::clear ()
593 template <
class RelaxationType,
typename VectorType>
594 template <
typename MatrixType2>
596 SmootherRelaxation<RelaxationType, VectorType>::initialize
598 const typename RelaxationType::AdditionalData &data)
605 for (
unsigned int i=
min; i<=
max; ++i)
606 (*
this)[i].initialize(m[i], data);
610 template <
class RelaxationType,
typename VectorType>
611 template <
typename MatrixType2,
class DATA>
613 SmootherRelaxation<RelaxationType, VectorType>::initialize
622 for (
unsigned int i=
min; i<=
max; ++i)
623 (*
this)[i].initialize(m[i], data[i]);
627 template <
class RelaxationType,
typename VectorType>
629 SmootherRelaxation<RelaxationType, VectorType>::smooth (
const unsigned int level,
631 const VectorType &rhs)
const 633 unsigned int maxlevel = this->max_level();
634 unsigned int steps2 = this->steps;
637 steps2 *= (1<<(maxlevel-level));
640 if (this->symmetric && (steps2 % 2 == 0))
643 deallog <<
'S' << level <<
' ';
645 for (
unsigned int i=0; i<steps2; ++i)
648 (*this)[level].Tstep(u, rhs);
650 (*
this)[level].step(u, rhs);
657 template <
class RelaxationType,
typename VectorType>
660 SmootherRelaxation<RelaxationType, VectorType>::
661 memory_consumption ()
const 673 template <
typename MatrixType,
class RelaxationType,
typename VectorType>
676 (
const unsigned int steps,
678 const bool symmetric,
686 template <
typename MatrixType,
class RelaxationType,
typename VectorType>
690 smoothers.clear_elements();
692 unsigned int i=matrices.min_level(),
693 max_level=matrices.max_level();
694 for (; i<=max_level; ++i)
699 template <
typename MatrixType,
class RelaxationType,
typename VectorType>
700 template <
typename MatrixType2>
704 const typename RelaxationType::AdditionalData &data)
709 matrices.resize(
min,
max);
710 smoothers.resize(
min,
max);
712 for (
unsigned int i=
min; i<=
max; ++i)
715 smoothers[i].initialize(m[i], data);
719 template <
typename MatrixType,
class RelaxationType,
typename VectorType>
720 template <
typename MatrixType2,
class DATA>
734 matrices.resize(
min,
max);
735 smoothers.resize(
min,
max);
737 for (
unsigned int i=
min; i<=
max; ++i)
740 smoothers[i].initialize(m[i], data[i]);
744 template <
typename MatrixType,
class RelaxationType,
typename VectorType>
745 template <
typename MatrixType2,
class DATA>
750 const unsigned int row,
751 const unsigned int col)
756 matrices.resize(
min,
max);
757 smoothers.resize(
min,
max);
759 for (
unsigned int i=
min; i<=
max; ++i)
761 matrices[i] = &(m[i].block(row, col));
762 smoothers[i].initialize(m[i].block(row, col), data);
766 template <
typename MatrixType,
class RelaxationType,
typename VectorType>
767 template <
typename MatrixType2,
class DATA>
772 const unsigned int row,
773 const unsigned int col)
783 matrices.resize(
min,
max);
784 smoothers.resize(
min,
max);
786 for (
unsigned int i=
min; i<=
max; ++i)
788 matrices[i] = &(m[i].block(row, col));
789 smoothers[i].initialize(m[i].block(row, col), data[i]);
794 template <
typename MatrixType,
class RelaxationType,
typename VectorType>
798 const VectorType &rhs)
const 800 unsigned int maxlevel = smoothers.max_level();
801 unsigned int steps2 = this->steps;
804 steps2 *= (1<<(maxlevel-level));
807 if (this->symmetric && (steps2 % 2 == 0))
810 deallog <<
'S' << level <<
' ';
812 for (
unsigned int i=0; i<steps2; ++i)
815 smoothers[level].Tstep(u, rhs);
817 smoothers[level].step(u, rhs);
825 template <
typename MatrixType,
class RelaxationType,
typename VectorType>
832 + matrices.memory_consumption()
833 + smoothers.memory_consumption()
834 + this->vector_memory.memory_consumption();
840 template <
typename MatrixType,
typename PreconditionerType,
typename VectorType>
843 (
const unsigned int steps,
845 const bool symmetric,
853 template <
typename MatrixType,
typename PreconditionerType,
typename VectorType>
857 smoothers.clear_elements();
859 unsigned int i=matrices.min_level(),
860 max_level=matrices.max_level();
861 for (; i<=max_level; ++i)
867 template <
typename MatrixType,
typename PreconditionerType,
typename VectorType>
868 template <
typename MatrixType2>
872 const typename PreconditionerType::AdditionalData &data)
877 matrices.resize(
min,
max);
878 smoothers.resize(
min,
max);
880 for (
unsigned int i=
min; i<=
max; ++i)
883 smoothers[i].initialize(m[i], data);
889 template <
typename MatrixType,
typename PreconditionerType,
typename VectorType>
890 template <
typename MatrixType2,
class DATA>
904 matrices.resize(
min,
max);
905 smoothers.resize(
min,
max);
907 for (
unsigned int i=
min; i<=
max; ++i)
910 smoothers[i].initialize(m[i], data[i]);
916 template <
typename MatrixType,
typename PreconditionerType,
typename VectorType>
917 template <
typename MatrixType2,
class DATA>
922 const unsigned int row,
923 const unsigned int col)
928 matrices.resize(
min,
max);
929 smoothers.resize(
min,
max);
931 for (
unsigned int i=
min; i<=
max; ++i)
933 matrices[i] = &(m[i].block(row, col));
934 smoothers[i].initialize(m[i].block(row, col), data);
940 template <
typename MatrixType,
typename PreconditionerType,
typename VectorType>
941 template <
typename MatrixType2,
class DATA>
946 const unsigned int row,
947 const unsigned int col)
957 matrices.resize(
min,
max);
958 smoothers.resize(
min,
max);
960 for (
unsigned int i=
min; i<=
max; ++i)
962 matrices[i] = &(m[i].block(row, col));
963 smoothers[i].initialize(m[i].block(row, col), data[i]);
969 template <
typename MatrixType,
typename PreconditionerType,
typename VectorType>
972 (
const unsigned int level,
974 const VectorType &rhs)
const 976 unsigned int maxlevel = matrices.max_level();
977 unsigned int steps2 = this->steps;
980 steps2 *= (1<<(maxlevel-level));
989 if (this->symmetric && (steps2 % 2 == 0))
992 deallog <<
'S' << level <<
' ';
994 for (
unsigned int i=0; i<steps2; ++i)
1000 if (i == 0 && u.all_zero())
1004 matrices[level].Tvmult(*r,u);
1005 r->sadd(-1.,1.,rhs);
1007 if (this->debug > 2)
1008 deallog <<
' ' << r->l2_norm() <<
' ';
1009 smoothers[level].Tvmult(*d, *r);
1010 if (this->debug > 1)
1011 deallog <<
' ' <<
d->l2_norm() <<
' ';
1015 if (this->debug > 0)
1017 if (i == 0 && u.all_zero())
1021 matrices[level].vmult(*r,u);
1024 if (this->debug > 2)
1025 deallog <<
' ' << r->l2_norm() <<
' ';
1026 smoothers[level].vmult(*d, *r);
1027 if (this->debug > 1)
1028 deallog <<
' ' <<
d->l2_norm() <<
' ';
1031 if (this->symmetric)
1034 if (this->debug > 0)
1035 deallog << std::endl;
1040 template <
typename MatrixType,
typename PreconditionerType,
typename VectorType>
1046 return sizeof(*this)
1047 + matrices.memory_consumption()
1048 + smoothers.memory_consumption()
1049 + this->vector_memory.memory_consumption();
1055 DEAL_II_NAMESPACE_CLOSE
std::size_t memory_consumption() const
unsigned int max_level() const
MGSmootherPrecondition(const unsigned int steps=1, const bool variable=false, const bool symmetric=false, const bool transpose=false)
void set_transpose(const bool)
virtual void smooth(const unsigned int level, VectorType &u, const VectorType &rhs) const
void set_symmetric(const bool)
SmootherRelaxation(const unsigned int steps=1, const bool variable=false, const bool symmetric=false, const bool transpose=false)
void set_variable(const bool)
MGSmootherRelaxation(const unsigned int steps=1, const bool variable=false, const bool symmetric=false, const bool transpose=false)
virtual void smooth(const unsigned int level, VectorType &u, const VectorType &rhs) const
std::size_t memory_consumption() const
MGLevelObject< RelaxationType > smoothers
unsigned int min_level() const
void set_debug(const unsigned int level)
MGLevelObject< PreconditionerType > smoothers
VectorizedArray< Number > min(const ::VectorizedArray< Number > &x, const ::VectorizedArray< Number > &y)
std::size_t memory_consumption() const
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
MGLevelObject< PointerMatrix< MatrixType, VectorType > > matrices
MGSmoother(const unsigned int steps=1, const bool variable=false, const bool symmetric=false, const bool transpose=false)
GrowingVectorMemory< VectorType > vector_memory
void initialize(const MGLevelObject< MatrixType2 > &matrices, const typename RelaxationType::AdditionalData &additional_data=typename RelaxationType::AdditionalData())
virtual void smooth(const unsigned int level, VectorType &u, const VectorType &rhs) const
std::size_t memory_consumption() const
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
void set_steps(const unsigned int)
MGLevelObject< PointerMatrix< MatrixType, VectorType > > matrices
void initialize(const MGLevelObject< MatrixType2 > &matrices, const typename RelaxationType::AdditionalData &additional_data=typename RelaxationType::AdditionalData())
VectorizedArray< Number > max(const ::VectorizedArray< Number > &x, const ::VectorizedArray< Number > &y)
void initialize(const MGLevelObject< MatrixType2 > &matrices, const typename PreconditionerType::AdditionalData &additional_data=typename PreconditionerType::AdditionalData())
virtual void smooth(const unsigned int level, VectorType &u, const VectorType &rhs) const