mirror of
https://github.com/paboyle/Grid.git
synced 2025-06-13 20:57:06 +01:00
multishift conjugate gradient added and a strong test: take a diagonal
but non-identity matrix l1 0 0 0 .... 0 l2 0 0 .... 0 0 l3 0 ... . . . . . . . . . And apply the multishift CG to it. Sum the poles and residues. Insist that this be the same as the exactly taken square root where l1,l2,l3 >= 0.
This commit is contained in:
@ -1,5 +1,5 @@
|
||||
|
||||
bin_PROGRAMS = Test_GaugeAction Test_cayley_cg Test_cayley_even_odd Test_contfrac_cg Test_contfrac_even_odd Test_cshift Test_cshift_red_black Test_dwf_cg_prec Test_dwf_cg_schur Test_dwf_cg_unprec Test_dwf_even_odd Test_gamma Test_main Test_nersc_io Test_remez Test_rng Test_rng_fixed Test_simd Test_stencil Test_wilson_cg_prec Test_wilson_cg_schur Test_wilson_cg_unprec Test_wilson_even_odd
|
||||
bin_PROGRAMS = Test_GaugeAction Test_cayley_cg Test_cayley_even_odd Test_contfrac_cg Test_contfrac_even_odd Test_cshift Test_cshift_red_black Test_dwf_cg_prec Test_dwf_cg_schur Test_dwf_cg_unprec Test_dwf_even_odd Test_gamma Test_main Test_multishift_sqrt Test_nersc_io Test_remez Test_rng Test_rng_fixed Test_simd Test_stencil Test_wilson_cg_prec Test_wilson_cg_schur Test_wilson_cg_unprec Test_wilson_even_odd
|
||||
|
||||
|
||||
Test_GaugeAction_SOURCES=Test_GaugeAction.cc
|
||||
@ -54,6 +54,10 @@ Test_main_SOURCES=Test_main.cc
|
||||
Test_main_LDADD=-lGrid
|
||||
|
||||
|
||||
Test_multishift_sqrt_SOURCES=Test_multishift_sqrt.cc
|
||||
Test_multishift_sqrt_LDADD=-lGrid
|
||||
|
||||
|
||||
Test_nersc_io_SOURCES=Test_nersc_io.cc
|
||||
Test_nersc_io_LDADD=-lGrid
|
||||
|
||||
|
110
tests/Test_multishift_sqrt.cc
Normal file
110
tests/Test_multishift_sqrt.cc
Normal file
@ -0,0 +1,110 @@
|
||||
#include <Grid.h>
|
||||
|
||||
using namespace std;
|
||||
using namespace Grid;
|
||||
using namespace Grid::QCD;
|
||||
|
||||
template<class Field> class DumbOperator : public LinearOperatorBase<Field> {
|
||||
public:
|
||||
LatticeComplex scale;
|
||||
LatticeComplex sqrtscale;
|
||||
DumbOperator(GridBase *grid)
|
||||
: scale(grid),
|
||||
sqrtscale(grid)
|
||||
{
|
||||
GridParallelRNG pRNG(grid);
|
||||
std::vector<int> seeds({5,6,7,8});
|
||||
pRNG.SeedFixedIntegers(seeds);
|
||||
|
||||
random(pRNG,sqrtscale);
|
||||
sqrtscale = sqrtscale * adj(sqrtscale);// force real pos def
|
||||
scale = sqrtscale * sqrtscale;
|
||||
}
|
||||
void Op (const Field &in, Field &out){
|
||||
out = scale * in;
|
||||
}
|
||||
void AdjOp (const Field &in, Field &out){
|
||||
out = scale * in;
|
||||
}
|
||||
void HermOpAndNorm(const Field &in, Field &out,double &n1,double &n2){
|
||||
ComplexD dot;
|
||||
|
||||
out = scale * in;
|
||||
|
||||
dot= innerProduct(in,out);
|
||||
n1=real(dot);
|
||||
|
||||
dot = innerProduct(out,out);
|
||||
n2=real(dot);
|
||||
}
|
||||
void ApplySqrt(const Field &in, Field &out){
|
||||
out = sqrtscale * in;
|
||||
}
|
||||
};
|
||||
|
||||
|
||||
int main (int argc, char ** argv)
|
||||
{
|
||||
Grid_init(&argc,&argv);
|
||||
|
||||
GridCartesian *grid = SpaceTimeGrid::makeFourDimGrid(GridDefaultLatt(),
|
||||
GridDefaultSimd(Nd,vComplexF::Nsimd()),
|
||||
GridDefaultMpi());
|
||||
|
||||
double lo=0.001;
|
||||
double hi=1.0;
|
||||
int precision=64;
|
||||
int degree=10;
|
||||
AlgRemez remez(lo,hi,precision);
|
||||
|
||||
////////////////////////////////////////
|
||||
// sqrt and inverse sqrt
|
||||
////////////////////////////////////////
|
||||
std::cout << "Generating degree "<<degree<<" for x^(1/2)"<<std::endl;
|
||||
remez.generateApprox(degree,1,2);
|
||||
|
||||
MultiShiftFunction Sqrt(remez,1.0e-6,false);
|
||||
|
||||
GridParallelRNG pRNG(grid);
|
||||
std::vector<int> seeds({1,2,3,4});
|
||||
pRNG.SeedFixedIntegers(seeds);
|
||||
|
||||
LatticeFermion src(grid); random(pRNG,src);
|
||||
LatticeFermion combined(grid);
|
||||
LatticeFermion reference(grid);
|
||||
LatticeFermion error(grid);
|
||||
std::vector<LatticeFermion> result(degree,grid);
|
||||
LatticeFermion summed(grid);
|
||||
|
||||
ConjugateGradientMultiShift<LatticeFermion> MSCG(10000,Sqrt);
|
||||
|
||||
DumbOperator<LatticeFermion> Diagonal(grid);
|
||||
|
||||
MSCG(Diagonal,src,result);
|
||||
|
||||
|
||||
// double a = norm;
|
||||
// for(int n=0;n<poles.size();n++){
|
||||
// a = a + residues[n]/(x+poles[n]);
|
||||
// }
|
||||
assert(Sqrt.order==degree);
|
||||
|
||||
combined = Sqrt.norm*src;
|
||||
for(int i=0;i<degree;i++){
|
||||
combined = combined + Sqrt.residues[i]*result[i];
|
||||
}
|
||||
|
||||
Diagonal.ApplySqrt(src,reference);
|
||||
|
||||
error = reference - combined;
|
||||
|
||||
std::cout << " Reference "<<norm2(reference)<<std::endl;
|
||||
std::cout << " combined "<<norm2(combined) <<std::endl;
|
||||
std::cout << " error "<<norm2(error) <<std::endl;
|
||||
|
||||
MSCG(Diagonal,src,summed);
|
||||
error = summed - combined;
|
||||
std::cout << " summed-combined "<<norm2(error) <<std::endl;
|
||||
|
||||
Grid_finalize();
|
||||
}
|
@ -4,43 +4,6 @@ using namespace std;
|
||||
using namespace Grid;
|
||||
using namespace Grid::QCD;
|
||||
|
||||
class MultiShiftFunction {
|
||||
public:
|
||||
std::vector<double> poles;
|
||||
std::vector<double> residues;
|
||||
double norm;
|
||||
double lo,hi;
|
||||
MultiShiftFunction(int n,double _lo,double _hi): poles(n), residues(n), lo(_lo), hi(_hi) {;};
|
||||
double approx(double x);
|
||||
void csv(std::ostream &out);
|
||||
void gnuplot(std::ostream &out);
|
||||
};
|
||||
double MultiShiftFunction::approx(double x)
|
||||
{
|
||||
double a = norm;
|
||||
for(int n=0;n<poles.size();n++){
|
||||
a = a + residues[n]/(x+poles[n]);
|
||||
}
|
||||
return a;
|
||||
}
|
||||
void MultiShiftFunction::gnuplot(std::ostream &out)
|
||||
{
|
||||
out<<"f(x) = "<<norm<<"";
|
||||
for(int n=0;n<poles.size();n++){
|
||||
out<<"+("<<residues[n]<<"/(x+"<<poles[n]<<"))";
|
||||
}
|
||||
out<<";"<<std::endl;
|
||||
}
|
||||
void MultiShiftFunction::csv(std::ostream &out)
|
||||
{
|
||||
for (double x=lo; x<hi; x*=1.05) {
|
||||
double f = approx(x);
|
||||
double r = sqrt(x);
|
||||
out<< x<<","<<r<<","<<f<<","<<r-f<<std::endl;
|
||||
}
|
||||
return;
|
||||
}
|
||||
|
||||
int main (int argc, char ** argv)
|
||||
{
|
||||
Grid_init(&argc,&argv);
|
||||
@ -56,22 +19,18 @@ int main (int argc, char ** argv)
|
||||
////////////////////////////////////////
|
||||
// sqrt and inverse sqrt
|
||||
////////////////////////////////////////
|
||||
MultiShiftFunction Sqrt(degree,lo,hi);
|
||||
MultiShiftFunction InvSqrt(degree,lo,hi);
|
||||
|
||||
MultiShiftFunction SqrtSqrt(degree,lo,hi);
|
||||
MultiShiftFunction InvSqrtSqrt(degree,lo,hi);
|
||||
|
||||
|
||||
std::cout << "Generating degree "<<degree<<" for x^(1/2)"<<std::endl;
|
||||
remez.generateApprox(degree,1,2);
|
||||
remez.getPFE (& Sqrt.residues[0],& Sqrt.poles[0],& Sqrt.norm);
|
||||
remez.getIPFE(&InvSqrt.residues[0],&InvSqrt.poles[0],&InvSqrt.norm);
|
||||
MultiShiftFunction Sqrt(remez,1.0,false);
|
||||
MultiShiftFunction InvSqrt(remez,1.0,true);
|
||||
|
||||
|
||||
std::cout << "Generating degree "<<degree<<" for x^(1/4)"<<std::endl;
|
||||
remez.generateApprox(degree,1,4);
|
||||
remez.getPFE (&SqrtSqrt.residues[0],&SqrtSqrt.poles[0],&SqrtSqrt.norm);
|
||||
remez.getIPFE(&InvSqrtSqrt.residues[0],&InvSqrtSqrt.poles[0],&InvSqrtSqrt.norm);
|
||||
MultiShiftFunction SqrtSqrt(remez,1.0,false);
|
||||
MultiShiftFunction InvSqrtSqrt(remez,1.0,true);
|
||||
|
||||
|
||||
ofstream gnuplot(std::string("Sqrt.gnu"),std::ios::out|std::ios::trunc);
|
||||
Sqrt.gnuplot(gnuplot);
|
||||
|
@ -20,7 +20,6 @@ int main (int argc, char ** argv)
|
||||
{
|
||||
Grid_init(&argc,&argv);
|
||||
|
||||
|
||||
std::vector<int> latt_size = GridDefaultLatt();
|
||||
std::vector<int> simd_layout = GridDefaultSimd(Nd,vComplexF::Nsimd());
|
||||
std::vector<int> mpi_layout = GridDefaultMpi();
|
||||
|
Reference in New Issue
Block a user