Skip to content

Commit 5c3a002

Browse files
authored
Fix a bug in the 2d mode of the IGF solver (#6199)
When there is only one node in the z-direction, the previous way of creating the domain for the IGF solver ends up with two nodes. That's not the desired behavior. Using the nodal BoxArray in MultiFab rho solves the issue.
1 parent 56f9b94 commit 5c3a002

File tree

3 files changed

+2
-6
lines changed

3 files changed

+2
-6
lines changed

Source/ablastr/fields/IntegratedGreenFunctionSolver.H

Lines changed: 0 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -123,14 +123,12 @@ namespace ablastr::fields
123123
* @param[in] rho the charge density amrex::MultiFab
124124
* @param[out] phi the electrostatic potential amrex::MultiFab
125125
* @param[in] cell_size an arreay of 3 reals dx dy dz
126-
* @param[in] ba amrex::BoxArray with the grid of a given level
127126
* @param[in] is_igf_2d_slices boolean to select between fully 3D Poisson solver and quasi-3D, i.e. one 2D Poisson solve on every z slice (default: false)
128127
*/
129128
void
130129
computePhiIGF (amrex::MultiFab const & rho,
131130
amrex::MultiFab & phi,
132131
std::array<amrex::Real, 3> const & cell_size,
133-
amrex::BoxArray const & ba,
134132
bool is_igf_2d_slices);
135133

136134
} // namespace ablastr::fields

Source/ablastr/fields/IntegratedGreenFunctionSolver.cpp

Lines changed: 1 addition & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -34,16 +34,14 @@ void
3434
computePhiIGF ( amrex::MultiFab const & rho,
3535
amrex::MultiFab & phi,
3636
std::array<amrex::Real, 3> const & cell_size,
37-
amrex::BoxArray const & ba,
3837
bool const is_igf_2d_slices)
3938
{
4039
using namespace amrex::literals;
4140

4241
BL_PROFILE("ablastr::fields::computePhiIGF");
4342

4443
// Define box that encompasses the full domain
45-
amrex::Box domain = ba.minimalBox();
46-
domain.surroundingNodes(); // get nodal points, since `phi` and `rho` are nodal
44+
amrex::Box domain = rho.boxArray().minimalBox();
4745
domain.grow( phi.nGrowVect() ); // include guard cells
4846

4947
int nprocs = amrex::ParallelDescriptor::NProcs();

Source/ablastr/fields/PoissonSolver.H

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -279,7 +279,7 @@ computePhi (
279279
if ( max_norm_b == 0 ) {
280280
phi[lev]->setVal(0);
281281
} else {
282-
computePhiIGF( *rho[lev], *phi[lev], dx_scaled, grids[lev], is_igf_2d);
282+
computePhiIGF( *rho[lev], *phi[lev], dx_scaled, is_igf_2d);
283283
}
284284
continue;
285285
}

0 commit comments

Comments
 (0)