5#ifndef AMREX_HYPRE_SOLVER_H_
6#define AMREX_HYPRE_SOLVER_H_
14#include "_hypre_utilities.h"
60 template <
class Marker,
class Filler>
69 std::string a_options_namespace =
"hypre");
71 template <class MF, std::enable_if_t<IsFabArray<MF>::value &&
72 std::is_same_v<
typename MF::value_type,
73 HYPRE_Real>,
int> = 0>
85 HYPRE_Real rel_tol, HYPRE_Real abs_tol,
int max_iter);
88 int getNumIters ()
const {
return m_hypre_ij->getNumIters(); }
92 return m_hypre_ij->getFinalResidualNorm();
96 HYPRE_IJMatrix
getA ()
const {
return m_hypre_ij->A(); }
98 HYPRE_IJVector
getb ()
const {
return m_hypre_ij->b(); }
100 HYPRE_IJVector
getx ()
const {
return m_hypre_ij->x(); }
107 template <
class Marker>
109 std::enable_if_t<IsCallable<Marker,int,int,int,int,int>::value>
111 std::enable_if_t<IsCallableR<bool,Marker,int,int,int,int,int>::value>
116 template <
typename AI>
127 template <
class Filler,
130 HYPRE_Int&, HYPRE_Int*,
131 HYPRE_Real*>::value,
int> FOO = 0>
140 template <class MF, std::enable_if_t<IsFabArray<MF>::value &&
141 std::is_same_v<
typename MF::value_type,
142 HYPRE_Real>,
int> = 0>
151 template <class MF, std::enable_if_t<IsFabArray<MF>::value &&
152 std::is_same_v<
typename MF::value_type,
153 HYPRE_Real>,
int> = 0>
166 std::string m_options_namespace;
182 HYPRE_Int m_nrows_proc;
184 std::unique_ptr<HypreIJIface> m_hypre_ij;
187 HYPRE_IJMatrix m_A =
nullptr;
188 HYPRE_IJVector m_b =
nullptr;
189 HYPRE_IJVector m_x =
nullptr;
193template <
class Marker,
class Filler>
202 std::string a_options_namespace)
203 : m_nvars (
int(a_index_type.size())),
204 m_index_type (a_index_type),
208 m_verbose (a_verbose),
209 m_options_namespace(std::move(a_options_namespace))
217 m_grids.resize(m_nvars);
218 m_local_id.resize(m_nvars);
219 m_global_id.resize(m_nvars);
220 m_nrows_grid.resize(m_nvars);
221 m_id_offset.resize(m_nvars);
223 for (
int ivar = 0; ivar < m_nvars; ++ivar) {
225 m_local_id [ivar].define(m_grids[ivar], m_dmap, 1, 0);
226 m_global_id [ivar].define(m_grids[ivar], m_dmap, 1, m_nghost);
227 m_nrows_grid[ivar].define(m_grids[0], m_dmap);
228 m_id_offset [ivar].define(m_grids[0], m_dmap);
229 nrows_max += m_grids[ivar].numPts();
231 m_global_id_vec.define(m_grids[0], m_dmap);
232 m_nrows.
define (m_grids[0], m_dmap);
234 "Need to configure Hypre with --enable-bigint");
236 m_owner_mask.resize(m_nvars);
237 for (
int ivar = 0; ivar < m_nvars; ++ivar) {
242 m_cell_offset.define(m_grids[0], m_dmap);
259 if (nrows_allprocs.
size() > 1) {
260 MPI_Allgather(&m_nrows_proc,
sizeof(HYPRE_Int), MPI_CHAR,
261 nrows_allprocs.data(),
sizeof(HYPRE_Int), MPI_CHAR, m_comm);
265 nrows_allprocs[0] = m_nrows_proc;
268 HYPRE_Int proc_begin = 0;
269 for (
int i = 0; i < myproc; ++i) {
270 proc_begin += nrows_allprocs[i];
273 HYPRE_Int proc_end = proc_begin;
275 for (
int ivar = 0; ivar < m_nvars; ++ivar) {
276 m_id_offset[ivar][mfi] = proc_end;
277 proc_end += m_nrows_grid[ivar][mfi];
287 using AtomicInt = std::conditional_t<
sizeof(HYPRE_Int) == 4,
288 HYPRE_Int,
unsigned long long>;
289 fill_global_id<AtomicInt>();
292 HYPRE_Int ilower = proc_begin;
293 HYPRE_Int iupper = proc_end-1;
294 m_hypre_ij = std::make_unique<HypreIJIface>(m_comm, ilower, iupper, m_verbose);
295 m_hypre_ij->parse_inputs(m_options_namespace);
298 m_A = m_hypre_ij->A();
299 m_b = m_hypre_ij->b();
300 m_x = m_hypre_ij->x();
306template <
class Marker>
308 std::enable_if_t<IsCallable<Marker,int,int,int,int,int>::value>
310 std::enable_if_t<IsCallableR<bool,Marker,int,int,int,int,int>::value>
319 int boxno = mfi.LocalIndex();
321 for (
int ivar = 0; ivar < m_nvars; ++ivar) {
325 m_cell_offset[mfi].resize(npts_tot);
327 int* p_cell_offset = m_cell_offset[mfi].data();
328 for (
int ivar = 0; ivar < m_nvars; ++ivar) {
330 auto const& lid = m_local_id[ivar].array(mfi);
331 auto const& owner = m_owner_mask[ivar]->const_array(mfi);
333 const auto npts =
static_cast<int>(bx.
numPts());
334 int npts_box = amrex::Scan::PrefixSum<int>(npts,
338 int id = (owner ( cell.
x,cell.
y,cell.
z ) &&
339 marker(boxno,cell.
x,cell.
y,cell.
z,ivar)) ? 1 : 0;
340 lid(cell.
x,cell.
y,cell.
z) = id;
346 if (lid(cell.
x,cell.
y,cell.
z)) {
347 lid(cell.
x,cell.
y,cell.
z) = ps;
348 p_cell_offset[ps] =
offset;
350 lid(cell.
x,cell.
y,cell.
z) = std::numeric_limits<int>::lowest();
354 m_nrows_grid[ivar][mfi] = npts_box;
355 npts_tot += npts_box;
356 p_cell_offset += npts_box;
358 m_cell_offset[mfi].resize(npts_tot);
367 int boxno = mfi.LocalIndex();
368 for (
int ivar = 0; ivar < m_nvars; ++ivar) {
370 auto const& lid = m_local_id[ivar].array(mfi);
371 auto const& owner = m_owner_mask[ivar]->const_array(mfi);
375 for (
int k = lo.z; k <= hi.z; ++k) {
376 for (
int j = lo.y; j <= hi.y; ++j) {
377 for (
int i = lo.x; i <= hi.x; ++i) {
378 if (owner(i,j,k) && marker(boxno,i,j,k,ivar)) {
381 lid(i,j,k) = std::numeric_limits<int>::lowest();
384 m_nrows_grid[ivar][mfi] = id;
392 for (
int ivar = 0; ivar < m_nvars; ++ivar) {
393 nrows += m_nrows_grid[ivar][mfi];
395 m_nrows[mfi] = nrows;
396 m_nrows_proc += nrows;
401template <
typename AI>
410 if constexpr (std::is_same_v<HYPRE_Int,AI>) {
411 for (
int ivar = 0; ivar < m_nvars; ++ivar) {
412 p_global_id.push_back(&(m_global_id[ivar]));
415 for (
int ivar = 0; ivar < m_nvars; ++ivar) {
416 global_id_raii.emplace_back(m_global_id[ivar].
boxArray(),
417 m_global_id[ivar].DistributionMap(),
418 1, m_global_id[ivar].nGrowVect());
419 p_global_id.push_back(&(global_id_raii[ivar]));
424#pragma omp parallel if (Gpu::notInLaunchRegion())
427 auto& rows_vec = m_global_id_vec[mfi];
428 rows_vec.resize(m_nrows[mfi]);
431 for (
int ivar = 0; ivar < m_nvars; ++ivar) {
432 HYPRE_Int
const os = m_id_offset[ivar][mfi];
433 Box bx = mfi.validbox();
435 Array4<AI> const& gid = p_global_id[ivar]->array(mfi);
436 auto const& lid = m_local_id[ivar].const_array(mfi);
437 HYPRE_Int* rows = rows_vec.data() + nrows;
438 nrows += m_nrows_grid[ivar][mfi];
441 if (lid.contains(i,j,k) && lid(i,j,k) >= 0) {
442 const auto id = lid(i,j,k) + os;
443 rows[lid(i,j,k)] = id;
444 gid(i,j,k) =
static_cast<AI
>(id);
446 gid(i,j,k) =
static_cast<AI
>
447 (std::numeric_limits<HYPRE_Int>::max());
453 for (
int ivar = 0; ivar < m_nvars; ++ivar) {
455 m_geom.periodicity());
456 p_global_id[ivar]->FillBoundary(m_geom.periodicity());
458 if constexpr (!std::is_same<HYPRE_Int, AI>()) {
459 auto const& dst = m_global_id[ivar].arrays();
460 auto const& src = p_global_id[ivar]->const_arrays();
464 dst[b](i,j,k) =
static_cast<HYPRE_Int
>(src[b](i,j,k));
479 auto* p_cols_tmp = cols_tmp.
data();
480 auto* p_mat_tmp = mat_tmp.
data();
481 auto const* p_cols = cols.
data();
482 auto const* p_mat = mat.
data();
484 Scan::PrefixSum<T>(N,
487 return static_cast<T
>(p_cols[i] >= 0);
491 if (p_cols[i] >= 0) {
492 p_cols_tmp[s] = p_cols[i];
493 p_mat_tmp[s] = p_mat[i];
497 std::swap(cols_tmp, cols);
498 std::swap(mat_tmp, mat);
505template <
class Filler,
507 Array4<HYPRE_Int const>
const*,
508 HYPRE_Int&, HYPRE_Int*,
509 HYPRE_Real*>::value,
int> FOO>
521 for (
MFIter mfi(m_local_id[0],mfitinfo); mfi.
isValid(); ++mfi)
523 int boxno = mfi.LocalIndex();
524 const HYPRE_Int nrows = m_nrows[mfi];
529 HYPRE_Int* ncols = ncols_vec.
data();
533 HYPRE_Int* cols = cols_vec.
data();
537 HYPRE_Real* mat = mat_vec.
data();
540 for (
int ivar = 0; ivar < m_nvars; ++ivar) {
541 gid_v[ivar] = m_global_id[ivar].const_array(mfi);
546 (gid_v.data(), gid_v.
size());
547 auto const* pgid = gid_buf.
data();
548 auto const* p_cell_offset = m_cell_offset[mfi].data();
550 for (
int ivar = 0; ivar < m_nvars; ++ivar) {
551 const HYPRE_Int nrows_var = m_nrows_grid[ivar][mfi];
554 ntot += Reduce::Sum<Long>(nrows_var,
558 filler(boxno, cell.
x, cell.
y, cell.
z, ivar, pgid,
563 p_cell_offset += nrows_var;
565 cols +=
Long(nrows_var)*MSS;
566 mat +=
Long(nrows_var)*MSS;
571 if (ntot >=
Long(std::numeric_limits<int>::max())) {
572 detail::pack_matrix_gpu<Long>(cols_tmp, mat_tmp, cols_vec, mat_vec);
574 detail::pack_matrix_gpu<int>(cols_tmp, mat_tmp, cols_vec, mat_vec);
577 auto* pgid = gid_v.data();
578 for (
int ivar = 0; ivar < m_nvars; ++ivar) {
579 if (m_nrows_grid[ivar][mfi] > 0) {
580 auto const& lid = m_local_id[ivar].const_array(mfi);
582 [=,&ncols,&cols,&mat] (
int i,
int j,
int k)
584 if (lid(i,j,k) >= 0) {
585 filler(boxno, i, j, k, ivar, pgid, *ncols, cols, mat);
595 const auto& rows_vec = m_global_id_vec[mfi];
596 HYPRE_Int
const* rows = rows_vec.data();
599 HYPRE_IJMatrixSetValues(m_A, nrows, ncols_vec.
data(), rows,
601 Gpu::hypreSynchronize();
604 HYPRE_IJMatrixAssemble(m_A);
608template <class MF, std::enable_if_t<IsFabArray<MF>::value &&
609 std::is_same_v<
typename MF::value_type,
610 HYPRE_Real>,
int> FOO>
614 HYPRE_Real rel_tol, HYPRE_Real abs_tol,
int max_iter)
620 HYPRE_IJVectorInitialize(m_b);
621 HYPRE_IJVectorInitialize(m_x);
623 load_vectors(a_soln, a_rhs);
625 HYPRE_IJVectorAssemble(m_x);
626 HYPRE_IJVectorAssemble(m_b);
628 m_hypre_ij->solve(rel_tol, abs_tol, max_iter);
630 get_solution(a_soln);
634template <class MF, std::enable_if_t<IsFabArray<MF>::value &&
635 std::is_same_v<
typename MF::value_type,
636 HYPRE_Real>,
int> FOO>
650 const HYPRE_Int nrows = m_nrows[mfi];
657 auto* xp = xvec.
data();
658 auto* bp = bvec.
data();
660 HYPRE_Int
const* rows = m_global_id_vec[mfi].data();
663 for (
int ivar = 0; ivar < m_nvars; ++ivar) {
664 if (m_nrows_grid[ivar][mfi] > 0) {
665 auto const& xfab = a_soln[ivar]->const_array(mfi);
666 auto const& bfab = a_rhs [ivar]->const_array(mfi);
667 auto const& lid = m_local_id[ivar].const_array(mfi);
669 HYPRE_Real* b = bp +
offset;
673 if (lid(i,j,k) >= 0) {
674 x[lid(i,j,k)] = xfab(i,j,k);
675 b[lid(i,j,k)] = bfab(i,j,k);
678 offset += m_nrows_grid[ivar][mfi];
682 Gpu::streamSynchronize();
683 HYPRE_IJVectorSetValues(m_x, nrows, rows, xp);
684 HYPRE_IJVectorSetValues(m_b, nrows, rows, bp);
685 Gpu::hypreSynchronize();
691template <class MF, std::enable_if_t<IsFabArray<MF>::value &&
692 std::is_same_v<
typename MF::value_type,
693 HYPRE_Real>,
int> FOO>
705 const HYPRE_Int nrows = m_nrows[mfi];
710 auto* xp = xvec.
data();
712 HYPRE_Int
const* rows = m_global_id_vec[mfi].data();
714 HYPRE_IJVectorGetValues(m_x, nrows, rows, xp);
715 Gpu::hypreSynchronize();
718 for (
int ivar = 0; ivar < m_nvars; ++ivar) {
719 if (m_nrows_grid[ivar][mfi] > 0) {
720 auto const& xfab = a_soln[ivar]->array(mfi);
721 auto const& lid = m_local_id[ivar].const_array(mfi);
726 if (lid(i,j,k) >= 0) {
727 xfab(i,j,k) =
x[lid(i,j,k)];
730 offset += m_nrows_grid[ivar][mfi];
733 Gpu::streamSynchronize();
737 for (
int ivar = 0; ivar < m_nvars; ++ivar) {
739 m_geom.periodicity());
#define BL_PROFILE(a)
Definition AMReX_BLProfiler.H:551
#define AMREX_ALWAYS_ASSERT_WITH_MESSAGE(EX, MSG)
Definition AMReX_BLassert.H:49
#define AMREX_ASSERT(EX)
Definition AMReX_BLassert.H:38
#define AMREX_GPU_DEVICE
Definition AMReX_GpuQualifiers.H:18
Array4< int const > offset
Definition AMReX_HypreMLABecLap.cpp:1139
A collection of Boxes stored in an Array.
Definition AMReX_BoxArray.H:568
__host__ __device__ BoxND & grow(int i) noexcept
Definition AMReX_Box.H:641
__host__ __device__ Long numPts() const noexcept
Return the number of points contained in the BoxND.
Definition AMReX_Box.H:356
__host__ __device__ BoxND & convert(IndexTypeND< dim > typ) noexcept
Convert the BoxND from the current type into the argument type. This may change the BoxND coordinates...
Definition AMReX_Box.H:974
__host__ __device__ IntVectND< dim > atOffset(Long offset) const noexcept
Given the offset, compute IntVectND<dim>
Definition AMReX_Box.H:1071
Calculates the distribution of FABs to MPI processes.
Definition AMReX_DistributionMapping.H:43
Rectangular problem domain geometry.
Definition AMReX_Geometry.H:74
Periodicity periodicity() const noexcept
Definition AMReX_Geometry.H:356
Definition AMReX_GpuBuffer.H:18
T const * data() const noexcept
Definition AMReX_GpuBuffer.H:45
Solve Ax = b using HYPRE's generic IJ matrix format where A is a sparse matrix specified using the co...
Definition AMReX_HypreSolver.H:31
HypreSolver(Vector< IndexType > const &a_index_type, IntVect const &a_nghost, Geometry const &a_geom, BoxArray const &a_grids, DistributionMapping const &a_dmap, Marker &&a_marker, Filler &&a_filler, int a_verbose=0, std::string a_options_namespace="hypre")
Definition AMReX_HypreSolver.H:194
int getNumIters() const
Number of iterations from the last solve().
Definition AMReX_HypreSolver.H:88
std::enable_if_t< IsCallable< Marker, int, int, int, int, int >::value > fill_local_id(Marker const &marker)
Assign local ids to each owned DOF by invoking marker.
Definition AMReX_HypreSolver.H:312
void get_solution(Vector< MF * > const &a_soln)
Copy IJ solution entries back into AMReX storage.
Definition AMReX_HypreSolver.H:695
HYPRE_IJVector getx() const
Access the IJ solution handle (non-owning).
Definition AMReX_HypreSolver.H:100
void load_vectors(Vector< MF * > const &a_soln, Vector< MF const * > const &a_rhs)
Copy AMReX RHS/initial guess data into the IJ vectors.
Definition AMReX_HypreSolver.H:638
HYPRE_IJMatrix getA() const
Access the assembled IJ matrix handle (non-owning).
Definition AMReX_HypreSolver.H:96
void fill_matrix(Filler const &filler)
Fill each CSR row using the supplied filler functor.
Definition AMReX_HypreSolver.H:511
void fill_global_id()
Convert the local ids to globally unique ids visible to HYPRE.
Definition AMReX_HypreSolver.H:403
void solve(Vector< MF * > const &a_soln, Vector< MF const * > const &a_rhs, HYPRE_Real rel_tol, HYPRE_Real abs_tol, int max_iter)
Solve Ax=b after the constructor assembled the IJ matrix.
Definition AMReX_HypreSolver.H:612
HYPRE_Real getFinalResidualNorm() const
Final residual norm from the last solve().
Definition AMReX_HypreSolver.H:91
HYPRE_IJVector getb() const
Access the assembled IJ RHS handle (non-owning).
Definition AMReX_HypreSolver.H:98
a one-thingy-per-box distributed object
Definition AMReX_LayoutData.H:13
void define(const BoxArray &a_grids, const DistributionMapping &a_dm)
Definition AMReX_LayoutData.H:25
Iterator for looping ever tiles and boxes of amrex::FabArray based containers.
Definition AMReX_MFIter.H:85
bool isValid() const noexcept
Is the iterator valid i.e. is it associated with a FAB?
Definition AMReX_MFIter.H:169
Dynamically allocated vector for trivially copyable data.
Definition AMReX_PODVector.H:308
size_type size() const noexcept
Definition AMReX_PODVector.H:648
void resize(size_type a_new_size, GrowthStrategy strategy=GrowthStrategy::Poisson)
Definition AMReX_PODVector.H:728
void clear() noexcept
Definition AMReX_PODVector.H:646
T * data() noexcept
Definition AMReX_PODVector.H:666
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
Long size() const noexcept
Definition AMReX_Vector.H:53
amrex_long Long
Definition AMReX_INT.H:30
__host__ __device__ Dim3 ubound(Array4< T > const &a) noexcept
Return the inclusive upper bounds of an Array4 in Dim3 form.
Definition AMReX_Array4.H:1331
__host__ __device__ Dim3 lbound(Array4< T > const &a) noexcept
Return the inclusive lower bounds of an Array4 in Dim3 form.
Definition AMReX_Array4.H:1317
void streamSynchronize() noexcept
Definition AMReX_GpuDevice.H:310
MPI_Comm CommunicatorSub() noexcept
sub-communicator for current frame
Definition AMReX_ParallelContext.H:70
int MyProcSub() noexcept
my sub-rank in current frame
Definition AMReX_ParallelContext.H:76
int NProcsSub() noexcept
number of ranks in current frame
Definition AMReX_ParallelContext.H:74
static constexpr struct amrex::Scan::Type::Exclusive exclusive
static constexpr RetSum noRetSum
Definition AMReX_Scan.H:33
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
__host__ __device__ BoxND< dim > convert(const BoxND< dim > &b, const IntVectND< dim > &typ) noexcept
Return a BoxND with different type.
Definition AMReX_Box.H:1558
std::enable_if_t< std::is_integral_v< T > > ParallelFor(TypeList< CTOs... > ctos, std::array< int, sizeof...(CTOs)> const &runtime_options, T N, F &&f)
Definition AMReX_CTOParallelForImpl.H:193
IntVect nGrowVect(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2856
void OverrideSync(FabArray< FAB > &fa, FabArray< IFAB > const &msk, const Periodicity &period)
Definition AMReX_FabArrayUtility.H:1442
std::unique_ptr< iMultiFab > OwnerMask(FabArrayBase const &mf, const Periodicity &period, const IntVect &ngrow)
Definition AMReX_iMultiFab.cpp:699
const int[]
Definition AMReX_BLProfiler.cpp:1664
__host__ __device__ void Loop(Dim3 lo, Dim3 hi, F const &f) noexcept
Definition AMReX_Loop.H:127
BoxArray const & boxArray(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2861
A multidimensional array accessor.
Definition AMReX_Array4.H:283
Definition AMReX_Dim3.H:12
int x
Definition AMReX_Dim3.H:12
int z
Definition AMReX_Dim3.H:12
int y
Definition AMReX_Dim3.H:12
Test if a given type T is callable with arguments of type Args...
Definition AMReX_TypeTraits.H:213
Definition AMReX_MFIter.H:20
MFItInfo & DisableDeviceSync() noexcept
Definition AMReX_MFIter.H:44
MFItInfo & UseDefaultStream() noexcept
Definition AMReX_MFIter.H:72