Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_HypreMLABecLap.H
Go to the documentation of this file.
1#ifndef AMREX_HYPRE_ML_ABECLAP_H_
2#define AMREX_HYPRE_ML_ABECLAP_H_
3#include <AMReX_Config.H>
4
5#include <AMReX_Geometry.H>
6#include <AMReX_LO_BCTYPES.H>
7#include <AMReX_MLMGBndry.H>
8#include <AMReX_MultiFab.H>
9#include <AMReX_MultiMask.H>
10
11#include <HYPRE_sstruct_ls.h>
12
13#include <limits>
14#include <utility>
15
23namespace amrex {
24
25enum struct HypreSolverID {
27};
28
29// single component only, cell centered only
30
35{
36public:
37 HypreMLABecLap (HypreMLABecLap const&) = delete;
41
52 Vector<BoxArray> a_grids,
54 HypreSolverID a_hypre_solver_id,
55 std::string a_parmparse_prefix = "hypre_mlabeclap");
56
58
60 void setVerbose (int v) { m_verbose = v; }
62 void setMaxIter (int v) { m_maxiter = v; }
64 void setIsSingular (bool v) { m_is_singular = v; }
65
78 void setup (Real a_ascalar, Real a_bscalar,
79 Vector<MultiFab const*> const& a_acoefs,
83 Vector<MultiFab const*> const& a_levelbcdata,
84 std::pair<MultiFab const*, IntVect> const& a_coarse_bc = {nullptr, IntVect(0)});
85
94 void solve (Vector<MultiFab*> const& a_sol, Vector<MultiFab const*> const& a_rhs,
95 Real a_reltol, Real a_abstol);
96
97 // update? updateDirichleBC? or updateCoeffs?
98
99 // public for cuda
100
107 void commBCoefs (int flev, Array<MultiFab const*,AMREX_SPACEDIM> const& a_bcoefs);
115 void commBCoefs_local (int flev, Array<MultiFab const*,AMREX_SPACEDIM> const& a_bcoefs,
116 Vector<FabArrayBase::CopyComTag> const& tags);
117
118private:
119
120 void addNonStencilEntriesToGraph ();
121
122 int m_verbose = 0;
123 int m_maxiter = 200;
124 bool m_is_singular = false;
125
126 Vector<Geometry> m_geom;
127 Vector<BoxArray> m_grids;
128 Vector<DistributionMapping> m_dmap;
129 std::string m_parmparse_prefix;
130 int m_nlevels = 0;
131 MPI_Comm m_comm = MPI_COMM_NULL;
132
133 Vector<IntVect> m_ref_ratio;
134
135 Real m_ascalar = std::numeric_limits<Real>::max();
136 Real m_bscalar = std::numeric_limits<Real>::max();
137 Array<LinOpBCType,AMREX_SPACEDIM> m_lobc{AMREX_D_DECL(LinOpBCType::bogus,
140 Array<LinOpBCType,AMREX_SPACEDIM> m_hibc{AMREX_D_DECL(LinOpBCType::bogus,
143
144 Vector<std::unique_ptr<MLMGBndry>> m_bndry;
145 Vector<std::unique_ptr<BndryRegister>> m_bndry_rhs;
146 Vector<iMultiFab> m_fine_masks;
147 Vector<iMultiFab> m_crse_masks;
148
149 // For coarse cells at coarse/fine interface. The vector is for AMR
150 // levels.
151 Vector<iMultiFab> m_c2f_offset_from; // offset for sparse coarse from-cells
152 Vector<LayoutData<int>> m_c2f_total_from; // # of coarse from-cells w/ c2f entries
153 Vector<iMultiFab> m_c2f_nentries; // # of non-stencil entries
154 Vector<iMultiFab> m_c2f_offset_to; // offset for sparse to-cells, including fine (and coarse in 3d) cells
155 Vector<LayoutData<int>> m_c2f_total_to; // total sum of non-stencil entries in a Box
156
157 // B coefficients at coarse/fine interface
158 Vector<Array<iMultiFab,AMREX_SPACEDIM>> m_offset_cf_bcoefs;
159 Vector<Array<LayoutData<std::unique_ptr<Gpu::DeviceVector<Real>>>,AMREX_SPACEDIM>> m_cf_bcoefs;
160
161#ifdef AMREX_USE_GPU
162 template <class T> using HostVector = Gpu::PinnedVector<T>;
163#else
164 template <class T> using HostVector = Vector<T>;
165#endif
166
167 // For fine cells at coarse/fine interface. The non-stencil entries are
168 // from fine to coarse. The outer vector is for AMR levels.
169 Vector<HostVector<int>> m_f2c_bno; // local box number
170 Vector<HostVector<IntVect>> m_f2c_cell; // fine cell
171 Vector<Vector<HYPRE_Int>> m_f2c_nentries; // # of non-stencil entries
172 Vector<HostVector<std::size_t>> m_f2c_offset; // offset into m_f2c_values
173 Vector<HostVector<Real>> m_f2c_values; // values for non-stencil entries
174
175 HYPRE_SStructGrid m_ss_grid = nullptr;
176 HYPRE_SStructStencil m_ss_stencil = nullptr;
177 HYPRE_SStructGraph m_ss_graph = nullptr;
178 HYPRE_SStructSolver m_ss_solver = nullptr;
179 HYPRE_SStructSolver m_ss_precond = nullptr;
180 HYPRE_SStructMatrix m_ss_A = nullptr;
181 HYPRE_SStructVector m_ss_x = nullptr;
182 HYPRE_SStructVector m_ss_b = nullptr;
183
184 HYPRE_Solver m_solver = nullptr;
185
186 HypreSolverID m_hypre_solver_id = HypreSolverID::BoomerAMG;
187 HYPRE_Int m_hypre_object_type = HYPRE_PARCSR;
188};
189
190}
191
192#endif
Constants describing linear-operator boundary condition types.
MultiFab-style container that stores Mask objects per grid.
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
Multilevel HYPRE solver that mirrors MLMG's ABecLaplacian structure.
Definition AMReX_HypreMLABecLap.H:35
void solve(Vector< MultiFab * > const &a_sol, Vector< MultiFab const * > const &a_rhs, Real a_reltol, Real a_abstol)
Solve the assembled system with HYPRE (single-component).
Definition AMReX_HypreMLABecLap.cpp:947
void setMaxIter(int v)
Override the maximum number of HYPRE iterations (v).
Definition AMReX_HypreMLABecLap.H:62
HypreMLABecLap(HypreMLABecLap const &)=delete
void commBCoefs_local(int flev, Array< MultiFab const *, 3 > const &a_bcoefs, Vector< FabArrayBase::CopyComTag > const &tags)
Local portion of commBCoefs used to issue device-friendly copies.
Definition AMReX_HypreMLABecLap.cpp:1499
void setVerbose(int v)
Set HYPRE's logging verbosity (v).
Definition AMReX_HypreMLABecLap.H:60
~HypreMLABecLap()
Definition AMReX_HypreMLABecLap.cpp:176
void commBCoefs(int flev, Array< MultiFab const *, 3 > const &a_bcoefs)
Communicate b coefficients needed by level flev at coarse/fine interfaces.
Definition AMReX_HypreMLABecLap.cpp:1177
void setIsSingular(bool v)
Flag whether the assembled operator is singular (v).
Definition AMReX_HypreMLABecLap.H:64
void setup(Real a_ascalar, Real a_bscalar, Vector< MultiFab const * > const &a_acoefs, Vector< Array< MultiFab const *, 3 > > const &a_bcoefs, Array< LinOpBCType, 3 > const &a_lobc, Array< LinOpBCType, 3 > const &a_hibc, Vector< MultiFab const * > const &a_levelbcdata, std::pair< MultiFab const *, IntVect > const &a_coarse_bc={nullptr, IntVect(0)})
Assemble the HYPRE operators/storage for the provided coefficients.
Definition AMReX_HypreMLABecLap.cpp:514
HypreMLABecLap(HypreMLABecLap &&)=delete
HypreMLABecLap & operator=(HypreMLABecLap const &)=delete
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
int MPI_Comm
Definition AMReX_ccse-mpi.H:51
static constexpr int MPI_COMM_NULL
Definition AMReX_ccse-mpi.H:59
Definition AMReX_Amr.cpp:49
IntVectND< 3 > IntVect
IntVect is an alias for amrex::IntVectND instantiated with AMREX_SPACEDIM.
Definition AMReX_BaseFwd.H:33
HypreSolverID
Definition AMReX_HypreMLABecLap.H:25