Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_OpenBC.H
Go to the documentation of this file.
1#ifndef AMREX_OPENBC_H_
2#define AMREX_OPENBC_H_
3#include <AMReX_Config.H>
4
5#include <AMReX_MLMG.H>
6#include <AMReX_MLPoisson.H>
7
8namespace amrex
9{
10
18namespace openbc {
19
21
22 static constexpr int M = 7; // highest order of moments
23 static constexpr int P = 3;
24
25 struct Moments
26 {
27 using array_type = GpuArray<Real,(M+2)*(M+1)/2>;
28 array_type mom;
29 Real x, y, z;
30 Orientation face;
31 };
32
33 struct MomTag
34 {
35 Array4<Real const> gp;
36 Box b2d;
37 Orientation face;
38 int offset;
39 };
40
42
43 std::ostream& operator<< (std::ostream& os, Moments const& mom);
44}
45
46#if defined(AMREX_USE_CUDA) || defined(AMREX_USE_HIP)
47template<>
48struct Gpu::SharedMemory<openbc::Moments::array_type>
49{
50 AMREX_GPU_DEVICE openbc::Moments::array_type* dataPtr () noexcept {
51 AMREX_HIP_OR_CUDA(HIP_DYNAMIC_SHARED(openbc::Moments::array_type,amrex_openbc_momarray);,
52 extern __shared__ openbc::Moments::array_type amrex_openbc_momarray[];)
53 return amrex_openbc_momarray;
54 }
55};
56#endif
57
70{
71public:
73 OpenBCSolver () = default;
74
83 OpenBCSolver (const Vector<Geometry>& a_geom,
84 const Vector<BoxArray>& a_grids,
85 const Vector<DistributionMapping>& a_dmap,
86 const LPInfo& a_info = LPInfo());
87
88 ~OpenBCSolver () = default;
89
90 OpenBCSolver (const OpenBCSolver&) = delete;
94
105 void define (const Vector<Geometry>& a_geom,
106 const Vector<BoxArray>& a_grids,
107 const Vector<DistributionMapping>& a_dmap,
108 const LPInfo& a_info = LPInfo());
109
111 void setVerbose (int v) noexcept;
113 void setBottomVerbose (int v) noexcept;
114
123 void useHypre (bool use_hypre) noexcept;
124
136 Real solve (const Vector<MultiFab*>& a_sol, const Vector<MultiFab const*>& a_rhs,
137 Real a_tol_rel, Real a_tol_abs);
138
139// public for cuda
157 void interpolate_potential (MultiFab& solg);
158
159private:
160
161#ifdef AMREX_USE_MPI
162 void bcast_moments (Gpu::DeviceVector<openbc::Moments>& moments);
163#endif
164
165 int m_verbose = 0;
166 int m_bottom_verbose = 0;
167 Vector<Geometry> m_geom;
168 Vector<BoxArray> m_grids;
170 LPInfo m_info;
171 std::unique_ptr<MLPoisson> m_poisson_1;
172 std::unique_ptr<MLPoisson> m_poisson_2;
173 std::unique_ptr<MLMG> m_mlmg_1;
174 std::unique_ptr<MLMG> m_mlmg_2;
175 BottomSolver m_bottom_solver_type = BottomSolver::bicgstab;
176
177 int m_coarsen_ratio = 0;
180#ifdef AMREX_USE_GPU
182 Gpu::PinnedVector<int> m_ngpublocks_h;
183 Gpu::DeviceVector<int> m_ngpublocks_d;
184 int m_nthreads_momtag;
185#endif
186
187 int m_nblocks_local = 0;
188 int m_nblocks = 0;
189#ifdef AMREX_USE_MPI
190 Vector<int> m_countvec;
191 Vector<int> m_offset;
192#endif
193
194 IntVect m_ngrowdomain;
195 MultiFab m_crse_grown_faces_phi;
196 MultiFab m_phind;
197 BoxArray m_bag;
198
199 Vector<IntVect> m_box_offset;
200 Vector<BoxArray> m_ba_all;
202 Vector<Geometry> m_geom_all;
203};
204
205}
206
207#endif
#define AMREX_HIP_OR_CUDA(a, b)
Definition AMReX_GpuControl.H:21
#define AMREX_GPU_DEVICE
Definition AMReX_GpuQualifiers.H:18
Array4< int const > offset
Definition AMReX_HypreMLABecLap.cpp:1139
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:568
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:40
Open Boundary Poisson Solver.
Definition AMReX_OpenBC.H:70
~OpenBCSolver()=default
OpenBCSolver()=default
Construct an empty solver; call define() before solving.
OpenBCSolver & operator=(const OpenBCSolver &)=delete
void interpolate_potential(MultiFab &solg)
Interpolate coarse potential values back onto the grown solution domain.
Definition AMReX_OpenBC.cpp:803
void setVerbose(int v) noexcept
Control verbosity (0 = silent) via v.
Definition AMReX_OpenBC.cpp:197
Real solve(const Vector< MultiFab * > &a_sol, const Vector< MultiFab const * > &a_rhs, Real a_tol_rel, Real a_tol_abs)
Solve the isolated Poisson problem.
Definition AMReX_OpenBC.cpp:218
void useHypre(bool use_hypre) noexcept
Toggle the Hypre bottom solver (requires Hypre-enabled build).
Definition AMReX_OpenBC.cpp:207
OpenBCSolver(const OpenBCSolver &)=delete
void compute_potential(Gpu::DeviceVector< openbc::Moments > const &moments)
Evaluate the coarse potential contribution from precomputed moments.
Definition AMReX_OpenBC.cpp:738
OpenBCSolver(OpenBCSolver &&)=delete
void compute_moments(Gpu::DeviceVector< openbc::Moments > &moments)
Build per-face multipole moments that summarize solution derivatives.
Definition AMReX_OpenBC.cpp:396
void setBottomVerbose(int v) noexcept
Control bottom solver verbosity via v.
Definition AMReX_OpenBC.cpp:202
void define(const Vector< Geometry > &a_geom, const Vector< BoxArray > &a_grids, const Vector< DistributionMapping > &a_dmap, const LPInfo &a_info=LPInfo())
Initialize/refresh the solver.
Definition AMReX_OpenBC.cpp:16
Dynamically allocated vector for trivially copyable data.
Definition AMReX_PODVector.H:308
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
std::array< T, N > Array
Definition AMReX_Array.H:26
std::ostream & operator<<(std::ostream &os, Moments const &mom)
Definition AMReX_OpenBC.cpp:880
Definition AMReX_Amr.cpp:49
BoxND< 3 > Box
Box is an alias for amrex::BoxND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:30
BottomSolver
Definition AMReX_MLLinOp.H:31
__device__ openbc::Moments::array_type * dataPtr() noexcept
Definition AMReX_OpenBC.H:50
Definition AMReX_GpuMemory.H:125
Definition AMReX_MLLinOp.H:36