1#ifndef AMREX_HYPRE_ML_ABECLAP_H_
2#define AMREX_HYPRE_ML_ABECLAP_H_
3#include <AMReX_Config.H>
11#include <HYPRE_sstruct_ls.h>
55 std::string a_parmparse_prefix =
"hypre_mlabeclap");
84 std::pair<MultiFab const*, IntVect>
const& a_coarse_bc = {
nullptr,
IntVect(0)});
94 void solve (Vector<MultiFab*>
const& a_sol, Vector<MultiFab const*>
const& a_rhs,
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);
120 void addNonStencilEntriesToGraph ();
124 bool m_is_singular =
false;
126 Vector<Geometry> m_geom;
127 Vector<BoxArray> m_grids;
128 Vector<DistributionMapping> m_dmap;
129 std::string m_parmparse_prefix;
133 Vector<IntVect> m_ref_ratio;
135 Real m_ascalar = std::numeric_limits<Real>::max();
136 Real m_bscalar = std::numeric_limits<Real>::max();
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;
151 Vector<iMultiFab> m_c2f_offset_from;
152 Vector<LayoutData<int>> m_c2f_total_from;
153 Vector<iMultiFab> m_c2f_nentries;
154 Vector<iMultiFab> m_c2f_offset_to;
155 Vector<LayoutData<int>> m_c2f_total_to;
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;
162 template <
class T>
using HostVector = Gpu::PinnedVector<T>;
164 template <
class T>
using HostVector = Vector<T>;
169 Vector<HostVector<int>> m_f2c_bno;
170 Vector<HostVector<IntVect>> m_f2c_cell;
171 Vector<Vector<HYPRE_Int>> m_f2c_nentries;
172 Vector<HostVector<std::size_t>> m_f2c_offset;
173 Vector<HostVector<Real>> m_f2c_values;
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;
184 HYPRE_Solver m_solver =
nullptr;
187 HYPRE_Int m_hypre_object_type = HYPRE_PARCSR;
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