7 #ifndef ABLASTR_POISSON_SOLVER_H
8 #define ABLASTR_POISSON_SOLVER_H
21 #if defined(AMREX_USE_EB) || defined(WARPX_DIM_RZ)
36 #include <AMReX_Config.H>
127 typename T_BoundaryHandler,
128 typename T_PostPhiCalculationFunctor = std::nullopt_t,
129 typename T_FArrayBoxFactory =
void
134 std::array<amrex::Real, 3>
const beta,
135 amrex::Real
const relative_tolerance,
136 amrex::Real absolute_tolerance,
142 T_BoundaryHandler
const boundary_handler,
143 bool const do_single_precision_comms =
false,
145 [[maybe_unused]] T_PostPhiCalculationFunctor post_phi_calculation = std::nullopt,
146 [[maybe_unused]] std::optional<amrex::Real const> current_time = std::nullopt,
150 using namespace amrex::literals;
152 if (!rel_ref_ratio.has_value()) {
154 "rel_ref_ratio must be set if mesh-refinement is used");
158 int const finest_level = rho.
size() - 1u;
161 amrex::Real max_norm_b = 0.0;
162 for (
int lev=0; lev<=finest_level; lev++) {
164 max_norm_b =
amrex::max(max_norm_b, rho[lev]->norm0());
168 const bool always_use_bnorm = (max_norm_b > 0);
169 if (!always_use_bnorm) {
170 if (absolute_tolerance == 0.0) absolute_tolerance = amrex::Real(1e-6);
172 "ElectrostaticSolver",
173 "Max norm of rho is 0",
179 for (
int lev=0; lev<=finest_level; lev++) {
182 #if defined(WARPX_DIM_1D_Z)
184 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
185 {{ beta[0], beta[2] }};
187 {{ beta[0], beta[1], beta[2] }};
190 #if !(defined(AMREX_USE_EB) || defined(WARPX_DIM_RZ))
193 {
AMREX_D_DECL(geom[lev].CellSize(0)/std::sqrt(1._rt-beta_solver[0]*beta_solver[0]),
194 geom[lev].CellSize(1)/std::sqrt(1._rt-beta_solver[1]*beta_solver[1]),
195 geom[lev].CellSize(2)/std::sqrt(1._rt-beta_solver[2]*beta_solver[2]))};
196 int max_semicoarsening_level = 0;
197 int semicoarsening_direction = -1;
198 const int min_dir = std::distance(dx_scaled.begin(),
199 std::min_element(dx_scaled.begin(),dx_scaled.end()));
200 const int max_dir = std::distance(dx_scaled.begin(),
201 std::max_element(dx_scaled.begin(),dx_scaled.end()));
202 if (dx_scaled[max_dir] > dx_scaled[min_dir]) {
203 semicoarsening_direction = max_dir;
204 max_semicoarsening_level =
static_cast<int>
205 (std::log2(dx_scaled[max_dir]/dx_scaled[min_dir]));
207 if (max_semicoarsening_level > 0) {
214 #if defined(AMREX_USE_EB) || defined(WARPX_DIM_RZ)
219 #if defined(AMREX_USE_EB)
220 , {eb_farray_box_factory.value()[lev]}
227 #if defined(WARPX_DIM_RZ)
228 linop.
setSigma({0._rt, 1._rt-beta_solver[1]*beta_solver[1]});
231 1._rt-beta_solver[0]*beta_solver[0],
232 1._rt-beta_solver[1]*beta_solver[1],
233 1._rt-beta_solver[2]*beta_solver[2])});
236 #if defined(AMREX_USE_EB)
239 if (boundary_handler.phi_EB_only_t) {
240 linop.setEBDirichlet(boundary_handler.potential_eb_t(current_time.value()));
243 linop.setEBDirichlet(boundary_handler.getPhiEB(current_time.value()));
254 linop.setDomainBC( boundary_handler.lobc, boundary_handler.hibc );
264 mlmg.
solve( {phi[lev]}, {rho[lev]},
265 relative_tolerance, absolute_tolerance );
272 if (lev < finest_level) {
278 const int ncomp = linop.getNComp();
293 do_single_precision_comms,
299 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
307 amrex::Box const b = mfi.tilebox(phi[lev + 1]->ixType().toIntVect());
315 if (post_phi_calculation.has_value())
316 post_phi_calculation.value()(mlmg, lev);
#define AMREX_FORCE_INLINE
#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< FArrayBox >::value_type const > array(const MFIter &mfi) const noexcept
AMREX_GPU_HOST_DEVICE static constexpr AMREX_FORCE_INLINE IntVect TheUnitVector() noexcept
void setSigma(Array< Real, AMREX_SPACEDIM > const &a_sigma) noexcept
void setVerbose(int v) noexcept
void setAlwaysUseBNorm(int flag) 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
void setBeta(Array< Real, AMREX_SPACEDIM > const &a_beta) noexcept
Long size() const 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 ep0
vacuum permittivity: dielectric permittivity of vacuum [F/m]
Definition: constant.H:46
Definition: PoissonSolver.H:59
void computePhi(amrex::Vector< amrex::MultiFab * > const &rho, amrex::Vector< amrex::MultiFab * > &phi, std::array< amrex::Real, 3 > const beta, 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_PostPhiCalculationFunctor post_phi_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: PoissonSolver.H:132
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(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:310
void ReduceRealMax(Vector< std::reference_wrapper< Real > > &&)
AMREX_GPU_HOST_DEVICE constexpr AMREX_FORCE_INLINE const T & max(const T &a, const T &b) noexcept
std::enable_if_t< std::is_integral< T >::value > ParallelFor(TypeList< CTOs... >, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE void mf_nodebilin_interp(int i, int, int, int n, Array4< Real > const &fine, int fcomp, Array4< Real const > const &crse, int ccomp, IntVect const &ratio) noexcept
bool TilingIfNotGPU() noexcept
i
Definition: check_interp_points_and_weights.py:174
value
Definition: updateAMReX.py:141
Definition: PoissonSolver.H:74
PoissonInterpCPtoFP(amrex::Array4< amrex::Real > const phi_fp_arr, amrex::Array4< amrex::Real const > const phi_cp_arr, amrex::IntVect const refratio)
Definition: PoissonSolver.H:75
amrex::Array4< amrex::Real > const m_phi_fp_arr
Definition: PoissonSolver.H:90
amrex::Array4< amrex::Real const > const m_phi_cp_arr
Definition: PoissonSolver.H:91
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void operator()(long i, long j, long k) const noexcept
Definition: PoissonSolver.H:84
amrex::IntVect const m_refratio
Definition: PoissonSolver.H:92
LPInfo & setSemicoarsening(bool x) noexcept
LPInfo & setMaxSemicoarseningLevel(int n) noexcept
LPInfo & setSemicoarseningDirection(int n) noexcept