15#ifndef dealii_mg_smoother_h
16#define dealii_mg_smoother_h
48template <
typename VectorType>
138template <
typename VectorType>
148 smooth(
const unsigned int level, VectorType &u,
const VectorType &rhs)
const;
183 template <
typename RelaxationType,
typename VectorType>
204 template <
typename MatrixType2>
207 const typename RelaxationType::AdditionalData &additional_data =
208 typename RelaxationType::AdditionalData());
217 template <
typename MatrixType2,
typename DataType>
234 const VectorType &rhs)
const override;
255 const VectorType &rhs)
const override;
294template <
typename MatrixType,
typename RelaxationType,
typename VectorType>
314 template <
typename MatrixType2>
317 const typename RelaxationType::AdditionalData &additional_data =
318 typename RelaxationType::AdditionalData());
328 template <
typename MatrixType2,
typename DataType>
342 template <
typename MatrixType2,
typename DataType>
345 const DataType &additional_data,
346 const unsigned int block_row,
347 const unsigned int block_col);
358 template <
typename MatrixType2,
typename DataType>
362 const unsigned int block_row,
363 const unsigned int block_col);
375 smooth(
const unsigned int level, VectorType &u,
const VectorType &rhs)
const;
393 apply(
const unsigned int level, VectorType &u,
const VectorType &rhs)
const;
444template <
typename MatrixType,
typename PreconditionerType,
typename VectorType>
465 template <
typename MatrixType2>
468 const typename PreconditionerType::AdditionalData &
469 additional_data =
typename PreconditionerType::AdditionalData());
478 template <
typename MatrixType2>
491 template <
typename MatrixType2,
typename DataType>
506 template <
typename MatrixType2,
typename DataType>
509 const DataType &additional_data,
510 const unsigned int block_row,
511 const unsigned int block_col);
523 template <
typename MatrixType2,
typename DataType>
527 const unsigned int block_row,
528 const unsigned int block_col);
542 const VectorType &rhs)
const override;
562 const VectorType &rhs)
const override;
590template <
typename VectorType>
594 const VectorType &)
const
597template <
typename VectorType>
604template <
typename VectorType>
607 const bool symmetric,
617template <
typename VectorType>
625template <
typename VectorType>
633template <
typename VectorType>
641template <
typename VectorType>
649template <
typename VectorType>
660 template <
typename RelaxationType,
typename VectorType>
662 const unsigned int steps,
664 const bool symmetric,
670 template <
typename RelaxationType,
typename VectorType>
672 SmootherRelaxation<RelaxationType, VectorType>::clear()
678 template <
typename RelaxationType,
typename VectorType>
679 template <
typename MatrixType2>
681 SmootherRelaxation<RelaxationType, VectorType>::initialize(
683 const typename RelaxationType::AdditionalData &data)
688 this->resize(min, max);
690 for (
unsigned int i = min; i <=
max; ++i)
695 template <
typename RelaxationType,
typename VectorType>
696 template <
typename MatrixType2,
typename DataType>
698 SmootherRelaxation<RelaxationType, VectorType>::initialize(
705 this->resize(min, max);
707 for (
unsigned int i = min; i <=
max; ++i)
712 template <
typename RelaxationType,
typename VectorType>
714 SmootherRelaxation<RelaxationType, VectorType>::smooth(
715 const unsigned int level,
717 const VectorType &rhs)
const
719 unsigned int maxlevel = this->max_level();
720 unsigned int steps2 = this->steps;
723 steps2 *= (1 << (maxlevel -
level));
726 if (this->symmetric && (steps2 % 2 == 0))
731 for (
unsigned int i = 0; i < steps2; ++i)
734 (*this)[
level].Tstep(u, rhs);
736 (*
this)[
level].step(u, rhs);
743 template <
typename RelaxationType,
typename VectorType>
745 SmootherRelaxation<RelaxationType, VectorType>::apply(
746 const unsigned int level,
748 const VectorType &rhs)
const
750 unsigned int maxlevel = this->max_level();
751 unsigned int steps2 = this->steps;
754 steps2 *= (1 << (maxlevel -
level));
757 if (this->symmetric && (steps2 % 2 == 0))
763 (*this)[
level].Tvmult(u, rhs);
765 (*
this)[
level].vmult(u, rhs);
768 for (
unsigned int i = 1; i < steps2; ++i)
771 (*this)[
level].Tstep(u, rhs);
773 (*
this)[
level].step(u, rhs);
780 template <
typename RelaxationType,
typename VectorType>
782 SmootherRelaxation<RelaxationType, VectorType>::memory_consumption()
const
793template <
typename MatrixType,
typename RelaxationType,
typename VectorType>
797 const bool symmetric,
804template <
typename MatrixType,
typename RelaxationType,
typename VectorType>
808 smoothers.clear_elements();
810 unsigned int i = matrices.min_level(), max_level = matrices.max_level();
811 for (; i <= max_level; ++i)
816template <
typename MatrixType,
typename RelaxationType,
typename VectorType>
817template <
typename MatrixType2>
821 const typename RelaxationType::AdditionalData &data)
826 matrices.resize(min, max);
827 smoothers.resize(min, max);
829 for (
unsigned int i = min; i <=
max; ++i)
841template <
typename MatrixType,
typename RelaxationType,
typename VectorType>
842template <
typename MatrixType2,
typename DataType>
854 matrices.resize(min, max);
855 smoothers.resize(min, max);
857 for (
unsigned int i = min; i <=
max; ++i)
869template <
typename MatrixType,
typename RelaxationType,
typename VectorType>
870template <
typename MatrixType2,
typename DataType>
874 const DataType &data,
875 const unsigned int row,
876 const unsigned int col)
881 matrices.resize(min, max);
882 smoothers.resize(min, max);
884 for (
unsigned int i = min; i <=
max; ++i)
890 m[i].block(row, col));
891 smoothers[i].initialize(m[i].block(row, col), data);
895template <
typename MatrixType,
typename RelaxationType,
typename VectorType>
896template <
typename MatrixType2,
typename DataType>
901 const unsigned int row,
902 const unsigned int col)
910 matrices.resize(min, max);
911 smoothers.resize(min, max);
913 for (
unsigned int i = min; i <=
max; ++i)
919 m[i].block(row, col));
920 smoothers[i].initialize(m[i].block(row, col), data[i]);
925template <
typename MatrixType,
typename RelaxationType,
typename VectorType>
928 const unsigned int level,
930 const VectorType &rhs)
const
932 unsigned int maxlevel = smoothers.max_level();
933 unsigned int steps2 = this->steps;
936 steps2 *= (1 << (maxlevel -
level));
939 if (this->symmetric && (steps2 % 2 == 0))
944 for (
unsigned int i = 0; i < steps2; ++i)
947 smoothers[
level].Tstep(u, rhs);
949 smoothers[
level].step(u, rhs);
956template <
typename MatrixType,
typename RelaxationType,
typename VectorType>
959 const unsigned int level,
961 const VectorType &rhs)
const
963 unsigned int maxlevel = smoothers.max_level();
964 unsigned int steps2 = this->steps;
967 steps2 *= (1 << (maxlevel -
level));
970 if (this->symmetric && (steps2 % 2 == 0))
976 smoothers[
level].Tvmult(u, rhs);
978 smoothers[
level].vmult(u, rhs);
981 for (
unsigned int i = 1; i < steps2; ++i)
984 smoothers[
level].Tstep(u, rhs);
986 smoothers[
level].step(u, rhs);
994template <
typename MatrixType,
typename RelaxationType,
typename VectorType>
999 return sizeof(*this) + matrices.memory_consumption() +
1000 smoothers.memory_consumption() +
1001 this->vector_memory.memory_consumption();
1007template <
typename MatrixType,
typename PreconditionerType,
typename VectorType>
1010 const bool variable,
1011 const bool symmetric,
1018template <
typename MatrixType,
typename PreconditionerType,
typename VectorType>
1022 smoothers.clear_elements();
1024 unsigned int i = matrices.min_level(), max_level = matrices.max_level();
1025 for (; i <= max_level; ++i)
1031template <
typename MatrixType,
typename PreconditionerType,
typename VectorType>
1032template <
typename MatrixType2>
1036 const typename PreconditionerType::AdditionalData &data)
1041 matrices.resize(min, max);
1042 smoothers.resize(min, max);
1044 for (
unsigned int i = min; i <=
max; ++i)
1058template <
typename MatrixType,
typename PreconditionerType,
typename VectorType>
1059template <
typename MatrixType2>
1067 matrices.resize(min, max);
1068 smoothers.resize(min, max);
1070 for (
unsigned int i = min; i <=
max; ++i)
1083template <
typename MatrixType,
typename PreconditionerType,
typename VectorType>
1084template <
typename MatrixType2,
typename DataType>
1096 matrices.resize(min, max);
1097 smoothers.resize(min, max);
1099 for (
unsigned int i = min; i <=
max; ++i)
1113template <
typename MatrixType,
typename PreconditionerType,
typename VectorType>
1114template <
typename MatrixType2,
typename DataType>
1118 const DataType &data,
1119 const unsigned int row,
1120 const unsigned int col)
1125 matrices.resize(min, max);
1126 smoothers.resize(min, max);
1128 for (
unsigned int i = min; i <=
max; ++i)
1130 matrices[i] = &(m[i].block(row, col));
1131 smoothers[i].initialize(m[i].block(row, col), data);
1137template <
typename MatrixType,
typename PreconditionerType,
typename VectorType>
1138template <
typename MatrixType2,
typename DataType>
1143 const unsigned int row,
1144 const unsigned int col)
1152 matrices.resize(min, max);
1153 smoothers.resize(min, max);
1155 for (
unsigned int i = min; i <=
max; ++i)
1157 matrices[i] = &(m[i].block(row, col));
1158 smoothers[i].initialize(m[i].block(row, col), data[i]);
1164template <
typename MatrixType,
typename PreconditionerType,
typename VectorType>
1167 const unsigned int level,
1169 const VectorType &rhs)
const
1171 unsigned int maxlevel = matrices.max_level();
1172 unsigned int steps2 = this->steps;
1175 steps2 *= (1 << (maxlevel -
level));
1184 if (this->symmetric && (steps2 % 2 == 0))
1186 if (this->debug > 0)
1189 for (
unsigned int i = 0; i < steps2; ++i)
1193 if (this->debug > 0)
1195 matrices[
level].Tvmult(*r, u);
1196 r->sadd(-1., 1., rhs);
1197 if (this->debug > 2)
1198 deallog <<
' ' << r->l2_norm() <<
' ';
1199 smoothers[
level].Tvmult(*d, *r);
1200 if (this->debug > 1)
1201 deallog <<
' ' <<
d->l2_norm() <<
' ';
1205 if (this->debug > 0)
1207 matrices[
level].vmult(*r, u);
1209 if (this->debug > 2)
1210 deallog <<
' ' << r->l2_norm() <<
' ';
1211 smoothers[
level].vmult(*d, *r);
1212 if (this->debug > 1)
1213 deallog <<
' ' <<
d->l2_norm() <<
' ';
1216 if (this->symmetric)
1219 if (this->debug > 0)
1225template <
typename MatrixType,
typename PreconditionerType,
typename VectorType>
1228 const unsigned int level,
1230 const VectorType &rhs)
const
1232 unsigned int maxlevel = matrices.max_level();
1233 unsigned int steps2 = this->steps;
1236 steps2 *= (1 << (maxlevel -
level));
1239 if (this->symmetric && (steps2 % 2 == 0))
1241 if (this->debug > 0)
1245 if (this->debug > 2)
1246 deallog <<
' ' << rhs.l2_norm() <<
' ';
1247 if (this->debug > 0)
1250 smoothers[
level].Tvmult(u, rhs);
1252 smoothers[
level].vmult(u, rhs);
1253 if (this->debug > 1)
1254 deallog <<
' ' << u.l2_norm() <<
' ';
1255 if (this->symmetric)
1267 for (
unsigned int i = 1; i < steps2; ++i)
1271 if (this->debug > 0)
1273 matrices[
level].Tvmult(*r, u);
1274 r->sadd(-1., 1., rhs);
1275 if (this->debug > 2)
1276 deallog <<
' ' << r->l2_norm() <<
' ';
1277 smoothers[
level].Tvmult(*d, *r);
1278 if (this->debug > 1)
1279 deallog <<
' ' <<
d->l2_norm() <<
' ';
1283 if (this->debug > 0)
1285 matrices[
level].vmult(*r, u);
1287 if (this->debug > 2)
1288 deallog <<
' ' << r->l2_norm() <<
' ';
1289 smoothers[
level].vmult(*d, *r);
1290 if (this->debug > 1)
1291 deallog <<
' ' <<
d->l2_norm() <<
' ';
1294 if (this->symmetric)
1297 if (this->debug > 0)
1303template <
typename MatrixType,
typename PreconditionerType,
typename VectorType>
1309 smoothers.memory_consumption() +
1310 this->vector_memory.memory_consumption();
std::size_t memory_consumption() const
std::size_t memory_consumption() const
unsigned int max_level() const
unsigned int min_level() const
virtual void smooth(const unsigned int level, VectorType &u, const VectorType &rhs) const
std::size_t memory_consumption() const
MGLevelObject< PreconditionerType > smoothers
void initialize(const MGLevelObject< MatrixType2 > &matrices, const MGLevelObject< DataType > &additional_data)
void initialize_matrices(const MGLevelObject< MatrixType2 > &matrices)
virtual void smooth(const unsigned int level, VectorType &u, const VectorType &rhs) const override
void initialize(const MGLevelObject< MatrixType2 > &matrices, const typename PreconditionerType::AdditionalData &additional_data=typename PreconditionerType::AdditionalData())
virtual void apply(const unsigned int level, VectorType &u, const VectorType &rhs) const override
void initialize(const MGLevelObject< MatrixType2 > &matrices, const DataType &additional_data, const unsigned int block_row, const unsigned int block_col)
MGLevelObject< LinearOperator< VectorType > > matrices
MGSmootherPrecondition(const unsigned int steps=1, const bool variable=false, const bool symmetric=false, const bool transpose=false)
void initialize(const MGLevelObject< MatrixType2 > &matrices, const MGLevelObject< DataType > &additional_data, const unsigned int block_row, const unsigned int block_col)
virtual void smooth(const unsigned int level, VectorType &u, const VectorType &rhs) const
MGLevelObject< LinearOperator< VectorType > > matrices
void initialize(const MGLevelObject< MatrixType2 > &matrices, const DataType &additional_data, const unsigned int block_row, const unsigned int block_col)
void initialize(const MGLevelObject< MatrixType2 > &matrices, const MGLevelObject< DataType > &additional_data)
void initialize(const MGLevelObject< MatrixType2 > &matrices, const MGLevelObject< DataType > &additional_data, const unsigned int block_row, const unsigned int block_col)
std::size_t memory_consumption() const
virtual void apply(const unsigned int level, VectorType &u, const VectorType &rhs) const
MGLevelObject< RelaxationType > smoothers
MGSmootherRelaxation(const unsigned int steps=1, const bool variable=false, const bool symmetric=false, const bool transpose=false)
void initialize(const MGLevelObject< MatrixType2 > &matrices, const typename RelaxationType::AdditionalData &additional_data=typename RelaxationType::AdditionalData())
GrowingVectorMemory< VectorType > vector_memory
MGSmoother(const unsigned int steps=1, const bool variable=false, const bool symmetric=false, const bool transpose=false)
void set_debug(const unsigned int level)
void set_steps(const unsigned int)
void set_symmetric(const bool)
void set_transpose(const bool)
void set_variable(const bool)
std::size_t memory_consumption() const
void initialize(const MGLevelObject< MatrixType2 > &matrices, const MGLevelObject< DataType > &additional_data)
SmootherRelaxation(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 override
virtual void apply(const unsigned int level, VectorType &u, const VectorType &rhs) const override
void initialize(const MGLevelObject< MatrixType2 > &matrices, const typename RelaxationType::AdditionalData &additional_data=typename RelaxationType::AdditionalData())
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
LinearOperator< Range, Domain, Payload > linear_operator(const OperatorExemplar &, const Matrix &)
@ symmetric
Matrix is symmetric.
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
T & get_underlying_value(T &p)
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)