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

GMRES/MLMG: Set the number of MG V-cycles per GMRES iteration #3875

Merged
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.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
27 changes: 21 additions & 6 deletions Src/LinearSolvers/AMReX_GMRES_MLMG.H
Expand Up @@ -91,14 +91,19 @@ public:

void precond (MF& lhs, MF const& rhs) const;

//! Control whether or not to use MLMG as preconditioner.
bool usePrecond (bool new_flag) { return std::exchange(m_use_precond, new_flag); }

//! Set the number of MLMG preconditioner iterations per GMRES iteration.
void setPrecondNumIters (int precond_niters) { m_precond_niters = precond_niters; }

private:
GM m_gmres;
MG* m_mlmg;
MLLinOpT<MF>* m_linop;
bool m_use_precond = true;
bool m_prop_zero = false;
int m_precond_niters = 1;
};

template <typename MF>
Expand Down Expand Up @@ -184,12 +189,22 @@ void GMRESMLMGT<MF>::precond (MF& lhs, MF const& rhs) const
if (m_use_precond) {
m_mlmg->prepareMGcycle();

LocalCopy(m_mlmg->res[0][0], rhs, 0, 0, nComp(rhs), IntVect(0));

m_mlmg->mgVcycle(0,0);

LocalCopy(lhs, m_mlmg->cor[0][0], 0, 0, nComp(rhs), IntVect(0));

for (int icycle = 0; icycle < m_precond_niters; ++icycle) {
if (icycle == 0) {
LocalCopy(m_mlmg->res[0][0], rhs, 0, 0, nComp(rhs), IntVect(0));
} else {
m_mlmg->computeResOfCorrection(0,0);
LocalCopy(m_mlmg->res[0][0], m_mlmg->rescor[0][0], 0, 0, nComp(rhs), IntVect(0));
}

m_mlmg->mgVcycle(0,0);

if (icycle == 0) {
LocalCopy(lhs, m_mlmg->cor[0][0], 0, 0, nComp(rhs), IntVect(0));
} else {
increment(lhs, m_mlmg->cor[0][0], RT(1));
}
}
} else {
LocalCopy(lhs, rhs, 0, 0, nComp(lhs), IntVect(0));
}
Expand Down
1 change: 1 addition & 0 deletions Tests/LinearSolvers/CurlCurl/MyTest.H
Expand Up @@ -32,6 +32,7 @@ private:

bool use_gmres = false;
bool gmres_use_precond = true;
int gmres_precond_niters = 1;

amrex::Geometry geom;
amrex::BoxArray grids;
Expand Down
2 changes: 2 additions & 0 deletions Tests/LinearSolvers/CurlCurl/MyTest.cpp
Expand Up @@ -62,6 +62,7 @@ MyTest::solve ()
{
GMRESMLMGT<V> gmsolver(mlmg);
gmsolver.usePrecond(gmres_use_precond);
gmsolver.setPrecondNumIters(gmres_precond_niters);

// This system has homogeneous BC unlike
// Tests/LinearSolvers/ABecLaplacian_C, so the setup is simpler.
Expand Down Expand Up @@ -106,6 +107,7 @@ MyTest::readParameters ()

pp.query("use_gmres", use_gmres);
pp.query("gmres_use_precond", gmres_use_precond);
pp.query("gmres_precond_niters", gmres_precond_niters);

pp.query("beta_factor", beta_factor);
pp.query("alpha", alpha);
Expand Down