7 #ifndef ABLASTR_VECTOR_POISSON_SOLVER_H
8 #define ABLASTR_VECTOR_POISSON_SOLVER_H
36 #include <AMReX_Config.H>
90 typename T_BoundaryHandler,
91 typename T_PostACalculationFunctor = std::nullopt_t,
92 typename T_FArrayBoxFactory =
void
97 amrex::Real
const relative_tolerance,
98 amrex::Real absolute_tolerance,
104 T_BoundaryHandler
const boundary_handler,
105 bool const do_single_precision_comms =
false,
107 [[maybe_unused]] T_PostACalculationFunctor post_A_calculation = std::nullopt,
108 [[maybe_unused]] std::optional<amrex::Real const> current_time = std::nullopt,
112 using namespace amrex::literals;
114 if (!rel_ref_ratio.has_value()) {
116 "rel_ref_ratio must be set if mesh-refinement is used");
120 int const finest_level = curr.
size() - 1u;
123 amrex::Real max_comp_J = 0.0;
124 for (
int lev=0; lev<=finest_level; lev++) {
125 for (
int adim=0; adim<3; adim++) {
127 max_comp_J =
amrex::max(max_comp_J, curr[lev][adim]->norm0());
132 bool always_use_bnorm = (max_comp_J > 0);
133 if (!always_use_bnorm) {
134 if (absolute_tolerance == 0.0) absolute_tolerance = amrex::Real(1e-6);
136 "MagnetostaticSolver",
137 "Max norm of J is 0",
145 for (
int lev=0; lev<=finest_level; lev++) {
147 {geom[lev]}, {grids[lev]}, {dmap[lev]}, info
148 #if defined(AMREX_USE_EB)
149 , {eb_farray_box_factory.value()[lev]}
153 {geom[lev]}, {grids[lev]}, {dmap[lev]}, info
154 #if defined(AMREX_USE_EB)
155 , {eb_farray_box_factory.value()[lev]}
159 {geom[lev]}, {grids[lev]}, {dmap[lev]}, info
160 #if defined(AMREX_USE_EB)
161 , {eb_farray_box_factory.value()[lev]}
168 for (
int adim=0; adim<3; adim++) {
175 #if defined(AMREX_USE_EB)
177 linop[adim]->setEBDirichlet(0_rt);
181 linop[adim]->setRZ(
true);
184 linop[adim]->setDomainBC( boundary_handler.lobc[adim], boundary_handler.hibc[adim] );
186 mlmg[adim] = std::make_unique<amrex::MLMG>(*linop[adim]);
188 mlmg[adim]->setVerbose(verbosity);
189 mlmg[adim]->setMaxIter(max_iters);
190 mlmg[adim]->setConvergenceNormType(always_use_bnorm ? amrex::MLMGNormType::bnorm : amrex::MLMGNormType::greater);
193 mlmg[adim]->solve( {A[lev][adim]}, {curr[lev][adim]},
194 relative_tolerance, absolute_tolerance );
198 A[lev][adim]->nGrowVect(),
200 geom[lev].periodicity(),
208 if (lev < finest_level) {
214 const int ncomp = linop[adim]->getNComp();
229 do_single_precision_comms,
235 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
244 amrex::Box const b = mfi.tilebox(A[lev + 1][adim]->ixType().toIntVect());
255 if (post_A_calculation.has_value())
256 post_A_calculation.value()(mlmg, lev);
#define AMREX_D_DECL(a, b, c)
#define ABLASTR_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition: TextMsg.H:77
static bool do_single_precision_comms
perform field communications in single precision
Definition: WarpX.H:358
BoxArray & coarsen(int refinement_ratio)
Array4< typename FabArray< FAB >::value_type const > array(const MFIter &mfi) const noexcept
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE constexpr IntVectND< dim > TheUnitVector() noexcept
AMREX_GPU_HOST_DEVICE static AMREX_FORCE_INLINE constexpr std::size_t size() noexcept
void setSigma(Array< Real, AMREX_SPACEDIM > const &a_sigma) noexcept
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void interp(int j, int k, int l, amrex::Array4< amrex::Real > const &fine, amrex::Array4< amrex::Real const > const &crse, const amrex::IntVect r_ratio, IntVect const &type) noexcept
Definition: Interpolate_K.H:9
static constexpr auto mu0
vacuum permeability: magnetic permeability of vacuum = 4.0e-7 * pi [H/m]
Definition: constant.H:48
Definition: PoissonSolver.H:60
void computeVectorPotential(amrex::Vector< amrex::Array< amrex::MultiFab *, 3 > > const &curr, amrex::Vector< amrex::Array< amrex::MultiFab *, 3 > > &A, amrex::Real const relative_tolerance, amrex::Real absolute_tolerance, int const max_iters, int const 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 const do_single_precision_comms=false, std::optional< amrex::Vector< amrex::IntVect > > rel_ref_ratio=std::nullopt, [[maybe_unused]] T_PostACalculationFunctor post_A_calculation=std::nullopt, [[maybe_unused]] std::optional< amrex::Real const > current_time=std::nullopt, [[maybe_unused]] std::optional< amrex::Vector< T_FArrayBoxFactory const * > > eb_farray_box_factory=std::nullopt)
Definition: VectorPoissonSolver.H:95
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:65
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:22
void WMRecordWarning(std::string topic, std::string text, WarnPriority priority=WarnPriority::medium)
Helper function to abbreviate the call to WarnManager::GetInstance().RecordWarning (recording a warni...
Definition: WarnManager.cpp:306
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)
IntVectND< AMREX_SPACEDIM > IntVect
bool TilingIfNotGPU() noexcept
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE constexpr const T & max(const T &a, const T &b) noexcept
value
Definition: updateAMReX.py:141
Definition: PoissonSolver.H:75