mirror of
https://github.com/paboyle/Grid.git
synced 2025-04-25 13:15:55 +01:00
Unified integrator and integrator algorithm into virtual class used as a policy for the
HMC.
This commit is contained in:
parent
eed889ea05
commit
29fd004d54
@ -7,16 +7,15 @@ template<class GaugeField>
|
|||||||
class Action {
|
class Action {
|
||||||
|
|
||||||
public:
|
public:
|
||||||
virtual void init (const GaugeField &U, GridParallelRNG& pRNG) = 0; //
|
// Boundary conditions? // Heatbath?
|
||||||
|
virtual void refresh(const GaugeField &U, GridParallelRNG& pRNG) = 0;// refresh pseudofermions
|
||||||
virtual RealD S (const GaugeField &U) = 0; // evaluate the action
|
virtual RealD S (const GaugeField &U) = 0; // evaluate the action
|
||||||
virtual void deriv(const GaugeField &U,GaugeField & dSdU ) = 0; // evaluate the action derivative
|
virtual void deriv(const GaugeField &U,GaugeField & dSdU ) = 0; // evaluate the action derivative
|
||||||
virtual void refresh(const GaugeField & ) {}; // Default to no-op for actions with no internal fields
|
|
||||||
// Boundary conditions?
|
|
||||||
// Heatbath?
|
|
||||||
virtual ~Action() {};
|
virtual ~Action() {};
|
||||||
};
|
};
|
||||||
|
|
||||||
// Could derive PseudoFermion action with a PF field, FermionField, and a Grid; implement refresh
|
// Could derive PseudoFermion action with a PF field, FermionField, and a Grid; implement refresh
|
||||||
|
/*
|
||||||
template<class GaugeField, class FermionField>
|
template<class GaugeField, class FermionField>
|
||||||
class PseudoFermionAction : public Action<GaugeField> {
|
class PseudoFermionAction : public Action<GaugeField> {
|
||||||
public:
|
public:
|
||||||
@ -32,5 +31,28 @@ class PseudoFermionAction : public Action<GaugeField> {
|
|||||||
};
|
};
|
||||||
|
|
||||||
};
|
};
|
||||||
|
*/
|
||||||
|
|
||||||
|
template<class GaugeField> struct ActionLevel{
|
||||||
|
public:
|
||||||
|
|
||||||
|
typedef Action<GaugeField>* ActPtr; // now force the same colours as the rest of the code
|
||||||
|
|
||||||
|
int multiplier;
|
||||||
|
|
||||||
|
std::vector<ActPtr> actions;
|
||||||
|
|
||||||
|
ActionLevel(int mul = 1) : multiplier(mul) {
|
||||||
|
assert (mul > 0);
|
||||||
|
};
|
||||||
|
|
||||||
|
void push_back(ActPtr ptr){
|
||||||
|
actions.push_back(ptr);
|
||||||
|
}
|
||||||
|
};
|
||||||
|
|
||||||
|
template<class GaugeField> using ActionSet = std::vector<ActionLevel< GaugeField > >;
|
||||||
|
|
||||||
|
|
||||||
}}
|
}}
|
||||||
#endif
|
#endif
|
||||||
|
@ -18,7 +18,7 @@ namespace Grid{
|
|||||||
public:
|
public:
|
||||||
WilsonGaugeAction(RealD b):beta(b){};
|
WilsonGaugeAction(RealD b):beta(b){};
|
||||||
|
|
||||||
virtual void init(const GaugeField &U, GridParallelRNG& pRNG) {};
|
virtual void refresh(const GaugeField &U, GridParallelRNG& pRNG) {}; // noop as no pseudoferms
|
||||||
|
|
||||||
virtual RealD S(const GaugeField &U) {
|
virtual RealD S(const GaugeField &U) {
|
||||||
RealD plaq = WilsonLoops<GaugeField>::avgPlaquette(U);
|
RealD plaq = WilsonLoops<GaugeField>::avgPlaquette(U);
|
||||||
|
@ -61,7 +61,7 @@ namespace Grid{
|
|||||||
PowerNegQuarter.Init(remez,param.tolerance,true);
|
PowerNegQuarter.Init(remez,param.tolerance,true);
|
||||||
};
|
};
|
||||||
|
|
||||||
virtual void init(const GaugeField &U, GridParallelRNG& pRNG) {
|
virtual void refresh(const GaugeField &U, GridParallelRNG& pRNG) {
|
||||||
|
|
||||||
// P(phi) = e^{- phi^dag (MpcdagMpc)^-1/2 phi}
|
// P(phi) = e^{- phi^dag (MpcdagMpc)^-1/2 phi}
|
||||||
// = e^{- phi^dag (MpcdagMpc)^-1/4 (MpcdagMpc)^-1/4 phi}
|
// = e^{- phi^dag (MpcdagMpc)^-1/4 (MpcdagMpc)^-1/4 phi}
|
||||||
|
@ -61,7 +61,7 @@ namespace Grid{
|
|||||||
PowerNegQuarter.Init(remez,param.tolerance,true);
|
PowerNegQuarter.Init(remez,param.tolerance,true);
|
||||||
};
|
};
|
||||||
|
|
||||||
virtual void init(const GaugeField &U, GridParallelRNG& pRNG) {
|
virtual void refresh(const GaugeField &U, GridParallelRNG& pRNG) {
|
||||||
|
|
||||||
// S_f = chi^dag* P(V^dag*V)/Q(V^dag*V)* N(M^dag*M)/D(M^dag*M)* P(V^dag*V)/Q(V^dag*V)* chi
|
// S_f = chi^dag* P(V^dag*V)/Q(V^dag*V)* N(M^dag*M)/D(M^dag*M)* P(V^dag*V)/Q(V^dag*V)* chi
|
||||||
//
|
//
|
||||||
|
@ -57,7 +57,7 @@ namespace Grid{
|
|||||||
PowerNegQuarter.Init(remez,param.tolerance,true);
|
PowerNegQuarter.Init(remez,param.tolerance,true);
|
||||||
};
|
};
|
||||||
|
|
||||||
virtual void init(const GaugeField &U, GridParallelRNG& pRNG) {
|
virtual void refresh(const GaugeField &U, GridParallelRNG& pRNG) {
|
||||||
|
|
||||||
// P(phi) = e^{- phi^dag (MdagM)^-1/2 phi}
|
// P(phi) = e^{- phi^dag (MdagM)^-1/2 phi}
|
||||||
// = e^{- phi^dag (MdagM)^-1/4 (MdagM)^-1/4 phi}
|
// = e^{- phi^dag (MdagM)^-1/4 (MdagM)^-1/4 phi}
|
||||||
|
@ -55,7 +55,7 @@ namespace Grid{
|
|||||||
PowerNegQuarter.Init(remez,param.tolerance,true);
|
PowerNegQuarter.Init(remez,param.tolerance,true);
|
||||||
};
|
};
|
||||||
|
|
||||||
virtual void init(const GaugeField &U, GridParallelRNG& pRNG) {
|
virtual void refresh(const GaugeField &U, GridParallelRNG& pRNG) {
|
||||||
|
|
||||||
// S_f = chi^dag* P(V^dag*V)/Q(V^dag*V)* N(M^dag*M)/D(M^dag*M)* P(V^dag*V)/Q(V^dag*V)* chi
|
// S_f = chi^dag* P(V^dag*V)/Q(V^dag*V)* N(M^dag*M)/D(M^dag*M)* P(V^dag*V)/Q(V^dag*V)* chi
|
||||||
//
|
//
|
||||||
|
@ -35,7 +35,7 @@ namespace Grid{
|
|||||||
//////////////////////////////////////////////////////////////////////////////////////
|
//////////////////////////////////////////////////////////////////////////////////////
|
||||||
// Push the gauge field in to the dops. Assume any BC's and smearing already applied
|
// Push the gauge field in to the dops. Assume any BC's and smearing already applied
|
||||||
//////////////////////////////////////////////////////////////////////////////////////
|
//////////////////////////////////////////////////////////////////////////////////////
|
||||||
virtual void init(const GaugeField &U, GridParallelRNG& pRNG) {
|
virtual void refresh(const GaugeField &U, GridParallelRNG& pRNG) {
|
||||||
|
|
||||||
// P(phi) = e^{- phi^dag (MdagM)^-1 phi}
|
// P(phi) = e^{- phi^dag (MdagM)^-1 phi}
|
||||||
// Phi = Mdag eta
|
// Phi = Mdag eta
|
||||||
|
@ -44,7 +44,7 @@ namespace Grid{
|
|||||||
//////////////////////////////////////////////////////////////////////////////////////
|
//////////////////////////////////////////////////////////////////////////////////////
|
||||||
// Push the gauge field in to the dops. Assume any BC's and smearing already applied
|
// Push the gauge field in to the dops. Assume any BC's and smearing already applied
|
||||||
//////////////////////////////////////////////////////////////////////////////////////
|
//////////////////////////////////////////////////////////////////////////////////////
|
||||||
virtual void init(const GaugeField &U, GridParallelRNG& pRNG) {
|
virtual void refresh(const GaugeField &U, GridParallelRNG& pRNG) {
|
||||||
|
|
||||||
// P(phi) = e^{- phi^dag (MpcdagMpc)^-1 phi}
|
// P(phi) = e^{- phi^dag (MpcdagMpc)^-1 phi}
|
||||||
// Phi = McpDag eta
|
// Phi = McpDag eta
|
||||||
|
@ -40,7 +40,7 @@ namespace Grid{
|
|||||||
conformable(_NumOp.GaugeRedBlackGrid(), _DenOp.GaugeRedBlackGrid());
|
conformable(_NumOp.GaugeRedBlackGrid(), _DenOp.GaugeRedBlackGrid());
|
||||||
};
|
};
|
||||||
|
|
||||||
virtual void init(const GaugeField &U, GridParallelRNG& pRNG) {
|
virtual void refresh(const GaugeField &U, GridParallelRNG& pRNG) {
|
||||||
|
|
||||||
// P(phi) = e^{- phi^dag Vpc (MpcdagMpc)^-1 Vpcdag phi}
|
// P(phi) = e^{- phi^dag Vpc (MpcdagMpc)^-1 Vpcdag phi}
|
||||||
//
|
//
|
||||||
|
@ -28,7 +28,7 @@ namespace Grid{
|
|||||||
OperatorFunction<FermionField> & AS
|
OperatorFunction<FermionField> & AS
|
||||||
) : NumOp(_NumOp), DenOp(_DenOp), DerivativeSolver(DS), ActionSolver(AS), Phi(_NumOp.FermionGrid()) {};
|
) : NumOp(_NumOp), DenOp(_DenOp), DerivativeSolver(DS), ActionSolver(AS), Phi(_NumOp.FermionGrid()) {};
|
||||||
|
|
||||||
virtual void init(const GaugeField &U, GridParallelRNG& pRNG) {
|
virtual void refresh(const GaugeField &U, GridParallelRNG& pRNG) {
|
||||||
|
|
||||||
// P(phi) = e^{- phi^dag V (MdagM)^-1 Vdag phi}
|
// P(phi) = e^{- phi^dag V (MdagM)^-1 Vdag phi}
|
||||||
//
|
//
|
||||||
|
@ -28,7 +28,7 @@ namespace Grid{
|
|||||||
};
|
};
|
||||||
|
|
||||||
// template <class GaugeField, class Integrator, class Smearer, class Boundary>
|
// template <class GaugeField, class Integrator, class Smearer, class Boundary>
|
||||||
template <class GaugeField, class Algorithm>
|
template <class GaugeField, class IntegratorType>
|
||||||
class HybridMonteCarlo {
|
class HybridMonteCarlo {
|
||||||
private:
|
private:
|
||||||
|
|
||||||
@ -36,7 +36,7 @@ namespace Grid{
|
|||||||
|
|
||||||
GridSerialRNG &sRNG; // Fixme: need a RNG management strategy.
|
GridSerialRNG &sRNG; // Fixme: need a RNG management strategy.
|
||||||
GridParallelRNG &pRNG; // Fixme: need a RNG management strategy.
|
GridParallelRNG &pRNG; // Fixme: need a RNG management strategy.
|
||||||
typedef Integrator<GaugeField,Algorithm> IntegratorType;
|
|
||||||
IntegratorType &TheIntegrator;
|
IntegratorType &TheIntegrator;
|
||||||
|
|
||||||
/////////////////////////////////////////////////////////
|
/////////////////////////////////////////////////////////
|
||||||
@ -69,7 +69,7 @@ namespace Grid{
|
|||||||
/////////////////////////////////////////////////////////
|
/////////////////////////////////////////////////////////
|
||||||
RealD evolve_step(GaugeField& U){
|
RealD evolve_step(GaugeField& U){
|
||||||
|
|
||||||
TheIntegrator.init(U); // set U and initialize P and phi's
|
TheIntegrator.refresh(U,pRNG); // set U and initialize P and phi's
|
||||||
|
|
||||||
RealD H0 = TheIntegrator.S(U); // initial state action
|
RealD H0 = TheIntegrator.S(U); // initial state action
|
||||||
|
|
||||||
|
@ -24,9 +24,9 @@ namespace Grid{
|
|||||||
RealD trajL; // trajectory length
|
RealD trajL; // trajectory length
|
||||||
RealD stepsize;
|
RealD stepsize;
|
||||||
|
|
||||||
IntegratorParameters(int Nexp_,
|
IntegratorParameters(int MDsteps_,
|
||||||
int MDsteps_,
|
RealD trajL_=1.0,
|
||||||
RealD trajL_):
|
int Nexp_=12):
|
||||||
Nexp(Nexp_),
|
Nexp(Nexp_),
|
||||||
MDsteps(MDsteps_),
|
MDsteps(MDsteps_),
|
||||||
trajL(trajL_),
|
trajL(trajL_),
|
||||||
@ -37,6 +37,24 @@ namespace Grid{
|
|||||||
|
|
||||||
};
|
};
|
||||||
|
|
||||||
|
/*! @brief Class for Molecular Dynamics management */
|
||||||
|
template<class GaugeField>
|
||||||
|
class Integrator {
|
||||||
|
|
||||||
|
protected:
|
||||||
|
|
||||||
|
typedef IntegratorParameters ParameterType;
|
||||||
|
|
||||||
|
IntegratorParameters Params;
|
||||||
|
|
||||||
|
const ActionSet<GaugeField> as;
|
||||||
|
|
||||||
|
int levels; //
|
||||||
|
double t_U; // Track time passing on each level and for U and for P
|
||||||
|
std::vector<double> t_P; //
|
||||||
|
|
||||||
|
GaugeField P;
|
||||||
|
|
||||||
// Should match any legal (SU(n)) gauge field
|
// Should match any legal (SU(n)) gauge field
|
||||||
// Need to use this template to match Ncol to pass to SU<N> class
|
// Need to use this template to match Ncol to pass to SU<N> class
|
||||||
template<int Ncol,class vec> void generate_momenta(Lattice< iVector< iScalar< iMatrix<vec,Ncol> >, Nd> > & P,GridParallelRNG& pRNG){
|
template<int Ncol,class vec> void generate_momenta(Lattice< iVector< iScalar< iMatrix<vec,Ncol> >, Nd> > & P,GridParallelRNG& pRNG){
|
||||||
@ -49,42 +67,6 @@ namespace Grid{
|
|||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
template<class GaugeField> struct ActionLevel{
|
|
||||||
public:
|
|
||||||
|
|
||||||
typedef Action<GaugeField>* ActPtr; // now force the same colours as the rest of the code
|
|
||||||
|
|
||||||
int multiplier;
|
|
||||||
|
|
||||||
std::vector<ActPtr> actions;
|
|
||||||
|
|
||||||
ActionLevel(int mul = 1) : multiplier(mul) {
|
|
||||||
assert (mul > 0);
|
|
||||||
};
|
|
||||||
|
|
||||||
void push_back(ActPtr ptr){
|
|
||||||
actions.push_back(ptr);
|
|
||||||
}
|
|
||||||
};
|
|
||||||
|
|
||||||
template<class GaugeField> using ActionSet = std::vector<ActionLevel< GaugeField > >;
|
|
||||||
|
|
||||||
/*! @brief Class for Molecular Dynamics management */
|
|
||||||
template<class GaugeField, class Algorithm >
|
|
||||||
class Integrator : public Algorithm {
|
|
||||||
|
|
||||||
private:
|
|
||||||
|
|
||||||
IntegratorParameters Params;
|
|
||||||
|
|
||||||
GridParallelRNG pRNG; // Store this somewhere more sensible and pass as reference
|
|
||||||
const ActionSet<GaugeField> as;
|
|
||||||
|
|
||||||
int levels; //
|
|
||||||
double t_U; // Track time passing on each level and for U and for P
|
|
||||||
std::vector<double> t_P; //
|
|
||||||
|
|
||||||
GaugeField P; // is a pointer really necessary?
|
|
||||||
|
|
||||||
//ObserverList observers; // not yet
|
//ObserverList observers; // not yet
|
||||||
// typedef std::vector<Observer*> ObserverList;
|
// typedef std::vector<Observer*> ObserverList;
|
||||||
@ -125,10 +107,14 @@ namespace Grid{
|
|||||||
|
|
||||||
}
|
}
|
||||||
|
|
||||||
|
/*
|
||||||
friend void Algorithm::step (GaugeField& U,
|
friend void Algorithm::step (GaugeField& U,
|
||||||
int level,
|
int level,
|
||||||
std::vector<int>& clock,
|
std::vector<int>& clock,
|
||||||
Integrator<GaugeField,Algorithm>* Integ);
|
Integrator<GaugeField,Algorithm>* Integ);
|
||||||
|
*/
|
||||||
|
|
||||||
|
virtual void step (GaugeField& U,int level, std::vector<int>& clock)=0;
|
||||||
|
|
||||||
public:
|
public:
|
||||||
|
|
||||||
@ -138,25 +124,21 @@ namespace Grid{
|
|||||||
Params(Par),
|
Params(Par),
|
||||||
as(Aset),
|
as(Aset),
|
||||||
P(grid),
|
P(grid),
|
||||||
pRNG(grid),
|
|
||||||
levels(Aset.size())
|
levels(Aset.size())
|
||||||
{
|
{
|
||||||
std::vector<int> seeds({1,2,3,4,5}); // Fixme; Pass it the RNG as a ref
|
|
||||||
pRNG.SeedFixedIntegers(seeds);
|
|
||||||
|
|
||||||
t_P.resize(levels,0.0);
|
t_P.resize(levels,0.0);
|
||||||
t_U=0.0;
|
t_U=0.0;
|
||||||
};
|
};
|
||||||
|
|
||||||
~Integrator(){}
|
virtual ~Integrator(){}
|
||||||
|
|
||||||
//Initialization of momenta and actions
|
//Initialization of momenta and actions
|
||||||
void init(GaugeField& U){
|
void refresh(GaugeField& U,GridParallelRNG &pRNG){
|
||||||
std::cout<<GridLogMessage<< "Integrator init\n";
|
std::cout<<GridLogMessage<< "Integrator refresh\n";
|
||||||
generate_momenta(P,pRNG);
|
generate_momenta(P,pRNG);
|
||||||
for(int level=0; level< as.size(); ++level){
|
for(int level=0; level< as.size(); ++level){
|
||||||
for(int actionID=0; actionID<as[level].actions.size(); ++actionID){
|
for(int actionID=0; actionID<as[level].actions.size(); ++actionID){
|
||||||
as[level].actions.at(actionID)->init(U, pRNG);
|
as[level].actions.at(actionID)->refresh(U, pRNG);
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
@ -197,12 +179,14 @@ namespace Grid{
|
|||||||
|
|
||||||
// All the clock stuff is removed if we pass first, last to the step down the way
|
// All the clock stuff is removed if we pass first, last to the step down the way
|
||||||
for(int step=0; step< Params.MDsteps; ++step){ // MD step
|
for(int step=0; step< Params.MDsteps; ++step){ // MD step
|
||||||
Algorithm::step(U,0,clock, (this));
|
int first_step = (step==0);
|
||||||
|
int last_step = (step==Params.MDsteps-1);
|
||||||
|
this->step(U,0,clock);
|
||||||
}
|
}
|
||||||
|
|
||||||
// Check the clocks all match
|
// Check the clocks all match
|
||||||
for(int level=0; level<as.size(); ++level){
|
for(int level=0; level<as.size(); ++level){
|
||||||
assert(fabs(t_U - t_P[level])<1.0e-4);
|
assert(fabs(t_U - t_P[level])<1.0e-6); // must be the same
|
||||||
std::cout<<GridLogMessage<<" times["<<level<<"]= "<<t_P[level]<< " " << t_U <<std::endl;
|
std::cout<<GridLogMessage<<" times["<<level<<"]= "<<t_P[level]<< " " << t_U <<std::endl;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
@ -63,30 +63,32 @@ namespace Grid{
|
|||||||
* P 1/2 P 1/2
|
* P 1/2 P 1/2
|
||||||
*/
|
*/
|
||||||
|
|
||||||
template<class GaugeField> class LeapFrog {
|
template<class GaugeField> class LeapFrog : public Integrator<GaugeField> {
|
||||||
public:
|
public:
|
||||||
|
|
||||||
typedef LeapFrog<GaugeField> Algorithm;
|
typedef LeapFrog<GaugeField> Algorithm;
|
||||||
|
|
||||||
void step (GaugeField& U,
|
LeapFrog(GridBase* grid,
|
||||||
int level, std::vector<int>& clock,
|
IntegratorParameters Par,
|
||||||
Integrator<GaugeField,Algorithm> * Integ){
|
ActionSet<GaugeField> & Aset): Integrator<GaugeField>(grid,Par,Aset) {};
|
||||||
|
|
||||||
|
|
||||||
|
void step (GaugeField& U, int level, std::vector<int>& clock){
|
||||||
|
|
||||||
|
int fl = this->as.size() -1;
|
||||||
// level : current level
|
// level : current level
|
||||||
// fl : final level
|
// fl : final level
|
||||||
// eps : current step size
|
// eps : current step size
|
||||||
|
|
||||||
int fl = Integ->as.size() -1;
|
|
||||||
|
|
||||||
// Get current level step size
|
// Get current level step size
|
||||||
int fin = 2*Integ->Params.MDsteps;
|
int fin = 2*this->Params.MDsteps;
|
||||||
for(int l=0; l<=level; ++l) fin*= Integ->as[l].multiplier;
|
for(int l=0; l<=level; ++l) fin*= this->as[l].multiplier;
|
||||||
fin = fin-1;
|
fin = fin-1;
|
||||||
|
|
||||||
double eps = Integ->Params.stepsize;
|
double eps = this->Params.stepsize;
|
||||||
for(int l=0; l<=level; ++l) eps/= Integ->as[l].multiplier;
|
for(int l=0; l<=level; ++l) eps/= this->as[l].multiplier;
|
||||||
|
|
||||||
int multiplier = Integ->as[level].multiplier;
|
int multiplier = this->as[level].multiplier;
|
||||||
for(int e=0; e<multiplier; ++e){
|
for(int e=0; e<multiplier; ++e){
|
||||||
|
|
||||||
int first_step,last_step;
|
int first_step,last_step;
|
||||||
@ -94,52 +96,52 @@ namespace Grid{
|
|||||||
first_step = (clock[level]==0);
|
first_step = (clock[level]==0);
|
||||||
|
|
||||||
if(first_step){ // initial half step
|
if(first_step){ // initial half step
|
||||||
Integ->update_P(U, level,eps/2.0); ++clock[level];
|
this->update_P(U, level,eps/2.0); ++clock[level];
|
||||||
}
|
}
|
||||||
|
|
||||||
if(level == fl){ // lowest level
|
if(level == fl){ // lowest level
|
||||||
Integ->update_U(U, eps);
|
this->update_U(U, eps);
|
||||||
}else{ // recursive function call
|
}else{ // recursive function call
|
||||||
step(U, level+1,clock, Integ);
|
this->step(U, level+1,clock);
|
||||||
}
|
}
|
||||||
|
|
||||||
last_step = (clock[level]==fin);
|
last_step = (clock[level]==fin);
|
||||||
int mm = last_step ? 1 : 2;
|
int mm = last_step ? 1 : 2;
|
||||||
Integ->update_P(U, level,mm*eps/2.0);
|
this->update_P(U, level,mm*eps/2.0);
|
||||||
clock[level]+=mm;
|
clock[level]+=mm;
|
||||||
|
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
};
|
};
|
||||||
|
|
||||||
template<class GaugeField> class MinimumNorm2 {
|
template<class GaugeField> class MinimumNorm2 : public Integrator<GaugeField> {
|
||||||
public:
|
|
||||||
typedef MinimumNorm2<GaugeField> Algorithm;
|
|
||||||
|
|
||||||
private:
|
private:
|
||||||
const double lambda = 0.1931833275037836;
|
const double lambda = 0.1931833275037836;
|
||||||
|
|
||||||
public:
|
public:
|
||||||
|
|
||||||
void step (GaugeField& U,
|
MinimumNorm2(GridBase* grid,
|
||||||
int level, std::vector<int>& clock,
|
IntegratorParameters Par,
|
||||||
Integrator<GaugeField,Algorithm>* Integ){
|
ActionSet<GaugeField> & Aset): Integrator<GaugeField>(grid,Par,Aset) {};
|
||||||
|
|
||||||
|
void step (GaugeField& U, int level, std::vector<int>& clock){
|
||||||
|
|
||||||
// level : current level
|
// level : current level
|
||||||
// fl : final level
|
// fl : final level
|
||||||
// eps : current step size
|
// eps : current step size
|
||||||
|
|
||||||
int fl = Integ->as.size() -1;
|
int fl = this->as.size() -1;
|
||||||
|
|
||||||
double eps = Integ->Params.stepsize;
|
double eps = this->Params.stepsize;
|
||||||
|
|
||||||
for(int l=0; l<=level; ++l) eps/= 2.0*Integ->as[l].multiplier;
|
for(int l=0; l<=level; ++l) eps/= 2.0*this->as[l].multiplier;
|
||||||
|
|
||||||
// which is final half step
|
// which is final half step
|
||||||
int fin = Integ->as[0].multiplier;
|
int fin = this->as[0].multiplier;
|
||||||
for(int l=1; l<=level; ++l) fin*= 2.0*Integ->as[l].multiplier;
|
for(int l=1; l<=level; ++l) fin*= 2.0*this->as[l].multiplier;
|
||||||
fin = 3*Integ->Params.MDsteps*fin -1;
|
fin = 3*this->Params.MDsteps*fin -1;
|
||||||
|
|
||||||
int multiplier = Integ->as[level].multiplier;
|
int multiplier = this->as[level].multiplier;
|
||||||
for(int e=0; e<multiplier; ++e){ // steps per step
|
for(int e=0; e<multiplier; ++e){ // steps per step
|
||||||
|
|
||||||
int first_step,last_step;
|
int first_step,last_step;
|
||||||
@ -147,26 +149,26 @@ namespace Grid{
|
|||||||
first_step = (clock[level]==0);
|
first_step = (clock[level]==0);
|
||||||
|
|
||||||
if(first_step){ // initial half step
|
if(first_step){ // initial half step
|
||||||
Integ->update_P(U,level,lambda*eps); ++clock[level];
|
this->update_P(U,level,lambda*eps); ++clock[level];
|
||||||
}
|
}
|
||||||
|
|
||||||
if(level == fl){ // lowest level
|
if(level == fl){ // lowest level
|
||||||
Integ->update_U(U,0.5*eps);
|
this->update_U(U,0.5*eps);
|
||||||
}else{ // recursive function call
|
}else{ // recursive function call
|
||||||
step(U,level+1,clock, Integ);
|
this->step(U,level+1,clock);
|
||||||
}
|
}
|
||||||
|
|
||||||
Integ->update_P(U,level,(1.0-2.0*lambda)*eps); ++clock[level];
|
this->update_P(U,level,(1.0-2.0*lambda)*eps); ++clock[level];
|
||||||
|
|
||||||
if(level == fl){ // lowest level
|
if(level == fl){ // lowest level
|
||||||
Integ->update_U(U,0.5*eps);
|
this->update_U(U,0.5*eps);
|
||||||
}else{ // recursive function call
|
}else{ // recursive function call
|
||||||
step(U,level+1,clock, Integ);
|
this->step(U,level+1,clock);
|
||||||
}
|
}
|
||||||
|
|
||||||
last_step = (clock[level]==fin);
|
last_step = (clock[level]==fin);
|
||||||
int mm = (last_step) ? 1 : 2;
|
int mm = (last_step) ? 1 : 2;
|
||||||
Integ->update_P(U,level,lambda*eps*mm); clock[level]+=mm;
|
this->update_P(U,level,lambda*eps*mm); clock[level]+=mm;
|
||||||
|
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
@ -47,14 +47,13 @@ int main (int argc, char ** argv)
|
|||||||
FullSet.push_back(Level1);
|
FullSet.push_back(Level1);
|
||||||
|
|
||||||
// Create integrator
|
// Create integrator
|
||||||
// typedef LeapFrog IntegratorAlgorithm;// change here to change the algorithm
|
typedef MinimumNorm2<LatticeGaugeField> IntegratorType;// change here to change the algorithm
|
||||||
typedef MinimumNorm2<LatticeGaugeField> IntegratorAlgorithm;// change here to change the algorithm
|
IntegratorParameters MDpar(20);
|
||||||
IntegratorParameters MDpar(12,20,1.0);
|
IntegratorType MDynamics(UGrid,MDpar, FullSet);
|
||||||
Integrator<LatticeGaugeField,IntegratorAlgorithm> MDynamics(UGrid,MDpar, FullSet);
|
|
||||||
|
|
||||||
// Create HMC
|
// Create HMC
|
||||||
HMCparameters HMCpar;
|
HMCparameters HMCpar;
|
||||||
HybridMonteCarlo<LatticeGaugeField,IntegratorAlgorithm> HMC(HMCpar, MDynamics,sRNG,pRNG);
|
HybridMonteCarlo<LatticeGaugeField,IntegratorType> HMC(HMCpar, MDynamics,sRNG,pRNG);
|
||||||
|
|
||||||
// Run it
|
// Run it
|
||||||
HMC.evolve(U);
|
HMC.evolve(U);
|
||||||
|
@ -52,8 +52,8 @@ int main (int argc, char ** argv)
|
|||||||
// typedef LeapFrog IntegratorAlgorithm;// change here to change the algorithm
|
// typedef LeapFrog IntegratorAlgorithm;// change here to change the algorithm
|
||||||
// IntegratorParameters MDpar(12,40,1.0);
|
// IntegratorParameters MDpar(12,40,1.0);
|
||||||
typedef MinimumNorm2<LatticeGaugeField> IntegratorAlgorithm;// change here to change the algorithm
|
typedef MinimumNorm2<LatticeGaugeField> IntegratorAlgorithm;// change here to change the algorithm
|
||||||
IntegratorParameters MDpar(12,10,1.0);
|
IntegratorParameters MDpar(10);
|
||||||
Integrator<LatticeGaugeField,IntegratorAlgorithm> MDynamics(&Fine,MDpar, FullSet);
|
IntegratorAlgorithm MDynamics(&Fine,MDpar, FullSet);
|
||||||
|
|
||||||
// Create HMC
|
// Create HMC
|
||||||
HMCparameters HMCpar;
|
HMCparameters HMCpar;
|
||||||
|
@ -51,8 +51,8 @@ int main (int argc, char ** argv)
|
|||||||
// Create integrator
|
// Create integrator
|
||||||
typedef MinimumNorm2<LatticeGaugeField> IntegratorAlgorithm;// change here to change the algorithm
|
typedef MinimumNorm2<LatticeGaugeField> IntegratorAlgorithm;// change here to change the algorithm
|
||||||
// typedef LeapFrog IntegratorAlgorithm;// change here to change the algorithm
|
// typedef LeapFrog IntegratorAlgorithm;// change here to change the algorithm
|
||||||
IntegratorParameters MDpar(12,20,1.0);
|
IntegratorParameters MDpar(20);
|
||||||
Integrator<LatticeGaugeField,IntegratorAlgorithm> MDynamics(&Fine,MDpar, FullSet);
|
IntegratorAlgorithm MDynamics(&Fine,MDpar, FullSet);
|
||||||
|
|
||||||
// Create HMC
|
// Create HMC
|
||||||
HMCparameters HMCpar;
|
HMCparameters HMCpar;
|
||||||
|
@ -52,8 +52,8 @@ int main (int argc, char ** argv)
|
|||||||
// Create integrator
|
// Create integrator
|
||||||
typedef MinimumNorm2<LatticeGaugeField> IntegratorAlgorithm;// change here to change the algorithm
|
typedef MinimumNorm2<LatticeGaugeField> IntegratorAlgorithm;// change here to change the algorithm
|
||||||
// typedef LeapFrog IntegratorAlgorithm;// change here to change the algorithm
|
// typedef LeapFrog IntegratorAlgorithm;// change here to change the algorithm
|
||||||
IntegratorParameters MDpar(12,20,1.0);
|
IntegratorParameters MDpar(20);
|
||||||
Integrator<LatticeGaugeField,IntegratorAlgorithm> MDynamics(&Fine,MDpar, FullSet);
|
IntegratorAlgorithm MDynamics(&Fine,MDpar, FullSet);
|
||||||
|
|
||||||
// Create HMC
|
// Create HMC
|
||||||
HMCparameters HMCpar;
|
HMCparameters HMCpar;
|
||||||
|
@ -47,8 +47,8 @@ int main (int argc, char ** argv)
|
|||||||
|
|
||||||
// Create integrator
|
// Create integrator
|
||||||
typedef MinimumNorm2<LatticeGaugeField> IntegratorAlgorithm;// change here to modify the algorithm
|
typedef MinimumNorm2<LatticeGaugeField> IntegratorAlgorithm;// change here to modify the algorithm
|
||||||
IntegratorParameters MDpar(12,20,1.0);
|
IntegratorParameters MDpar(20);
|
||||||
Integrator<LatticeGaugeField,IntegratorAlgorithm> MDynamics(&Fine,MDpar, FullSet);
|
IntegratorAlgorithm MDynamics(&Fine,MDpar, FullSet);
|
||||||
|
|
||||||
// Create HMC
|
// Create HMC
|
||||||
HMCparameters HMCpar;
|
HMCparameters HMCpar;
|
||||||
|
@ -51,8 +51,8 @@ int main (int argc, char ** argv)
|
|||||||
// Create integrator
|
// Create integrator
|
||||||
typedef MinimumNorm2<LatticeGaugeField> IntegratorAlgorithm;// change here to change the algorithm
|
typedef MinimumNorm2<LatticeGaugeField> IntegratorAlgorithm;// change here to change the algorithm
|
||||||
// typedef LeapFrog IntegratorAlgorithm;// change here to change the algorithm
|
// typedef LeapFrog IntegratorAlgorithm;// change here to change the algorithm
|
||||||
IntegratorParameters MDpar(12,20,1.0);
|
IntegratorParameters MDpar(20);
|
||||||
Integrator<LatticeGaugeField,IntegratorAlgorithm> MDynamics(&Fine,MDpar, FullSet);
|
IntegratorAlgorithm MDynamics(&Fine,MDpar, FullSet);
|
||||||
|
|
||||||
// Create HMC
|
// Create HMC
|
||||||
HMCparameters HMCpar;
|
HMCparameters HMCpar;
|
||||||
|
@ -52,8 +52,8 @@ int main (int argc, char ** argv)
|
|||||||
// Create integrator
|
// Create integrator
|
||||||
typedef MinimumNorm2<LatticeGaugeField> IntegratorAlgorithm;// change here to change the algorithm
|
typedef MinimumNorm2<LatticeGaugeField> IntegratorAlgorithm;// change here to change the algorithm
|
||||||
|
|
||||||
IntegratorParameters MDpar(12,20,1.0);
|
IntegratorParameters MDpar(20);
|
||||||
Integrator<LatticeGaugeField,IntegratorAlgorithm> MDynamics(&Fine,MDpar, FullSet);
|
IntegratorAlgorithm MDynamics(&Fine,MDpar, FullSet);
|
||||||
|
|
||||||
// Create HMC
|
// Create HMC
|
||||||
HMCparameters HMCpar;
|
HMCparameters HMCpar;
|
||||||
|
@ -54,8 +54,8 @@ int main (int argc, char ** argv)
|
|||||||
// Create integrator
|
// Create integrator
|
||||||
typedef MinimumNorm2<LatticeGaugeField> IntegratorAlgorithm;// change here to change the algorithm
|
typedef MinimumNorm2<LatticeGaugeField> IntegratorAlgorithm;// change here to change the algorithm
|
||||||
// typedef LeapFrog IntegratorAlgorithm;// change here to change the algorithm
|
// typedef LeapFrog IntegratorAlgorithm;// change here to change the algorithm
|
||||||
IntegratorParameters MDpar(12,20,1.0);
|
IntegratorParameters MDpar(20);
|
||||||
Integrator<LatticeGaugeField,IntegratorAlgorithm> MDynamics(&Fine,MDpar, FullSet);
|
IntegratorAlgorithm MDynamics(&Fine,MDpar, FullSet);
|
||||||
|
|
||||||
// Create HMC
|
// Create HMC
|
||||||
HMCparameters HMCpar;
|
HMCparameters HMCpar;
|
||||||
|
@ -52,8 +52,8 @@ int main (int argc, char ** argv)
|
|||||||
// Create integrator
|
// Create integrator
|
||||||
typedef MinimumNorm2<LatticeGaugeField> IntegratorAlgorithm;// change here to change the algorithm
|
typedef MinimumNorm2<LatticeGaugeField> IntegratorAlgorithm;// change here to change the algorithm
|
||||||
|
|
||||||
IntegratorParameters MDpar(12,20,1.0);
|
IntegratorParameters MDpar(20);
|
||||||
Integrator<LatticeGaugeField,IntegratorAlgorithm> MDynamics(&Fine,MDpar, FullSet);
|
IntegratorAlgorithm MDynamics(&Fine,MDpar, FullSet);
|
||||||
|
|
||||||
// Create HMC
|
// Create HMC
|
||||||
HMCparameters HMCpar;
|
HMCparameters HMCpar;
|
||||||
|
@ -54,8 +54,8 @@ int main (int argc, char ** argv)
|
|||||||
// Create integrator
|
// Create integrator
|
||||||
typedef MinimumNorm2<LatticeGaugeField> IntegratorAlgorithm;// change here to change the algorithm
|
typedef MinimumNorm2<LatticeGaugeField> IntegratorAlgorithm;// change here to change the algorithm
|
||||||
// typedef LeapFrog IntegratorAlgorithm;// change here to change the algorithm
|
// typedef LeapFrog IntegratorAlgorithm;// change here to change the algorithm
|
||||||
IntegratorParameters MDpar(12,20,1.0);
|
IntegratorParameters MDpar(20);
|
||||||
Integrator<LatticeGaugeField,IntegratorAlgorithm> MDynamics(&Fine,MDpar, FullSet);
|
IntegratorAlgorithm MDynamics(&Fine,MDpar, FullSet);
|
||||||
|
|
||||||
// Create HMC
|
// Create HMC
|
||||||
HMCparameters HMCpar;
|
HMCparameters HMCpar;
|
||||||
|
Loading…
x
Reference in New Issue
Block a user