Skip to content
Open
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
5 changes: 5 additions & 0 deletions Docs/source/usage/parameters.rst
Original file line number Diff line number Diff line change
Expand Up @@ -391,6 +391,11 @@ Overall simulation parameters
non-zero value is specified by the user via
:pp:param:`warpx.self_fields_absolute_tolerance`).

* ``gmres``: Poisson's equation is solved using a standard GMRES solver.
The GMRES algorithm is the classic right-preconditioned variation presented in the original GMRES
paper in `Y. Saad et al. (SIAM J. Sci. Stat. Comput., 7, 3, 1986) <https://doi.org/10.1137/0907058>`__.
The preconditioner is 1 cycle of the AMReX-based MLMG solvers.

* ``fft``: Poisson's equation is solved using an Integrated Green Function method (which requires FFT calculations).
See these references for more details :cite:t:`param-QiangPhysRevSTAB2006`, :cite:t:`param-QiangPhysRevSTAB2006err`.
It only works in 3D and it requires the compilation flag ``-DWarpX_FFT=ON``.
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -81,7 +81,7 @@ void PoissonBoundaryHandler::DefinePhiBCs (const amrex::Geometry& geom)
amrex::ignore_unused(geom);
#endif
for (int idim=dim_start; idim<AMREX_SPACEDIM; idim++){
if (WarpX::poisson_solver_id == PoissonSolverAlgo::Multigrid){
if (WarpX::poisson_solver_id == PoissonSolverAlgo::Multigrid || WarpX::poisson_solver_id == PoissonSolverAlgo::GMRES){
if ( WarpX::field_boundary_lo[idim] == FieldBoundaryType::Periodic
&& WarpX::field_boundary_hi[idim] == FieldBoundaryType::Periodic ) {
lobc[idim] = LinOpBCType::Periodic;
Expand Down
2 changes: 1 addition & 1 deletion Source/Initialization/DivCleaner/ProjectionDivCleaner.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -401,7 +401,7 @@ WarpX::ProjectionCleanDivB() {
|| WarpX::electromagnetic_solver_id == ElectromagneticSolverAlgo::HybridPIC
|| ( (WarpX::electrostatic_solver_id == ElectrostaticSolverAlgo::LabFrame
|| WarpX::electrostatic_solver_id == ElectrostaticSolverAlgo::LabFrameElectroMagnetostatic)
&& WarpX::poisson_solver_id == PoissonSolverAlgo::Multigrid))
&& (WarpX::poisson_solver_id == PoissonSolverAlgo::Multigrid || WarpX::poisson_solver_id == PoissonSolverAlgo::GMRES)))
#if defined(WARPX_DIM_RZ)
&& WarpX::grid_type == GridType::Staggered
#endif
Expand Down
3 changes: 3 additions & 0 deletions Source/Initialization/WarpXInitData.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -591,6 +591,9 @@ WarpX::PrintMainPICparameters ()
else if(poisson_solver_id == PoissonSolverAlgo::Multigrid){
amrex::Print() << "Poisson solver: | multigrid" << "\n";
}
else if(poisson_solver_id == PoissonSolverAlgo::GMRES){
amrex::Print() << "Poisson solver: | gmres" << "\n";
}
}

amrex::Print() << "-------------------------------------------------------------------------------\n";
Expand Down
1 change: 1 addition & 0 deletions Source/Utils/WarpXAlgorithmSelection.H
Original file line number Diff line number Diff line change
Expand Up @@ -67,6 +67,7 @@ AMREX_ENUM(ElectrostaticSolverAlgo,

AMREX_ENUM(PoissonSolverAlgo,
Multigrid,
GMRES,
IntegratedGreenFunction,
fft = IntegratedGreenFunction,
Default = Multigrid);
Expand Down
2 changes: 1 addition & 1 deletion Source/WarpX.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -1147,7 +1147,7 @@ WarpX::ReadParameters ()
|| WarpX::electromagnetic_solver_id == ElectromagneticSolverAlgo::HybridPIC
|| ( (WarpX::electrostatic_solver_id == ElectrostaticSolverAlgo::LabFrame
|| WarpX::electrostatic_solver_id == ElectrostaticSolverAlgo::LabFrameElectroMagnetostatic)
&& WarpX::poisson_solver_id == PoissonSolverAlgo::Multigrid)))
&& (WarpX::poisson_solver_id == PoissonSolverAlgo::Multigrid || WarpX::poisson_solver_id == PoissonSolverAlgo::GMRES))))
{
m_do_initial_div_cleaning = true;
}
Expand Down
56 changes: 50 additions & 6 deletions Source/ablastr/fields/PoissonSolver.H
Original file line number Diff line number Diff line change
Expand Up @@ -51,6 +51,8 @@
#ifdef AMREX_USE_EB
# include <AMReX_EBFabFactory.H>
#endif
#include <AMReX_GMRES_MLMG.H>
#include "WarpX.H"

