54class MultigridAlgorithm
57 MultigridAlgorithm(dealii::MGLevelObject<std::shared_ptr<MatrixType>>
const & matrix,
60 dealii::MGLevelObject<std::shared_ptr<SmootherType>>
const & smoother,
61 MPI_Comm
const & comm,
62 unsigned int const n_cycles = 1)
63 : minlevel(matrix.min_level()),
64 maxlevel(matrix.max_level()),
65 defect(minlevel, maxlevel),
66 solution(minlevel, maxlevel),
67 t(minlevel, maxlevel),
68 matrix(&matrix,
typeid(*this).name()),
69 coarse(&coarse,
typeid(*this).name()),
71 smoother(&smoother,
typeid(*this).name()),
75 AssertThrow(n_cycles == 1, dealii::ExcNotImplemented());
77 for(
unsigned int level = minlevel; level <= maxlevel; ++level)
79 matrix[level]->initialize_dof_vector(solution[level]);
80 defect[level] = solution[level];
81 t[level] = solution[level];
84 timer_tree = std::make_shared<TimerTree>();
87 template<
class OtherVectorType>
89 vmult(OtherVectorType & dst, OtherVectorType
const & src)
const
95 for(
unsigned int i = minlevel; i < maxlevel; i++)
99 defect[maxlevel].copy_locally_owned_data_from(src);
101 v_cycle(maxlevel,
false);
103 dst.copy_locally_owned_data_from(solution[maxlevel]);
106 timer_tree->insert({
"Multigrid"}, timer.wall_time());
110 template<
class OtherVectorType>
112 solve(OtherVectorType & dst, OtherVectorType
const & src)
const
114 defect[maxlevel].copy_locally_owned_data_from(src);
116 solution[maxlevel].copy_locally_owned_data_from(dst);
119 (*matrix)[maxlevel]->initialize_dof_vector(residual);
122 double const norm_r_0 = calculate_residual(residual);
123 double norm_r = norm_r_0;
125 int const max_iter = 1000;
126 double const abstol = 1.e-12;
127 double const reltol = 1.e-6;
130 bool converged = norm_r_0 < abstol;
133 for(
unsigned int i = minlevel; i < maxlevel; i++)
138 v_cycle(maxlevel,
true);
141 norm_r = calculate_residual(residual);
142 std::cout <<
"Norm of residual = " << norm_r << std::endl;
143 converged = (norm_r < abstol or norm_r / norm_r_0 < reltol or n_iter >= max_iter);
148 dst.copy_locally_owned_data_from(solution[maxlevel]);
153 template<
class OtherVectorType>
155 calculate_residual(OtherVectorType & residual)
const
157 (*matrix)[maxlevel]->vmult(residual, solution[maxlevel]);
158 residual.sadd(-1.0, 1.0, defect[maxlevel]);
159 return residual.l2_norm();
162 std::shared_ptr<TimerTree>
173 v_cycle(
unsigned int const level,
bool const multigrid_is_a_solver)
const
180 if(level == minlevel)
186 (*coarse)(level, solution[level], defect[level]);
189 timer_tree->insert({
"Multigrid",
"level " + std::to_string(level)}, timer.wall_time());
199 if(multigrid_is_a_solver)
203 (*smoother)[level]->step(solution[level], defect[level]);
211 (*smoother)[level]->vmult(solution[level], defect[level]);
215 (*matrix)[level]->vmult_interface_down(t[level], solution[level]);
216 t[level].sadd(-1.0, 1.0, defect[level]);
217 transfer.restrict_and_add(level, defect[level - 1], t[level]);
220 timer_tree->insert({
"Multigrid",
"level " + std::to_string(level)}, timer.wall_time());
224 v_cycle(level - 1,
false);
231 transfer.prolongate_and_add(level, solution[level], solution[level - 1]);
234 (*smoother)[level]->step(solution[level], defect[level]);
237 timer_tree->insert({
"Multigrid",
"level " + std::to_string(level)}, timer.wall_time());
245 unsigned int minlevel;
250 unsigned int maxlevel;
256 mutable dealii::MGLevelObject<VectorType> defect;
261 mutable dealii::MGLevelObject<VectorType> solution;
266 mutable dealii::MGLevelObject<VectorType> t;
271 dealii::SmartPointer<dealii::MGLevelObject<std::shared_ptr<MatrixType>>
const> matrix;
276 dealii::SmartPointer<dealii::MGCoarseGridBase<VectorType>
const> coarse;
286 dealii::SmartPointer<dealii::MGLevelObject<std::shared_ptr<SmootherType>>
const> smoother;
288 MPI_Comm
const mpi_comm;
290 unsigned int const n_cycles;
292 std::shared_ptr<TimerTree> timer_tree;