mirror of
				https://github.com/paboyle/Grid.git
				synced 2025-11-03 21:44:33 +00:00 
			
		
		
		
	Added a mixed precision multishift algorithm for which the matrix multiplies are performed in single precision but the search directions are accumulated in double precision.
A reliable update step is performed at a tunable frequency to correct the residual. A final mixed-prec single-shift solve is performed on each pole to perform cleanup if necessary. A test is provided to demonstrate the algorithm.
This commit is contained in:
		@@ -54,6 +54,7 @@ NAMESPACE_CHECK(BiCGSTAB);
 | 
			
		||||
#include <Grid/algorithms/iterative/SchurRedBlack.h>
 | 
			
		||||
#include <Grid/algorithms/iterative/ConjugateGradientMultiShift.h>
 | 
			
		||||
#include <Grid/algorithms/iterative/ConjugateGradientMixedPrec.h>
 | 
			
		||||
#include <Grid/algorithms/iterative/ConjugateGradientMultiShiftMixedPrec.h>
 | 
			
		||||
#include <Grid/algorithms/iterative/BiCGSTABMixedPrec.h>
 | 
			
		||||
#include <Grid/algorithms/iterative/BlockConjugateGradient.h>
 | 
			
		||||
#include <Grid/algorithms/iterative/ConjugateGradientReliableUpdate.h>
 | 
			
		||||
 
 | 
			
		||||
							
								
								
									
										410
									
								
								Grid/algorithms/iterative/ConjugateGradientMultiShiftMixedPrec.h
									
									
									
									
									
										Normal file
									
								
							
							
						
						
									
										410
									
								
								Grid/algorithms/iterative/ConjugateGradientMultiShiftMixedPrec.h
									
									
									
									
									
										Normal file
									
								
							@@ -0,0 +1,410 @@
 | 
			
		||||
/*************************************************************************************
 | 
			
		||||
 | 
			
		||||
    Grid physics library, www.github.com/paboyle/Grid 
 | 
			
		||||
 | 
			
		||||
    Source file: ./lib/algorithms/iterative/ConjugateGradientMultiShift.h
 | 
			
		||||
 | 
			
		||||
    Copyright (C) 2015
 | 
			
		||||
 | 
			
		||||
Author: Azusa Yamaguchi <ayamaguc@staffmail.ed.ac.uk>
 | 
			
		||||
Author: Peter Boyle <paboyle@ph.ed.ac.uk>
 | 
			
		||||
Author: Christopher Kelly <ckelly@bnl.gov>
 | 
			
		||||
 | 
			
		||||
    This program is free software; you can redistribute it and/or modify
 | 
			
		||||
    it under the terms of the GNU General Public License as published by
 | 
			
		||||
    the Free Software Foundation; either version 2 of the License, or
 | 
			
		||||
    (at your option) any later version.
 | 
			
		||||
 | 
			
		||||
    This program is distributed in the hope that it will be useful,
 | 
			
		||||
    but WITHOUT ANY WARRANTY; without even the implied warranty of
 | 
			
		||||
    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 | 
			
		||||
    GNU General Public License for more details.
 | 
			
		||||
 | 
			
		||||
    You should have received a copy of the GNU General Public License along
 | 
			
		||||
    with this program; if not, write to the Free Software Foundation, Inc.,
 | 
			
		||||
    51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
 | 
			
		||||
 | 
			
		||||
    See the full license in the file "LICENSE" in the top level distribution directory
 | 
			
		||||
*************************************************************************************/
 | 
			
		||||
/*  END LEGAL */
 | 
			
		||||
#ifndef GRID_CONJUGATE_GRADIENT_MULTI_SHIFT_MIXEDPREC_H
 | 
			
		||||
#define GRID_CONJUGATE_GRADIENT_MULTI_SHIFT_MIXEDPREC_H
 | 
			
		||||
 | 
			
		||||
NAMESPACE_BEGIN(Grid);
 | 
			
		||||
 | 
			
		||||
//CK 2020: A variant of the multi-shift conjugate gradient with the matrix multiplication in single precision. 
 | 
			
		||||
//The residual is stored in single precision, but the search directions and solution are stored in double precision. 
 | 
			
		||||
//Every update_freq iterations the residual is corrected in double precision. 
 | 
			
		||||
    
 | 
			
		||||
//For safety the a final regular CG is applied to clean up if necessary
 | 
			
		||||
 | 
			
		||||
//Linop to add shift to input linop, used in cleanup CG
 | 
			
		||||