#include <algorithm>
#include <array>
Expand Down Expand Up @@ -416,8 +418,15 @@ computePhi (
rho[lev]->OverrideSync(geom[lev].periodicity());

// Solve Poisson equation at lev
mlmg.solve( {phi[lev]}, {rho[lev]},
relative_tolerance, absolute_tolerance );
if (WarpX::poisson_solver_id == PoissonSolverAlgo::Multigrid) {
mlmg.solve( {phi[lev]}, {rho[lev]},
relative_tolerance, absolute_tolerance );
} else {
amrex::GMRESMLMG gmsolve(mlmg);
gmsolve.setMaxIters(max_iters);
gmsolve.solve( *phi[lev], *rho[lev], relative_tolerance, absolute_tolerance );
Comment on lines +425 to +427
Copy link
Copy Markdown
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@RemiLehe commented: can we also handle max_iters here? Otherwise it is silently ignored if user-set

linop->postSolve({phi[lev]});
Copy link
Copy Markdown
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

What is the need for linop->postSolve() here ? Why is it not needed when using mlmg?

Copy link
Copy Markdown
Contributor Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

MLMG calls postSolve within, which sets the value of phi in covered cells to the Dirichlet value. GMRES initially left phi in covered cells at zero; I added this call so the plotfiles would match up in the covered cells.

}

const amrex::IntVect& refratio = rel_ref_ratio.value()[lev];
const int ncomp = linop->getNComp();
Expand All @@ -437,12 +446,47 @@ computePhi (
ng);
}

// Run additional operations, such as calculation of the E field for embedded boundaries
if constexpr (!std::is_same_v<T_PostPhiCalculationFunctor, std::nullopt_t>) {
if (post_phi_calculation.has_value()) {
post_phi_calculation.value()(mlmg, lev);
if (WarpX::poisson_solver_id == PoissonSolverAlgo::Multigrid) {
// Run additional operations, such as calculation of the E field for embedded boundaries
if constexpr (!std::is_same_v<T_PostPhiCalculationFunctor, std::nullopt_t>) {
if (post_phi_calculation.has_value()) {
post_phi_calculation.value()(mlmg, lev);
}
}
}

if (WarpX::poisson_solver_id == PoissonSolverAlgo::GMRES && eb_enabled) {

// create an alieas to the WarpX Efield
auto & warpx = WarpX::GetInstance();

amrex::Vector< amrex::Array<amrex::MultiFab *, AMREX_SPACEDIM> > e_field;
e_field.push_back(
#if defined(WARPX_DIM_1D_Z)
amrex::Array<amrex::MultiFab*, 1>{
warpx.m_fields.get(warpx::fields::FieldType::Efield_fp, Direction{2}, lev)
}
#elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
amrex::Array<amrex::MultiFab*, 2>{
warpx.m_fields.get(warpx::fields::FieldType::Efield_fp, Direction{0}, lev),
warpx.m_fields.get(warpx::fields::FieldType::Efield_fp, Direction{2}, lev)
}
#elif defined(WARPX_DIM_3D)
amrex::Array<amrex::MultiFab *, 3>{
warpx.m_fields.get(warpx::fields::FieldType::Efield_fp, Direction{0}, lev),
warpx.m_fields.get(warpx::fields::FieldType::Efield_fp, Direction{1}, lev),
warpx.m_fields.get(warpx::fields::FieldType::Efield_fp, Direction{2}, lev)
}
#endif
Comment on lines +465 to +480
Copy link
Copy Markdown
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@dpgrote Is electrostatic mode in WarpX not setup for RCYLINDER and RSPHERE geometries?

Copy link
Copy Markdown
Member

@dpgrote dpgrote Apr 16, 2026

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

@JustinRayAngus Those geometries are not supported in electrostatic. Also, this block of code is for EB which is also not supported.

);

// compute E = grad(phi)
linop->compGrad(lev, e_field[lev], *phi[lev], amrex::MLLinOp::Location::CellCenter);
for(int dir=0; dir<AMREX_SPACEDIM; ++dir) {
e_field[lev][dir]->mult(-1.);
}
}

rho[lev]->mult(-ablastr::constant::SI::epsilon_0); // Multiply rho by epsilon again

} // loop over lev(els)
Expand Down
Loading