7 #ifndef ABLASTR_POISSON_SOLVER_H
8 #define ABLASTR_POISSON_SOLVER_H
22 #if defined(AMREX_USE_EB) || defined(WARPX_DIM_RZ)
37 #include <AMReX_Config.H>
128 typename T_BoundaryHandler,
129 typename T_PostPhiCalculationFunctor = std::nullopt_t,
130 typename T_FArrayBoxFactory =
void
135 std::array<amrex::Real, 3>
const beta,
136 amrex::Real
const relative_tolerance,
137 amrex::Real absolute_tolerance,
143 T_BoundaryHandler
const boundary_handler,
144 bool const do_single_precision_comms =
false,
146 [[maybe_unused]] T_PostPhiCalculationFunctor post_phi_calculation = std::nullopt,
147 [[maybe_unused]] std::optional<amrex::Real const> current_time = std::nullopt,
151 using namespace amrex::literals;
153 if (!rel_ref_ratio.has_value()) {
155 "rel_ref_ratio must be set if mesh-refinement is used");
159 int const finest_level = rho.
size() - 1u;
162 amrex::Real max_norm_b = 0.0;
163 for (
int lev=0; lev<=finest_level; lev++) {
165 max_norm_b =
amrex::max(max_norm_b, rho[lev]->norm0());
169 bool always_use_bnorm = (max_norm_b > 0);
170 if (!always_use_bnorm) {
171 if (absolute_tolerance == 0.0) absolute_tolerance = amrex::Real(1e-6);
173 "ElectrostaticSolver",
174 "Max norm of rho is 0",
180 for (
int lev=0; lev<=finest_level; lev++) {
183 #if defined(WARPX_DIM_1D_Z)
185 #elif defined(WARPX_DIM_XZ) || defined(WARPX_DIM_RZ)
186 {{ beta[0], beta[2] }};
188 {{ beta[0], beta[1], beta[2] }};
191 #if !(defined(AMREX_USE_EB) || defined(WARPX_DIM_RZ))
194 {
AMREX_D_DECL(geom[lev].CellSize(0)/std::sqrt(1._rt-beta_solver[0]*beta_solver[0]),
195 geom[lev].CellSize(1)/std::sqrt(1._rt-beta_solver[1]*beta_solver[1]),
196 geom[lev].CellSize(2)/std::sqrt(1._rt-beta_solver[2]*beta_solver[2]))};
197 int max_semicoarsening_level = 0;
198 int semicoarsening_direction = -1;
199 int min_dir = std::distance(dx_scaled.begin(),
200 std::min_element(dx_scaled.begin(),dx_scaled.end()));
201 int max_dir = std::distance(dx_scaled.begin(),
202 std::max_element(dx_scaled.begin(),dx_scaled.end()));
203 if (dx_scaled[max_dir] > dx_scaled[min_dir]) {
204 semicoarsening_direction = max_dir;
205 max_semicoarsening_level =
static_cast<int>
206 (std::log2(dx_scaled[max_dir]/dx_scaled[min_dir]));
208 if (max_semicoarsening_level > 0) {
215 #if defined(AMREX_USE_EB) || defined(WARPX_DIM_RZ)
220 #if defined(AMREX_USE_EB)
221 , {eb_farray_box_factory.value()[lev]}
228 #if defined(WARPX_DIM_RZ)
229 linop.
setSigma({0._rt, 1._rt-beta_solver[1]*beta_solver[1]});
232 1._rt-beta_solver[0]*beta_solver[0],
233 1._rt-beta_solver[1]*beta_solver[1],
234 1._rt-beta_solver[2]*beta_solver[2])});
237 #if defined(AMREX_USE_EB)
240 if (boundary_handler.phi_EB_only_t) {
241 linop.setEBDirichlet(boundary_handler.potential_eb_t(current_time.value()));
244 linop.setEBDirichlet(boundary_handler.getPhiEB(current_time.value()));
255 linop.setDomainBC( boundary_handler.lobc, boundary_handler.hibc );
262 mlmg.
setConvergenceNormType(always_use_bnorm ? amrex::MLMGNormType::bnorm : amrex::MLMGNormType::greater);
265 mlmg.
solve( {phi[lev]}, {rho[lev]},
266 relative_tolerance, absolute_tolerance );
273 if (lev < finest_level) {
279 const int ncomp = linop.getNComp();
294 do_single_precision_comms,
300 #pragma omp parallel if (amrex::Gpu::notInLaunchRegion())
308 amrex::Box const b = mfi.tilebox(phi[lev + 1]->ixType().toIntVect());
316 if (post_phi_calculation.has_value())
317 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:77
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
void setSigma(Array< Real, AMREX_SPACEDIM > const &a_sigma) noexcept
void setVerbose(int v) 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
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:60
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:133
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
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
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE constexpr const T & max(const T &a, const T &b) noexcept
i
Definition: check_interp_points_and_weights.py:174
value
Definition: updateAMReX.py:141
Definition: PoissonSolver.H:75
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:76
amrex::Array4< amrex::Real > const m_phi_fp_arr
Definition: PoissonSolver.H:91
amrex::Array4< amrex::Real const > const m_phi_cp_arr
Definition: PoissonSolver.H:92
AMREX_GPU_DEVICE AMREX_FORCE_INLINE void operator()(long i, long j, long k) const noexcept
Definition: PoissonSolver.H:85
amrex::IntVect const m_refratio
Definition: PoissonSolver.H:93
LPInfo & setSemicoarsening(bool x) noexcept
LPInfo & setSemicoarseningDirection(int n) noexcept
LPInfo & setMaxSemicoarseningLevel(int n) noexcept