TBCI Numerical high perf. C++ Library 2.8.0
cgs.h
Go to the documentation of this file.
1
4/* $Id: cgs.h,v 1.11.2.15 2019/05/28 11:13:02 garloff Exp $ */
5
6#ifndef TBCI_SOLVER_CGS_H
7#define TBCI_SOLVER_CGS_H
8
9#include "tbci/basics.h"
10
11#ifdef SOLVER_LOG
12# include <fstream>
13# include <iomanip>
14#endif
15
17
18/* NOTE: Template parameter T and T_tol should really be double (or
19 * similar: float, long double are OK), as fabs() is used, which
20 * returns double. This could be changed to abs(), but will then cause
21 * some trouble for complex numbers.
22 */
23
24
57
58INST3(template <typename T, Matrix<T>, Vector<T> > class NN friend \
59int CGS (const Matrix<T>&, Vector<T>&, const Vector<T>&,\
60 const Preconditioner_Sig<T, Matrix<T> >&, unsigned int&, double&, const unsigned);)
61INST3(template <typename T, BdMatrix<T>, Vector<T> > class NN friend \
62int CGS (const BdMatrix<T>&, Vector<T>&, const Vector<T>&,\
63 const Preconditioner_Sig<T, BdMatrix<T> >&, unsigned int&, double&, const unsigned);)
64INST3(template <typename T, Symm_BdMatrix<T>, Vector<T> > class NN friend \
65int CGS (const Symm_BdMatrix<T>&, Vector<T>&, const Vector<T>&,\
66 const Preconditioner_Sig<T, Symm_BdMatrix<T> >&, unsigned int&, double&, const unsigned);)
67
68
69template < typename T, typename SysMatrix, typename SysVector >
70int CGS(const SysMatrix &A, SysVector &x, const SysVector &b,
71 const Preconditioner_Sig<T, SysMatrix> &M, unsigned int &max_iter,
72 double &tol, const unsigned off = 0)
73{
74 const unsigned int dim = A.rows();
75 SysVector p(dim), phat(dim), q(dim), qhat(dim), vhat(dim), u(dim), uhat(dim);
76 typename SysVector::value_type rho_1, rho_2((typename SysVector::value_type)1), alpha, beta;
77
78 double residsqr;
79 double tolsqr = TBCI__ sqr (tol);
80 double normbsqr = b.fabssqr();
81 if (normbsqr < 1e-32)
82 normbsqr = 1e32;
83 else
84 normbsqr = 1.0 / normbsqr;
85
86 SysVector r(b - A*x);
87 SysVector rtilde(r);
88 //M.update (A);
89
90 if ((residsqr = r.fabssqr() * normbsqr) <= tolsqr)
91 {
92 tol = MATH__ sqrt(residsqr);
93 max_iter = 0;
94 return 0;
95 }
96#ifdef SOLVER_LOG
97 STD__ ofstream cnvg;
98 if (off) cnvg.open ("cgs_cnvg.gnu", STD__ ios::app);
99 else cnvg.open ("cgs_cnvg.gnu");
100#endif
101
102 for (unsigned int i = 1; i <= max_iter; i++)
103 {
104#ifdef SOLVER_LOG
105 cnvg << i+off << "\t" << MATH__ sqrt(residsqr) << "\n";
106#endif
107 rho_1 = dot(rtilde, r);
108 if (rho_1 == (typename SysVector::value_type)0)
109 {
110 tol = MATH__ sqrt (residsqr);
111 return 2;
112 }
113 if (i == 1)
114 {
115 u = r;
116 p = u;
117 } else
118 {
119 beta = rho_1 / rho_2;
120 u = r + beta * q;
121 p = u + beta * (q + beta * p);
122 }
123 phat = M.solve(p);
124 vhat = A*phat;
125 alpha = rho_1 / dot(rtilde, vhat);
126 q = u - alpha * vhat;
127 uhat = M.solve(u + q);
128 x += alpha * uhat;
129 qhat = A * uhat;
130 r -= alpha * qhat;
131 rho_2 = rho_1;
132 if ((residsqr = r.fabssqr() * normbsqr) < tolsqr)
133 {
134 tol = MATH__ sqrt(residsqr);
135 max_iter = i;
136#ifdef SOLVER_LOG
137 cnvg << i+off << "\t" << MATH__ sqrt(residsqr) << STD__ endl; cnvg.close ();
138#endif
139 return 0;
140 }
141 }
142
143 tol = MATH__ sqrt(residsqr);
144 return 1;
145}
146
148
149#endif /* TBCI_SOLVER_CGS_H */
const Vector< T > const Vector< T > const Vector< T > & p
Definition LM_fit.h:98
const Vector< T > const Vector< T > & x
Definition LM_fit.h:97
int i
Definition LM_fit.h:71
#define STD__
Definition basics.h:338
#define INST3(x, y, z)
Definition basics.h:240
#define NAMESPACE_END
Definition basics.h:323
#define NAMESPACE_TBCI
Definition basics.h:317
#define TBCI__
Definition basics.h:332
#define MATH__
Definition basics.h:339
#define T
Definition bdmatlib.cc:20
SysMatrix
Definition bicgstab.h:90
NAMESPACE_TBCI Vector< T > class NN friend int CGS(const BdMatrix< T > &, Vector< T > &, const Vector< T > &, const Preconditioner_Sig< T, BdMatrix< T > > &, unsigned int &, double &, const unsigned)
The class BdMatrix is an implementation to store and do operations on sparse Matrices with a band str...
Abstract base class for all Preconditioners.
Definition precond.h:41
virtual TVector< T > solve(const Vector< T > &) const =0
Matrix class with optimized Matrix-Vector multiplication for symmetrical Matrices.
cplx< T > dot(const cplx< T > &a, const cplx< T > &b)
Definition cplx.h:300
cplx< T > sqr(const cplx< T > &c)
Definition cplx.h:449
NAMESPACE_END NAMESPACE_CSTD TBCI__ cplx< T > sqrt(const TBCI__ cplx< T > &z)
Definition cplx.h:751
F_TMatrix< T > b
Definition f_matrix.h:736