WarpX
Loading...
Searching...
No Matches
VectorPoissonSolver.H
Go to the documentation of this file.
1/* Copyright 2022-2024 S. Eric Clark (Helion Energy, formerly LLNL)
2 *
3 * This file is part of WarpX.
4 *
5 * License: BSD-3-Clause-LBNL
6 */
7#ifndef ABLASTR_VECTOR_POISSON_SOLVER_H
8#define ABLASTR_VECTOR_POISSON_SOLVER_H
9
10#include <ablastr/constant.H>
15
16#include <AMReX_Array.H>
17#include <AMReX_Array4.H>
18#include <AMReX_BLassert.H>
19#include <AMReX_Box.H>
20#include <AMReX_BoxArray.H>
21#include <AMReX_Config.H>
23#include <AMReX_FArrayBox.H>
24#include <AMReX_FabArray.H>
25#include <AMReX_Geometry.H>
26#include <AMReX_GpuControl.H>
27#include <AMReX_GpuLaunch.H>
28#include <AMReX_GpuQualifiers.H>
29#include <AMReX_IndexType.H>
30#include <AMReX_IntVect.H>
31#include <AMReX_LO_BCTYPES.H>
32#include <AMReX_MFIter.H>
33#include <AMReX_MFInterp_C.H>
35#include <AMReX_MLLinOp.H>
36#include <AMReX_MLMG.H>
37#include <AMReX_MultiFab.H>
38#include <AMReX_Parser.H>
39#include <AMReX_REAL.H>
40#include <AMReX_SPACE.H>
41#include <AMReX_Vector.H>
42#ifdef AMREX_USE_EB
43# include <AMReX_EBFabFactory.H>
44#endif
45
46#include <array>
47#include <optional>
48
49namespace ablastr::fields {
50
81template<
82 typename T_BoundaryHandler,
83 typename T_PostACalculationFunctor = std::nullopt_t,
84 typename T_FArrayBoxFactory = void
85>
86void
89 amrex::Real relative_tolerance,
90 amrex::Real absolute_tolerance,
91 int max_iters,
92 int verbosity,
96 T_BoundaryHandler const boundary_handler,
97 bool eb_enabled = false,
98 bool do_single_precision_comms = false,
99 std::optional<amrex::Vector<amrex::IntVect> > rel_ref_ratio = std::nullopt,
100 [[maybe_unused]] T_PostACalculationFunctor post_A_calculation = std::nullopt,
101 [[maybe_unused]] std::optional<amrex::Real const> current_time = std::nullopt, // only used for EB
102 [[maybe_unused]] std::optional<amrex::Vector<T_FArrayBoxFactory const *> > eb_farray_box_factory = std::nullopt // only used for EB
103)
104{
105 using namespace amrex::literals;
106
107 if (!rel_ref_ratio.has_value()) {
108 ABLASTR_ALWAYS_ASSERT_WITH_MESSAGE(curr.size() == 1u,
109 "rel_ref_ratio must be set if mesh-refinement is used");
110 rel_ref_ratio = amrex::Vector<amrex::IntVect>{{amrex::IntVect(AMREX_D_DECL(1, 1, 1))}};
111 }
112
113 auto const finest_level = static_cast<int>(curr.size()) - 1;
114
115 // scale J appropriately; also determine if current is zero everywhere
116 amrex::Real max_comp_J = 0.0;
117 for (int lev=0; lev<=finest_level; lev++) {
118 for (int adim=0; adim<3; adim++) {
119 curr[lev][adim]->mult(-1._rt*ablastr::constant::SI::mu0); // Unscaled below
120 max_comp_J = amrex::max(max_comp_J, curr[lev][adim]->norm0());
121 }
122 }
124
125 const bool always_use_bnorm = (max_comp_J > 0);
126 if (!always_use_bnorm) {
127 if (absolute_tolerance == 0.0) { absolute_tolerance = amrex::Real(1e-6); }
129 "MagnetostaticSolver",
130 "Max norm of J is 0",
132 );
133 }
134
135 // Loop over dimensions of A to solve each component individually
136 for (int lev=0; lev<=finest_level; lev++) {
137 amrex::LPInfo info;
138
139#ifdef WARPX_DIM_RZ
140 constexpr bool is_rz = true;
141#else
142 constexpr bool is_rz = false;
143#endif
144
146 {AMREX_D_DECL(geom[lev].CellSize(0),
147 geom[lev].CellSize(1),
148 geom[lev].CellSize(2))};
149
150
151 if (!eb_enabled && !is_rz) {
152 // Determine whether to use semi-coarsening
153 int max_semicoarsening_level = 0;
154 int semicoarsening_direction = -1;
155 const auto min_dir = static_cast<int>(std::distance(dx.begin(),
156 std::min_element(dx.begin(), dx.end())));
157 const auto max_dir = static_cast<int>(std::distance(dx.begin(),
158 std::max_element(dx.begin(), dx.end())));
159 if (dx[max_dir] > dx[min_dir]) {
160 semicoarsening_direction = max_dir;
161 max_semicoarsening_level = static_cast<int>(std::log2(dx[max_dir] / dx[min_dir]));
162 }
163 if (max_semicoarsening_level > 0) {
164 info.setSemicoarsening(true);
165 info.setMaxSemicoarseningLevel(max_semicoarsening_level);
166 info.setSemicoarseningDirection(semicoarsening_direction);
167 }
168 }
169
170 amrex::MLEBNodeFDLaplacian linopx, linopy, linopz;
171 if (eb_enabled) {
172#ifdef AMREX_USE_EB
173 linopx.define({geom[lev]}, {grids[lev]}, {dmap[lev]}, info, {eb_farray_box_factory.value()[lev]});
174 linopy.define({geom[lev]}, {grids[lev]}, {dmap[lev]}, info, {eb_farray_box_factory.value()[lev]});
175 linopz.define({geom[lev]}, {grids[lev]}, {dmap[lev]}, info, {eb_farray_box_factory.value()[lev]});
176#endif
177 } else {
178 linopx.define({geom[lev]}, {grids[lev]}, {dmap[lev]}, info);
179 linopy.define({geom[lev]}, {grids[lev]}, {dmap[lev]}, info);
180 linopz.define({geom[lev]}, {grids[lev]}, {dmap[lev]}, info);
181 }
182
183 amrex::Array<amrex::MLEBNodeFDLaplacian*,3> linop = {&linopx,&linopy,&linopz};
185
186 for (int adim=0; adim<3; adim++) {
187 // Solve the Poisson equation
188 // This is solving the self fields using the magnetostatic solver in the lab frame
189
190 // Note: this assumes that beta is zero
191 linop[adim]->setSigma({AMREX_D_DECL(1._rt, 1._rt, 1._rt)});
192
193 // Set Homogeneous Dirichlet Boundary on EB
194#if defined(AMREX_USE_EB)
195 if (eb_enabled) { linop[adim]->setEBDirichlet(0_rt); }
196#endif
197
198#ifdef WARPX_DIM_RZ
199 linop[adim]->setRZ(true);
200#endif
201
202 linop[adim]->setDomainBC( boundary_handler.lobc[adim], boundary_handler.hibc[adim] );
203
204 mlmg[adim] = std::make_unique<amrex::MLMG>(*linop[adim]);
205
206 mlmg[adim]->setVerbose(verbosity);
207 mlmg[adim]->setMaxIter(max_iters);
208 mlmg[adim]->setConvergenceNormType(always_use_bnorm ? amrex::MLMGNormType::bnorm : amrex::MLMGNormType::greater);
209
210 // Solve Poisson equation at lev
211 mlmg[adim]->solve( {A[lev][adim]}, {curr[lev][adim]},
212 relative_tolerance, absolute_tolerance );
213
214 // Synchronize the ghost cells, do halo exchange
216 *A[lev][adim], A[lev][adim]->nGrowVect(),
217 do_single_precision_comms,
218 geom[lev].periodicity(), false);
219
220 // needed for solving the levels by levels:
221 // - coarser level is initial guess for finer level
222 // - coarser level provides boundary values for finer level patch
223 // Interpolation from phi[lev] to phi[lev+1]
224 // (This provides both the boundary conditions and initial guess for phi[lev+1])
225 if (lev < finest_level) {
226
227 // Allocate A_cp for lev+1
228 amrex::BoxArray ba = A[lev+1][adim]->boxArray();
229 const amrex::IntVect& refratio = rel_ref_ratio.value()[lev];
230 ba.coarsen(refratio);
231 const int ncomp = linop[adim]->getNComp();
232 amrex::MultiFab A_cp(ba, A[lev+1][adim]->DistributionMap(), ncomp, 1);
233
234 // Copy from A[lev] to A_cp (in parallel)
236 const amrex::Periodicity& crse_period = geom[lev].periodicity();
237
239 A_cp,
240 *A[lev][adim],
241 0,
242 0,
243 1,
244 ng,
245 ng,
246 do_single_precision_comms,
247 crse_period
248 );
249
250 // Local interpolation from A_cp to A[lev+1]
251#ifdef AMREX_USE_OMP
252#pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
253#endif
254 for (amrex::MFIter mfi(*A[lev+1][adim],amrex::TilingIfNotGPU()); mfi.isValid(); ++mfi)
255 {
256 amrex::Array4<amrex::Real> const& A_fp_arr = A[lev+1][adim]->array(mfi);
257 amrex::Array4<amrex::Real> const& A_cp_arr = A_cp.array(mfi);
258
259 details::PoissonInterpCPtoFP const interp(A_fp_arr, A_cp_arr, refratio);
260
261 amrex::Box const b = mfi.tilebox(A[lev + 1][adim]->ixType().toIntVect());
262 amrex::ParallelFor(b, interp);
263 }
264
265 }
266
267 // Unscale current
268 curr[lev][adim]->mult(-1._rt/ablastr::constant::SI::mu0);
269 } // Loop over adim
270 // Run additional operations, such as calculation of the B fields for embedded boundaries
271 if constexpr (!std::is_same_v<T_PostACalculationFunctor, std::nullopt_t>) {
272 if (post_A_calculation.has_value()) {
273 post_A_calculation.value()(mlmg, lev);
274 }
275 }
276 } // loop over lev(els)
277}
278} // namepace Magnetostatic
279#endif
#define AMREX_D_DECL(a, b, c)
#define ABLASTR_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition TextMsg.H:75
BoxArray & coarsen(int refinement_ratio)
Array4< typename FabArray< FAB >::value_type const > array(const MFIter &mfi) const noexcept
__host__ static __device__ constexpr IntVectND< dim > TheUnitVector() noexcept
void define(const Vector< Geometry > &a_geom, const Vector< BoxArray > &a_grids, const Vector< DistributionMapping > &a_dmap, const LPInfo &a_info, const Vector< EBFArrayBoxFactory const * > &a_factory)
constexpr auto mu0
vacuum permeability: magnetic permeability of vacuum = 4.0e-7 * pi [H/m]
Definition constant.H:155
Definition EffectivePotentialPoissonSolver.H:63
void computeVectorPotential(amrex::Vector< amrex::Array< amrex::MultiFab *, 3 > > const &curr, amrex::Vector< amrex::Array< amrex::MultiFab *, 3 > > &A, 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, T_BoundaryHandler const boundary_handler, bool eb_enabled=false, bool do_single_precision_comms=false, std::optional< amrex::Vector< amrex::IntVect > > rel_ref_ratio=std::nullopt, T_PostACalculationFunctor post_A_calculation=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 VectorPoissonSolver.H:87
void FillBoundary(amrex::MultiFab &mf, amrex::IntVect ng, bool do_single_precision_comms, const amrex::Periodicity &period, std::optional< bool > nodal_sync)
Definition Communication.cpp:71
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
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)
BoxND< 3 > Box
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