1
0
mirror of https://github.com/paboyle/Grid.git synced 2025-06-12 20:27:06 +01:00

Getting closer to having a wilson solver... introducing a first and untested

cut at Conjugate gradient. Also copied in Remez, Zolotarev, Chebyshev from
Mike Clark, Tony Kennedy and my BFM package respectively since we know we will
need these. I wanted the structure of

algorithms/approx
algorithms/iterative

etc.. to start taking shape.
This commit is contained in:
Peter Boyle
2015-05-18 07:47:05 +01:00
parent 8e99e4671f
commit d0e4673a3f
22 changed files with 1798 additions and 157 deletions

View File

@ -1,19 +1,25 @@
#ifndef GRID_CONJUGATE_GRADIENT_H
#define GRID_CONJUGATE_GRADIENT_H
#include<Grid.h>
#include<algorithms/LinearOperator.h>
namespace Grid {
template<class Field> class ConjugateGradient : public IterativeProcess<Field> {
public:
/////////////////////////////////////////////////////////////
// Base classes for iterative processes based on operators
// single input vec, single output vec.
/////////////////////////////////////////////////////////////
ConjugateGradient(RealD tol,Integer maxit): IterativeProces(tol,maxit) {};
template<class Field>
class ConjugateGradient : public OperatorFunction<Field> {
public:
RealD Tolerance;
Integer MaxIterations;
void operator() (HermitianOperatorBase<Field> *Linop,const Field &src, Field &psi){
ConjugateGradient(RealD tol,Integer maxit) : Tolerance(tol), MaxIterations(maxit) {
std::cout << Tolerance<<std::endl;
};
void operator() (HermitianOperatorBase<Field> &Linop,const Field &src, Field &psi){
RealD residual = Tolerance;
RealD cp,c,a,d,b,ssq,qq,b_pred;
Field p(src);
@ -32,24 +38,24 @@ namespace Grid {
cp =a;
ssq=norm2(src);
std::cout <<setprecision(4)<< "ConjugateGradient: guess "<<guess<<std::endl;
std::cout <<setprecision(4)<< "ConjugateGradient: src "<<ssq <<std::endl;
std::cout <<setprecision(4)<< "ConjugateGradient: mp "<<d <<std::endl;
std::cout <<setprecision(4)<< "ConjugateGradient: mmp "<<b <<std::endl;
std::cout <<setprecision(4)<< "ConjugateGradient: r "<<cp <<std::endl;
std::cout <<setprecision(4)<< "ConjugateGradient: p "<<a <<std::endl;
std::cout <<std::setprecision(4)<< "ConjugateGradient: guess "<<guess<<std::endl;
std::cout <<std::setprecision(4)<< "ConjugateGradient: src "<<ssq <<std::endl;
std::cout <<std::setprecision(4)<< "ConjugateGradient: mp "<<d <<std::endl;
std::cout <<std::setprecision(4)<< "ConjugateGradient: mmp "<<b <<std::endl;
std::cout <<std::setprecision(4)<< "ConjugateGradient: r "<<cp <<std::endl;
std::cout <<std::setprecision(4)<< "ConjugateGradient: p "<<a <<std::endl;
RealD rsq = residual* residual*ssq;
RealD rsq = Tolerance* Tolerance*ssq;
//Check if guess is really REALLY good :)
if ( cp <= rsq ) {
return 0;
return;
}
std::cout << setprecision(4)<< "ConjugateGradient: k=0 residual "<<cp<<" rsq"<<rsq<<std::endl;
std::cout << std::setprecision(4)<< "ConjugateGradient: k=0 residual "<<cp<<" rsq"<<rsq<<std::endl;
int k;
for (k=1;k<=max_iter;k++){
for (k=1;k<=MaxIterations;k++){
c=cp;
@ -62,10 +68,10 @@ namespace Grid {
b = cp/c;
// Fuse these loops ; should be really easy
// Update psi; New (conjugate/M-orthogonal) search direction
psi= a*p+psi;
p = p*b+r;
std::cout << "Iteration " <<k<<" residual "<<cp<< " target"<< rsq<<std::endl;
// Stopping condition
if ( cp <= rsq ) {
@ -78,10 +84,11 @@ namespace Grid {
RealD resnorm = sqrt(norm2(p));
RealD true_residual = resnorm/srcnorm;
std::cout<<"ConjugateGradient: true residual is "<<true_residual<<" sol "<<psinorm<<" src "<<srcnorm<<std::endl;
std::cout<<"ConjugateGradient: target residual was "<<residual<<std::endl;
return k;
std::cout<<"ConjugateGradient: target residual was "<<Tolerance<<std::endl;
}
}
std::cout<<"ConjugateGradient did NOT converge"<<std::endl;
assert(0);
}
};
}

View File

@ -0,0 +1,34 @@
#ifndef GRID_NORMAL_EQUATIONS_H
#define GRID_NORMAL_EQUATIONS_H
namespace Grid {
///////////////////////////////////////////////////////////////////////////////////////////////////////
// Take a matrix and form a Red Black solver calling a Herm solver
// Use of RB info prevents making SchurRedBlackSolve conform to standard interface
///////////////////////////////////////////////////////////////////////////////////////////////////////
template<class Field> class NormalEquations : public OperatorFunction<Field>{
private:
SparseMatrixBase<Field> & _Matrix;
OperatorFunction<Field> & _HermitianSolver;
public:
/////////////////////////////////////////////////////
// Wrap the usual normal equations Schur trick
/////////////////////////////////////////////////////
NormalEquations(SparseMatrixBase<Field> &Matrix, OperatorFunction<Field> &HermitianSolver)
: _Matrix(Matrix), _HermitianSolver(HermitianSolver) {};
void operator() (const Field &in, Field &out){
Field src(in._grid);
_Matrix.Mdag(in,src);
_HermitianSolver(src,out); // Mdag M out = Mdag in
}
};
}
#endif

View File

@ -0,0 +1,106 @@
#ifndef GRID_SCHUR_RED_BLACK_H
#define GRID_SCHUR_RED_BLACK_H
/*
* Red black Schur decomposition
*
* M = (Mee Meo) = (1 0 ) (Mee 0 ) (1 Mee^{-1} Meo)
* (Moe Moo) (Moe Mee^-1 1 ) (0 Moo-Moe Mee^-1 Meo) (0 1 )
* = L D U
*
* L^-1 = (1 0 )
* (-MoeMee^{-1} 1 )
* L^{dag} = ( 1 Mee^{-dag} Moe^{dag} )
* ( 0 1 )
* L^{-d} = ( 1 -Mee^{-dag} Moe^{dag} )
* ( 0 1 )
*
* U^-1 = (1 -Mee^{-1} Meo)
* (0 1 )
* U^{dag} = ( 1 0)
* (Meo^dag Mee^{-dag} 1)
* U^{-dag} = ( 1 0)
* (-Meo^dag Mee^{-dag} 1)
***********************
* M psi = eta
***********************
*Odd
* i) (D_oo)^{\dag} D_oo psi_o = (D_oo)^\dag L^{-1} eta_o
* eta_o' = D_oo (eta_o - Moe Mee^{-1} eta_e)
*Even
* ii) Mee psi_e + Meo psi_o = src_e
*
* => sol_e = M_ee^-1 * ( src_e - Meo sol_o )...
*
*/
namespace Grid {
///////////////////////////////////////////////////////////////////////////////////////////////////////
// Take a matrix and form a Red Black solver calling a Herm solver
// Use of RB info prevents making SchurRedBlackSolve conform to standard interface
///////////////////////////////////////////////////////////////////////////////////////////////////////
template<class Field> class SchurRedBlackSolve : public OperatorFunction<Field>{
private:
SparseMatrixBase<Field> & _Matrix;
OperatorFunction<Field> & _HermitianRBSolver;
int CBfactorise;
public:
/////////////////////////////////////////////////////
// Wrap the usual normal equations Schur trick
/////////////////////////////////////////////////////
SchurRedBlackSolve(SparseMatrixBase<Field> &Matrix, OperatorFunction<Field> &HermitianRBSolver)
: _Matrix(Matrix), _HermitianRBSolver(HermitianRBSolver) {
CBfactorise=0;
};
void operator() (const Field &in, Field &out){
// FIXME CGdiagonalMee not implemented virtual function
// FIXME need to make eo grid from full grid.
// FIXME use CBfactorise to control schur decomp
const int Even=0;
const int Odd =1;
// Make a cartesianRedBlack from full Grid
GridRedBlackCartesian grid(in._grid);
Field src_e(&grid);
Field src_o(&grid);
Field sol_e(&grid);
Field sol_o(&grid);
Field tmp(&grid);
Field Mtmp(&grid);
pickCheckerboard(Even,src_e,in);
pickCheckerboard(Odd ,src_o,in);
/////////////////////////////////////////////////////
// src_o = Mdag * (source_o - Moe MeeInv source_e)
/////////////////////////////////////////////////////
_Matrix.MooeeInv(src_e,tmp); // MooeeInv(source[Even],tmp,DaggerNo,Even);
_Matrix.Meooe (tmp,Mtmp); // Meo (tmp,src,Odd,DaggerNo);
tmp=src_o-Mtmp; // axpy (tmp,src,source[Odd],-1.0);
_Matrix.MpcDag(tmp,src_o); // Mprec(tmp,src,Mtmp,DaggerYes);
//////////////////////////////////////////////////////////////
// Call the red-black solver
//////////////////////////////////////////////////////////////
_HermitianRBSolver(src_o,sol_o); // CGNE_prec_MdagM(solution[Odd],src);
///////////////////////////////////////////////////
// sol_e = M_ee^-1 * ( src_e - Meo sol_o )...
///////////////////////////////////////////////////
_Matrix.Meooe(sol_o,tmp); // Meo(solution[Odd],tmp,Even,DaggerNo);
src_e = src_e-tmp; // axpy(src,tmp,source[Even],-1.0);
_Matrix.MooeeInv(src_e,sol_e); // MooeeInv(src,solution[Even],DaggerNo,Even);
setCheckerboard(out,sol_e);
setCheckerboard(out,sol_o);
}
};
}
#endif

View File

@ -0,0 +1,15 @@
- ConjugateGradientMultiShift
- MCR
- Potentially Useful Boost libraries
- MultiArray
- Aligned allocator; memory pool
- Remez -- Mike or Boost?
- Multiprecision
- quaternians
- Tokenize
- Serialization
- Regex
- Proto (ET)
- uBlas