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
35DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
36# include <ml_MultiLevelPreconditioner.h>
37DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
41#include <exadg/solvers_and_preconditioners/multigrid/multigrid_parameters.h>
42#include <exadg/solvers_and_preconditioners/preconditioners/preconditioner_base.h>
43#include <exadg/solvers_and_preconditioners/solvers/iterative_solvers_dealii_wrapper.h>
44#include <exadg/solvers_and_preconditioners/utilities/linear_algebra_utilities.h>
45#include <exadg/utilities/print_functions.h>
49template<
int dim,
int spacedim>
50std::unique_ptr<MPI_Comm, void (*)(MPI_Comm *)>
51create_subcommunicator(dealii::DoFHandler<dim, spacedim>
const & dof_handler)
53 unsigned int n_locally_owned_cells = 0;
54 for(
auto const & cell : dof_handler.active_cell_iterators())
55 if(cell->is_locally_owned())
56 ++n_locally_owned_cells;
58 MPI_Comm
const mpi_comm = dof_handler.get_mpi_communicator();
66 if(dealii::Utilities::MPI::min(n_locally_owned_cells, mpi_comm) == 0)
68 std::unique_ptr<MPI_Comm, void (*)(MPI_Comm *)> subcommunicator(
new MPI_Comm,
73 MPI_Comm_split(mpi_comm,
74 n_locally_owned_cells > 0,
75 dealii::Utilities::MPI::this_mpi_process(mpi_comm),
76 subcommunicator.get());
78 return subcommunicator;
82 std::unique_ptr<MPI_Comm, void (*)(MPI_Comm *)> communicator(
new MPI_Comm, [](MPI_Comm * comm) {
85 *communicator = mpi_comm;
91#ifdef DEAL_II_WITH_TRILINOS
92inline Teuchos::ParameterList
93get_ML_parameter_list(dealii::TrilinosWrappers::PreconditionAMG::AdditionalData
const & ml_data,
96 Teuchos::ParameterList parameter_list;
99 if(ml_data.elliptic ==
true)
101 ML_Epetra::SetDefaults(
"SA", parameter_list);
105 ML_Epetra::SetDefaults(
"NSSA", parameter_list);
106 parameter_list.set(
"aggregation: block scaling",
true);
108 parameter_list.set(
"aggregation: type",
"Uncoupled");
109 parameter_list.set(
"smoother: type", ml_data.smoother_type);
110 parameter_list.set(
"coarse: type", ml_data.coarse_type);
114# if DEAL_II_TRILINOS_VERSION_GTE(12, 4, 0)
115 parameter_list.set(
"initialize random seed",
true);
118 parameter_list.set(
"smoother: sweeps",
static_cast<int>(ml_data.smoother_sweeps));
119 parameter_list.set(
"cycle applications",
static_cast<int>(ml_data.n_cycles));
120 if(ml_data.w_cycle ==
true)
122 parameter_list.set(
"prec type",
"MGW");
126 parameter_list.set(
"prec type",
"MGV");
129 parameter_list.set(
"smoother: Chebyshev alpha", 10.);
130 parameter_list.set(
"smoother: ifpack overlap",
static_cast<int>(ml_data.smoother_overlap));
131 parameter_list.set(
"aggregation: threshold", ml_data.aggregation_threshold);
135 parameter_list.set(
"coarse: max size", 2000);
138 parameter_list.set(
"repartition: enable", 1);
139 parameter_list.set(
"repartition: max min ratio", 1.3);
140 parameter_list.set(
"repartition: min per proc", 300);
141 parameter_list.set(
"repartition: partitioner",
"Zoltan");
142 parameter_list.set(
"repartition: Zoltan dimensions", dimension);
144 if(ml_data.output_details)
146 parameter_list.set(
"ML output", 10);
150 parameter_list.set(
"ML output", 0);
153 return parameter_list;
156template<
typename Operator>
160 typedef dealii::LinearAlgebra::distributed::Vector<double> VectorType;
162 typedef dealii::TrilinosWrappers::PreconditionAMG::AdditionalData MLData;
166 dealii::TrilinosWrappers::SparseMatrix system_matrix;
169 dealii::TrilinosWrappers::PreconditionAMG amg;
172 PreconditionerML(Operator
const & op,
bool const initialize, MLData ml_data = MLData())
173 : pde_operator(op), ml_data(ml_data)
176 pde_operator.init_system_matrix(system_matrix,
177 op.get_matrix_free().get_dof_handler().get_mpi_communicator());
186 vmult(VectorType & dst, VectorType
const & src)
const override
192 apply_krylov_solver_with_amg_preconditioner(VectorType & dst,
193 VectorType
const & src,
194 MultigridCoarseGridSolver
const & solver_type,
195 SolverData
const & solver_data)
const
197 dealii::ReductionControl solver_control(solver_data.max_iter,
199 solver_data.rel_tol);
201 if(solver_type == MultigridCoarseGridSolver::CG)
203 dealii::SolverCG<VectorType> solver(solver_control);
204 solver.solve(system_matrix, dst, src, *
this);
206 else if(solver_type == MultigridCoarseGridSolver::GMRES)
208 typename dealii::SolverGMRES<VectorType>::AdditionalData gmres_data;
209 gmres_data.max_n_tmp_vectors = solver_data.max_krylov_size;
210 gmres_data.right_preconditioning =
true;
212 dealii::SolverGMRES<VectorType> solver(solver_control, gmres_data);
213 solver.solve(system_matrix, dst, src, *
this);
217 AssertThrow(
false, dealii::ExcMessage(
"Not implemented."));
225 system_matrix *= 0.0;
228 pde_operator.calculate_system_matrix(system_matrix);
231 unsigned int const dimension = pde_operator.get_matrix_free().dimension;
232 Teuchos::ParameterList parameter_list = get_ML_parameter_list(ml_data, dimension);
237 std::vector<std::vector<bool>> constant_modes;
238 std::vector<std::vector<double>> constant_modes_values;
239 pde_operator.get_constant_modes(constant_modes, constant_modes_values);
240 MPI_Comm
const & mpi_comm = system_matrix.get_mpi_communicator();
241 bool const some_constant_mode_set =
242 dealii::Utilities::MPI::logical_or(constant_modes.size() > 0, mpi_comm);
243 if(not some_constant_mode_set)
245 bool const some_constant_mode_values_set =
246 dealii::Utilities::MPI::logical_or(constant_modes_values.size() > 0, mpi_comm);
247 AssertThrow(some_constant_mode_values_set > 0,
249 "Neither `constant_modes` nor `constant_modes_values` were provided. "
250 "AMG setup requires near null space basis vectors."));
253 if(constant_modes_values.size() > 0 and constant_modes_values[0].size() > 0)
255 ml_data.constant_modes_values = constant_modes_values;
261 if(constant_modes.size() > 0 and constant_modes[0].size() > 0)
263 ml_data.constant_modes = constant_modes;
269 std::unique_ptr<Epetra_MultiVector> ptr_operator_modes;
270 ml_data.set_operator_null_space(parameter_list,
272 system_matrix.trilinos_matrix());
275 amg.initialize(system_matrix, parameter_list);
277 this->update_needed =
false;
282 Operator
const & pde_operator;
288#ifdef DEAL_II_WITH_PETSC
292template<
typename Operator,
typename Number>
296 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
298 typedef dealii::PETScWrappers::PreconditionBoomerAMG::AdditionalData BoomerData;
302 std::unique_ptr<MPI_Comm, void (*)(MPI_Comm *)> subcommunicator;
306 dealii::PETScWrappers::MPI::SparseMatrix system_matrix;
309 dealii::PETScWrappers::PreconditionBoomerAMG amg;
311 PreconditionerBoomerAMG(Operator
const & op,
312 bool const initialize,
313 BoomerData boomer_data = BoomerData())
315 create_subcommunicator(op.get_matrix_free().get_dof_handler(op.get_dof_index()))),
317 boomer_data(boomer_data)
320 pde_operator.init_system_matrix(system_matrix, *subcommunicator);
328 ~PreconditionerBoomerAMG()
330 if(system_matrix.m() > 0)
332 PetscErrorCode ierr = VecDestroy(&petsc_vector_dst);
333 AssertThrow(ierr == 0, dealii::ExcPETScError(ierr));
334 ierr = VecDestroy(&petsc_vector_src);
335 AssertThrow(ierr == 0, dealii::ExcPETScError(ierr));
340 vmult(VectorType & dst, VectorType
const & src)
const override
342 if(system_matrix.m() > 0)
343 apply_petsc_operation(dst,
347 [&](dealii::PETScWrappers::VectorBase & petsc_dst,
348 dealii::PETScWrappers::VectorBase
const & petsc_src) {
349 amg.vmult(petsc_dst, petsc_src);
354 apply_krylov_solver_with_amg_preconditioner(VectorType & dst,
355 VectorType
const & src,
356 MultigridCoarseGridSolver
const & solver_type,
357 SolverData
const & solver_data)
const
359 apply_petsc_operation(dst,
361 system_matrix.get_mpi_communicator(),
362 [&](dealii::PETScWrappers::VectorBase & petsc_dst,
363 dealii::PETScWrappers::VectorBase
const & petsc_src) {
364 dealii::ReductionControl solver_control(solver_data.max_iter,
366 solver_data.rel_tol);
368 if(solver_type == MultigridCoarseGridSolver::CG)
370 dealii::PETScWrappers::SolverCG solver(solver_control);
371 solver.solve(system_matrix, petsc_dst, petsc_src, amg);
373 else if(solver_type == MultigridCoarseGridSolver::GMRES)
375 dealii::PETScWrappers::SolverGMRES solver(solver_control);
376 solver.solve(system_matrix, petsc_dst, petsc_src, amg);
380 AssertThrow(false, dealii::ExcMessage(
"Not implemented."));
392 if(system_matrix.m() > 0)
395 calculate_preconditioner();
397 this->update_needed =
false;
402 calculate_preconditioner()
405 if(system_matrix.m() > 0)
407 pde_operator.calculate_system_matrix(system_matrix);
409 amg.initialize(system_matrix, boomer_data);
412 dealii::LinearAlgebra::distributed::Vector<typename Operator::value_type> vector;
413 pde_operator.initialize_dof_vector(vector);
414 VecCreateMPI(system_matrix.get_mpi_communicator(),
415 vector.get_partitioner()->locally_owned_size(),
418 VecCreateMPI(system_matrix.get_mpi_communicator(),
419 vector.get_partitioner()->locally_owned_size(),
426 Operator
const & pde_operator;
428 BoomerData boomer_data;
431 mutable Vec petsc_vector_src;
432 mutable Vec petsc_vector_dst;
439template<
typename Operator,
typename Number>
440class PreconditionerAMG :
public PreconditionerBase<Number>
443 typedef typename PreconditionerBase<Number>::VectorType VectorType;
446 PreconditionerAMG(Operator
const & pde_operator,
bool const initialize,
AMGData const & data)
452 if(data.amg_type == AMGType::BoomerAMG)
454#ifdef DEAL_II_WITH_PETSC
455 preconditioner_boomer =
456 std::make_shared<PreconditionerBoomerAMG<Operator, Number>>(pde_operator,
460 AssertThrow(
false, dealii::ExcMessage(
"deal.II is not compiled with PETSc!"));
463 else if(data.amg_type == AMGType::ML)
465#ifdef DEAL_II_WITH_TRILINOS
467 std::make_shared<PreconditionerML<Operator>>(pde_operator, initialize, data.ml_data);
469 AssertThrow(
false, dealii::ExcMessage(
"deal.II is not compiled with Trilinos!"));
474 AssertThrow(
false, dealii::ExcNotImplemented());
479 vmult(VectorType & dst, VectorType
const & src)
const final
481 if(data.amg_type == AMGType::BoomerAMG)
483#ifdef DEAL_II_WITH_PETSC
484 preconditioner_boomer->vmult(dst, src);
486 AssertThrow(
false, dealii::ExcMessage(
"deal.II is not compiled with PETSc!"));
489 else if(data.amg_type == AMGType::ML)
491#ifdef DEAL_II_WITH_TRILINOS
495 [&](dealii::LinearAlgebra::distributed::Vector<double> & dst_double,
496 dealii::LinearAlgebra::distributed::Vector<double>
const & src_double) {
497 preconditioner_ml->vmult(dst_double, src_double);
500 AssertThrow(
false, dealii::ExcMessage(
"deal.II is not compiled with Trilinos!"));
505 AssertThrow(
false, dealii::ExcNotImplemented());
510 apply_krylov_solver_with_amg_preconditioner(VectorType & dst,
511 VectorType
const & src,
512 MultigridCoarseGridSolver
const & solver_type,
515 if(data.amg_type == AMGType::BoomerAMG)
517#ifdef DEAL_II_WITH_PETSC
518 std::shared_ptr<PreconditionerBoomerAMG<Operator, Number>> preconditioner =
519 std::dynamic_pointer_cast<PreconditionerBoomerAMG<Operator, Number>>(preconditioner_boomer);
521 preconditioner->apply_krylov_solver_with_amg_preconditioner(dst,
526 AssertThrow(
false, dealii::ExcMessage(
"deal.II is not compiled with PETSc!"));
529 else if(data.amg_type == AMGType::ML)
531#ifdef DEAL_II_WITH_TRILINOS
532 std::shared_ptr<PreconditionerML<Operator>> preconditioner =
533 std::dynamic_pointer_cast<PreconditionerML<Operator>>(preconditioner_ml);
538 [&](dealii::LinearAlgebra::distributed::Vector<double> & dst_double,
539 dealii::LinearAlgebra::distributed::Vector<double>
const & src_double) {
540 preconditioner->apply_krylov_solver_with_amg_preconditioner(dst_double,
546 AssertThrow(
false, dealii::ExcMessage(
"deal.II is not compiled with Trilinos!"));
551 AssertThrow(
false, dealii::ExcNotImplemented());
558 if(data.amg_type == AMGType::BoomerAMG)
560#ifdef DEAL_II_WITH_PETSC
561 preconditioner_boomer->update();
563 AssertThrow(
false, dealii::ExcMessage(
"deal.II is not compiled with PETSc!"));
566 else if(data.amg_type == AMGType::ML)
568#ifdef DEAL_II_WITH_TRILINOS
569 preconditioner_ml->update();
571 AssertThrow(
false, dealii::ExcMessage(
"deal.II is not compiled with Trilinos!"));
576 AssertThrow(
false, dealii::ExcNotImplemented());
579 this->update_needed =
false;
585 std::shared_ptr<PreconditionerBase<Number>> preconditioner_boomer;
587 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:92