Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
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>
10
11#include <algorithm>
12
19namespace amrex
20{
21
25class HypreABecLap3 final
26 : public Hypre
27{
28public:
29
39 HypreABecLap3 (const BoxArray& grids, const DistributionMapping& dmap,
40 const Geometry& geom_, MPI_Comm comm_,
41 const iMultiFab* overset_mask = nullptr);
42 ~HypreABecLap3 () override;
43
44 HypreABecLap3 (HypreABecLap3 const&) = delete;
48
60 void solve (MultiFab& soln, const MultiFab& rhs, Real rel_tol, Real abs_tol,
61 int max_iter, const BndryData& bndry, int max_bndry_order) override;
62
63#ifdef AMREX_USE_EB
65 void setEBDirichlet (MultiFab const* beb) { m_eb_b_coeffs = beb; }
66#endif
67
68private :
69 std::unique_ptr<HypreIJIface> hypre_ij;
70
71 // Non-owning references to HYPRE matrix, rhs, and solution data
72 HYPRE_IJMatrix A = nullptr;
73 HYPRE_IJVector b = nullptr;
74 HYPRE_IJVector x = nullptr;
75
76 LayoutData<HYPRE_Int> ncells_grid;
78 FabArray<BaseFab<HYPRE_Int> > cell_id_vec;
79
80 MultiFab const* m_eb_b_coeffs = nullptr;
81
82 iMultiFab const* m_overset_mask = nullptr;
83
84public: // for CUDA
86 void prepareSolver ();
88 void getSolution (MultiFab& soln);
90 void loadVectors (MultiFab& soln, const MultiFab& rhs);
91};
92
93}
94
95#endif
A BndryData stores and manipulates boundary data information on each side of each box in a BoxArray.
Definition AMReX_BndryData.H:46
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:568
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
An Array of FortranArrayBox(FAB)-like Objects.
Definition AMReX_FabArray.H:349
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:74
HYPRE IJ solver for ABecLaplacian (optionally with overset masks).
Definition AMReX_HypreABecLap3.H:27
void setEBDirichlet(MultiFab const *beb)
Provide b coefficients stored in beb for embedded-boundary faces.
Definition AMReX_HypreABecLap3.H:65
void loadVectors(MultiFab &soln, const MultiFab &rhs)
Load RHS rhs and initial guess soln into HYPRE vectors.
Definition AMReX_HypreABecLap3.cpp:515
void getSolution(MultiFab &soln)
Copy HYPRE solution back into AMReX storage (soln).
Definition AMReX_HypreABecLap3.cpp:61
HypreABecLap3 & operator=(HypreABecLap3 const &)=delete
void prepareSolver()
Build IJ structures and configure the HYPRE solver/preconditioner.
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
Solve with the IJ wrapper (supports overset/EB contributions).
Definition AMReX_HypreABecLap3.cpp:30
~HypreABecLap3() override
HypreABecLap3(HypreABecLap3 &&)=delete
HypreABecLap3(HypreABecLap3 const &)=delete
Lightweight base class shared by HYPRE-backed linear solvers, HypreABecLap, HypreABecLap2,...
Definition AMReX_Hypre.H:33
a one-thingy-per-box distributed object
Definition AMReX_LayoutData.H:13
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:40
A Collection of IArrayBoxes.
Definition AMReX_iMultiFab.H:34
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
int MPI_Comm
Definition AMReX_ccse-mpi.H:51
Definition AMReX_Amr.cpp:49