namespace ConjugateGradientMultiShiftMixedPrecSupport{
 | 
			
		||||
template<typename Field>
 | 
			
		||||
class ShiftedLinop: public LinearOperatorBase<Field>{
 | 
			
		||||
public:
 | 
			
		||||
  LinearOperatorBase<Field> &linop_base;
 | 
			
		||||
  RealD shift;
 | 
			
		||||
 | 
			
		||||
  ShiftedLinop(LinearOperatorBase<Field> &_linop_base, RealD _shift): linop_base(_linop_base), shift(_shift){}
 | 
			
		||||
 | 
			
		||||
  void OpDiag (const Field &in, Field &out){ assert(0); }
 | 
			
		||||
  void OpDir  (const Field &in, Field &out,int dir,int disp){ assert(0); }
 | 
			
		||||
  void OpDirAll  (const Field &in, std::vector<Field> &out){ assert(0); }
 | 
			
		||||
  
 | 
			
		||||
  void Op     (const Field &in, Field &out){ assert(0); }
 | 
			
		||||
  void AdjOp  (const Field &in, Field &out){ assert(0); }
 | 
			
		||||
 | 
			
		||||
  void HermOp(const Field &in, Field &out){
 | 
			
		||||
    linop_base.HermOp(in, out);
 | 
			
		||||
    axpy(out, shift, in, out);
 | 
			
		||||
  }    
 | 
			
		||||
 | 
			
		||||
  void HermOpAndNorm(const Field &in, Field &out,RealD &n1,RealD &n2){
 | 
			
		||||
    HermOp(in,out);
 | 
			
		||||
    ComplexD dot = innerProduct(in,out);
 | 
			
		||||
    n1=real(dot);
 | 
			
		||||
    n2=norm2(out);
 | 
			
		||||
  }
 | 
			
		||||
};
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
template<class FieldD, class FieldF,
 | 
			
		||||
	 typename std::enable_if< getPrecision<FieldD>::value == 2, int>::type = 0,
 | 
			
		||||
	 typename std::enable_if< getPrecision<FieldF>::value == 1, int>::type = 0> 
 | 
			
		||||
class ConjugateGradientMultiShiftMixedPrec : public OperatorMultiFunction<FieldD>,
 | 
			
		||||
					     public OperatorFunction<FieldD>
 | 
			
