Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Allow to iterate intermediate GMG levels #1211

Merged
merged 4 commits into from
Aug 1, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
6 changes: 6 additions & 0 deletions CHANGELOG.md
Original file line number Diff line number Diff line change
Expand Up @@ -3,6 +3,12 @@
All notable changes to the Lethe project will be documented in this file.
The format is based on [Keep a Changelog](http://keepachangelog.com/).

## [Master] - 2024-08-01

### Added

- MINOR The possibility to specify an intermediate level as coarse grid solver for the gcmg preconditioner was added as a parameter. It allows to perform several v-cycles at the level chosen. [#1211](https://github.com/chaos-polymtl/lethe/pull/1211)

## [Master] - 2024-07-31

### Fixed
Expand Down
4 changes: 4 additions & 0 deletions doc/source/parameters/cfd/linear_solver_control.rst
Original file line number Diff line number Diff line change
Expand Up @@ -208,6 +208,7 @@ Different parameters for the main components of the two geometric multigrid algo
set mg verbosity = quiet
set mg min level = -1
set mg level min cells = -1
set mg int level = -1
set mg enable hessians in jacobian = true

# Relaxation smoother parameters
Expand Down Expand Up @@ -259,6 +260,9 @@ Different parameters for the main components of the two geometric multigrid algo
.. tip::
Evaluating terms involving the hessian is expensive. Therefore, one can turn on or off those terms in the mg level operators to improve performance by setting ``mg enable hessians in jacobian`` to ``false``. This is useful for certain problems and must be used carefully.

.. tip::
The ``mg int level`` option only works for ``gcmg`` preconditioner. It allows to choose an intermediate level as coarse grid solver and to perform several multigrid v-cycles there.

In addition, Lethe supports `p-multigrid` through the ``gcmg`` preconditioner. It can be used by specifying two additional parameters:

.. code-block:: text
Expand Down
3 changes: 3 additions & 0 deletions include/core/parameters.h
Original file line number Diff line number Diff line change
Expand Up @@ -1183,6 +1183,9 @@ namespace Parameters
/// MG minimum number of cells per level
int mg_level_min_cells;

/// MG intermediate level
int mg_int_level;

/// MG enable hessians in jacobian
bool mg_enable_hessians_jacobian;

Expand Down
21 changes: 20 additions & 1 deletion include/solvers/mf_navier_stokes.h
Original file line number Diff line number Diff line change
Expand Up @@ -147,6 +147,9 @@ class MFNavierStokesPreconditionGMG
/// Max level of the multigrid hierarchy
unsigned int maxlevel;

/// Intermediate level of the multigrid hierarchy
unsigned int intlevel;

/// Triangulations for the global coarsening case
std::vector<std::shared_ptr<const Triangulation<dim>>>
coarse_grid_triangulations;
Expand Down Expand Up @@ -206,14 +209,30 @@ class MFNavierStokesPreconditionGMG
/// Multigrid wrapper for the coarse grid solver
std::shared_ptr<MGCoarseGridBase<VectorType>> mg_coarse;

/// Solver control for the coarse grid solver (intermediate level)
std::shared_ptr<SolverControl> coarse_grid_solver_control_intermediate;

/// Multigrid wrapper for the coarse grid solver (intermediate level)
std::shared_ptr<MGCoarseGridBase<VectorType>> mg_coarse_intermediate;

/// GMRES as coarse grid solver (intermediate level)
std::shared_ptr<SolverGMRES<VectorType>> coarse_grid_solver_intermediate;

/// Multigrid method (intermediate level)
std::shared_ptr<Multigrid<VectorType>> mg_intermediate;

/// Global coarsening multigrid preconditioner object (intermediate level)
std::shared_ptr<PreconditionMG<dim, VectorType, GCTransferType>>
gc_multigrid_preconditioner_intermediate;

/// Multigrid method
std::shared_ptr<Multigrid<VectorType>> mg;

/// Local smoothing multigrid preconditioner object
std::shared_ptr<PreconditionMG<dim, VectorType, LSTransferType>>
ls_multigrid_preconditioner;

/// Global coarsening multigrid preconiditoner object
/// Global coarsening multigrid preconditioner object
std::shared_ptr<PreconditionMG<dim, VectorType, GCTransferType>>
gc_multigrid_preconditioner;

Expand Down
6 changes: 6 additions & 0 deletions source/core/parameters.cc
Original file line number Diff line number Diff line change
Expand Up @@ -2501,6 +2501,11 @@ namespace Parameters
Patterns::Integer(),
"mg minimum number of cells for coarse level");

prm.declare_entry("mg int level",
"-1",
Patterns::Integer(),
"mg int level");

prm.declare_entry(
"mg enable hessians in jacobian",
"true",
Expand Down Expand Up @@ -2685,6 +2690,7 @@ namespace Parameters

mg_min_level = prm.get_integer("mg min level");
mg_level_min_cells = prm.get_integer("mg level min cells");
mg_int_level = prm.get_integer("mg int level");
mg_enable_hessians_jacobian =
prm.get_bool("mg enable hessians in jacobian");
Assert(enable_hessians_jacobian || !mg_enable_hessians_jacobian,
Expand Down
76 changes: 70 additions & 6 deletions source/solvers/mf_navier_stokes.cc
Original file line number Diff line number Diff line change
Expand Up @@ -800,6 +800,10 @@ MFNavierStokesPreconditionGMG<dim>::MFNavierStokesPreconditionGMG(
this->simulation_parameters.linear_solver.at(PhysicsID::fluid_dynamics)
.mg_min_level;

int mg_int_level =
this->simulation_parameters.linear_solver.at(PhysicsID::fluid_dynamics)
.mg_int_level;

AssertThrow(
(mg_min_level + 1) <=
static_cast<int>(this->coarse_grid_triangulations.size()),
Expand Down Expand Up @@ -912,6 +916,8 @@ MFNavierStokesPreconditionGMG<dim>::MFNavierStokesPreconditionGMG(
this->minlevel = 0;
this->maxlevel = levels.size() - 1;

this->intlevel = (mg_int_level == -1) ? this->minlevel : mg_int_level;

// Local object for constraints of the different levels
MGLevelObject<AffineConstraints<typename VectorType::value_type>>
constraints;
Expand Down Expand Up @@ -1522,7 +1528,8 @@ MFNavierStokesPreconditionGMG<dim>::initialize(
*this->mg_transfer_ls,
*this->mg_smoother,
*this->mg_smoother,
this->minlevel);
this->minlevel,
this->maxlevel);

if (this->dof_handler.get_triangulation().has_hanging_nodes())
this->mg->set_edge_in_matrix(*this->mg_interface_matrix_in);
Expand All @@ -1537,12 +1544,53 @@ MFNavierStokesPreconditionGMG<dim>::initialize(
.preconditioner ==
Parameters::LinearSolver::PreconditionerType::gcmg)
{
if (this->minlevel != this->intlevel)
{
// Create main MG object
this->mg_intermediate =
std::make_shared<Multigrid<VectorType>>(*this->mg_matrix,
*this->mg_coarse,
*this->mg_transfer_gc,
*this->mg_smoother,
*this->mg_smoother,
this->minlevel,
this->intlevel);

// Create MG preconditioner
this->gc_multigrid_preconditioner_intermediate =
std::make_shared<PreconditionMG<dim, VectorType, GCTransferType>>(
this->dof_handler, *this->mg_intermediate, *this->mg_transfer_gc);

// TODO: introduce parameters
lpsaavedra marked this conversation as resolved.
Show resolved Hide resolved
this->coarse_grid_solver_control_intermediate =
std::make_shared<ReductionControl>(1000, 1e-20, 1e-2, false, false);

this->coarse_grid_solver_intermediate =
std::make_shared<SolverGMRES<VectorType>>(
*this->coarse_grid_solver_control_intermediate);

this->mg_coarse_intermediate =
std::make_shared<MGCoarseGridIterativeSolver<
VectorType,
SolverGMRES<VectorType>,
OperatorType,
PreconditionMG<dim, VectorType, GCTransferType>>>(
*this->coarse_grid_solver_intermediate,
*this->mg_operators[this->intlevel],
*this->gc_multigrid_preconditioner_intermediate);
}

// Create main MG object
this->mg = std::make_shared<Multigrid<VectorType>>(*this->mg_matrix,
*this->mg_coarse,
*this->mg_transfer_gc,
*this->mg_smoother,
*this->mg_smoother);
this->mg = std::make_shared<Multigrid<VectorType>>(
*this->mg_matrix,
(this->minlevel != this->intlevel) ? (*this->mg_coarse_intermediate) :
(*this->mg_coarse),
*this->mg_transfer_gc,
*this->mg_smoother,
*this->mg_smoother,
this->intlevel,
this->maxlevel,
Multigrid<VectorType>::Cycle::v_cycle);

// Create MG preconditioner
this->gc_multigrid_preconditioner =
Expand Down Expand Up @@ -1576,6 +1624,22 @@ MFNavierStokesPreconditionGMG<dim>::initialize(
this->mg->connect_post_smoother_step(
create_mg_timer_function("5_post_smoother_step"));

if (mg_intermediate)
{
this->mg_intermediate->connect_pre_smoother_step(
create_mg_timer_function("0_pre_smoother_step"));
this->mg_intermediate->connect_residual_step(
create_mg_timer_function("1_residual_step"));
this->mg_intermediate->connect_restriction(
create_mg_timer_function("2_restriction"));
this->mg_intermediate->connect_coarse_solve(create_mg_timer_function(""));
this->mg_intermediate->connect_prolongation(
create_mg_timer_function("3_prolongation"));
this->mg_intermediate->connect_edge_prolongation(
create_mg_timer_function("4_edge_prolongation"));
this->mg_intermediate->connect_post_smoother_step(
create_mg_timer_function("5_post_smoother_step"));
}

const auto create_mg_precon_timer_function = [&](const std::string &label) {
return [label, this](const bool flag) {
Expand Down
Loading