Block-Structured AMR Software Framework
Loading...
Searching...
No Matches
AMReX_GMRES_MV.H
Go to the documentation of this file.
1#ifndef AMREX_GMRES_MV_H_
2#define AMREX_GMRES_MV_H_
3
4#include <AMReX_GMRES.H>
5#include <AMReX_Algebra.H>
6#include <AMReX_Smoother_MV.H>
7#include <functional>
8
9namespace amrex {
10
17template <typename T>
19{
20public:
21 using RT = Real; // double or float
25 using PC = std::function<void(VEC&,VEC const&)>;
26
32 GMRES_MV (MAT const* a_mat);
33
41 void setPrecond (PC a_pc) { m_pc = std::move(a_pc); }
42
51 void solve (VEC& a_sol, VEC const& a_rhs, T a_tol_rel, T a_tol_abs);
52
54 void setVerbose (int v) { m_gmres.setVerbose(v); }
55
57 GM& getGMRES () { return m_gmres; }
58
60 [[nodiscard]] VEC makeVecRHS () const;
61
63 [[nodiscard]] VEC makeVecLHS () const;
64
66 static T norm2 (VEC const& vec);
67
69 static void scale (VEC& vec, T scale_factor);
70
72 static T dotProduct (VEC const& vec1, VEC const& vec2);
73
79 static void setToZero (VEC& lhs);
80
87 static void assign (VEC& lhs, VEC const& rhs);
88
96 static void increment (VEC& lhs, VEC const& rhs, T a);
97
107 static void linComb (VEC& lhs, T a, VEC const& rhs_a, T b, VEC const& rhs_b);
108
115 void apply (VEC& lhs, VEC& rhs) const;
116
118 void precond (VEC& lhs, VEC const& rhs) const;
119
120private:
121 GM m_gmres;
122 MAT const* m_mat = nullptr;
123 PC m_pc;
124};
125
126template <typename T>
128 : m_mat(a_mat)
129{
130 m_gmres.define(*this);
131}
132
133template <typename T>
134void GMRES_MV<T>::solve (VEC& a_sol, VEC const& a_rhs, T a_tol_rel, T a_tol_abs)
135{
136 m_gmres.solve(a_sol, a_rhs, a_tol_rel, a_tol_abs);
137}
138
139template <typename T>
141{
142 return VEC(m_mat->partition());
143}
144
145template <typename T>
147{
148 return VEC(m_mat->partition());
149}
150
151template <typename T>
153{
154 return vec.norm2();
155}
156
157template <typename T>
158void GMRES_MV<T>::scale (VEC& vec, T scale_factor)
159{
160 vec.scaleAsync(scale_factor);
161}
162
163template <typename T>
164T GMRES_MV<T>::dotProduct (VEC const& vec1, VEC const& vec2)
165{
166 return amrex::Dot(vec1,vec2);
167}
168
169template <typename T>
171{
172 lhs.setValAsync(0);
173}
174
175template <typename T>
176void GMRES_MV<T>::assign (VEC& lhs, VEC const& rhs)
177{
178 lhs.copyAsync(rhs);
179}
180
181template <typename T>
182void GMRES_MV<T>::increment (VEC& lhs, VEC const& rhs, T a)
183{
184 amrex::Axpy(lhs, a, rhs);
185}
186
187template <typename T>
188void GMRES_MV<T>::linComb (VEC& lhs, T a, VEC const& rhs_a, T b, VEC const& rhs_b)
189{
190 amrex::LinComb(lhs, a, rhs_a, b, rhs_b);
191}
192
193template <typename T>
194void GMRES_MV<T>::apply (VEC& lhs, VEC& rhs) const
195{
196 amrex::SpMV(lhs, *m_mat, rhs);
197}
198
199template <typename T>
200void GMRES_MV<T>::precond (VEC& lhs, VEC const& rhs) const
201{
202 if (m_pc) {
203 m_pc(lhs, rhs);
204 } else {
205 lhs.copyAsync(rhs);
206 }
207}
208
209}
210#endif
Distributed dense vector that mirrors the layout of an AlgPartition.
Definition AMReX_AlgVector.H:29
void setValAsync(T val)
Definition AMReX_AlgVector.H:282
T norm2(bool local=false) const
Return the 2-norm.
Definition AMReX_AlgVector.H:384
void copyAsync(AlgVector< T, Allocator > const &rhs)
Definition AMReX_AlgVector.H:297
void scaleAsync(T scale_factor)
Definition AMReX_AlgVector.H:331
Definition AMReX_GMRES_MV.H:19
GM & getGMRES()
Access the underlying GMRES object for additional tuning.
Definition AMReX_GMRES_MV.H:57
static void assign(VEC &lhs, VEC const &rhs)
Copy rhs into lhs.
Definition AMReX_GMRES_MV.H:176
std::function< void(VEC &, VEC const &)> PC
Definition AMReX_GMRES_MV.H:25
static T norm2(VEC const &vec)
Euclidean norm of vec.
Definition AMReX_GMRES_MV.H:152
void setVerbose(int v)
Set verbosity level v for the underlying GMRES solver.
Definition AMReX_GMRES_MV.H:54
VEC makeVecLHS() const
Return another AlgVector with the same partition for LHS storage.
Definition AMReX_GMRES_MV.H:146
AlgVector< T > VEC
Definition AMReX_GMRES_MV.H:22
void apply(VEC &lhs, VEC &rhs) const
Apply the sparse operator to rhs and store the result in lhs.
Definition AMReX_GMRES_MV.H:194
amrex::GMRES< VEC, GMRES_MV< T > > GM
Definition AMReX_GMRES_MV.H:24
static void scale(VEC &vec, T scale_factor)
Scale vec in place by scale_factor.
Definition AMReX_GMRES_MV.H:158
void solve(VEC &a_sol, VEC const &a_rhs, T a_tol_rel, T a_tol_abs)
Solve the linear system.
Definition AMReX_GMRES_MV.H:134
static void linComb(VEC &lhs, T a, VEC const &rhs_a, T b, VEC const &rhs_b)
Form the linear combination .
Definition AMReX_GMRES_MV.H:188
Real RT
Definition AMReX_GMRES_MV.H:21
SpMatrix< T > MAT
Definition AMReX_GMRES_MV.H:23
static void setToZero(VEC &lhs)
Reset lhs to zero.
Definition AMReX_GMRES_MV.H:170
static void increment(VEC &lhs, VEC const &rhs, T a)
Accumulate .
Definition AMReX_GMRES_MV.H:182
void precond(VEC &lhs, VEC const &rhs) const
Apply the optional preconditioner (or copy if none is provided).
Definition AMReX_GMRES_MV.H:200
GMRES_MV(MAT const *a_mat)
Bind GMRES to a sparse matrix described by a_mat.
Definition AMReX_GMRES_MV.H:127
static T dotProduct(VEC const &vec1, VEC const &vec2)
Dot product between vec1 and vec2.
Definition AMReX_GMRES_MV.H:164
VEC makeVecRHS() const
Return an AlgVector that mirrors the matrix partition for RHS storage.
Definition AMReX_GMRES_MV.H:140
void setPrecond(PC a_pc)
Supply an optional right-preconditioner functor.
Definition AMReX_GMRES_MV.H:41
void define(M &linop)
Bind the solver to a linear operator.
Definition AMReX_GMRES.H:208
void setVerbose(int v)
Set verbosity level v (0 = silent).
Definition AMReX_GMRES.H:117
Distributed CSR matrix that manages storage and GPU-friendly partitions.
Definition AMReX_SpMatrix.H:61
amrex_real Real
Floating Point Type for Fields.
Definition AMReX_REAL.H:79
Definition AMReX_Amr.cpp:49
void LinComb(MF &dst, typename MF::value_type a, MF const &src_a, int acomp, typename MF::value_type b, MF const &src_b, int bcomp, int dcomp, int ncomp, IntVect const &nghost)
dst = a*src_a + b*src_b
Definition AMReX_FabArrayUtility.H:2009
void SpMV(Long nrows, Long ncols, T *__restrict__ py, CsrView< T const > const &A, T const *__restrict__ px)
Perform y = A * x using CSR data (GPU/CPU aware).
Definition AMReX_SpMV.H:28
FAB::value_type Dot(FabArray< FAB > const &x, int xcomp, FabArray< FAB > const &y, int ycomp, int ncomp, IntVect const &nghost, bool local=false)
Compute dot products of two FabArrays.
Definition AMReX_FabArrayUtility.H:1673
void Axpy(AlgVector< T, Allocator > &y, T a, AlgVector< T, Allocator > const &x)
y = ax + y. For GPU builds this is asynchronous with respect to the host.
Definition AMReX_AlgVecUtil.H:183