mirror of
https://github.com/paboyle/Grid.git
synced 2025-06-13 20:57:06 +01:00
Global edit on HMC sector -- making GaugeField a template parameter and
preparing to pass integrator, smearing, bc's as policy classes to hmc. Propose to unify "integrator" and integrator algorithm in a base/derived way to override step. Want to read through ForceGradient to ensure that abstraction covers the force gradient case.
This commit is contained in:
@ -62,7 +62,6 @@ namespace QCD {
|
||||
template<typename vtype> using iGparitySpinColourVector = iVector<iVector<iVector<vtype, Nc>, Nhs>, Ngp >;
|
||||
template<typename vtype> using iGparityHalfSpinColourVector = iVector<iVector<iVector<vtype, Nc>, Nhs>, Ngp >;
|
||||
|
||||
|
||||
// Spin matrix
|
||||
typedef iSpinMatrix<Complex > SpinMatrix;
|
||||
typedef iSpinMatrix<ComplexF > SpinMatrixF;
|
||||
@ -154,7 +153,7 @@ namespace QCD {
|
||||
typedef iHalfSpinColourVector<vComplexD> vHalfSpinColourVectorD;
|
||||
|
||||
// singlets
|
||||
typedef iSinglet<Complex > TComplex; // FIXME This is painful. Tensor singlet complex type.
|
||||
typedef iSinglet<Complex > TComplex; // FIXME This is painful. Tensor singlet complex type.
|
||||
typedef iSinglet<ComplexF> TComplexF; // FIXME This is painful. Tensor singlet complex type.
|
||||
typedef iSinglet<ComplexD> TComplexD; // FIXME This is painful. Tensor singlet complex type.
|
||||
|
||||
@ -162,7 +161,7 @@ namespace QCD {
|
||||
typedef iSinglet<vComplexF> vTComplexF; // what if we don't know the tensor structure
|
||||
typedef iSinglet<vComplexD> vTComplexD; // what if we don't know the tensor structure
|
||||
|
||||
typedef iSinglet<Real > TReal; // Shouldn't need these; can I make it work without?
|
||||
typedef iSinglet<Real > TReal; // Shouldn't need these; can I make it work without?
|
||||
typedef iSinglet<RealF> TRealF; // Shouldn't need these; can I make it work without?
|
||||
typedef iSinglet<RealD> TRealD; // Shouldn't need these; can I make it work without?
|
||||
|
||||
@ -251,6 +250,8 @@ namespace QCD {
|
||||
typedef LatticeDoubleStoredColourMatrixF LatticeDoubledGaugeFieldF;
|
||||
typedef LatticeDoubleStoredColourMatrixD LatticeDoubledGaugeFieldD;
|
||||
|
||||
template<class GF> using LorentzScalar = Lattice<iScalar<typename GF::vector_object::element> >;
|
||||
|
||||
// Uhgg... typing this hurt ;)
|
||||
// (my keyboard got burning hot when I typed this, must be the anti-Fermion)
|
||||
typedef Lattice<vColourVector> LatticeStaggeredFermion;
|
||||
|
@ -16,10 +16,6 @@
|
||||
#include <qcd/action/ActionBase.h>
|
||||
#include <qcd/action/ActionParams.h>
|
||||
|
||||
////////////////////////////////////////////
|
||||
// Gauge Actions
|
||||
////////////////////////////////////////////
|
||||
#include <qcd/action/gauge/WilsonGaugeAction.h>
|
||||
|
||||
|
||||
////////////////////////////////////////////
|
||||
@ -31,6 +27,17 @@
|
||||
#include <qcd/action/fermion/WilsonKernels.h> //used by all wilson type fermions
|
||||
|
||||
|
||||
////////////////////////////////////////////
|
||||
// Gauge Actions
|
||||
////////////////////////////////////////////
|
||||
#include <qcd/action/gauge/WilsonGaugeAction.h>
|
||||
namespace Grid {
|
||||
namespace QCD {
|
||||
typedef WilsonGaugeAction<LatticeGaugeField> WilsonGaugeActionR;
|
||||
typedef WilsonGaugeAction<LatticeGaugeFieldF> WilsonGaugeActionF;
|
||||
typedef WilsonGaugeAction<LatticeGaugeFieldD> WilsonGaugeActionD;
|
||||
}}
|
||||
|
||||
////////////////////////////////////////////////////////////////////////////////////////////////////
|
||||
// Explicit explicit template instantiation is still required in the .cc files
|
||||
//
|
||||
|
@ -7,8 +7,12 @@ namespace Grid{
|
||||
////////////////////////////////////////////////////////////////////////
|
||||
// Wilson Gauge Action .. should I template the Nc etc..
|
||||
////////////////////////////////////////////////////////////////////////
|
||||
template<class GaugeField, class MatrixField>
|
||||
class WilsonGaugeAction : public Action<GaugeField> {
|
||||
template<class GaugeField>
|
||||
class WilsonGaugeAction : public Action<GaugeField> {
|
||||
public:
|
||||
|
||||
typedef LorentzScalar<GaugeField> GaugeLinkField;
|
||||
|
||||
private:
|
||||
RealD beta;
|
||||
public:
|
||||
@ -17,23 +21,26 @@ namespace Grid{
|
||||
virtual void init(const GaugeField &U, GridParallelRNG& pRNG) {};
|
||||
|
||||
virtual RealD S(const GaugeField &U) {
|
||||
RealD plaq = WilsonLoops<MatrixField,GaugeField>::avgPlaquette(U);
|
||||
RealD plaq = WilsonLoops<GaugeField>::avgPlaquette(U);
|
||||
std::cout<<GridLogMessage << "Plaq : "<<plaq << "\n";
|
||||
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) {
|
||||
//not optimal implementation FIXME
|
||||
//extend Ta to include Lorentz indexes
|
||||
RealD factor = 0.5*beta/RealD(Nc);
|
||||
MatrixField dSdU_mu(U._grid);
|
||||
MatrixField Umu(U._grid);
|
||||
GaugeLinkField Umu(U._grid);
|
||||
GaugeLinkField dSdU_mu(U._grid);
|
||||
for (int mu=0; mu < Nd; mu++){
|
||||
|
||||
Umu = PeekIndex<LorentzIndex>(U,mu);
|
||||
|
||||
// Staple in direction mu
|
||||
WilsonLoops<MatrixField,GaugeField>::Staple(dSdU_mu,U,mu);
|
||||
WilsonLoops<GaugeField>::Staple(dSdU_mu,U,mu);
|
||||
dSdU_mu = Ta(Umu*adj(dSdU_mu))*factor;
|
||||
pokeLorentz(dSdU, dSdU_mu, mu);
|
||||
}
|
||||
|
@ -15,6 +15,7 @@
|
||||
namespace Grid{
|
||||
namespace QCD{
|
||||
|
||||
|
||||
struct HMCparameters{
|
||||
Integer Nsweeps; /* @brief Number of sweeps in this run */
|
||||
Integer TotalSweeps; /* @brief If provided, the total number of sweeps */
|
||||
@ -26,14 +27,17 @@ namespace Grid{
|
||||
HMCparameters();
|
||||
};
|
||||
|
||||
template <class Algorithm>
|
||||
class HybridMonteCarlo{
|
||||
// template <class GaugeField, class Integrator, class Smearer, class Boundary>
|
||||
template <class GaugeField, class Algorithm>
|
||||
class HybridMonteCarlo {
|
||||
private:
|
||||
|
||||
const HMCparameters Params;
|
||||
|
||||
GridSerialRNG sRNG; // Fixme: need a RNG management strategy.
|
||||
|
||||
Integrator<Algorithm>& MD;
|
||||
|
||||
GridSerialRNG &sRNG; // Fixme: need a RNG management strategy.
|
||||
GridParallelRNG &pRNG; // Fixme: need a RNG management strategy.
|
||||
typedef Integrator<GaugeField,Algorithm> IntegratorType;
|
||||
IntegratorType &TheIntegrator;
|
||||
|
||||
/////////////////////////////////////////////////////////
|
||||
// Metropolis step
|
||||
@ -63,16 +67,20 @@ namespace Grid{
|
||||
/////////////////////////////////////////////////////////
|
||||
// Evolution
|
||||
/////////////////////////////////////////////////////////
|
||||
RealD evolve_step(LatticeGaugeField& U){
|
||||
MD.init(U); // set U and initialize P and phi's
|
||||
RealD evolve_step(GaugeField& U){
|
||||
|
||||
TheIntegrator.init(U); // set U and initialize P and phi's
|
||||
|
||||
RealD H0 = TheIntegrator.S(U); // initial state action
|
||||
|
||||
RealD H0 = MD.S(U); // initial state action
|
||||
std::cout<<GridLogMessage<<"Total H before = "<< H0 << "\n";
|
||||
|
||||
MD.integrate(U);
|
||||
TheIntegrator.integrate(U);
|
||||
|
||||
RealD H1 = MD.S(U); // updated state action
|
||||
RealD H1 = TheIntegrator.S(U); // updated state action
|
||||
|
||||
std::cout<<GridLogMessage<<"Total H after = "<< H1 << "\n";
|
||||
|
||||
return (H1-H0);
|
||||
}
|
||||
|
||||
@ -81,16 +89,18 @@ namespace Grid{
|
||||
/////////////////////////////////////////
|
||||
// Constructor
|
||||
/////////////////////////////////////////
|
||||
HybridMonteCarlo(HMCparameters Pms, Integrator<Algorithm>& MolDyn): Params(Pms),MD(MolDyn) {
|
||||
|
||||
//FIXME... initialize RNGs also with seed ; RNG management strategy
|
||||
sRNG.SeedRandomDevice();
|
||||
|
||||
HybridMonteCarlo(HMCparameters Pms, IntegratorType &_Int, GridSerialRNG &_sRNG, GridParallelRNG &_pRNG ) :
|
||||
Params(Pms),
|
||||
TheIntegrator(_Int),
|
||||
sRNG(_sRNG),
|
||||
pRNG(_pRNG)
|
||||
{
|
||||
}
|
||||
~HybridMonteCarlo(){};
|
||||
|
||||
|
||||
void evolve(LatticeGaugeField& Uin){
|
||||
void evolve(GaugeField& Uin){
|
||||
|
||||
Real DeltaH;
|
||||
|
||||
// Thermalizations
|
||||
@ -102,7 +112,7 @@ namespace Grid{
|
||||
}
|
||||
|
||||
// Actual updates (evolve a copy Ucopy then copy back eventually)
|
||||
LatticeGaugeField Ucopy(Uin._grid);
|
||||
GaugeField Ucopy(Uin._grid);
|
||||
for(int iter=Params.StartingConfig; iter < Params.Nsweeps+Params.StartingConfig; ++iter){
|
||||
std::cout<<GridLogMessage << "-- # Sweep = "<< iter << "\n";
|
||||
|
||||
|
@ -1,28 +0,0 @@
|
||||
/*!
|
||||
@file Integrator_base.cc
|
||||
@brief utilities for MD including funcs to generate initial HMC momentum
|
||||
*/
|
||||
|
||||
#include <Grid.h>
|
||||
|
||||
namespace Grid{
|
||||
namespace QCD{
|
||||
|
||||
void MDutils::generate_momenta(LatticeLorentzColourMatrix& P,GridParallelRNG& pRNG){
|
||||
// for future support of different groups
|
||||
MDutils::generate_momenta_su3(P, pRNG);
|
||||
}
|
||||
|
||||
void MDutils::generate_momenta_su3(LatticeLorentzColourMatrix& P,GridParallelRNG& pRNG){
|
||||
LatticeColourMatrix Pmu(P._grid);
|
||||
Pmu = zero;
|
||||
for(int mu=0;mu<Nd;mu++){
|
||||
SU3::GaussianLieAlgebraMatrix(pRNG, Pmu);
|
||||
PokeIndex<LorentzIndex>(P, Pmu, mu);
|
||||
}
|
||||
|
||||
}
|
||||
|
||||
|
||||
}
|
||||
}
|
@ -10,27 +10,15 @@
|
||||
#ifndef INTEGRATOR_INCLUDED
|
||||
#define INTEGRATOR_INCLUDED
|
||||
|
||||
class Observer;
|
||||
//class Observer;
|
||||
|
||||
#include <memory>
|
||||
|
||||
namespace Grid{
|
||||
namespace QCD{
|
||||
|
||||
typedef Action<LatticeGaugeField>* ActPtr; // now force the same colours as the rest of the code
|
||||
struct ActionLevel{
|
||||
int multiplier;
|
||||
public:
|
||||
std::vector<ActPtr> actions;
|
||||
explicit ActionLevel(int mul = 1):multiplier(mul){assert (mul > 0);};
|
||||
void push_back(ActPtr ptr){
|
||||
actions.push_back(ptr);
|
||||
}
|
||||
};
|
||||
typedef std::vector<ActionLevel> ActionSet;
|
||||
typedef std::vector<Observer*> ObserverList;
|
||||
|
||||
struct IntegratorParameters{
|
||||
|
||||
int Nexp;
|
||||
int MDsteps; // number of outer steps
|
||||
RealD trajL; // trajectory length
|
||||
@ -39,94 +27,133 @@ namespace Grid{
|
||||
IntegratorParameters(int Nexp_,
|
||||
int MDsteps_,
|
||||
RealD trajL_):
|
||||
Nexp(Nexp_),MDsteps(MDsteps_),trajL(trajL_),stepsize(trajL/MDsteps){};
|
||||
Nexp(Nexp_),
|
||||
MDsteps(MDsteps_),
|
||||
trajL(trajL_),
|
||||
stepsize(trajL/MDsteps)
|
||||
{
|
||||
// empty body constructor
|
||||
};
|
||||
|
||||
};
|
||||
|
||||
namespace MDutils{
|
||||
void generate_momenta(LatticeGaugeField&,GridParallelRNG&);
|
||||
void generate_momenta_su3(LatticeGaugeField&,GridParallelRNG&);
|
||||
// Should match any legal (SU(n)) gauge field
|
||||
// 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){
|
||||
typedef Lattice< iScalar< iScalar< iMatrix<vec,Ncol> > > > GaugeLinkField;
|
||||
GaugeLinkField Pmu(P._grid);
|
||||
Pmu = zero;
|
||||
for(int mu=0;mu<Nd;mu++){
|
||||
SU<Ncol>::GaussianLieAlgebraMatrix(pRNG, Pmu);
|
||||
PokeIndex<LorentzIndex>(P, Pmu, mu);
|
||||
}
|
||||
}
|
||||
|
||||
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 IntegratorAlgorithm >
|
||||
class Integrator{
|
||||
template<class GaugeField, class Algorithm >
|
||||
class Integrator : public Algorithm {
|
||||
|
||||
private:
|
||||
|
||||
int levels;
|
||||
std::vector<double> t_P;
|
||||
double t_U;
|
||||
|
||||
IntegratorParameters Params;
|
||||
const ActionSet as;
|
||||
std::unique_ptr<LatticeGaugeField> P;
|
||||
GridParallelRNG pRNG;
|
||||
|
||||
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
|
||||
|
||||
IntegratorAlgorithm TheIntegrator;
|
||||
// typedef std::vector<Observer*> ObserverList;
|
||||
// void register_observers();
|
||||
// void notify_observers();
|
||||
|
||||
void register_observers();
|
||||
|
||||
void notify_observers();
|
||||
|
||||
void update_P(LatticeGaugeField&U, int level,double ep){
|
||||
void update_P(GaugeField&U, int level,double ep){
|
||||
t_P[level]+=ep;
|
||||
for(int a=0; a<as[level].actions.size(); ++a){
|
||||
LatticeGaugeField force(U._grid);
|
||||
GaugeField force(U._grid);
|
||||
as[level].actions.at(a)->deriv(U,force);
|
||||
*P -= force*ep;
|
||||
P = P - force*ep;
|
||||
}
|
||||
|
||||
std::cout<<GridLogMessage;
|
||||
for(int l=0; l<level;++l) std::cout<<" ";
|
||||
std::cout<<"["<<level<<"] P " << " dt "<< ep <<" : t_P "<< t_P[level] <<std::endl;
|
||||
|
||||
}
|
||||
|
||||
void update_U(LatticeGaugeField&U, double ep){
|
||||
void update_U(GaugeField&U, double ep){
|
||||
//rewrite exponential to deal automatically with the lorentz index?
|
||||
LatticeColourMatrix Umu(U._grid);
|
||||
LatticeColourMatrix Pmu(U._grid);
|
||||
// GaugeLinkField Umu(U._grid);
|
||||
// GaugeLinkField Pmu(U._grid);
|
||||
for (int mu = 0; mu < Nd; mu++){
|
||||
Umu=PeekIndex<LorentzIndex>(U, mu);
|
||||
Pmu=PeekIndex<LorentzIndex>(*P, mu);
|
||||
auto Umu=PeekIndex<LorentzIndex>(U, mu);
|
||||
auto Pmu=PeekIndex<LorentzIndex>(P, mu);
|
||||
Umu = expMat(Pmu, ep, Params.Nexp)*Umu;
|
||||
PokeIndex<LorentzIndex>(U, Umu, mu);
|
||||
}
|
||||
|
||||
t_U+=ep;
|
||||
int fl = levels-1;
|
||||
std::cout<<GridLogMessage<<" ";
|
||||
for(int l=0; l<fl;++l) std::cout<<" ";
|
||||
std::cout<<"["<<fl<<"] U " << " dt "<< ep <<" : t_U "<< t_U <<std::endl;
|
||||
|
||||
}
|
||||
|
||||
friend void IntegratorAlgorithm::step (LatticeGaugeField& U,
|
||||
int level, std::vector<int>& clock,
|
||||
Integrator<IntegratorAlgorithm>* Integ);
|
||||
friend void Algorithm::step (GaugeField& U,
|
||||
int level,
|
||||
std::vector<int>& clock,
|
||||
Integrator<GaugeField,Algorithm>* Integ);
|
||||
|
||||
public:
|
||||
|
||||
Integrator(GridBase* grid,
|
||||
IntegratorParameters Par,
|
||||
ActionSet& Aset):
|
||||
ActionSet<GaugeField> & Aset):
|
||||
Params(Par),
|
||||
as(Aset),
|
||||
P(new LatticeGaugeField(grid)),
|
||||
P(grid),
|
||||
pRNG(grid),
|
||||
levels(Aset.size())
|
||||
{
|
||||
|
||||
std::vector<int> seeds({1,2,3,4,5});
|
||||
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_U=0.0;
|
||||
|
||||
};
|
||||
|
||||
~Integrator(){}
|
||||
|
||||
//Initialization of momenta and actions
|
||||
void init(LatticeGaugeField& U){
|
||||
void init(GaugeField& U){
|
||||
std::cout<<GridLogMessage<< "Integrator init\n";
|
||||
MDutils::generate_momenta(*P,pRNG);
|
||||
generate_momenta(P,pRNG);
|
||||
for(int level=0; level< as.size(); ++level){
|
||||
for(int actionID=0; actionID<as[level].actions.size(); ++actionID){
|
||||
as[level].actions.at(actionID)->init(U, pRNG);
|
||||
@ -135,12 +162,12 @@ namespace Grid{
|
||||
}
|
||||
|
||||
// Calculate action
|
||||
RealD S(LatticeGaugeField& U){
|
||||
LatticeComplex Hloc(U._grid);
|
||||
Hloc = zero;
|
||||
RealD S(GaugeField& U){
|
||||
|
||||
LatticeComplex Hloc(U._grid); Hloc = zero;
|
||||
// Momenta
|
||||
for (int mu=0; mu <Nd; mu++){
|
||||
LatticeColourMatrix Pmu = peekLorentz(*P, mu);
|
||||
auto Pmu = PeekIndex<LorentzIndex>(P, mu);
|
||||
Hloc -= trace(Pmu*Pmu);
|
||||
}
|
||||
Complex Hsum = sum(Hloc);
|
||||
@ -162,16 +189,23 @@ namespace Grid{
|
||||
return H;
|
||||
}
|
||||
|
||||
void integrate(LatticeGaugeField& U){
|
||||
void integrate(GaugeField& U){
|
||||
|
||||
std::vector<int> clock;
|
||||
|
||||
clock.resize(as.size(),0);
|
||||
|
||||
// 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
|
||||
TheIntegrator.step(U,0,clock, (this));
|
||||
Algorithm::step(U,0,clock, (this));
|
||||
}
|
||||
|
||||
// Check the clocks all match
|
||||
for(int level=0; level<as.size(); ++level){
|
||||
assert(fabs(t_U - t_P[level])<1.0e-3);
|
||||
assert(fabs(t_U - t_P[level])<1.0e-4);
|
||||
std::cout<<GridLogMessage<<" times["<<level<<"]= "<<t_P[level]<< " " << t_U <<std::endl;
|
||||
}
|
||||
|
||||
}
|
||||
};
|
||||
|
||||
|
@ -12,111 +12,13 @@
|
||||
namespace Grid{
|
||||
namespace QCD{
|
||||
|
||||
/*
|
||||
Chroma: Recursive min norm
|
||||
00132 Real dtau = traj_length / Real(n_steps);
|
||||
00133 Real lambda_dt = dtau*lambda;
|
||||
00134 Real dtauby2 = dtau / Real(2);
|
||||
00135 Real one_minus_2lambda_dt = (Real(1)-Real(2)*lambda)*dtau;
|
||||
00136 Real two_lambda_dt = lambda_dt*Real(2);
|
||||
00137
|
||||
00138 // Its sts so:
|
||||
00139 expSdt(s, lambda_dt);
|
||||
00140 for(int i=0; i < n_steps-1; i++) { // N-1 full steps
|
||||
00141 // Roll the exp(lambda_dt T) here and start
|
||||
00142 // Next iter into one
|
||||
00143 subIntegrator(s, dtauby2); <--- either leapU or next integrator
|
||||
00144 expSdt(s, one_minus_2lambda_dt);
|
||||
00145 subIntegrator(s, dtauby2); <--- either leapU or next integrator
|
||||
00146 expSdt(s, two_lambda_dt);
|
||||
00147 }
|
||||
00148 // Last step, can't roll the first and last exp(lambda_dt T)
|
||||
00149 // together.
|
||||
00150 subIntegrator(s, dtauby2);
|
||||
00151 expSdt(s, one_minus_2lambda_dt);
|
||||
00152 subIntegrator(s, dtauby2);
|
||||
00153 expSdt(s, lambda_dt);
|
||||
|
||||
*
|
||||
*/
|
||||
|
||||
|
||||
class MinimumNorm2{
|
||||
const double lambda = 0.1931833275037836;
|
||||
|
||||
public:
|
||||
void step (LatticeLorentzColourMatrix& U,
|
||||
int level, std::vector<int>& clock,
|
||||
Integrator<MinimumNorm2>* Integ){
|
||||
|
||||
// level : current level
|
||||
// fl : final level
|
||||
// eps : current step size
|
||||
|
||||
int fl = Integ->as.size() -1;
|
||||
|
||||
double eps = Integ->Params.stepsize;
|
||||
|
||||
for(int l=0; l<=level; ++l) eps/= 2.0*Integ->as[l].multiplier;
|
||||
|
||||
// which is final half step
|
||||
int fin = Integ->as[0].multiplier;
|
||||
for(int l=1; l<=level; ++l) fin*= 2.0*Integ->as[l].multiplier;
|
||||
fin = 3*Integ->Params.MDsteps*fin -1;
|
||||
|
||||
int multiplier = Integ->as[level].multiplier;
|
||||
for(int e=0; e<multiplier; ++e){ // steps per step
|
||||
|
||||
if(clock[level] == 0){ // initial half step
|
||||
Integ->update_P(U,level,lambda*eps); ++clock[level];
|
||||
}
|
||||
|
||||
if(level == fl){ // lowest level
|
||||
Integ->update_U(U,0.5*eps);
|
||||
}else{ // recursive function call
|
||||
step(U,level+1,clock, Integ);
|
||||
}
|
||||
|
||||
Integ->update_P(U,level,(1.0-2.0*lambda)*eps); ++clock[level];
|
||||
|
||||
if(level == fl){ // lowest level
|
||||
Integ->update_U(U,0.5*eps);
|
||||
}else{ // recursive function call
|
||||
step(U,level+1,clock, Integ);
|
||||
}
|
||||
|
||||
// Handle the final half step
|
||||
std::cout << GridLogMessage << " P[" <<level<<"] clock = "<<clock[level]<<"/"<<fin<<std::endl;
|
||||
int mm = (clock[level]==fin) ? 1 : 2;
|
||||
Integ->update_P(U,level,lambda*eps*mm); clock[level]+=mm;
|
||||
|
||||
}
|
||||
|
||||
}
|
||||
|
||||
};
|
||||
|
||||
|
||||
|
||||
/*
|
||||
|
||||
Chroma: Recursive leapfrog
|
||||
|
||||
|
||||
00124 Real dtau = traj_length / n_steps;
|
||||
00125 Real dtauby2 = dtau/2;
|
||||
00126
|
||||
00127 // Its sts so:
|
||||
00128 expSdt(s, dtauby2); // First half step
|
||||
00129 for(int i=0; i < n_steps-1; i++) { // N-1 full steps
|
||||
00130 subIntegrator(s, dtau); <--- either leapU or next integrator
|
||||
00131 expSdt(s, dtau);
|
||||
00132 }
|
||||
00133 subIntegrator(s, dtau); // Last Full Step
|
||||
00134 expSdt(s, dtauby2); // Last Half Step
|
||||
|
||||
/* PAB:
|
||||
*
|
||||
* Recursive leapfrog; explanation of nested stepping
|
||||
*
|
||||
* Nested 1:4; units in dt for top level integrator
|
||||
* CHROMA GUIDO
|
||||
*
|
||||
* CHROMA IroIro
|
||||
* 0 1 0
|
||||
* P 1/2 P 1/2
|
||||
* P 1/16 P1/16
|
||||
@ -129,7 +31,7 @@ Chroma: Recursive leapfrog
|
||||
* U 1/8 U1/8
|
||||
* P 1/16 P1/8
|
||||
* P 1 P 1
|
||||
* P 1/16 * skipped
|
||||
* P 1/16 * skipped --- avoids revaluating force
|
||||
* U 1/8 U1/8
|
||||
* P 1/8 P1/8
|
||||
* U 1/8 U1/8
|
||||
@ -159,13 +61,16 @@ Chroma: Recursive leapfrog
|
||||
* U 1/8 U1/8
|
||||
* P 1/16 P1/16
|
||||
* P 1/2 P 1/2
|
||||
* Total
|
||||
*/
|
||||
class LeapFrog{
|
||||
|
||||
template<class GaugeField> class LeapFrog {
|
||||
public:
|
||||
void step (LatticeLorentzColourMatrix& U,
|
||||
|
||||
typedef LeapFrog<GaugeField> Algorithm;
|
||||
|
||||
void step (GaugeField& U,
|
||||
int level, std::vector<int>& clock,
|
||||
Integrator<LeapFrog>* Integ){
|
||||
Integrator<GaugeField,Algorithm> * Integ){
|
||||
|
||||
// level : current level
|
||||
// fl : final level
|
||||
@ -186,13 +91,7 @@ Chroma: Recursive leapfrog
|
||||
|
||||
int first_step,last_step;
|
||||
|
||||
if ( level==0 ) {
|
||||
first_step = (clock[level]==0);
|
||||
last_step = (clock[level]==fin);
|
||||
} else {
|
||||
first_step = (e==0);
|
||||
last_step = (e==multiplier-1);
|
||||
}
|
||||
first_step = (clock[level]==0);
|
||||
|
||||
if(first_step){ // initial half step
|
||||
Integ->update_P(U, level,eps/2.0); ++clock[level];
|
||||
@ -204,8 +103,7 @@ Chroma: Recursive leapfrog
|
||||
step(U, level+1,clock, Integ);
|
||||
}
|
||||
|
||||
// Handle the final half step
|
||||
std::cout << GridLogMessage << " P[" <<level<<"] clock = "<<clock[level]<<"/"<<fin<<std::endl;
|
||||
last_step = (clock[level]==fin);
|
||||
int mm = last_step ? 1 : 2;
|
||||
Integ->update_P(U, level,mm*eps/2.0);
|
||||
clock[level]+=mm;
|
||||
@ -214,6 +112,65 @@ Chroma: Recursive leapfrog
|
||||
}
|
||||
};
|
||||
|
||||
template<class GaugeField> class MinimumNorm2 {
|
||||
public:
|
||||
typedef MinimumNorm2<GaugeField> Algorithm;
|
||||
|
||||
private:
|
||||
const double lambda = 0.1931833275037836;
|
||||
public:
|
||||
|
||||
void step (GaugeField& U,
|
||||
int level, std::vector<int>& clock,
|
||||
Integrator<GaugeField,Algorithm>* Integ){
|
||||
|
||||
// level : current level
|
||||
// fl : final level
|
||||
// eps : current step size
|
||||
|
||||
int fl = Integ->as.size() -1;
|
||||
|
||||
double eps = Integ->Params.stepsize;
|
||||
|
||||
for(int l=0; l<=level; ++l) eps/= 2.0*Integ->as[l].multiplier;
|
||||
|
||||
// which is final half step
|
||||
int fin = Integ->as[0].multiplier;
|
||||
for(int l=1; l<=level; ++l) fin*= 2.0*Integ->as[l].multiplier;
|
||||
fin = 3*Integ->Params.MDsteps*fin -1;
|
||||
|
||||
int multiplier = Integ->as[level].multiplier;
|
||||
for(int e=0; e<multiplier; ++e){ // steps per step
|
||||
|
||||
int first_step,last_step;
|
||||
|
||||
first_step = (clock[level]==0);
|
||||
|
||||
if(first_step){ // initial half step
|
||||
Integ->update_P(U,level,lambda*eps); ++clock[level];
|
||||
}
|
||||
|
||||
if(level == fl){ // lowest level
|
||||
Integ->update_U(U,0.5*eps);
|
||||
}else{ // recursive function call
|
||||
step(U,level+1,clock, Integ);
|
||||
}
|
||||
|
||||
Integ->update_P(U,level,(1.0-2.0*lambda)*eps); ++clock[level];
|
||||
|
||||
if(level == fl){ // lowest level
|
||||
Integ->update_U(U,0.5*eps);
|
||||
}else{ // recursive function call
|
||||
step(U,level+1,clock, Integ);
|
||||
}
|
||||
|
||||
last_step = (clock[level]==fin);
|
||||
int mm = (last_step) ? 1 : 2;
|
||||
Integ->update_P(U,level,lambda*eps*mm); clock[level]+=mm;
|
||||
|
||||
}
|
||||
}
|
||||
};
|
||||
|
||||
}
|
||||
}
|
||||
|
@ -4,9 +4,12 @@ namespace Grid {
|
||||
namespace QCD {
|
||||
|
||||
// Common wilson loop observables
|
||||
template<class GaugeMat,class GaugeLorentz>
|
||||
template<class GaugeLorentz>
|
||||
class WilsonLoops {
|
||||
public:
|
||||
|
||||
typedef LorentzScalar<GaugeLorentz> GaugeMat;
|
||||
|
||||
//////////////////////////////////////////////////
|
||||
// directed plaquette oriented in mu,nu plane
|
||||
//////////////////////////////////////////////////
|
||||
@ -141,10 +144,10 @@ void siteRectangle(GaugeMat &plaq,const std::vector<GaugeMat> &U, const int mu,
|
||||
};
|
||||
|
||||
|
||||
typedef WilsonLoops<LatticeColourMatrix,LatticeGaugeField> ColourWilsonLoops;
|
||||
typedef WilsonLoops<LatticeColourMatrix,LatticeGaugeField> U1WilsonLoops;
|
||||
typedef WilsonLoops<LatticeColourMatrix,LatticeGaugeField> SU2WilsonLoops;
|
||||
typedef WilsonLoops<LatticeColourMatrix,LatticeGaugeField> SU3WilsonLoops;
|
||||
typedef WilsonLoops<LatticeGaugeField> ColourWilsonLoops;
|
||||
typedef WilsonLoops<LatticeGaugeField> U1WilsonLoops;
|
||||
typedef WilsonLoops<LatticeGaugeField> SU2WilsonLoops;
|
||||
typedef WilsonLoops<LatticeGaugeField> SU3WilsonLoops;
|
||||
|
||||
}}
|
||||
|
||||
|
Reference in New Issue
Block a user