Block-Structured AMR Software Framework
AMReX_HypreABecLap3.H
Go to the documentation of this file.
1 #ifndef AMREX_HYPREABECLAP3_H_
2 #define AMREX_HYPREABECLAP3_H_
3 #include <AMReX_Config.H>
4 
5 #include <AMReX_Hypre.H>
6 
7 #include <AMReX_iMultiFab.H>
8 #include <AMReX_LayoutData.H>
9 #include <AMReX_HypreIJIface.H>
10 
11 #include <algorithm>
12 
13 namespace amrex
14 {
15 
16 class HypreABecLap3 final
17  : public Hypre
18 {
19 public:
20 
21  HypreABecLap3 (const BoxArray& grids, const DistributionMapping& dmap,
22  const Geometry& geom_, MPI_Comm comm_,
23  const iMultiFab* overset_mask = nullptr);
24  ~HypreABecLap3 () override;
25 
26  HypreABecLap3 (HypreABecLap3 const&) = delete;
30 
31  void solve (MultiFab& soln, const MultiFab& rhs, Real rel_tol, Real abs_tol,
32  int max_iter, const BndryData& bndry, int max_bndry_order) override;
33 
34 #ifdef AMREX_USE_EB
35  void setEBDirichlet (MultiFab const* beb) { m_eb_b_coeffs = beb; }
36 #endif
37 
38 private :
39  std::unique_ptr<HypreIJIface> hypre_ij;
40 
41  // Non-owning references to hypre matrix, rhs, and solution data
42  HYPRE_IJMatrix A = nullptr;
43  HYPRE_IJVector b = nullptr;
44  HYPRE_IJVector x = nullptr;
45 
49 
50  MultiFab const* m_eb_b_coeffs = nullptr;
51 
52  iMultiFab const* m_overset_mask = nullptr;
53 
54 public: // for CUDA
55  void prepareSolver ();
56  void getSolution (MultiFab& soln);
57  void loadVectors (MultiFab& soln, const MultiFab& rhs);
58 };
59 
60 }
61 
62 #endif
int MPI_Comm
Definition: AMReX_ccse-mpi.H:47
A BndryData stores and manipulates boundary data information on each side of each box in a BoxArray.
Definition: AMReX_BndryData.H:41
A collection of Boxes stored in an Array.
Definition: AMReX_BoxArray.H:549
Calculates the distribution of FABs to MPI processes.
Definition: AMReX_DistributionMapping.H:41
An Array of FortranArrayBox(FAB)-like Objects.
Definition: AMReX_FabArray.H:344
Rectangular problem domain geometry.
Definition: AMReX_Geometry.H:73
Definition: AMReX_HypreABecLap3.H:18
FabArray< BaseFab< HYPRE_Int > > cell_id_vec
Definition: AMReX_HypreABecLap3.H:48
void loadVectors(MultiFab &soln, const MultiFab &rhs)
Definition: AMReX_HypreABecLap3.cpp:511
LayoutData< HYPRE_Int > ncells_grid
Definition: AMReX_HypreABecLap3.H:46
HYPRE_IJMatrix A
Definition: AMReX_HypreABecLap3.H:42
void getSolution(MultiFab &soln)
Definition: AMReX_HypreABecLap3.cpp:61
FabArray< BaseFab< HYPRE_Int > > cell_id
Definition: AMReX_HypreABecLap3.H:47
std::unique_ptr< HypreIJIface > hypre_ij
Definition: AMReX_HypreABecLap3.H:39
iMultiFab const * m_overset_mask
Definition: AMReX_HypreABecLap3.H:52
void prepareSolver()
Definition: AMReX_HypreABecLap3.cpp:88
void solve(MultiFab &soln, const MultiFab &rhs, Real rel_tol, Real abs_tol, int max_iter, const BndryData &bndry, int max_bndry_order) override
Definition: AMReX_HypreABecLap3.cpp:30
MultiFab const * m_eb_b_coeffs
Definition: AMReX_HypreABecLap3.H:50
HypreABecLap3 & operator=(HypreABecLap3 const &)=delete
~HypreABecLap3() override
HYPRE_IJVector b
Definition: AMReX_HypreABecLap3.H:43
HypreABecLap3(HypreABecLap3 &&)=delete
HYPRE_IJVector x
Definition: AMReX_HypreABecLap3.H:44
HypreABecLap3(HypreABecLap3 const &)=delete
HypreABecLap3(const BoxArray &grids, const DistributionMapping &dmap, const Geometry &geom_, MPI_Comm comm_, const iMultiFab *overset_mask=nullptr)
Definition: AMReX_HypreABecLap3.cpp:20
Definition: AMReX_Hypre.H:18
A collection (stored as an array) of FArrayBox objects.
Definition: AMReX_MultiFab.H:38
Definition: AMReX_iMultiFab.H:32
Definition: AMReX_Amr.cpp:49