CHROMA
invcg1_array.h
Go to the documentation of this file.
1 // -*- C++ -*-
2 /*! \file
3  * \brief Conjugate-Gradient algorithm for a generic Linear Operator
4  */
5 
6 #ifndef __invcg1_array__
7 #define __invcg1_array__
8 
9 #include "linearop.h"
10 
11 namespace Chroma {
12 
13 //! Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator
14 /*! \ingroup invert
15  * This subroutine uses the Conjugate Gradient (CG) algorithm to find
16  * the solution of the set of linear equations
17  *
18  * Chi = A . Psi
19  *
20  * where A is hermitian
21  *
22  * Algorithm:
23 
24  * Psi[0] := initial guess; Linear interpolation (argument)
25  * r[0] := Chi - A . Psi[0] ; Initial residual
26  * p[1] := r[0] ; Initial direction
27  * IF |r[0]| <= RsdCG |Chi| THEN RETURN; Converged?
28  * FOR k FROM 1 TO MaxCG DO CG iterations
29  * a[k] := |r[k-1]|**2 / <A p[k],p[k]> ;
30  * Psi[k] += a[k] p[k] ; New solution std::vector
31  * r[k] -= a[k] A. p[k] ; New residual
32  * IF |r[k]| <= RsdCG |Chi| THEN RETURN; Converged?
33  * b[k+1] := |r[k]|**2 / |r[k-1]|**2 ;
34  * p[k+1] := r[k] + b[k+1] p[k]; New direction
35  *
36  * Arguments:
37  *
38  * \param M Linear Operator (Read)
39  * \param chi Source (Read)
40  * \param psi Solution (Modify)
41  * \param RsdCG CG residual accuracy (Read)
42  * \param MaxCG Maximum CG iterations (Read)
43  * \param n_count Number of CG iteration (Write)
44  *
45  * Local Variables:
46  *
47  * p Direction std::vector
48  * r Residual std::vector
49  * cp | r[k] |**2
50  * c | r[k-1] |**2
51  * k CG iteration counter
52  * a a[k]
53  * b b[k+1]
54  * d < p[k], A.p[k] >
55  * Mp Temporary for M.p
56  *
57  * Subroutines:
58  * +
59  *
60  *
61  * Operations:
62  *
63  * 2 A + 2 Nc Ns + N_Count ( 2 A + 10 Nc Ns )
64  */
65 template<typename T>
66 void InvCG1(const LinearOperatorArray<T>& A,
67  const multi1d<T>& chi,
68  multi1d<T>& psi,
69  const Real& RsdCG,
70  int MaxCG,
71  int& n_count);
72 
73 } // end namespace Chroma
74 
75 #endif
SystemSolverResults_t InvCG1(const LinearOperator< LatticeFermion > &A, const LatticeFermion &chi, LatticeFermion &psi, const Real &RsdCG, int MaxCG, int MinCG)
Conjugate-Gradient (CGNE) algorithm for a generic Linear Operator.
Definition: invcg1.cc:215
Linear Operators.
Asqtad Staggered-Dirac operator.
Definition: klein_gord.cc:10
const WilsonTypeFermAct< multi1d< LatticeFermion > > Handle< const ConnectState > const multi1d< Real > enum InvType invType const multi1d< Real > & RsdCG
Definition: pbg5p_w.cc:30
int n_count
Definition: invbicg.cc:78
const WilsonTypeFermAct< multi1d< LatticeFermion > > Handle< const ConnectState > const multi1d< Real > enum InvType invType const multi1d< Real > int MaxCG
Definition: pbg5p_w.cc:32
multi1d< LatticeFermion > chi(Ncb)
LatticeFermion psi
Definition: mespbg5p_w.cc:35
A(A, psi, r, Ncb, PLUS)