52template<
class ContainerType>
69 LGMRES(
const ContainerType&
copyable,
unsigned max_inner,
unsigned max_outer,
unsigned max_restarts):
73 m_maxRestarts( max_restarts),
74 m_inner_m( max_inner),
75 m_outer_k( max_outer),
76 m_krylovDimension( max_inner+max_outer)
78 if( m_inner_m < m_outer_k)
79 std::cerr <<
"WARNING (LGMRES): max_inner is smaller than the restart dimension max_outer. Did you swap the constructor parameters?\n";
81 m_H.assign( m_krylovDimension+1, std::vector<value_type>( m_krylovDimension, 0));
84 m_givens.assign( m_krylovDimension+1, {0,0});
86 m_s.assign(m_krylovDimension+1,0);
89 m_outer_w.assign(m_outer_k,
copyable);
90 m_outer_Az.assign(m_outer_k,
copyable);
91 m_V.assign(m_krylovDimension+1,
copyable);
99 template<
class ...Params>
103 *
this =
LGMRES( std::forward<Params>( ps)...);
107 void set_max(
unsigned new_Restarts) {m_maxRestarts = new_Restarts;}
110 unsigned get_max()
const {
return m_maxRestarts;}
113 m_throw_on_fail = throw_on_fail;
117 const ContainerType&
copyable()
const{
return m_tmp;}
141 template<
class MatrixType0,
class ContainerType0,
class ContainerType1,
class MatrixType1,
class ContainerType2 >
142 unsigned solve( MatrixType0&& A, ContainerType0& x,
const ContainerType1& b, MatrixType1&& P,
const ContainerType2& W,
value_type eps = 1e-12,
value_type nrmb_correction = 1);
154 template <
class Preconditioner,
class ContainerType0>
155 void Update(Preconditioner&& P, ContainerType &dx, ContainerType0 &x,
156 unsigned dimension,
const std::vector<std::vector<value_type>> &H,
157 std::vector<value_type> &s,
const std::vector<const ContainerType*> &W);
158 std::vector<std::array<value_type,2>> m_givens;
159 std::vector<std::vector<value_type>> m_H, m_HH;
160 ContainerType m_tmp, m_dx, m_residual;
161 std::vector<ContainerType> m_V, m_outer_w, m_outer_Az;
162 std::vector<value_type> m_s;
163 unsigned m_maxRestarts, m_inner_m, m_outer_k, m_krylovDimension;
164 bool m_converged =
true, m_throw_on_fail =
true;
168template<
class ContainerType>
169template <
class Preconditioner,
class ContainerType0>
170void LGMRES<ContainerType>::Update(Preconditioner&& P, ContainerType &dx,
172 unsigned dimension,
const std::vector<std::vector<value_type>> &H,
173 std::vector<value_type> &s,
const std::vector<const ContainerType*> &W)
177 for (
int lupe = dimension; lupe >= 0; --lupe)
179 s[lupe] = s[lupe]/H[lupe][lupe];
180 for (
int innerLupe = lupe - 1; innerLupe >= 0; --innerLupe)
183 s[innerLupe] = DG_FMA( -s[lupe],H[innerLupe][lupe], s[innerLupe]);
194template<
class ContainerType>
195template<
class Matrix,
class ContainerType0,
class ContainerType1,
class Preconditioner,
class ContainerType2>
196unsigned LGMRES< ContainerType>::solve( Matrix&& A, ContainerType0& x,
const ContainerType1& b, Preconditioner&& P,
const ContainerType2& S, value_type eps, value_type nrmb_correction)
204 value_type tol = eps*(nrmb + nrmb_correction);
212 unsigned restartCycle = 0;
213 unsigned counter = 0;
216 std::vector<ContainerType const*> m_W, m_Vptr;
217 m_W.assign(m_krylovDimension,
nullptr);
218 m_Vptr.assign(m_krylovDimension+1,
nullptr);
219 for(
unsigned i=0; i<m_krylovDimension+1; i++)
233 for(
unsigned lupe=1;lupe<=m_krylovDimension;++lupe)
237 for(
unsigned iteration=0;iteration<m_krylovDimension;++iteration)
239 unsigned outer_w_count = std::min(restartCycle,m_outer_k);
240 if(iteration < m_krylovDimension-outer_w_count){
241 m_W[iteration] = &m_V[iteration];
242 dg::blas2::gemv(std::forward<Preconditioner>(P),*m_W[iteration],m_tmp);
245 }
else if( iteration < m_krylovDimension){
246 unsigned w_idx = iteration - (m_krylovDimension - outer_w_count);
247 m_W[iteration] = &m_outer_w[w_idx];
253 for(
unsigned row=0;row<=iteration;++row)
255 m_HH[row][iteration] = m_H[row][iteration]
260 m_HH[iteration+1][iteration] = m_H[iteration+1][iteration]
271 for (
unsigned row = 0; row < iteration; row++)
273 tmp = m_givens[row][0]*m_H[row][iteration] +
274 m_givens[row][1]*m_H[row+1][iteration];
275 m_H[row+1][iteration] = -m_givens[row][1]*m_H[row][iteration]
276 + m_givens[row][0]*m_H[row+1][iteration];
277 m_H[row][iteration] = tmp;
281 if(m_H[iteration+1][iteration] == 0.0)
284 m_givens[iteration][0] = 1.0;
285 m_givens[iteration][1] = 0.0;
287 else if (fabs(m_H[iteration+1][iteration]) > fabs(m_H[iteration][iteration]))
292 tmp = m_H[iteration][iteration]/m_H[iteration+1][iteration];
293 m_givens[iteration][1] = 1.0/sqrt(1.0+tmp*tmp);
294 m_givens[iteration][0] = tmp*m_givens[iteration][1];
301 tmp = m_H[iteration+1][iteration]/m_H[iteration][iteration];
302 m_givens[iteration][0] = 1.0/sqrt(1.0+tmp*tmp);
303 m_givens[iteration][1] = tmp*m_givens[iteration][0];
306 tmp = m_givens[iteration][0]*m_H[iteration][iteration] +
307 m_givens[iteration][1]*m_H[iteration+1][iteration];
308 m_H[iteration+1][iteration] = -m_givens[iteration][1]*m_H[iteration][iteration] +
309 m_givens[iteration][0]*m_H[iteration+1][iteration];
310 m_H[iteration][iteration] = tmp;
312 tmp = m_givens[iteration][0]*m_s[iteration] + m_givens[iteration][1]*m_s[iteration+1];
313 m_s[iteration+1] = -m_givens[iteration][1]*m_s[iteration] + m_givens[iteration][1]*m_s[iteration+1];
314 m_s[iteration] = tmp;
316 rho = fabs(m_s[iteration+1]);
319 Update(std::forward<Preconditioner>(P),m_dx,x,iteration,m_H,m_s,m_W);
323 Update(std::forward<Preconditioner>(P),m_dx,x,m_krylovDimension-1,m_H,m_s,m_W);
326 std::rotate(m_outer_w.rbegin(),m_outer_w.rbegin()+1,m_outer_w.rend());
327 std::rotate(m_outer_Az.rbegin(),m_outer_Az.rbegin()+1,m_outer_Az.rend());
333 std::vector<value_type> coeffs( m_krylovDimension+1, 0.);
334 for(
unsigned i=0; i<m_krylovDimension+1; i++)
337 for(
unsigned k=0; k<m_krylovDimension; k++)
338 coeffs[i] = DG_FMA( m_HH[i][k],m_s[k], coeffs[i]);
345 }
while( (restartCycle < m_maxRestarts) && (rho > tol));
351 <<
"After "<<counter<<
" LGMRES iterations");
Functor class for the right preconditioned LGMRES method to solve .
Definition lgmres.h:54
bool converged() const
If last call to solve converged or not.
Definition lgmres.h:149
void set_throw_on_fail(bool throw_on_fail)
Set or unset a throw on failure-to-converge.
Definition lgmres.h:112
LGMRES(const ContainerType ©able, unsigned max_inner, unsigned max_outer, unsigned max_restarts)
Allocate memory for the preconditioned LGMRES method.
Definition lgmres.h:69
const ContainerType & copyable() const
Return an object of same size as the object used for construction.
Definition lgmres.h:117
LGMRES()=default
Allocate nothing, Call construct method before usage.
void construct(Params &&...ps)
Perfect forward parameters to one of the constructors.
Definition lgmres.h:100
unsigned get_max() const
Get the current maximum number of restarts.
Definition lgmres.h:110
void set_max(unsigned new_Restarts)
Set the number of restarts.
Definition lgmres.h:107
ContainerType container_type
Definition lgmres.h:56
unsigned solve(MatrixType0 &&A, ContainerType0 &x, const ContainerType1 &b, MatrixType1 &&P, const ContainerType2 &W, value_type eps=1e-12, value_type nrmb_correction=1)
Solve using a right preconditioned LGMRES method.
dg::get_value_type< ContainerType > value_type
Definition lgmres.h:57
#define _ping_
Definition exceptions.h:12
void copy(const ContainerTypeIn &source, ContainerTypeOut &target)
Definition blas1.h:243
void axpby(value_type alpha, const ContainerType1 &x, value_type1 beta, ContainerType &y)
Definition blas1.h:306
void scal(ContainerType &x, value_type alpha)
Definition blas1.h:263
void gemv(get_value_type< ContainerType1 > alpha, MatrixType &&M, const ContainerType1 &x, get_value_type< ContainerType1 > beta, ContainerType2 &y)
Alias for blas2::symv ;.
Definition blas2.h:339
auto dot(const ContainerType1 &x, const MatrixType &m, const ContainerType2 &y)
; Binary reproducible general dot product
Definition blas2.h:94
auto asDenseMatrix(const std::vector< const ContainerType * > &in)
Lightweight DenseMatrix for dg::blas2::gemv.
Definition densematrix.h:75
typename TensorTraits< std::decay_t< Vector > >::value_type get_value_type
Definition tensor_traits.h:45
This is the namespace for all functions and classes defined and used by the discontinuous Galerkin li...
Indicate failure to converge.
Definition exceptions.h:126