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