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
16namespace amrex {
17
18enum struct HypreSolverID {
20};
21
22// single component only, cell centered only
23
25{
26public:
27 HypreMLABecLap (HypreMLABecLap const&) = delete;
31
33 Vector<BoxArray> a_grids,
35 HypreSolverID a_hypre_solver_id,
36 std::string a_parmparse_prefix = "hypre_mlabeclap");
37
39
40 void setVerbose (int v) { m_verbose = v; }
41 void setMaxIter (int v) { m_maxiter = v; }
42 void setIsSingular (bool v) { m_is_singular = v; }
43
44 void setup (Real a_ascalar, Real a_bscalar,
45 Vector<MultiFab const*> const& a_acoefs,
49 Vector<MultiFab const*> const& a_levelbcdata,
50 std::pair<MultiFab const*, IntVect> const& a_coarse_bc = {nullptr, IntVect(0)});
51
52 void solve (Vector<MultiFab*> const& a_sol, Vector<MultiFab const*> const& a_rhs,
53 Real a_reltol, Real a_abstol);
54
55 // update? updateDirichleBC? or updateCoeffs?
56
57 // public for cuda
58
59 void commBCoefs (int flev, Array<MultiFab const*,AMREX_SPACEDIM> const& a_bcoefs);
60 void commBCoefs_local (int flev, Array<MultiFab const*,AMREX_SPACEDIM> const& a_bcoefs,
61 Vector<FabArrayBase::CopyComTag> const& tags);
62
63private:
64
65 void addNonStencilEntriesToGraph ();
66
67 int m_verbose = 0;
68 int m_maxiter = 200;
69 bool m_is_singular = false;
70
71 Vector<Geometry> m_geom;
72 Vector<BoxArray> m_grids;
73 Vector<DistributionMapping> m_dmap;
74 std::string m_parmparse_prefix;
75 int m_nlevels = 0;
76 MPI_Comm m_comm = MPI_COMM_NULL;
77
78 Vector<IntVect> m_ref_ratio;
79
80 Real m_ascalar = std::numeric_limits<Real>::max();
81 Real m_bscalar = std::numeric_limits<Real>::max();
82 Array<LinOpBCType,AMREX_SPACEDIM> m_lobc{AMREX_D_DECL(LinOpBCType::bogus,
85 Array<LinOpBCType,AMREX_SPACEDIM> m_hibc{AMREX_D_DECL(LinOpBCType::bogus,
88
89 Vector<std::unique_ptr<MLMGBndry>> m_bndry;
90 Vector<std::unique_ptr<BndryRegister>> m_bndry_rhs;
91 Vector<iMultiFab> m_fine_masks;
92 Vector<iMultiFab> m_crse_masks;
93
94 // For coarse cells at coarse/fine interface. The vector is for AMR
95 // levels.
96 Vector<iMultiFab> m_c2f_offset_from; // offset for sparse coarse from-cells
97 Vector<LayoutData<int>> m_c2f_total_from; // # of coarse from-cells w/ c2f entries
98 Vector<iMultiFab> m_c2f_nentries; // # of non-stencil entries
99 Vector<iMultiFab> m_c2f_offset_to; // offset for sparse to-cells, including fine (and coarse in 3d) cells
100 Vector<LayoutData<int>> m_c2f_total_to; // total sum of non-stencil entries in a Box
101
102 // B coefficients at coarse/fine interface
103 Vector<Array<iMultiFab,AMREX_SPACEDIM>> m_offset_cf_bcoefs;
104 Vector<Array<LayoutData<std::unique_ptr<Gpu::DeviceVector<Real>>>,AMREX_SPACEDIM>> m_cf_bcoefs;
105
106#ifdef AMREX_USE_GPU
107 template <class T> using HostVector = Gpu::PinnedVector<T>;
108#else
109 template <class T> using HostVector = Vector<T>;
110#endif
111
112 // For fine cells at coarse/fine interface. The non-stencil entries are
113 // from fine to coarse. The outer vector is for AMR levels.
114 Vector<HostVector<int>> m_f2c_bno; // local box number
115 Vector<HostVector<IntVect>> m_f2c_cell; // fine cell
116 Vector<Vector<HYPRE_Int>> m_f2c_nentries; // # of non-stencil entries
117 Vector<HostVector<std::size_t>> m_f2c_offset; // offset into m_f2c_values
118 Vector<HostVector<Real>> m_f2c_values; // values for non-stencil entries
119
120 HYPRE_SStructGrid m_ss_grid = nullptr;
121 HYPRE_SStructStencil m_ss_stencil = nullptr;
122 HYPRE_SStructGraph m_ss_graph = nullptr;
123 HYPRE_SStructSolver m_ss_solver = nullptr;
124 HYPRE_SStructSolver m_ss_precond = nullptr;
125 HYPRE_SStructMatrix m_ss_A = nullptr;
126 HYPRE_SStructVector m_ss_x = nullptr;
127 HYPRE_SStructVector m_ss_b = nullptr;
128
129 HYPRE_Solver m_solver = nullptr;
130
131 HypreSolverID m_hypre_solver_id = HypreSolverID::BoomerAMG;
132 HYPRE_Int m_hypre_object_type = HYPRE_PARCSR;
133};
134
135}
136
137#endif
#define AMREX_D_DECL(a, b, c)
Definition AMReX_SPACE.H:171
Definition AMReX_HypreMLABecLap.H:25
void solve(Vector< MultiFab * > const &a_sol, Vector< MultiFab const * > const &a_rhs, Real a_reltol, Real a_abstol)
Definition AMReX_HypreMLABecLap.cpp:907
void setMaxIter(int v)
Definition AMReX_HypreMLABecLap.H:41
HypreMLABecLap(HypreMLABecLap const &)=delete
void commBCoefs_local(int flev, Array< MultiFab const *, 3 > const &a_bcoefs, Vector< FabArrayBase::CopyComTag > const &tags)
Definition AMReX_HypreMLABecLap.cpp:1440
void setVerbose(int v)
Definition AMReX_HypreMLABecLap.H:40
~HypreMLABecLap()
Definition AMReX_HypreMLABecLap.cpp:167
void commBCoefs(int flev, Array< MultiFab const *, 3 > const &a_bcoefs)
Definition AMReX_HypreMLABecLap.cpp:1127
void setIsSingular(bool v)
Definition AMReX_HypreMLABecLap.H:42
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)})
Definition AMReX_HypreMLABecLap.cpp:487
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:18