WarpX
Loading...
Searching...
No Matches
PoissonSolver.H
Go to the documentation of this file.
1/* Copyright 2019-2022 Axel Huebl, Remi Lehe
2 *
3 * This file is part of WarpX.
4 *
5 * License: BSD-3-Clause-LBNL
6 */
7#ifndef ABLASTR_POISSON_SOLVER_H
8#define ABLASTR_POISSON_SOLVER_H
9
10#include <ablastr/constant.H>
12#include <ablastr/utils/Enums.H>
19
20#if defined(ABLASTR_USE_FFT) && defined(WARPX_DIM_3D)
22#endif
23
24#include <AMReX_Array.H>
25#include <AMReX_Array4.H>
26#include <AMReX_BLassert.H>
27#include <AMReX_Box.H>
28#include <AMReX_BoxArray.H>
29#include <AMReX_Config.H>
31#include <AMReX_FArrayBox.H>
32#include <AMReX_FabArray.H>
33#include <AMReX_Geometry.H>
34#include <AMReX_GpuControl.H>
35#include <AMReX_GpuLaunch.H>
36#include <AMReX_GpuQualifiers.H>
37#include <AMReX_IndexType.H>
38#include <AMReX_IntVect.H>
39#include <AMReX_LO_BCTYPES.H>
40#include <AMReX_MFIter.H>
41#include <AMReX_MFInterp_C.H>
42#include <AMReX_MLMG.H>
43#include <AMReX_MLLinOp.H>
45#include <AMReX_MultiFab.H>
46#include <AMReX_Parser.H>
47#include <AMReX_REAL.H>
48#include <AMReX_SPACE.H>
49#include <AMReX_Vector.H>
51#ifdef AMREX_USE_EB
52# include <AMReX_EBFabFactory.H>
53#endif
54
55#include <algorithm>
56#include <array>
57#include <optional>
58#include <stdexcept>
59
60
61namespace ablastr::fields {
62
70inline amrex::Real getMaxNormRho (
72 int finest_level,
73 amrex::Real & absolute_tolerance)
74{
75 amrex::Real max_norm_rho = 0.0;
76 for (int lev=0; lev<=finest_level; lev++) {
77 max_norm_rho = amrex::max(max_norm_rho, rho[lev]->norm0());
78 }
80
81 if (max_norm_rho == 0) {
82 if (absolute_tolerance == 0.0) { absolute_tolerance = amrex::Real(1e-6); }
84 "ElectrostaticSolver",
85 "Max norm of rho is 0",
87 );
88 }
89 return max_norm_rho;
90}
91
108 amrex::MultiFab const* phi_lev,
109 amrex::MultiFab* phi_lev_plus_one,
110 amrex::Geometry const & geom_lev,
111 bool do_single_precision_comms,
112 const amrex::IntVect& refratio,
113 const int ncomp,
114 const int ng)
115{
116 using namespace amrex::literals;
117
118 // Allocate phi_cp for lev+1
119 amrex::BoxArray ba = phi_lev_plus_one->boxArray();
120 ba.coarsen(refratio);
121 amrex::MultiFab phi_cp(ba, phi_lev_plus_one->DistributionMap(), ncomp, ng);
122 if (ng > 0) {
123 // Set all values outside the domain to zero
124 phi_cp.setDomainBndry(0.0_rt, geom_lev);
125 }
126
127 // Copy from phi[lev] to phi_cp (in parallel)
128 const amrex::Periodicity& crse_period = geom_lev.periodicity();
129
131 phi_cp,
132 *phi_lev,
133 0,
134 0,
135 1,
137 amrex::IntVect(ng),
138 do_single_precision_comms,
139 crse_period
140 );
141
142 // Local interpolation from phi_cp to phi[lev+1]
143#ifdef AMREX_USE_OMP
144#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
145#endif
146 for (amrex::MFIter mfi(*phi_lev_plus_one, amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi) {
147 amrex::Array4<amrex::Real> const phi_fp_arr = phi_lev_plus_one->array(mfi);
148 amrex::Array4<amrex::Real const> const phi_cp_arr = phi_cp.array(mfi);
149
150 details::PoissonInterpCPtoFP const interp(phi_fp_arr, phi_cp_arr, refratio);
151
152 amrex::Box const& b = mfi.growntilebox(ng);
153 amrex::ParallelFor(b, interp);
154 }
155}
156
191template<
192 typename T_PostPhiCalculationFunctor = std::nullopt_t,
193 typename T_BoundaryHandler = std::nullopt_t,
194 typename T_FArrayBoxFactory = void
195>
196void
200 std::array<amrex::Real, 3> const beta,
201 amrex::Real relative_tolerance,
202 amrex::Real absolute_tolerance,
203 int max_iters,
204 int verbosity,
208 utils::enums::GridType grid_type,
209 bool is_solver_igf_on_lev0,
210 [[maybe_unused]] bool const is_igf_2d,
211 bool eb_enabled = false,
212 bool do_single_precision_comms = false,
213 std::optional<amrex::Vector<amrex::IntVect> > rel_ref_ratio = std::nullopt,
214 [[maybe_unused]] T_PostPhiCalculationFunctor post_phi_calculation = std::nullopt,
215 [[maybe_unused]] T_BoundaryHandler const boundary_handler = std::nullopt,
216 [[maybe_unused]] std::optional<amrex::Real const> current_time = std::nullopt, // only used for EB
217 [[maybe_unused]] std::optional<amrex::Vector<T_FArrayBoxFactory const *> > eb_farray_box_factory = std::nullopt // only used for EB
218)
219{
220 using namespace amrex::literals;
221
222 ABLASTR_PROFILE("computePhi");
223
224 auto const finest_level = static_cast<int>(rho.size() - 1);
225
226 if (!rel_ref_ratio.has_value()) {
227 ABLASTR_ALWAYS_ASSERT_WITH_MESSAGE(finest_level == 0u,
228 "rel_ref_ratio must be set if mesh-refinement is used");
229 rel_ref_ratio = amrex::Vector<amrex::IntVect>{{amrex::IntVect(AMREX_D_DECL(1, 1, 1))}};
230 }
231
232#if !defined(AMREX_USE_EB)
234 "Embedded boundary solve requested but not compiled in");
235#endif
236
237#if !defined(ABLASTR_USE_FFT)
238 ABLASTR_ALWAYS_ASSERT_WITH_MESSAGE( !is_solver_igf_on_lev0,
239 "Must compile with FFT support to use the IGF solver!");
240#endif
241
242#if !defined(WARPX_DIM_3D)
243 ABLASTR_ALWAYS_ASSERT_WITH_MESSAGE( !is_solver_igf_on_lev0,
244 "The FFT Poisson solver is currently only implemented for 3D!");
245#endif
246
247 // Set the value of beta
249#if defined(WARPX_DIM_1D_Z)
250 {{ beta[2] }}; // beta_x and beta_z
251#elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
252 {{ beta[0], beta[2] }}; // beta_x and beta_z
253#elif defined(WARPX_DIM_RCYLINDER) || defined(WARPX_DIM_RSPHERE)
254 {{ beta[0] }}; // beta_x
255#else
256 {{ beta[0], beta[1], beta[2] }};
257#endif
258
260 std::all_of(beta_solver.begin(), beta_solver.end(),
261 [](const auto b){return (b<1.0_rt);}),
262 "Components of beta_solver must be < 1.");
263
264 // determine if rho is zero everywhere
265 const amrex::Real max_norm_b = getMaxNormRho(
266 amrex::GetVecOfConstPtrs(rho), finest_level, absolute_tolerance);
267
268 amrex::LPInfo info;
269
270 for (int lev=0; lev<=finest_level; lev++) {
272 {AMREX_D_DECL(geom[lev].CellSize(0)/std::sqrt(1._rt-beta_solver[0]*beta_solver[0]),
273 geom[lev].CellSize(1)/std::sqrt(1._rt-beta_solver[1]*beta_solver[1]),
274 geom[lev].CellSize(2)/std::sqrt(1._rt-beta_solver[2]*beta_solver[2]))};
275
276#if (defined(ABLASTR_USE_FFT) && defined(WARPX_DIM_3D))
277 // Use the Integrated Green Function solver (FFT) on the coarsest level if it was selected
278 if(is_solver_igf_on_lev0 && lev==0){
279 if ( max_norm_b == 0 ) {
280 phi[lev]->setVal(0);
281 } else {
282 computePhiIGF( *rho[lev], *phi[lev], dx_scaled, is_igf_2d);
283 }
284 continue;
285 }
286#endif
287 // Use the Multigrid (MLMG) solver if selected or on refined patches
288 // but first scale rho appropriately
289 rho[lev]->mult(-1._rt / ablastr::constant::SI::epsilon_0);
290
291#ifdef WARPX_DIM_RZ
292 constexpr bool is_rz = true;
293#else
294 constexpr bool is_rz = false;
295#endif
296
297 if (!eb_enabled && !is_rz) {
298 // Determine whether to use semi-coarsening
299 int max_semicoarsening_level = 0;
300 int semicoarsening_direction = -1;
301 const auto min_dir = static_cast<int>(std::distance(dx_scaled.begin(),
302 std::min_element(dx_scaled.begin(), dx_scaled.end())));
303 const auto max_dir = static_cast<int>(std::distance(dx_scaled.begin(),
304 std::max_element(dx_scaled.begin(), dx_scaled.end())));
305 if (dx_scaled[max_dir] > dx_scaled[min_dir]) {
306 semicoarsening_direction = max_dir;
307 max_semicoarsening_level = static_cast<int>(std::log2(dx_scaled[max_dir] / dx_scaled[min_dir]));
308 }
309 if (max_semicoarsening_level > 0) {
310 info.setSemicoarsening(true);
311 info.setMaxSemicoarseningLevel(max_semicoarsening_level);
312 info.setSemicoarseningDirection(semicoarsening_direction);
313 }
314 }
315
316 std::unique_ptr<amrex::MLNodeLinOp> linop;
317 if (eb_enabled || is_rz) {
318 // In the presence of EB or RZ: the solver assumes that the beam is
319 // propagating along one of the axes of the grid, i.e. that only *one*
320 // of the components of `beta` is non-negligible.
321 auto linop_nodelap = std::make_unique<amrex::MLEBNodeFDLaplacian>();
322 if (eb_enabled) {
323#if defined(AMREX_USE_EB)
324 if constexpr(std::is_same_v<void, T_FArrayBoxFactory>) {
325 throw std::runtime_error("EB requested by eb_farray_box_factory not provided!");
326 } else {
327 linop_nodelap->define(
331 info,
332 amrex::Vector<amrex::EBFArrayBoxFactory const*>{eb_farray_box_factory.value()[lev]}
333 );
334 }
335#endif
336 }
337 else {
338 // TODO: rather use MLNodeTensorLaplacian (for RZ w/o EB) here? Semi-Coarsening would be nice here
339 linop_nodelap->define(
343 info
344 );
345 }
346
347 // Note: this assumes that the beam is propagating along
348 // one of the axes of the grid, i.e. that only *one* of the
349 // components of `beta` is non-negligible. // we use this
350#if defined(WARPX_DIM_RZ)
351 linop_nodelap->setRZ(true);
352 linop_nodelap->setSigma({0._rt, 1._rt-beta_solver[1]*beta_solver[1]});
353#else
354 linop_nodelap->setSigma({AMREX_D_DECL(
355 1._rt-beta_solver[0]*beta_solver[0],
356 1._rt-beta_solver[1]*beta_solver[1],
357 1._rt-beta_solver[2]*beta_solver[2])});
358#endif
359#if defined(AMREX_USE_EB)
360 if (eb_enabled) {
361 if constexpr (!std::is_same_v<T_BoundaryHandler, std::nullopt_t>) {
362 // if the EB potential only depends on time, the potential can be passed
363 // as a float instead of a callable
364 if (boundary_handler.phi_EB_only_t) {
365 linop_nodelap->setEBDirichlet(boundary_handler.potential_eb_t(current_time.value()));
366 } else {
367 linop_nodelap->setEBDirichlet(boundary_handler.getPhiEB(current_time.value()));
368 }
369 } else
370 {
371 ABLASTR_ALWAYS_ASSERT_WITH_MESSAGE( !is_solver_igf_on_lev0,
372 "EB Poisson solver enabled but no 'boundary_handler' passed!");
373 }
374 }
375#endif
376 linop = std::move(linop_nodelap);
377 } else {
378 // In the absence of EB and RZ: use a more generic solver
379 // that can handle beams propagating in any direction
380 auto linop_tenslap = std::make_unique<amrex::MLNodeTensorLaplacian>(
384 info
385 );
386 linop_tenslap->setBeta(beta_solver); // for the non-axis-aligned solver
387 linop = std::move(linop_tenslap);
388 }
389
390 // Level 0 domain boundary
391 if constexpr (std::is_same_v<T_BoundaryHandler, std::nullopt_t>) {
393 amrex::LinOpBCType::Dirichlet,
394 amrex::LinOpBCType::Dirichlet,
395 amrex::LinOpBCType::Dirichlet
396 )};
398 linop->setDomainBC(lobc, hibc);
399 } else {
400 linop->setDomainBC(boundary_handler.lobc, boundary_handler.hibc);
401 }
402
403 // Solve the Poisson equation
404 amrex::MLMG mlmg(*linop); // actual solver defined here
405 mlmg.setVerbose(verbosity);
406 mlmg.setMaxIter(max_iters);
407 mlmg.setConvergenceNormType((max_norm_b > 0) ? amrex::MLMGNormType::bnorm : amrex::MLMGNormType::greater);
408
409 const int ng = int(grid_type == utils::enums::GridType::Collocated); // ghost cells
410 if (ng) {
411 // In this case, computeE needs to use ghost nodes data. So we
412 // ask MLMG to fill BC for us after it solves the problem.
413 mlmg.setFinalFillBC(true);
414 }
415
416 // Solve Poisson equation at lev
417 mlmg.solve( {phi[lev]}, {rho[lev]},
418 relative_tolerance, absolute_tolerance );
419
420 const amrex::IntVect& refratio = rel_ref_ratio.value()[lev];
421 const int ncomp = linop->getNComp();
422
423 // needed for solving the levels by levels:
424 // - coarser level is initial guess for finer level
425 // - coarser level provides boundary values for finer level patch
426 // Interpolation from phi[lev] to phi[lev+1]
427 // (This provides both the boundary conditions and initial guess for phi[lev+1])
428 if (lev < finest_level) {
430 phi[lev+1],
431 geom[lev],
432 do_single_precision_comms,
433 refratio,
434 ncomp,
435 ng);
436 }
437
438 // Run additional operations, such as calculation of the E field for embedded boundaries
439 if constexpr (!std::is_same_v<T_PostPhiCalculationFunctor, std::nullopt_t>) {
440 if (post_phi_calculation.has_value()) {
441 post_phi_calculation.value()(mlmg, lev);
442 }
443 }
444 rho[lev]->mult(-ablastr::constant::SI::epsilon_0); // Multiply rho by epsilon again
445
446 } // loop over lev(els)
447} // computePhi
448} // namespace ablastr::fields
449
450#endif // ABLASTR_POISSON_SOLVER_H
#define AMREX_D_DECL(a, b, c)
#define ABLASTR_PROFILE(fname)
Definition ProfilerWrapper.H:13
#define ABLASTR_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition TextMsg.H:75
BoxArray & coarsen(int refinement_ratio)
const DistributionMapping & DistributionMap() const noexcept
const BoxArray & boxArray() const noexcept
void setDomainBndry(value_type val, const Geometry &geom)
Array4< typename FabArray< FAB >::value_type const > array(const MFIter &mfi) const noexcept
Periodicity periodicity() const noexcept
void setVerbose(int v) noexcept
void setFinalFillBC(int flag) noexcept
void setConvergenceNormType(MLMGNormType norm) noexcept
RT solve(const Vector< AMF * > &a_sol, const Vector< AMF const * > &a_rhs, RT a_tol_rel, RT a_tol_abs, const char *checkpoint_file=nullptr)
void setMaxIter(int n) noexcept
Long size() const noexcept
constexpr auto epsilon_0
vacuum permittivity: dielectric permittivity of vacuum [F/m]
Definition constant.H:152
Definition EffectivePotentialPoissonSolver.H:63
void interpolatePhiBetweenLevels(amrex::MultiFab const *phi_lev, amrex::MultiFab *phi_lev_plus_one, amrex::Geometry const &geom_lev, bool do_single_precision_comms, const amrex::IntVect &refratio, const int ncomp, const int ng)
Definition PoissonSolver.H:107
void computePhiIGF(amrex::MultiFab const &rho, amrex::MultiFab &phi, std::array< amrex::Real, 3 > const &cell_size, bool const is_igf_2d_slices)
Compute the electrostatic potential using the Integrated Green Function method as in http://dx....
Definition IntegratedGreenFunctionSolver.cpp:34
amrex::Real getMaxNormRho(ablastr::fields::ConstMultiLevelScalarField const &rho, int finest_level, amrex::Real &absolute_tolerance)
Definition PoissonSolver.H:70
void computePhi(ablastr::fields::MultiLevelScalarField const &rho, ablastr::fields::MultiLevelScalarField const &phi, std::array< amrex::Real, 3 > const beta, amrex::Real relative_tolerance, amrex::Real absolute_tolerance, int max_iters, int verbosity, amrex::Vector< amrex::Geometry > const &geom, amrex::Vector< amrex::DistributionMapping > const &dmap, amrex::Vector< amrex::BoxArray > const &grids, utils::enums::GridType grid_type, bool is_solver_igf_on_lev0, bool const is_igf_2d, bool eb_enabled=false, bool do_single_precision_comms=false, std::optional< amrex::Vector< amrex::IntVect > > rel_ref_ratio=std::nullopt, T_PostPhiCalculationFunctor post_phi_calculation=std::nullopt, T_BoundaryHandler const boundary_handler=std::nullopt, std::optional< amrex::Real const > current_time=std::nullopt, std::optional< amrex::Vector< T_FArrayBoxFactory const * > > eb_farray_box_factory=std::nullopt)
Definition PoissonSolver.H:197
amrex::Vector< ScalarField > MultiLevelScalarField
Definition MultiFabRegister.H:200
amrex::Vector< ConstScalarField > ConstMultiLevelScalarField
Definition MultiFabRegister.H:204
void ParallelCopy(amrex::MultiFab &dst, const amrex::MultiFab &src, int src_comp, int dst_comp, int num_comp, const amrex::IntVect &src_nghost, const amrex::IntVect &dst_nghost, bool do_single_precision_comms, const amrex::Periodicity &period, amrex::FabArrayBase::CpOp op)
Definition Communication.cpp:28
GridType
Definition Enums.H:23
@ Collocated
Definition Enums.H:23
void WMRecordWarning(const std::string &topic, const std::string &text, const WarnPriority &priority=WarnPriority::medium)
Helper function to abbreviate the call to WarnManager::GetInstance().RecordWarning (recording a warni...
Definition WarnManager.cpp:320
@ low
Definition WarnManager.H:31
void ReduceRealMax(Vector< std::reference_wrapper< Real > > const &)
std::enable_if_t< std::is_integral_v< T > > ParallelFor(TypeList< CTOs... > ctos, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
Vector< const T * > GetVecOfConstPtrs(const Vector< T > &a)
BoxND< 3 > Box
MLMGT< MultiFab > MLMG
IntVectND< 3 > IntVect
__host__ __device__ constexpr const T & max(const T &a, const T &b) noexcept
bool TilingIfNotGPU() noexcept
std::array< T, N > Array
LPInfo & setSemicoarsening(bool x) noexcept
LPInfo & setSemicoarseningDirection(int n) noexcept
LPInfo & setMaxSemicoarseningLevel(int n) noexcept