22#ifndef EXADG_SOLVERS_AND_PRECONDITIONERS_PRECONDITIONERS_PRECONDITIONER_AMG_H_
23#define EXADG_SOLVERS_AND_PRECONDITIONERS_PRECONDITIONERS_PRECONDITIONER_AMG_H_
26#include <deal.II/lac/la_parallel_vector.h>
27#include <deal.II/lac/petsc_precondition.h>
28#include <deal.II/lac/petsc_solver.h>
29#include <deal.II/lac/petsc_sparse_matrix.h>
30#include <deal.II/lac/trilinos_precondition.h>
31#include <deal.II/lac/trilinos_sparse_matrix.h>
34#ifdef DEAL_II_WITH_TRILINOS
35# include <ml_MultiLevelPreconditioner.h>
39#include <exadg/solvers_and_preconditioners/multigrid/multigrid_parameters.h>
40#include <exadg/solvers_and_preconditioners/preconditioners/preconditioner_base.h>
41#include <exadg/solvers_and_preconditioners/solvers/iterative_solvers_dealii_wrapper.h>
42#include <exadg/solvers_and_preconditioners/utilities/linear_algebra_utilities.h>
43#include <exadg/utilities/print_functions.h>
47template<
int dim,
int spacedim>
48std::unique_ptr<MPI_Comm, void (*)(MPI_Comm *)>
49create_subcommunicator(dealii::DoFHandler<dim, spacedim>
const & dof_handler)
51 unsigned int n_locally_owned_cells = 0;
52 for(
auto const & cell : dof_handler.active_cell_iterators())
53 if(cell->is_locally_owned())
54 ++n_locally_owned_cells;
56 MPI_Comm
const mpi_comm = dof_handler.get_mpi_communicator();
64 if(dealii::Utilities::MPI::min(n_locally_owned_cells, mpi_comm) == 0)
66 std::unique_ptr<MPI_Comm, void (*)(MPI_Comm *)> subcommunicator(
new MPI_Comm,
71 MPI_Comm_split(mpi_comm,
72 n_locally_owned_cells > 0,
73 dealii::Utilities::MPI::this_mpi_process(mpi_comm),
74 subcommunicator.get());
76 return subcommunicator;
80 std::unique_ptr<MPI_Comm, void (*)(MPI_Comm *)> communicator(
new MPI_Comm, [](MPI_Comm * comm) {
83 *communicator = mpi_comm;
89#ifdef DEAL_II_WITH_TRILINOS
90inline Teuchos::ParameterList
91get_ML_parameter_list(dealii::TrilinosWrappers::PreconditionAMG::AdditionalData
const & ml_data,
94 Teuchos::ParameterList parameter_list;
97 if(ml_data.elliptic ==
true)
99 ML_Epetra::SetDefaults(
"SA", parameter_list);
103 ML_Epetra::SetDefaults(
"NSSA", parameter_list);
104 parameter_list.set(
"aggregation: block scaling",
true);
106 parameter_list.set(
"aggregation: type",
"Uncoupled");
107 parameter_list.set(
"smoother: type", ml_data.smoother_type);
108 parameter_list.set(
"coarse: type", ml_data.coarse_type);
112# if DEAL_II_TRILINOS_VERSION_GTE(12, 4, 0)
113 parameter_list.set(
"initialize random seed",
true);
116 parameter_list.set(
"smoother: sweeps",
static_cast<int>(ml_data.smoother_sweeps));
117 parameter_list.set(
"cycle applications",
static_cast<int>(ml_data.n_cycles));
118 if(ml_data.w_cycle ==
true)
120 parameter_list.set(
"prec type",
"MGW");
124 parameter_list.set(
"prec type",
"MGV");
127 parameter_list.set(
"smoother: Chebyshev alpha", 10.);
128 parameter_list.set(
"smoother: ifpack overlap",
static_cast<int>(ml_data.smoother_overlap));
129 parameter_list.set(
"aggregation: threshold", ml_data.aggregation_threshold);
133 parameter_list.set(
"coarse: max size", 2000);
136 parameter_list.set(
"repartition: enable", 1);
137 parameter_list.set(
"repartition: max min ratio", 1.3);
138 parameter_list.set(
"repartition: min per proc", 300);
139 parameter_list.set(
"repartition: partitioner",
"Zoltan");
140 parameter_list.set(
"repartition: Zoltan dimensions", dimension);
142 if(ml_data.output_details)
144 parameter_list.set(
"ML output", 10);
148 parameter_list.set(
"ML output", 0);
151 return parameter_list;
154template<
typename Operator>
158 typedef dealii::LinearAlgebra::distributed::Vector<double> VectorType;
160 typedef dealii::TrilinosWrappers::PreconditionAMG::AdditionalData MLData;
164 dealii::TrilinosWrappers::SparseMatrix system_matrix;
167 dealii::TrilinosWrappers::PreconditionAMG amg;
170 PreconditionerML(Operator
const & op,
bool const initialize, MLData ml_data = MLData())
171 : pde_operator(op), ml_data(ml_data)
174 pde_operator.init_system_matrix(system_matrix,
175 op.get_matrix_free().get_dof_handler().get_mpi_communicator());
184 vmult(VectorType & dst, VectorType
const & src)
const override
190 apply_krylov_solver_with_amg_preconditioner(VectorType & dst,
191 VectorType
const & src,
192 MultigridCoarseGridSolver
const & solver_type,
193 SolverData
const & solver_data)
const
195 dealii::ReductionControl solver_control(solver_data.max_iter,
197 solver_data.rel_tol);
199 if(solver_type == MultigridCoarseGridSolver::CG)
201 dealii::SolverCG<VectorType> solver(solver_control);
202 solver.solve(system_matrix, dst, src, *
this);
204 else if(solver_type == MultigridCoarseGridSolver::GMRES)
206 typename dealii::SolverGMRES<VectorType>::AdditionalData gmres_data;
207 gmres_data.max_n_tmp_vectors = solver_data.max_krylov_size;
208 gmres_data.right_preconditioning =
true;
210 dealii::SolverGMRES<VectorType> solver(solver_control, gmres_data);
211 solver.solve(system_matrix, dst, src, *
this);
215 AssertThrow(
false, dealii::ExcMessage(
"Not implemented."));
223 system_matrix *= 0.0;
226 pde_operator.calculate_system_matrix(system_matrix);
229 unsigned int const dimension = pde_operator.get_matrix_free().dimension;
230 Teuchos::ParameterList parameter_list = get_ML_parameter_list(ml_data, dimension);
235 std::vector<std::vector<bool>> constant_modes;
236 std::vector<std::vector<double>> constant_modes_values;
237 pde_operator.get_constant_modes(constant_modes, constant_modes_values);
238 MPI_Comm
const & mpi_comm = system_matrix.get_mpi_communicator();
239 bool const some_constant_mode_set =
240 dealii::Utilities::MPI::logical_or(constant_modes.size() > 0, mpi_comm);
241 if(not some_constant_mode_set)
243 bool const some_constant_mode_values_set =
244 dealii::Utilities::MPI::logical_or(constant_modes_values.size() > 0, mpi_comm);
245 AssertThrow(some_constant_mode_values_set > 0,
247 "Neither `constant_modes` nor `constant_modes_values` were provided. "
248 "AMG setup requires near null space basis vectors."));
251 if(constant_modes_values.size() > 0 and constant_modes_values[0].size() > 0)
253 ml_data.constant_modes_values = constant_modes_values;
259 if(constant_modes.size() > 0 and constant_modes[0].size() > 0)
261 ml_data.constant_modes = constant_modes;
267 std::unique_ptr<Epetra_MultiVector> ptr_operator_modes;
268 ml_data.set_operator_null_space(parameter_list,
270 system_matrix.trilinos_matrix());
273 amg.initialize(system_matrix, parameter_list);
275 this->update_needed =
false;
280 Operator
const & pde_operator;
286#ifdef DEAL_II_WITH_PETSC
290template<
typename Operator,
typename Number>
294 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
296 typedef dealii::PETScWrappers::PreconditionBoomerAMG::AdditionalData BoomerData;
300 std::unique_ptr<MPI_Comm, void (*)(MPI_Comm *)> subcommunicator;
304 dealii::PETScWrappers::MPI::SparseMatrix system_matrix;
307 dealii::PETScWrappers::PreconditionBoomerAMG amg;
309 PreconditionerBoomerAMG(Operator
const & op,
310 bool const initialize,
311 BoomerData boomer_data = BoomerData())
313 create_subcommunicator(op.get_matrix_free().get_dof_handler(op.get_dof_index()))),
315 boomer_data(boomer_data)
318 pde_operator.init_system_matrix(system_matrix, *subcommunicator);
326 ~PreconditionerBoomerAMG()
328 if(system_matrix.m() > 0)
330 PetscErrorCode ierr = VecDestroy(&petsc_vector_dst);
331 AssertThrow(ierr == 0, dealii::ExcPETScError(ierr));
332 ierr = VecDestroy(&petsc_vector_src);
333 AssertThrow(ierr == 0, dealii::ExcPETScError(ierr));
338 vmult(VectorType & dst, VectorType
const & src)
const override
340 if(system_matrix.m() > 0)
341 apply_petsc_operation(dst,
345 [&](dealii::PETScWrappers::VectorBase & petsc_dst,
346 dealii::PETScWrappers::VectorBase
const & petsc_src) {
347 amg.vmult(petsc_dst, petsc_src);
352 apply_krylov_solver_with_amg_preconditioner(VectorType & dst,
353 VectorType
const & src,
354 MultigridCoarseGridSolver
const & solver_type,
355 SolverData
const & solver_data)
const
357 apply_petsc_operation(dst,
359 system_matrix.get_mpi_communicator(),
360 [&](dealii::PETScWrappers::VectorBase & petsc_dst,
361 dealii::PETScWrappers::VectorBase
const & petsc_src) {
362 dealii::ReductionControl solver_control(solver_data.max_iter,
364 solver_data.rel_tol);
366 if(solver_type == MultigridCoarseGridSolver::CG)
368 dealii::PETScWrappers::SolverCG solver(solver_control);
369 solver.solve(system_matrix, petsc_dst, petsc_src, amg);
371 else if(solver_type == MultigridCoarseGridSolver::GMRES)
373 dealii::PETScWrappers::SolverGMRES solver(solver_control);
374 solver.solve(system_matrix, petsc_dst, petsc_src, amg);
378 AssertThrow(false, dealii::ExcMessage(
"Not implemented."));
390 if(system_matrix.m() > 0)
393 calculate_preconditioner();
395 this->update_needed =
false;
400 calculate_preconditioner()
403 if(system_matrix.m() > 0)
405 pde_operator.calculate_system_matrix(system_matrix);
407 amg.initialize(system_matrix, boomer_data);
410 dealii::LinearAlgebra::distributed::Vector<typename Operator::value_type> vector;
411 pde_operator.initialize_dof_vector(vector);
412 VecCreateMPI(system_matrix.get_mpi_communicator(),
413 vector.get_partitioner()->locally_owned_size(),
416 VecCreateMPI(system_matrix.get_mpi_communicator(),
417 vector.get_partitioner()->locally_owned_size(),
424 Operator
const & pde_operator;
426 BoomerData boomer_data;
429 mutable Vec petsc_vector_src;
430 mutable Vec petsc_vector_dst;
437template<
typename Operator,
typename Number>
438class PreconditionerAMG :
public PreconditionerBase<Number>
441 typedef typename PreconditionerBase<Number>::VectorType VectorType;
444 PreconditionerAMG(Operator
const & pde_operator,
bool const initialize,
AMGData const & data)
450 if(data.amg_type == AMGType::BoomerAMG)
452#ifdef DEAL_II_WITH_PETSC
453 preconditioner_boomer =
454 std::make_shared<PreconditionerBoomerAMG<Operator, Number>>(pde_operator,
458 AssertThrow(
false, dealii::ExcMessage(
"deal.II is not compiled with PETSc!"));
461 else if(data.amg_type == AMGType::ML)
463#ifdef DEAL_II_WITH_TRILINOS
465 std::make_shared<PreconditionerML<Operator>>(pde_operator, initialize, data.ml_data);
467 AssertThrow(
false, dealii::ExcMessage(
"deal.II is not compiled with Trilinos!"));
472 AssertThrow(
false, dealii::ExcNotImplemented());
477 vmult(VectorType & dst, VectorType
const & src)
const final
479 if(data.amg_type == AMGType::BoomerAMG)
481#ifdef DEAL_II_WITH_PETSC
482 preconditioner_boomer->vmult(dst, src);
484 AssertThrow(
false, dealii::ExcMessage(
"deal.II is not compiled with PETSc!"));
487 else if(data.amg_type == AMGType::ML)
489#ifdef DEAL_II_WITH_TRILINOS
493 [&](dealii::LinearAlgebra::distributed::Vector<double> & dst_double,
494 dealii::LinearAlgebra::distributed::Vector<double>
const & src_double) {
495 preconditioner_ml->vmult(dst_double, src_double);
498 AssertThrow(
false, dealii::ExcMessage(
"deal.II is not compiled with Trilinos!"));
503 AssertThrow(
false, dealii::ExcNotImplemented());
508 apply_krylov_solver_with_amg_preconditioner(VectorType & dst,
509 VectorType
const & src,
510 MultigridCoarseGridSolver
const & solver_type,
513 if(data.amg_type == AMGType::BoomerAMG)
515#ifdef DEAL_II_WITH_PETSC
516 std::shared_ptr<PreconditionerBoomerAMG<Operator, Number>> preconditioner =
517 std::dynamic_pointer_cast<PreconditionerBoomerAMG<Operator, Number>>(preconditioner_boomer);
519 preconditioner->apply_krylov_solver_with_amg_preconditioner(dst,
524 AssertThrow(
false, dealii::ExcMessage(
"deal.II is not compiled with PETSc!"));
527 else if(data.amg_type == AMGType::ML)
529#ifdef DEAL_II_WITH_TRILINOS
530 std::shared_ptr<PreconditionerML<Operator>> preconditioner =
531 std::dynamic_pointer_cast<PreconditionerML<Operator>>(preconditioner_ml);
536 [&](dealii::LinearAlgebra::distributed::Vector<double> & dst_double,
537 dealii::LinearAlgebra::distributed::Vector<double>
const & src_double) {
538 preconditioner->apply_krylov_solver_with_amg_preconditioner(dst_double,
544 AssertThrow(
false, dealii::ExcMessage(
"deal.II is not compiled with Trilinos!"));
549 AssertThrow(
false, dealii::ExcNotImplemented());
556 if(data.amg_type == AMGType::BoomerAMG)
558#ifdef DEAL_II_WITH_PETSC
559 preconditioner_boomer->update();
561 AssertThrow(
false, dealii::ExcMessage(
"deal.II is not compiled with PETSc!"));
564 else if(data.amg_type == AMGType::ML)
566#ifdef DEAL_II_WITH_TRILINOS
567 preconditioner_ml->update();
569 AssertThrow(
false, dealii::ExcMessage(
"deal.II is not compiled with Trilinos!"));
574 AssertThrow(
false, dealii::ExcNotImplemented());
577 this->update_needed =
false;
583 std::shared_ptr<PreconditionerBase<Number>> preconditioner_boomer;
585 std::shared_ptr<PreconditionerBase<double>> preconditioner_ml;
Definition preconditioner_base.h:35
void apply_function_in_double_precision(dealii::LinearAlgebra::distributed::Vector< Number > &dst, dealii::LinearAlgebra::distributed::Vector< Number > const &src, std::function< void(dealii::LinearAlgebra::distributed::Vector< double > &, dealii::LinearAlgebra::distributed::Vector< double > const &)> operation)
Definition linear_algebra_utilities.h:141
Definition multigrid_parameters.h:98
Definition solver_data.h:34