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

Two flavour HMC for Wilson/Wilson is conserving energy.

Still to check plaq and <e(-dH)>, but nevertheless this is
progress
This commit is contained in:
Peter Boyle
2015-07-29 17:53:39 +09:00
parent 4cc2ef84d3
commit 4fe110bd07
20 changed files with 259 additions and 76 deletions

View File

@ -79,4 +79,10 @@
///////////////////////////////////////////////////////////////////////////////
#include <qcd/action/fermion/g5HermitianLinop.h>
////////////////////////////////////////
// Pseudo fermion combinations
////////////////////////////////////////
#include <qcd/action/pseudofermion/TwoFlavour.h>
#endif

View File

@ -56,6 +56,10 @@ namespace Grid {
virtual void Mdiag (const FermionField &in, FermionField &out) { Mooee(in,out);}; // Same as Mooee applied to both CB's
virtual void Mdir (const FermionField &in, FermionField &out,int dir,int disp)=0; // case by case Wilson, Clover, Cayley, ContFrac, PartFrac
///////////////////////////////////////////////
// Updates gauge field during HMC
///////////////////////////////////////////////
virtual void ImportGauge(const GaugeField & _U);
};

View File

@ -24,6 +24,10 @@ WilsonFermion::WilsonFermion(LatticeGaugeField &_Umu,
{
// Allocate the required comms buffer
comm_buf.resize(Stencil._unified_buffer_size); // this is always big enough to contain EO
ImportGauge(_Umu);
}
void WilsonFermion::ImportGauge(const LatticeGaugeField &_Umu)
{
DoubleStore(Umu,_Umu);
pickCheckerboard(Even,UmuEven,Umu);
pickCheckerboard(Odd ,UmuOdd,Umu);

View File

@ -136,6 +136,7 @@ namespace Grid {
WilsonFermion(LatticeGaugeField &_Umu,GridCartesian &Fgrid,GridRedBlackCartesian &Hgrid,RealD _mass);
// DoubleStore
virtual void ImportGauge(const LatticeGaugeField &_Umu);
void DoubleStore(LatticeDoubledGaugeField &Uds,const LatticeGaugeField &Umu);
///////////////////////////////////////////////////////////////

View File

@ -65,7 +65,10 @@ namespace QCD {
// Allocate the required comms buffer
comm_buf.resize(Stencil._unified_buffer_size); // this is always big enough to contain EO
ImportGauge(_Umu);
}
void WilsonFermion5D::ImportGauge(const LatticeGaugeField &_Umu)
{
DoubleStore(Umu,_Umu);
pickCheckerboard(Even,UmuEven,Umu);
pickCheckerboard(Odd ,UmuOdd,Umu);

View File

@ -94,6 +94,7 @@ namespace Grid {
double _M5);
// DoubleStore
virtual void ImportGauge(const LatticeGaugeField &_Umu);
void DoubleStore(LatticeDoubledGaugeField &Uds,const LatticeGaugeField &Umu);
///////////////////////////////////////////////////////////////

View File

@ -19,8 +19,10 @@ namespace Grid{
virtual RealD S(const GaugeField &U) {
RealD plaq = WilsonLoops<MatrixField,GaugeField>::avgPlaquette(U);
std::cout<<GridLogMessage << "Plaq : "<<plaq << "\n";
double vol = U._grid->gSites();
return beta*(1.0 -plaq)*(Nd*(Nd-1.0))*vol*0.5;
RealD vol = U._grid->gSites();
RealD action=beta*(1.0 -plaq)*(Nd*(Nd-1.0))*vol*0.5;
std::cout << GridLogMessage << "WilsonGauge action "<<action<<std::endl;
return action;
};
virtual void deriv(const GaugeField &U,GaugeField & dSdU) {

View File

@ -95,18 +95,18 @@ namespace Grid{
////////////////////////////////////////////////////////////////////////
// Two flavour pseudofermion action for any dop
////////////////////////////////////////////////////////////////////////
template<class GaugeField,class MatrixField,class FermionField,class FermionOperator>
template<class GaugeField,class MatrixField,class FermionField>
class TwoFlavourPseudoFermionAction : public Action<GaugeField> {
private:
FermionOperator<FermionField,GaugeField> & FermOp;// the basic operator
OperatorFunction<FermionField> &DerivativeSolver;
OperatorFunction<FermionField> &ActionSolver;
GridBase *Grid;
GridBase &Grid;
FermionField Phi; // the pseudo fermion field for this trajectory
@ -114,10 +114,11 @@ namespace Grid{
/////////////////////////////////////////////////
// Pass in required objects.
/////////////////////////////////////////////////
TwoFlavourPseudoFermionAction(FermionOperator &Op,
OperatorFunction<FermionField> & DS,
OperatorFunction<FermionField> & AS
) : FermOp(Op), DerivativeSolver(DS), ActionSolver(AS) {
TwoFlavourPseudoFermionAction(FermionOperator<FermionField,GaugeField> &Op,
OperatorFunction<FermionField> & DS,
OperatorFunction<FermionField> & AS,
GridBase &_Grid
) : FermOp(Op), DerivativeSolver(DS), ActionSolver(AS), Phi(&_Grid), Grid(_Grid) {
};
//////////////////////////////////////////////////////////////////////////////////////
@ -125,9 +126,28 @@ namespace Grid{
//////////////////////////////////////////////////////////////////////////////////////
virtual void init(const GaugeField &U, GridParallelRNG& pRNG) {
// width? Must check
gaussian(Phi,pRNG);
// P(phi) = e^{- phi^dag (MdagM)^-1 phi}
// Phi = Mdag eta
// P(eta) = e^{- eta^dag eta}
//
// e^{x^2/2 sig^2} => sig^2 = 0.5.
//
// So eta should be of width sig = 1/sqrt(2).
// and must multiply by 0.707....
//
// Chroma has this scale factor: two_flavor_monomial_w.h
// IroIro: does not use this scale. It is absorbed by a change of vars
// in the Phi integral, and thus is only an irrelevant prefactor for the partition function.
//
RealD scale = std::sqrt(0.5);
FermionField eta(&Grid);
gaussian(pRNG,eta);
FermOp.Mdag(eta,Phi);
Phi=Phi*scale;
};
//////////////////////////////////////////////////////
@ -135,38 +155,49 @@ namespace Grid{
//////////////////////////////////////////////////////
virtual RealD S(const GaugeField &U) {
FermionField X(Grid);
FermionField Y(Grid);
MdagMLinearOperator<FermionOperator<FermionField,GaugeField>,FermionField> MdagMOp(FermOp);
FermOp.ImportGauge(U);
ActionSolver(MdagMop,Phi,X);
FermionField X(&Grid);
FermionField Y(&Grid);
MdagMLinearOperator<FermionOperator<FermionField,GaugeField> ,FermionField> MdagMOp(FermOp);
X=zero;
ActionSolver(MdagMOp,Phi,X);
MdagMOp.Op(X,Y);
RealD action = norm2(Y);
std::cout << GridLogMessage << "Pseudofermion action "<<action<<std::endl;
return action;
};
//////////////////////////////////////////////////////
// dS/du = - phi^dag (Mdag M)^-1 [ Mdag dM + dMdag M ] (Mdag M)^-1 phi
// = - phi^dag M^-1 dM (MdagM)^-1 phi - phi^dag (MdagM)^-1 dMdag dM (Mdag)^-1 phi
//
// = - Ydag dM X - Xdag dMdag Y
//
//////////////////////////////////////////////////////
virtual void deriv(const GaugeField &U,GaugeField & dSdU) {
FermionField X(Grid);
FermionField Y(Grid);
GaugeField tmp(Grid);
FermOp.ImportGauge(U);
MdagMLinearOperator<FermionOperator<FermionField,GaugeField>,FermionField> MdagMOp(FermOp);
FermionField X(&Grid);
FermionField Y(&Grid);
GaugeField tmp(&Grid);
DerivativeSolver(MdagMop,Phi,X);
MdagMLinearOperator<FermionOperator<FermionField,GaugeField> ,FermionField> MdagMOp(FermOp);
X=zero;
DerivativeSolver(MdagMOp,Phi,X);
MdagMOp.Op(X,Y);
// Our conventions really make this UdSdU; We do not differentiate wrt Udag here.
// So must take dSdU - adj(dSdU) and left multiply by mom to get dS/dt.
FermOp.MDeriv(tmp , Y, X,DaggerNo ); dSdU=tmp;
FermOp.MDeriv(tmp , X, Y,DaggerYes); dSdU=-UdSdU-tmp;
FermOp.MDeriv(tmp , X, Y,DaggerYes); dSdU=dSdU+tmp;
dSdU = Ta(dSdU);
};

View File

@ -7,8 +7,8 @@ namespace Grid{
// FIXME fill this constructor now just default values
////////////////////////////// Default values
Nsweeps = 100;
TotalSweeps = 20;
Nsweeps = 200;
TotalSweeps = 220;
ThermalizationSteps = 20;
StartingConfig = 0;
SaveInterval = 1;

View File

@ -73,6 +73,7 @@ namespace Grid{
RealD H1 = MD.S(U); // updated state action
std::cout<<GridLogMessage<<"Total H after = "<< H1 << "\n";
std::cout<<GridLogMessage<<"DeltaH is "<< H1-H0 << "\n";
return (H1-H0);
}

View File

@ -63,7 +63,6 @@ namespace Grid{
}
}
void update_U(LatticeLorentzColourMatrix&U, double ep){
//rewrite exponential to deal automatically with the lorentz index?
LatticeColourMatrix Umu(U._grid);
@ -77,7 +76,6 @@ namespace Grid{
}
friend void IntegratorAlgorithm::step (LatticeLorentzColourMatrix& U,
int level, std::vector<int>& clock,
@ -92,9 +90,9 @@ namespace Grid{
~Integrator(){}
//Initialization of momenta and actions
void init(LatticeLorentzColourMatrix& U){
std::cout<<GridLogMessage<< "Integrator init\n";
MDutils::generate_momenta(*P,pRNG);
@ -105,7 +103,6 @@ namespace Grid{
}
}
// Calculate action
RealD S(LatticeLorentzColourMatrix& U){
LatticeComplex Hloc(U._grid);
@ -119,12 +116,14 @@ namespace Grid{
RealD H = Hsum.real();
std::cout<<GridLogMessage << "H_p = "<< H << "\n";
std::cout<<GridLogMessage << "Momentum action H_p = "<< H << "\n";
// Actions
for(int level=0; level<as.size(); ++level)
for(int actionID=0; actionID<as.at(level).size(); ++actionID)
H += as[level].at(actionID)->S(U);
std::cout<<GridLogMessage << "Total action H = "<< H << "\n";
return H;
}

View File

@ -32,8 +32,7 @@ namespace Grid{
int fin = Integ->Nrel[0];
for(int l=1; l<=level; ++l) fin*= 2.0*Integ->Nrel[l];
fin = 3*Integ->Params.MDsteps*fin -1;
for(int e=0; e<Integ->Nrel[level]; ++e){
if(clock[level] == 0){ // initial half step
@ -45,7 +44,6 @@ namespace Grid{
if(level == fl){ // lowest level
Integ->update_U(U,0.5*eps);
for(int l=0; l<level;++l) std::cout<<GridLogMessage<<" ";
std::cout<<GridLogMessage<<"U "<< (clock[level]+1) <<std::endl;
}else{ // recursive function call
@ -81,9 +79,7 @@ namespace Grid{
std::cout<<GridLogMessage<<"P "<< clock[level] <<std::endl;
}
}
}
};
@ -93,6 +89,7 @@ namespace Grid{
void step (LatticeLorentzColourMatrix& U,
int level, std::vector<int>& clock,
Integrator<LeapFrog>* Integ){
// level : current level
// fl : final level
// eps : current step size
@ -115,6 +112,7 @@ namespace Grid{
for(int l=0; l<level;++l) std::cout<<GridLogMessage<<" ";
std::cout<<GridLogMessage<<"P "<< 0.5*clock[level] <<std::endl;
}
if(level == fl){ // lowest level
Integ->update_U(U, eps);
for(int l=0; l<level;++l) std::cout<<GridLogMessage<<" ";
@ -122,24 +120,21 @@ namespace Grid{
}else{ // recursive function call
step(U, level+1,clock, Integ);
}
if(clock[level] == fin){ // final half step
Integ->update_P(U, level,eps/2.0);
++clock[level];
for(int l=0; l<level;++l) std::cout<<GridLogMessage<<" ";
std::cout<<GridLogMessage<<"P "<< 0.5*clock[level] <<std::endl;
}else{ // bulk step
Integ->update_P(U, level,eps);
clock[level]+=2;
for(int l=0; l<level;++l) std::cout<<GridLogMessage<<" ";
std::cout<<GridLogMessage<<"P "<< 0.5*clock[level] <<std::endl;
}
}
}
};

View File

@ -65,6 +65,18 @@ namespace Grid{
for(int a=0; a<as[level].size(); ++a){
LatticeLorentzColourMatrix force(U._grid);
as[level].at(a)->deriv(U,force);
Complex dSdt=0.0;
for(int mu=0;mu<Nd;mu++){
LatticeColourMatrix forcemu(U._grid);
LatticeColourMatrix mommu(U._grid);
forcemu=PeekIndex<LorentzIndex>(force,mu);
mommu=PeekIndex<LorentzIndex>(*P,mu);
dSdt += sum(trace(forcemu*(*P)));
}
std::cout << GridLogMessage << " action "<<level<<","<<a<<" dSdt "<< dSdt << " dt "<<ep <<std::endl;
*P -= force*ep;
}
}

View File

@ -554,9 +554,7 @@ Note that in step D setting B ~ X - A and using B in place of A in step E will g
for(int a=0;a<generators();a++){
gaussian(pRNG,ca);
generator(a,ta);
la=toComplex(ca)*ci*ta;
out += la;
}