		||||
{
 | 
			
		||||
public:                                                
 | 
			
		||||
 | 
			
		||||
  using OperatorFunction<FieldD>::operator();
 | 
			
		||||
 | 
			
		||||
  RealD   Tolerance;
 | 
			
		||||
  Integer MaxIterations;
 | 
			
		||||
  Integer IterationsToComplete; //Number of iterations the CG took to finish. Filled in upon completion
 | 
			
		||||
  std::vector<int> IterationsToCompleteShift;  // Iterations for this shift
 | 
			
		||||
  int verbose;
 | 
			
		||||
  MultiShiftFunction shifts;
 | 
			
		||||
  std::vector<RealD> TrueResidualShift;
 | 
			
		||||
 | 
			
		||||
  int ReliableUpdateFreq; //number of iterations between reliable updates
 | 
			
		||||
 | 
			
		||||
  GridBase* SinglePrecGrid; //Grid for single-precision fields
 | 
			
		||||
  LinearOperatorBase<FieldF> &Linop_f; //single precision
 | 
			
		||||
 | 
			
		||||
  ConjugateGradientMultiShiftMixedPrec(Integer maxit, MultiShiftFunction &_shifts,
 | 
			
		||||
				       GridBase* _SinglePrecGrid, LinearOperatorBase<FieldF> &_Linop_f,
 | 
			
		||||
				       int _ReliableUpdateFreq
 | 
			
		||||
				       ) : 
 | 
			
		||||
    MaxIterations(maxit),  shifts(_shifts), SinglePrecGrid(_SinglePrecGrid), Linop_f(_Linop_f), ReliableUpdateFreq(_ReliableUpdateFreq)
 | 
			
		||||
  { 
 | 
			
		||||
    verbose=1;
 | 
			
		||||
    IterationsToCompleteShift.resize(_shifts.order);
 | 
			
		||||
    TrueResidualShift.resize(_shifts.order);
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  void operator() (LinearOperatorBase<FieldD> &Linop, const FieldD &src, FieldD &psi)
 | 
			
		||||
  {
 | 
			
		||||
    GridBase *grid = src.Grid();
 | 
			
		||||
    int nshift = shifts.order;
 | 
			
		||||
    std::vector<FieldD> results(nshift,grid);
 | 
			
		||||
    (*this)(Linop,src,results,psi);
 | 
			
		||||
  }
 | 
			
		||||
  void operator() (LinearOperatorBase<FieldD> &Linop, const FieldD &src, std::vector<FieldD> &results, FieldD &psi)
 | 
			
		||||
  {
 | 
			
		||||
    int nshift = shifts.order;
 | 
			
		||||
 | 
			
		||||
    (*this)(Linop,src,results);
 | 
			
		||||
  
 | 
			
		||||
    psi = shifts.norm*src;
 | 
			
		||||
    for(int i=0;i<nshift;i++){
 | 
			
		||||
      psi = psi + shifts.residues[i]*results[i];
 | 
			
		||||
    }
 | 
			
		||||
 | 
			
		||||
    return;
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  void operator() (LinearOperatorBase<FieldD> &Linop_d, const FieldD &src_d, std::vector<FieldD> &psi_d)
 | 
			
		||||
  { 
 | 
			
		||||
    GridBase *DoublePrecGrid = src_d.Grid();
 | 
			
		||||
    precisionChangeWorkspace wk_f_from_d(SinglePrecGrid, DoublePrecGrid);
 | 
			
		||||
    precisionChangeWorkspace wk_d_from_f(DoublePrecGrid, SinglePrecGrid);
 | 
			
		||||
 | 
			
		||||
    ////////////////////////////////////////////////////////////////////////
 | 
			
		||||
    // Convenience references to the info stored in "MultiShiftFunction"
 | 
			
		||||
    ////////////////////////////////////////////////////////////////////////
 | 
			
		||||
    int nshift = shifts.order;
 | 
			
		||||
 | 
			
		||||
    std::vector<RealD> &mass(shifts.poles); // Make references to array in "shifts"
 | 
			
		||||
    std::vector<RealD> &mresidual(shifts.tolerances);
 | 
			
		||||
    std::vector<RealD> alpha(nshift,1.0);
 | 
			
		||||
 | 
			
		||||
    //Double precision search directions
 | 
			
		||||
    FieldD p_d(DoublePrecGrid);
 | 
			
		||||
    std::vector<FieldD> ps_d(nshift, DoublePrecGrid);// Search directions (double precision)
 | 
			
		||||
 | 
			
		||||
    FieldD tmp_d(DoublePrecGrid);
 | 
			
		||||
    FieldD r_d(DoublePrecGrid);
 | 
			
		||||
    FieldD mmp_d(DoublePrecGrid);
 | 
			
		||||
 | 
			
		||||
    assert(psi_d.size()==nshift);
 | 
			
		||||
    assert(mass.size()==nshift);
 | 
			
		||||
    assert(mresidual.size()==nshift);
 | 
			
		||||
  
 | 
			
		||||
    // dynamic sized arrays on stack; 2d is a pain with vector
 | 
			
		||||
    RealD  bs[nshift];
 | 
			
		||||
    RealD  rsq[nshift];
 | 
			
		||||
    RealD  z[nshift][2];
 | 
			
		||||
    int     converged[nshift];
 | 
			
		||||
  
 | 
			
		||||
    const int       primary =0;
 | 
			
		||||
  
 | 
			
		||||
    //Primary shift fields CG iteration
 | 
			
		||||
    RealD a,b,c,d;
 | 
			
		||||
    RealD cp,bp,qq; //prev
 | 
			
		||||
  
 | 
			
		||||
    // Matrix mult fields
 | 
			
		||||
    FieldF r_f(SinglePrecGrid);
 | 
			
		||||
    FieldF p_f(SinglePrecGrid);
 | 
			
		||||
    FieldF tmp_f(SinglePrecGrid);
 | 
			
		||||
    FieldF mmp_f(SinglePrecGrid);
 | 
			
		||||
    FieldF src_f(SinglePrecGrid);
 | 
			
		||||
    precisionChange(src_f, src_d, wk_f_from_d);
 | 
			
		||||
 | 
			
		||||
    // Check lightest mass
 | 
			
		||||
    for(int s=0;s<nshift;s++){
 | 
			
		||||
      assert( mass[s]>= mass[primary] );
 | 
			
		||||
      converged[s]=0;
 | 
			
		||||
    }
 | 
			
		||||
  
 | 
			
		||||
    // Wire guess to zero
 | 
			
		||||
    // Residuals "r" are src
 | 
			
		||||
    // First search direction "p" is also src
 | 
			
		||||
    cp = norm2(src_d);
 | 
			
		||||
 | 
			
		||||
    // Handle trivial case of zero src.
 | 
			
		||||
    if( cp == 0. ){
 | 
			
		||||
      for(int s=0;s<nshift;s++){
 | 
			
		||||
	psi_d[s] = Zero();
 | 
			
		||||
	IterationsToCompleteShift[s] = 1;
 | 
			
		||||
	TrueResidualShift[s] = 0.;
 | 
			
		||||
      }
 | 
			
		||||
      return;
 | 
			
		||||
    }
 | 
			
		||||
 | 
			
		||||
    for(int s=0;s<nshift;s++){
 | 
			
		||||
      rsq[s] = cp * mresidual[s] * mresidual[s];
 | 
			
		||||
      std::cout<<GridLogMessage<<"ConjugateGradientMultiShiftMixedPrec: shift "<< s <<" target resid "<<rsq[s]<<std::endl;
 | 
			
		||||
      ps_d[s] = src_d;
 | 
			
		||||
    }
 | 
			
		||||
    // r and p for primary
 | 
			
		||||
    r_f=src_f; //residual maintained in single
 | 
			
		||||
    p_f=src_f;
 | 
			
		||||
    p_d = src_d; //primary copy --- make this a reference to ps_d to save axpys
 | 
			
		||||
  
 | 
			
		||||
    //MdagM+m[0]
 | 
			
		||||
    Linop_f.HermOpAndNorm(p_f,mmp_f,d,qq); // mmp = MdagM p        d=real(dot(p, mmp)),  qq=norm2(mmp)
 | 
			
		||||
    axpy(mmp_f,mass[0],p_f,mmp_f);
 | 
			
		||||
    RealD rn = norm2(p_f);
 | 
			
		||||
    d += rn*mass[0];
 | 
			
		||||
 | 
			
		||||
    b = -cp /d;
 | 
			
		||||
  
 | 
			
		||||
    // Set up the various shift variables
 | 
			
		||||
    int       iz=0;
 | 
			
		||||
    z[0][1-iz] = 1.0;
 | 
			
		||||
    z[0][iz]   = 1.0;
 | 
			
		||||
    bs[0]      = b;
 | 
			
		||||
    for(int s=1;s<nshift;s++){
 | 
			
		||||
      z[s][1-iz] = 1.0;
 | 
			
		||||
      z[s][iz]   = 1.0/( 1.0 - b*(mass[s]-mass[0]));
 | 
			
		||||
      bs[s]      = b*z[s][iz]; 
 | 
			
		||||
    }
 | 
			
		||||
  
 | 
			
		||||
    // r += b[0] A.p[0]
 | 
			
		||||
    // c= norm(r)
 | 
			
		||||
    c=axpy_norm(r_f,b,mmp_f,r_f);
 | 
			
		||||
  
 | 
			
		||||
    for(int s=0;s<nshift;s++) {
 | 
			
		||||
      axpby(psi_d[s],0.,-bs[s]*alpha[s],src_d,src_d);
 | 
			
		||||
    }
 | 
			
		||||
  
 | 
			
		||||
    ///////////////////////////////////////
 | 
			
		||||
    // Timers
 | 
			
		||||
    ///////////////////////////////////////
 | 
			
		||||
    GridStopWatch AXPYTimer, ShiftTimer, QRTimer, MatrixTimer, SolverTimer, PrecChangeTimer, CleanupTimer;
 | 
			
		||||
 | 
			
		||||
    SolverTimer.Start();
 | 
			
		||||
  
 | 
			
		||||
    // Iteration loop
 | 
			
		||||
    int k;
 | 
			
		||||
  
 | 
			
		||||
    for (k=1;k<=MaxIterations;k++){    
 | 
			
		||||
      a = c /cp;
 | 
			
		||||
 | 
			
		||||
      //Update double precision search direction by residual
 | 
			
		||||
      PrecChangeTimer.Start();
 | 
			
		||||
      precisionChange(r_d, r_f, wk_d_from_f);
 | 
			
		||||
      PrecChangeTimer.Stop();
 | 
			
		||||
 | 
			
		||||
      AXPYTimer.Start();
 | 
			
		||||
      axpy(p_d,a,p_d,r_d); 
 | 
			
		||||
 | 
			
		||||
      for(int s=0;s<nshift;s++){
 | 
			
		||||
	if ( ! converged[s] ) { 
 | 
			
		||||
	  if (s==0){
 | 
			
		||||
	    axpy(ps_d[s],a,ps_d[s],r_d);
 | 
			
		||||
	  } else{
 | 
			
		||||
	    RealD as =a *z[s][iz]*bs[s] /(z[s][1-iz]*b);
 | 
			
		||||
	    axpby(ps_d[s],z[s][iz],as,r_d,ps_d[s]);
 | 
			
		||||
	  }
 | 
			
		||||
	}
 | 
			
		||||
      }
 | 
			
		||||
      AXPYTimer.Stop();
 | 
			
		||||
 | 
			
		||||
      PrecChangeTimer.Start();
 | 
			
		||||
      precisionChange(p_f, p_d, wk_f_from_d); //get back single prec search direction for linop
 | 
			
		||||
      PrecChangeTimer.Stop();
 | 
			
		||||
 | 
			
		||||
      cp=c;
 | 
			
		||||
      MatrixTimer.Start();  
 | 
			
		||||
      Linop_f.HermOp(p_f,mmp_f); 
 | 
			
		||||
      d=real(innerProduct(p_f,mmp_f));    
 | 
			
		||||
      MatrixTimer.Stop();  
 | 
			
		||||
 | 
			
		||||
      AXPYTimer.Start();
 | 
			
		||||
      axpy(mmp_f,mass[0],p_f,mmp_f);
 | 
			
		||||
      AXPYTimer.Stop();
 | 
			
		||||
      RealD rn = norm2(p_f);
 | 
			
		||||
      d += rn*mass[0];
 | 
			
		||||
    
 | 
			
		||||
      bp=b;
 | 
			
		||||
      b=-cp/d;
 | 
			
		||||
    
 | 
			
		||||
      // Toggle the recurrence history
 | 
			
		||||
      bs[0] = b;
 | 
			
		||||
      iz = 1-iz;
 | 
			
		||||
      ShiftTimer.Start();
 | 
			
		||||
      for(int s=1;s<nshift;s++){
 | 
			
		||||
	if((!converged[s])){
 | 
			
		||||
	  RealD z0 = z[s][1-iz];
 | 
			
		||||
	  RealD z1 = z[s][iz];
 | 
			
		||||
	  z[s][iz] = z0*z1*bp
 | 
			
		||||
	    / (b*a*(z1-z0) + z1*bp*(1- (mass[s]-mass[0])*b)); 
 | 
			
		||||
	  bs[s] = b*z[s][iz]/z0; // NB sign  rel to Mike
 | 
			
		||||
	}
 | 
			
		||||
      }
 | 
			
		||||
      ShiftTimer.Stop();
 | 
			
		||||
 | 
			
		||||
      //Update double precision solutions
 | 
			
		||||
      AXPYTimer.Start();
 | 
			
		||||
      for(int s=0;s<nshift;s++){
 | 
			
		||||
	int ss = s;
 | 
			
		||||
	if( (!converged[s]) ) { 
 | 
			
		||||
	  axpy(psi_d[ss],-bs[s]*alpha[s],ps_d[s],psi_d[ss]);
 | 
			
		||||
	}
 | 
			
		||||
      }
 | 
			
		||||
 | 
			
		||||
      //Perform reliable update if necessary; otherwise update residual from single-prec mmp
 | 
			
		||||
      RealD c_f = axpy_norm(r_f,b,mmp_f,r_f);
 | 
			
		||||
      AXPYTimer.Stop();
 | 
			
		||||
 | 
			
		||||
      c = c_f;
 | 
			
		||||
 | 
			
		||||
      if(k % ReliableUpdateFreq == 0){
 | 
			
		||||
	//Replace r with true residual
 | 
			
		||||
	MatrixTimer.Start();  
 | 
			
		||||
	Linop_d.HermOp(psi_d[0],mmp_d); 
 | 
			
		||||
	MatrixTimer.Stop();  
 | 
			
		||||
 | 
			
		||||
	AXPYTimer.Start();
 | 
			
		||||
	axpy(mmp_d,mass[0],psi_d[0],mmp_d);
 | 
			
		||||
 | 
			
		||||
	RealD c_d = axpy_norm(r_d, -1.0, mmp_d, src_d);
 | 
			
		||||
	AXPYTimer.Stop();
 | 
			
		||||
 | 
			
		||||
	std::cout<<GridLogMessage<<"ConjugateGradientMultiShiftMixedPrec k="<<k<< ", replaced |r|^2 = "<<c_f <<" with |r|^2 = "<<c_d<<std::endl;
 | 
			
		||||
	
 | 
			
		||||
	PrecChangeTimer.Start();
 | 
			
		||||
	precisionChange(r_f, r_d, wk_f_from_d);
 | 
			
		||||
	PrecChangeTimer.Stop();
 | 
			
		||||
	c = c_d;
 | 
			
		||||
      }
 | 
			
		||||
    
 | 
			
		||||
      // Convergence checks
 | 
			
		||||
      int all_converged = 1;
 | 
			
		||||
      for(int s=0;s<nshift;s++){
 | 
			
		||||
      
 | 
			
		||||
	if ( (!converged[s]) ){
 | 
			
		||||
	  IterationsToCompleteShift[s] = k;
 | 
			
		||||
	
 | 
			
		||||
	  RealD css  = c * z[s][iz]* z[s][iz];
 | 
			
		||||
	
 | 
			
		||||
	  if(css<rsq[s]){
 | 
			
		||||
	    if ( ! converged[s] )
 | 
			
		||||
	      std::cout<<GridLogMessage<<"ConjugateGradientMultiShiftMixedPrec k="<<k<<" Shift "<<s<<" has converged"<<std::endl;
 | 
			
		||||
	    converged[s]=1;
 | 
			
		||||
	  } else {
 | 
			
		||||
	    all_converged=0;
 | 
			
		||||
	  }
 | 
			
		||||
 | 
			
		||||
	}
 | 
			
		||||
      }
 | 
			
		||||
 | 
			
		||||
      if ( all_converged ){
 | 
			
		||||
 | 
			
		||||
	SolverTimer.Stop();
 | 
			
		||||
	std::cout<<GridLogMessage<< "ConjugateGradientMultiShiftMixedPrec: All shifts have converged iteration "<<k<<std::endl;
 | 
			
		||||
	std::cout<<GridLogMessage<< "ConjugateGradientMultiShiftMixedPrec: Checking solutions"<<std::endl;
 | 
			
		||||
      
 | 
			
		||||
	// Check answers 
 | 
			
		||||
	for(int s=0; s < nshift; s++) { 
 | 
			
		||||
	  Linop_d.HermOpAndNorm(psi_d[s],mmp_d,d,qq);
 | 
			
		||||
	  axpy(tmp_d,mass[s],psi_d[s],mmp_d);
 | 
			
		||||
	  axpy(r_d,-alpha[s],src_d,tmp_d);
 | 
			
		||||
	  RealD rn = norm2(r_d);
 | 
			
		||||
	  RealD cn = norm2(src_d);
 | 
			
		||||
	  TrueResidualShift[s] = std::sqrt(rn/cn);
 | 
			
		||||
	  std::cout<<GridLogMessage<<"ConjugateGradientMultiShiftMixedPrec: shift["<<s<<"] true residual "<< TrueResidualShift[s] << " target " << mresidual[s] << std::endl;
 | 
			
		||||
 | 
			
		||||
	  //If we have not reached the desired tolerance, do a (mixed precision) CG cleanup
 | 
			
		||||
	  if(rn >= rsq[s]){
 | 
			
		||||
	    CleanupTimer.Start();
 | 
			
		||||
	    std::cout<<GridLogMessage<<"ConjugateGradientMultiShiftMixedPrec: performing cleanup step for shift " << s << std::endl;
 | 
			
		||||
 | 
			
		||||
	    //Setup linear operators for final cleanup
 | 
			
		||||
	    ConjugateGradientMultiShiftMixedPrecSupport::ShiftedLinop<FieldD> Linop_shift_d(Linop_d, mass[s]);
 | 
			
		||||
	    ConjugateGradientMultiShiftMixedPrecSupport::ShiftedLinop<FieldF> Linop_shift_f(Linop_f, mass[s]);
 | 
			
		||||
					       
 | 
			
		||||
	    MixedPrecisionConjugateGradient<FieldD,FieldF> cg(mresidual[s], MaxIterations, MaxIterations, SinglePrecGrid, Linop_shift_f, Linop_shift_d); 
 | 
			
		||||
	    cg(src_d, psi_d[s]);
 | 
			
		||||
	    
 | 
			
		||||
	    TrueResidualShift[s] = cg.TrueResidual;
 | 
			
		||||
	    CleanupTimer.Stop();
 | 
			
		||||
	  }
 | 
			
		||||
	}
 | 
			
		||||
 | 
			
		||||
	std::cout << GridLogMessage << "ConjugateGradientMultiShiftMixedPrec: Time Breakdown for body"<<std::endl;
 | 
			
		||||
	std::cout << GridLogMessage << "\tElapsed    " << SolverTimer.Elapsed()     <<std::endl;
 | 
			
		||||
	std::cout << GridLogMessage << "\tAXPY    " << AXPYTimer.Elapsed()     <<std::endl;
 | 
			
		||||
	std::cout << GridLogMessage << "\tMatrix    " << MatrixTimer.Elapsed()     <<std::endl;
 | 
			
		||||
	std::cout << GridLogMessage << "\tShift    " << ShiftTimer.Elapsed()     <<std::endl;
 | 
			
		||||
	std::cout << GridLogMessage << "\tPrecision Change " << PrecChangeTimer.Elapsed()     <<std::endl;
 | 
			
		||||
	std::cout << GridLogMessage << "\tFinal Cleanup " << CleanupTimer.Elapsed()     <<std::endl;
 | 
			
		||||
 | 
			
		||||
	IterationsToComplete = k;	
 | 
			
		||||
 | 
			
		||||
	return;
 | 
			
		||||
      }
 | 
			
		||||
 | 
			
		||||
   
 | 
			
		||||
    }
 | 
			
		||||
    // ugly hack
 | 
			
		||||
    std::cout<<GridLogMessage<<"CG multi shift did not converge"<<std::endl;
 | 
			
		||||
    //  assert(0);
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
};
 | 
			
		||||
NAMESPACE_END(Grid);
 | 
			
		||||
#endif
 | 
			
		||||
							
								
								
									
										93
									
								
								tests/solver/Test_dwf_multishift_mixedprec.cc
									
									
									
									
									
										Normal file
									
								
							
							
						
						
									
										93
									
								
								tests/solver/Test_dwf_multishift_mixedprec.cc
									
									
									
									
									
										Normal file
									
								
							@@ -0,0 +1,93 @@
 | 
			
		||||
    /*************************************************************************************
 | 
			
		||||
 | 
			
		||||
    Grid physics library, www.github.com/paboyle/Grid 
 | 
			
		||||
 | 
			
		||||
    Source file: ./tests/Test_dwf_multishift_mixedprec.cc
 | 
			
		||||
 | 
			
		||||
    Copyright (C) 2015
 | 
			
		||||
 | 
			
		||||
Author: Christopher Kelly <ckelly@bnl.gov>
 | 
			
		||||
 | 
			
		||||
    This program is free software; you can redistribute it and/or modify
 | 
			
		||||
    it under the terms of the GNU General Public License as published by
 | 
			
		||||
    the Free Software Foundation; either version 2 of the License, or
 | 
			
		||||
    (at your option) any later version.
 | 
			
		||||
 | 
			
		||||
    This program is distributed in the hope that it will be useful,
 | 
			
		||||
    but WITHOUT ANY WARRANTY; without even the implied warranty of
 | 
			
		||||
    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 | 
			
		||||
    GNU General Public License for more details.
 | 
			
		||||
 | 
			
		||||
    You should have received a copy of the GNU General Public License along
 | 
			
		||||
    with this program; if not, write to the Free Software Foundation, Inc.,
 | 
			
		||||
    51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
 | 
			
		||||
 | 
			
		||||
    See the full license in the file "LICENSE" in the top level distribution directory
 | 
			
		||||
    *************************************************************************************/
 | 
			
		||||
    /*  END LEGAL */
 | 
			
		||||
#include <Grid/Grid.h>
 | 
			
		||||
 | 
			
		||||
using namespace Grid;
 | 
			
		||||
 | 
			
		||||
int main (int argc, char ** argv)
 | 
			
		||||
{
 | 
			
		||||
  Grid_init(&argc, &argv);
 | 
			
		||||
 | 
			
		||||
  const int Ls = 16;
 | 
			
		||||
  GridCartesian* UGrid_d = SpaceTimeGrid::makeFourDimGrid(GridDefaultLatt(), GridDefaultSimd(Nd, vComplexD::Nsimd()), GridDefaultMpi());
 | 
			
		||||
  GridRedBlackCartesian* UrbGrid_d = SpaceTimeGrid::makeFourDimRedBlackGrid(UGrid_d);
 | 
			
		||||
  GridCartesian* FGrid_d = SpaceTimeGrid::makeFiveDimGrid(Ls, UGrid_d);
 | 
			
		||||
  GridRedBlackCartesian* FrbGrid_d = SpaceTimeGrid::makeFiveDimRedBlackGrid(Ls, UGrid_d);
 | 
			
		||||
 | 
			
		||||
  GridCartesian* UGrid_f = SpaceTimeGrid::makeFourDimGrid(GridDefaultLatt(), GridDefaultSimd(Nd, vComplexF::Nsimd()), GridDefaultMpi());
 | 
			
		||||
  GridRedBlackCartesian* UrbGrid_f = SpaceTimeGrid::makeFourDimRedBlackGrid(UGrid_f);
 | 
			
		||||
  GridCartesian* FGrid_f = SpaceTimeGrid::makeFiveDimGrid(Ls, UGrid_f);
 | 
			
		||||
  GridRedBlackCartesian* FrbGrid_f = SpaceTimeGrid::makeFiveDimRedBlackGrid(Ls, UGrid_f);
 | 
			
		||||
 | 
			
		||||
  std::vector<int> seeds4({1, 2, 3, 4});
 | 
			
		||||
  std::vector<int> seeds5({5, 6, 7, 8});
 | 
			
		||||
  GridParallelRNG RNG5(FGrid_d);
 | 
			
		||||
  RNG5.SeedFixedIntegers(seeds5);
 | 
			
		||||
  GridParallelRNG RNG4(UGrid_d);
 | 
			
		||||
  RNG4.SeedFixedIntegers(seeds4);
 | 
			
		||||
 | 
			
		||||
  LatticeFermionD src_d(FGrid_d);
 | 
			
		||||
  random(RNG5, src_d);
 | 
			
		||||
 | 
			
		||||
  LatticeGaugeFieldD Umu_d(UGrid_d);
 | 
			
		||||
  SU<Nc>::HotConfiguration(RNG4, Umu_d);
 | 
			
		||||
 | 
			
		||||
  LatticeGaugeFieldF Umu_f(UGrid_f);
 | 
			
		||||
  precisionChange(Umu_f, Umu_d);
 | 
			
		||||
 | 
			
		||||
  std::cout << GridLogMessage << "Lattice dimensions: " << GridDefaultLatt() << "   Ls: " << Ls << std::endl;
 | 
			
		||||
 | 
			
		||||
  RealD mass = 0.01;
 | 
			
		||||
  RealD M5 = 1.8;
 | 
			
		||||
  DomainWallFermionD Ddwf_d(Umu_d, *FGrid_d, *FrbGrid_d, *UGrid_d, *UrbGrid_d, mass, M5);
 | 
			
		||||
  DomainWallFermionF Ddwf_f(Umu_f, *FGrid_f, *FrbGrid_f, *UGrid_f, *UrbGrid_f, mass, M5);
 | 
			
		||||
 | 
			
		||||
  LatticeFermionD src_o_d(FrbGrid_d);
 | 
			
		||||
  pickCheckerboard(Odd, src_o_d, src_d);
 | 
			
		||||
 | 
			
		||||
  SchurDiagMooeeOperator<DomainWallFermionD, LatticeFermionD> HermOpEO_d(Ddwf_d);
 | 
			
		||||
  SchurDiagMooeeOperator<DomainWallFermionF, LatticeFermionF> HermOpEO_f(Ddwf_f);
 | 
			
		||||
 | 
			
		||||
  AlgRemez remez(1e-4, 64, 50);
 | 
			
		||||
  int order = 15;
 | 
			
		||||
  remez.generateApprox(order, 1, 2); //sqrt
 | 
			
		||||
 | 
			
		||||
  MultiShiftFunction shifts(remez, 1e-10, false);
 | 
			
		||||
 | 
			
		||||
  int relup_freq = 50;
 | 
			
		||||
  double t1=usecond();
 | 
			
		||||
  ConjugateGradientMultiShiftMixedPrec<LatticeFermionD,LatticeFermionF> mcg(10000, shifts, FrbGrid_f, HermOpEO_f, relup_freq);
 | 
			
		||||
 | 
			
		||||
  std::vector<LatticeFermionD> results_o_d(order, FrbGrid_d);
 | 
			
		||||
  mcg(HermOpEO_d, src_o_d, results_o_d);
 | 
			
		||||
  double t2=usecond();
 | 
			
		||||
 | 
			
		||||
  std::cout<<GridLogMessage << "Test: Total usec    =   "<< (t2-t1)<<std::endl;
 | 
			
		||||
 | 
			
		||||
  Grid_finalize();
 | 
			
		||||
}
 | 
			
		||||
		Reference in New Issue
	
	Block a user