mirror of
				https://github.com/paboyle/Grid.git
				synced 2025-11-03 21:44:33 +00:00 
			
		
		
		
	Integrator works now
This commit is contained in:
		@@ -138,7 +138,15 @@ class HMCWrapperTemplate: public HMCRunnerBase<ReaderClass> {
 | 
			
		||||
 | 
			
		||||
    // Can move this outside?
 | 
			
		||||
    typedef IntegratorType<SmearingPolicy> TheIntegrator;
 | 
			
		||||
    TheIntegrator MDynamics(UGrid, Parameters.MD, TheAction, Smearing);
 | 
			
		||||
    // Metric
 | 
			
		||||
    //TrivialMetric<typename Implementation::Field> Mtr;
 | 
			
		||||
    ConjugateGradient<LatticeGaugeField> CG(1.0e-8,10000);
 | 
			
		||||
    LaplacianParams LapPar(0.0001, 1.0, 1000, 1e-8, 12, 64);
 | 
			
		||||
    RealD Kappa = 0.9;
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
    LaplacianAdjointField<PeriodicGimplR> Laplacian(UGrid, CG, LapPar, Kappa);
 | 
			
		||||
    TheIntegrator MDynamics(UGrid, Parameters.MD, TheAction, Smearing, Laplacian);
 | 
			
		||||
 | 
			
		||||
    if (Parameters.StartingType == "HotStart") {
 | 
			
		||||
      // Hot start
 | 
			
		||||
 
 | 
			
		||||
@@ -202,7 +202,7 @@ class HybridMonteCarlo {
 | 
			
		||||
    RealD H0 = TheIntegrator.S(U);  // initial state action
 | 
			
		||||
 | 
			
		||||
    std::streamsize current_precision = std::cout.precision();
 | 
			
		||||
    std::cout.precision(17);
 | 
			
		||||
    std::cout.precision(15);
 | 
			
		||||
    std::cout << GridLogMessage << "Total H before trajectory = " << H0 << "\n";
 | 
			
		||||
    std::cout.precision(current_precision);
 | 
			
		||||
 | 
			
		||||
@@ -210,7 +210,19 @@ class HybridMonteCarlo {
 | 
			
		||||
 | 
			
		||||
    RealD H1 = TheIntegrator.S(U);  // updated state action
 | 
			
		||||
 | 
			
		||||
    std::cout.precision(17);
 | 
			
		||||
    ///////////////////////////////////////////////////////////
 | 
			
		||||
    if(0){
 | 
			
		||||
      std::cout << "------------------------- Reversibility test" << std::endl;
 | 
			
		||||
      TheIntegrator.reverse_momenta();
 | 
			
		||||
      TheIntegrator.integrate(U);
 | 
			
		||||
 | 
			
		||||
      H1 = TheIntegrator.S(U);  // updated state action
 | 
			
		||||
      std::cout << "--------------------------------------------" << std::endl;
 | 
			
		||||
    }
 | 
			
		||||
    ///////////////////////////////////////////////////////////
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
    std::cout.precision(15);
 | 
			
		||||
    std::cout << GridLogMessage << "Total H after trajectory  = " << H1
 | 
			
		||||
	      << "  dH = " << H1 - H0 << "\n";
 | 
			
		||||
    std::cout.precision(current_precision);
 | 
			
		||||
 
 | 
			
		||||
@@ -78,7 +78,7 @@ class Integrator {
 | 
			
		||||
  std::vector<double> t_P;  
 | 
			
		||||
 | 
			
		||||
  //MomentaField P;
 | 
			
		||||
  GeneralisedMomenta<FieldImplementation, TrivialMetric<MomentaField>> P;
 | 
			
		||||
  GeneralisedMomenta<FieldImplementation > P;
 | 
			
		||||
  SmearingPolicy& Smearer;
 | 
			
		||||
  RepresentationPolicy Representations;
 | 
			
		||||
  IntegratorParameters Params;
 | 
			
		||||
@@ -125,9 +125,15 @@ class Integrator {
 | 
			
		||||
      force = FieldImplementation::projectForce(force); // Ta for gauge fields
 | 
			
		||||
      Real force_abs = std::sqrt(norm2(force)/U._grid->gSites());
 | 
			
		||||
      std::cout << GridLogIntegrator << "Force average: " << force_abs << std::endl;
 | 
			
		||||
      Mom -= force * ep;
 | 
			
		||||
      Mom -= force * ep; 
 | 
			
		||||
    }
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
    MomentaField MomDer(P.Mom._grid);
 | 
			
		||||
    P.M.ImportGauge(U);
 | 
			
		||||
    P.DerivativeU(P.Mom, MomDer);
 | 
			
		||||
    Mom -= MomDer * ep;
 | 
			
		||||
 | 
			
		||||
    // Force from the other representations
 | 
			
		||||
    as[level].apply(update_P_hireps, Representations, Mom, U, ep);
 | 
			
		||||
  }
 | 
			
		||||
@@ -141,7 +147,7 @@ class Integrator {
 | 
			
		||||
    MomentaField Msum(P.Mom._grid);
 | 
			
		||||
    Msum = zero;
 | 
			
		||||
    for (int a = 0; a < as[level].actions.size(); ++a) {
 | 
			
		||||
      // Compute the force
 | 
			
		||||
      // Compute the force terms for the lagrangian part
 | 
			
		||||
      // We need to compute the derivative of the actions
 | 
			
		||||
      // only once
 | 
			
		||||
      Field force(U._grid);
 | 
			
		||||
@@ -153,7 +159,7 @@ class Integrator {
 | 
			
		||||
      if (as[level].actions.at(a)->is_smeared) Smearer.smeared_force(force);
 | 
			
		||||
      force = FieldImplementation::projectForce(force);  // Ta for gauge fields
 | 
			
		||||
      Real force_abs = std::sqrt(norm2(force) / U._grid->gSites());
 | 
			
		||||
      std::cout << GridLogIntegrator << "Force average: " << force_abs
 | 
			
		||||
      std::cout << GridLogIntegrator << "|Force| site average: " << force_abs
 | 
			
		||||
                << std::endl;
 | 
			
		||||
      Msum += force;
 | 
			
		||||
    }
 | 
			
		||||
@@ -162,21 +168,44 @@ class Integrator {
 | 
			
		||||
    MomentaField OldMom = P.Mom;
 | 
			
		||||
    double threshold = 1e-6;
 | 
			
		||||
    P.M.ImportGauge(U);
 | 
			
		||||
    MomentaField MomDer(P.Mom._grid);
 | 
			
		||||
    MomentaField MomDer1(P.Mom._grid);
 | 
			
		||||
    MomDer1 = zero;
 | 
			
		||||
    MomentaField diff(P.Mom._grid);
 | 
			
		||||
 | 
			
		||||
    // be careful here, we need the first step
 | 
			
		||||
    // in every trajectory
 | 
			
		||||
    static int call = 0;
 | 
			
		||||
    if (call == 1)
 | 
			
		||||
      P.DerivativeU(P.Mom, MomDer1);
 | 
			
		||||
 | 
			
		||||
    call = 1;
 | 
			
		||||
 | 
			
		||||
    // Here run recursively
 | 
			
		||||
    int counter = 1;
 | 
			
		||||
    RealD RelativeError;
 | 
			
		||||
    do {
 | 
			
		||||
      MomentaField MomDer(P.Mom._grid);
 | 
			
		||||
      MomentaField X(P.Mom._grid);
 | 
			
		||||
      OldMom = NewMom;
 | 
			
		||||
      std::cout << GridLogIntegrator << "UpdateP implicit step "<< counter << std::endl;
 | 
			
		||||
 | 
			
		||||
      // Compute the derivative of the kinetic term
 | 
			
		||||
      // with respect to the gauge field
 | 
			
		||||
      P.DerivativeU(NewMom, MomDer);
 | 
			
		||||
      NewMom = P.Mom - ep * (MomDer + Msum);
 | 
			
		||||
      Real force_abs = std::sqrt(norm2(MomDer) / U._grid->gSites());
 | 
			
		||||
      std::cout << GridLogIntegrator << "|Force| laplacian site average: " << force_abs
 | 
			
		||||
                << std::endl;
 | 
			
		||||
 | 
			
		||||
    } while (norm2(NewMom - OldMom) > threshold);
 | 
			
		||||
      NewMom = P.Mom - ep* 0.5 * (2.0*Msum + MomDer + MomDer1);
 | 
			
		||||
      diff = NewMom - OldMom;
 | 
			
		||||
      counter++;
 | 
			
		||||
      RelativeError = std::sqrt(norm2(diff))/std::sqrt(norm2(NewMom));
 | 
			
		||||
      std::cout << GridLogIntegrator << "UpdateP RelativeError: " << RelativeError << std::endl;
 | 
			
		||||
      OldMom = NewMom;
 | 
			
		||||
    } while (RelativeError > threshold);
 | 
			
		||||
 | 
			
		||||
    P.Mom = NewMom;
 | 
			
		||||
 | 
			
		||||
    // update the auxiliary fields momenta
 | 
			
		||||
    // todo
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
@@ -204,19 +233,44 @@ class Integrator {
 | 
			
		||||
    int fl = levels - 1;
 | 
			
		||||
    std::cout << GridLogIntegrator << "   " << "[" << fl << "] U " << " dt " << ep << " : t_U " << t_U << std::endl;
 | 
			
		||||
 | 
			
		||||
    Real threshold = 1e-6;
 | 
			
		||||
    P.M.ImportGauge(U);
 | 
			
		||||
    MomentaField Mom1(P.Mom._grid);
 | 
			
		||||
    MomentaField Mom2(P.Mom._grid);
 | 
			
		||||
    RealD RelativeError;
 | 
			
		||||
    Field diff(U._grid);
 | 
			
		||||
    Real threshold = 1e-6;
 | 
			
		||||
    int counter = 1;
 | 
			
		||||
    int MaxCounter = 1000;
 | 
			
		||||
 | 
			
		||||
    Field OldU = U;
 | 
			
		||||
    Field NewU = U;
 | 
			
		||||
 | 
			
		||||
    P.M.ImportGauge(U);
 | 
			
		||||
    P.DerivativeP(Mom1); // first term in the derivative 
 | 
			
		||||
 | 
			
		||||
    do {
 | 
			
		||||
      OldU = NewU; // some redundancy to be eliminated
 | 
			
		||||
      std::cout << GridLogIntegrator << "UpdateU implicit step "<< counter << std::endl;
 | 
			
		||||
      
 | 
			
		||||
      P.DerivativeP(Mom2); // second term in the derivative, on the updated U
 | 
			
		||||
      FieldImplementation::update_field(Mom1 + Mom2, NewU, ep);
 | 
			
		||||
      MomentaField sum = (Mom1 + Mom2);
 | 
			
		||||
      //std::cout << GridLogMessage << "sum Norm " << norm2(sum) << std::endl;
 | 
			
		||||
 | 
			
		||||
      for (int mu = 0; mu < Nd; mu++) {
 | 
			
		||||
        auto Umu = PeekIndex<LorentzIndex>(U, mu);
 | 
			
		||||
        auto Pmu = PeekIndex<LorentzIndex>(sum, mu);
 | 
			
		||||
        Umu = expMat(Pmu, ep * 0.5, 12) * Umu;
 | 
			
		||||
        PokeIndex<LorentzIndex>(NewU, ProjectOnGroup(Umu), mu);
 | 
			
		||||
      }
 | 
			
		||||
 | 
			
		||||
      diff = NewU - OldU;
 | 
			
		||||
      RelativeError = std::sqrt(norm2(diff))/std::sqrt(norm2(NewU));
 | 
			
		||||
      std::cout << GridLogIntegrator << "UpdateU RelativeError: " << RelativeError << std::endl;
 | 
			
		||||
      
 | 
			
		||||
      P.M.ImportGauge(NewU);
 | 
			
		||||
    } while (norm2(NewU - OldU) > threshold);
 | 
			
		||||
      OldU = NewU; // some redundancy to be eliminated
 | 
			
		||||
      counter++;
 | 
			
		||||
    } while (RelativeError > threshold && counter < MaxCounter);
 | 
			
		||||
 | 
			
		||||
    U = NewU;
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
@@ -225,10 +279,10 @@ class Integrator {
 | 
			
		||||
 public:
 | 
			
		||||
  Integrator(GridBase* grid, IntegratorParameters Par,
 | 
			
		||||
             ActionSet<Field, RepresentationPolicy>& Aset,
 | 
			
		||||
             SmearingPolicy& Sm)
 | 
			
		||||
             SmearingPolicy& Sm, Metric<MomentaField>& M)
 | 
			
		||||
      : Params(Par),
 | 
			
		||||
        as(Aset),
 | 
			
		||||
        P(grid),
 | 
			
		||||
        P(grid, M),
 | 
			
		||||
        levels(Aset.size()),
 | 
			
		||||
        Smearer(Sm),
 | 
			
		||||
        Representations(grid) {
 | 
			
		||||
@@ -260,6 +314,10 @@ class Integrator {
 | 
			
		||||
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  void reverse_momenta(){
 | 
			
		||||
    P.Mom *= 1.0;
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  // to be used by the actionlevel class to iterate
 | 
			
		||||
  // over the representations
 | 
			
		||||
  struct _refresh {
 | 
			
		||||
@@ -278,7 +336,10 @@ class Integrator {
 | 
			
		||||
  void refresh(Field& U, GridParallelRNG& pRNG) {
 | 
			
		||||
    assert(P.Mom._grid == U._grid);
 | 
			
		||||
    std::cout << GridLogIntegrator << "Integrator refresh\n";
 | 
			
		||||
 | 
			
		||||
    //FieldImplementation::generate_momenta(P, pRNG);
 | 
			
		||||
 | 
			
		||||
    P.M.ImportGauge(U);
 | 
			
		||||
    P.MomentaDistribution(pRNG);
 | 
			
		||||
 | 
			
		||||
    // Update the smeared fields, can be implemented as observer
 | 
			
		||||
@@ -325,7 +386,9 @@ class Integrator {
 | 
			
		||||
  // Calculate action
 | 
			
		||||
  RealD S(Field& U) {  // here also U not used
 | 
			
		||||
 | 
			
		||||
    RealD H = - FieldImplementation::FieldSquareNorm(P.Mom); // - trace (P*P)
 | 
			
		||||
    //RealD H = - FieldImplementation::FieldSquareNorm(P.Mom); // - trace (P*P)
 | 
			
		||||
    P.M.ImportGauge(U);
 | 
			
		||||
    RealD H = - P.MomentaAction();
 | 
			
		||||
    RealD Hterm;
 | 
			
		||||
    std::cout << GridLogMessage << "Momentum action H_p = " << H << "\n";
 | 
			
		||||
 | 
			
		||||
@@ -369,6 +432,7 @@ class Integrator {
 | 
			
		||||
 | 
			
		||||
    // and that we indeed got to the end of the trajectory
 | 
			
		||||
    assert(fabs(t_U - Params.trajL) < 1.0e-6);
 | 
			
		||||
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
 
 | 
			
		||||
@@ -294,7 +294,7 @@ class ForceGradient : public Integrator<FieldImplementation, SmearingPolicy,
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
// correct
 | 
			
		||||
template <class FieldImplementation, class SmearingPolicy,
 | 
			
		||||
          class RepresentationPolicy =
 | 
			
		||||
              Representations<FundamentalRepresentation> >
 | 
			
		||||
@@ -311,9 +311,9 @@ class ImplicitLeapFrog : public Integrator<FieldImplementation, SmearingPolicy,
 | 
			
		||||
  std::string integrator_name(){return "ImplicitLeapFrog";}
 | 
			
		||||
 | 
			
		||||
  ImplicitLeapFrog(GridBase* grid, IntegratorParameters Par,
 | 
			
		||||
           ActionSet<Field, RepresentationPolicy>& Aset, SmearingPolicy& Sm)
 | 
			
		||||
           ActionSet<Field, RepresentationPolicy>& Aset, SmearingPolicy& Sm, Metric<Field>& M)
 | 
			
		||||
      : Integrator<FieldImplementation, SmearingPolicy, RepresentationPolicy>(
 | 
			
		||||
            grid, Par, Aset, Sm){};
 | 
			
		||||
            grid, Par, Aset, Sm, M){};
 | 
			
		||||
 | 
			
		||||
  void step(Field& U, int level, int _first, int _last) {
 | 
			
		||||
    int fl = this->as.size() - 1;
 | 
			
		||||
@@ -335,19 +335,89 @@ class ImplicitLeapFrog : public Integrator<FieldImplementation, SmearingPolicy,
 | 
			
		||||
      }
 | 
			
		||||
 | 
			
		||||
      if (level == fl) {  // lowest level
 | 
			
		||||
        this->implicit_update_U(U, eps / 2.0);
 | 
			
		||||
        this->implicit_update_U(U, eps);
 | 
			
		||||
      } else {  // recursive function call
 | 
			
		||||
        this->step(U, level + 1, first_step, last_step);
 | 
			
		||||
      }
 | 
			
		||||
 | 
			
		||||
      int mm = last_step ? 1 : 2;
 | 
			
		||||
      this->update_P(U, level, mm * eps / 2.0);
 | 
			
		||||
      //int mm = last_step ? 1 : 2;
 | 
			
		||||
      if (last_step){
 | 
			
		||||
        this->update_P(U, level, eps / 2.0);
 | 
			
		||||
      } else {
 | 
			
		||||
      this->implicit_update_P(U, level, eps);
 | 
			
		||||
      }
 | 
			
		||||
    }
 | 
			
		||||
  }
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
// This is not completely tested
 | 
			
		||||
template <class FieldImplementation, class SmearingPolicy,
 | 
			
		||||
          class RepresentationPolicy =
 | 
			
		||||
              Representations<FundamentalRepresentation> >
 | 
			
		||||
class ImplicitMinimumNorm2 : public Integrator<FieldImplementation, SmearingPolicy,
 | 
			
		||||
                                       RepresentationPolicy> {
 | 
			
		||||
 private:
 | 
			
		||||
  const RealD lambda = 0.1931833275037836;
 | 
			
		||||
 | 
			
		||||
 public:
 | 
			
		||||
  INHERIT_FIELD_TYPES(FieldImplementation);
 | 
			
		||||
 | 
			
		||||
  ImplicitMinimumNorm2(GridBase* grid, IntegratorParameters Par,
 | 
			
		||||
               ActionSet<Field, RepresentationPolicy>& Aset, SmearingPolicy& Sm, Metric<Field>& M)
 | 
			
		||||
      : Integrator<FieldImplementation, SmearingPolicy, RepresentationPolicy>(
 | 
			
		||||
            grid, Par, Aset, Sm, M){};
 | 
			
		||||
 | 
			
		||||
  std::string integrator_name(){return "ImplicitMininumNorm2";}
 | 
			
		||||
 | 
			
		||||
  void step(Field& U, int level, int _first, int _last) {
 | 
			
		||||
    // level  : current level
 | 
			
		||||
    // fl     : final level
 | 
			
		||||
    // eps    : current step size
 | 
			
		||||
 | 
			
		||||
    int fl = this->as.size() - 1;
 | 
			
		||||
 | 
			
		||||
    RealD eps = this->Params.trajL/this->Params.MDsteps * 2.0;
 | 
			
		||||
    for (int l = 0; l <= level; ++l) eps /= 2.0 * this->as[l].multiplier;
 | 
			
		||||
 | 
			
		||||
    // Nesting:  2xupdate_U of size eps/2
 | 
			
		||||
    // Next level is eps/2/multiplier
 | 
			
		||||
 | 
			
		||||
    int multiplier = this->as[level].multiplier;
 | 
			
		||||
    for (int e = 0; e < multiplier; ++e) {  // steps per step
 | 
			
		||||
 | 
			
		||||
      int first_step = _first && (e == 0);
 | 
			
		||||
      int last_step = _last && (e == multiplier - 1);
 | 
			
		||||
 | 
			
		||||
      if (first_step) {  // initial half step
 | 
			
		||||
        this->implicit_update_P(U, level, lambda * eps);
 | 
			
		||||
      }
 | 
			
		||||
 | 
			
		||||
      if (level == fl) {  // lowest level
 | 
			
		||||
        this->implicit_update_U(U, 0.5 * eps);
 | 
			
		||||
      } else {  // recursive function call
 | 
			
		||||
        this->step(U, level + 1, first_step, 0);
 | 
			
		||||
      }
 | 
			
		||||
 | 
			
		||||
      this->implicit_update_P(U, level, (1.0 - 2.0 * lambda) * eps);
 | 
			
		||||
 | 
			
		||||
      if (level == fl) {  // lowest level
 | 
			
		||||
        this->implicit_update_U(U, 0.5 * eps);
 | 
			
		||||
      } else {  // recursive function call
 | 
			
		||||
        this->step(U, level + 1, 0, last_step);
 | 
			
		||||
      }
 | 
			
		||||
 | 
			
		||||
      //int mm = (last_step) ? 1 : 2;
 | 
			
		||||
      //this->update_P(U, level, lambda * eps * mm);
 | 
			
		||||
 | 
			
		||||
      if (last_step) {
 | 
			
		||||
        this->update_P(U, level, eps * lambda);
 | 
			
		||||
      } else {
 | 
			
		||||
        this->implicit_update_P(U, level, lambda * eps*2.0);
 | 
			
		||||
      }
 | 
			
		||||
    }
 | 
			
		||||
  }
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
 
 | 
			
		||||
@@ -77,7 +77,8 @@ template <class Impl>
 | 
			
		||||
class LaplacianAdjointField: public Metric<typename Impl::Field> {
 | 
			
		||||
  OperatorFunction<typename Impl::Field> &Solver;
 | 
			
		||||
  LaplacianParams param;
 | 
			
		||||
  MultiShiftFunction PowerNegHalf;    
 | 
			
		||||
  MultiShiftFunction PowerHalf;    
 | 
			
		||||
  MultiShiftFunction PowerInvHalf;    
 | 
			
		||||
 | 
			
		||||
 public:
 | 
			
		||||
  INHERIT_GIMPL_TYPES(Impl);
 | 
			
		||||
@@ -87,10 +88,15 @@ class LaplacianAdjointField: public Metric<typename Impl::Field> {
 | 
			
		||||
        AlgRemez remez(param.lo,param.hi,param.precision);
 | 
			
		||||
        std::cout<<GridLogMessage << "Generating degree "<<param.degree<<" for x^(1/2)"<<std::endl;
 | 
			
		||||
        remez.generateApprox(param.degree,1,2);
 | 
			
		||||
        PowerNegHalf.Init(remez,param.tolerance,true);
 | 
			
		||||
        PowerHalf.Init(remez,param.tolerance,false);
 | 
			
		||||
        PowerInvHalf.Init(remez,param.tolerance,true);
 | 
			
		||||
        
 | 
			
		||||
 | 
			
		||||
      };
 | 
			
		||||
 | 
			
		||||
  void Mdir(const GaugeField&, GaugeField&, int, int){ assert(0);}
 | 
			
		||||
  void Mdiag(const GaugeField&, GaugeField&){ assert(0);}
 | 
			
		||||
 | 
			
		||||
  void ImportGauge(const GaugeField& _U) {
 | 
			
		||||
    for (int mu = 0; mu < Nd; mu++) {
 | 
			
		||||
      U[mu] = PeekIndex<LorentzIndex>(_U, mu);
 | 
			
		||||
@@ -98,6 +104,11 @@ class LaplacianAdjointField: public Metric<typename Impl::Field> {
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  void M(const GaugeField& in, GaugeField& out) {
 | 
			
		||||
    // in is an antihermitian matrix
 | 
			
		||||
    // test
 | 
			
		||||
    //GaugeField herm = in + adj(in);
 | 
			
		||||
    //std::cout << "AHermiticity: " << norm2(herm) << std::endl;
 | 
			
		||||
 | 
			
		||||
    GaugeLinkField tmp(in._grid);
 | 
			
		||||
    GaugeLinkField tmp2(in._grid);
 | 
			
		||||
    GaugeLinkField sum(in._grid);
 | 
			
		||||
@@ -116,17 +127,21 @@ class LaplacianAdjointField: public Metric<typename Impl::Field> {
 | 
			
		||||
    }
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  void MDeriv(const GaugeField& in, GaugeField& der, bool dag) {
 | 
			
		||||
  void MDeriv(const GaugeField& in, GaugeField& der) {
 | 
			
		||||
    // in is anti-hermitian
 | 
			
		||||
    RealD factor = -kappa / (double(4 * Nd));
 | 
			
		||||
    for (int mu = 0; mu < Nd; mu++) {
 | 
			
		||||
      GaugeLinkField in_mu = PeekIndex<LorentzIndex>(in, mu);
 | 
			
		||||
    
 | 
			
		||||
    for (int mu = 0; mu < Nd; mu++){
 | 
			
		||||
      GaugeLinkField der_mu(der._grid);
 | 
			
		||||
      if (!dag)
 | 
			
		||||
        der_mu =
 | 
			
		||||
            factor * Cshift(in_mu, mu, +1) * adj(U[mu]) + adj(U[mu]) * in_mu;
 | 
			
		||||
      else
 | 
			
		||||
        der_mu = factor * U[mu] * Cshift(in_mu, mu, +1) + in_mu * U[mu];
 | 
			
		||||
    }
 | 
			
		||||
      der_mu = zero;
 | 
			
		||||
      for (int nu = 0; nu < Nd; nu++){
 | 
			
		||||
        GaugeLinkField in_nu = PeekIndex<LorentzIndex>(in, nu);
 | 
			
		||||
        der_mu += U[mu] * Cshift(in_nu, mu, 1) * adj(U[mu]) * in_nu;
 | 
			
		||||
      }
 | 
			
		||||
      // the minus sign comes by using the in_nu instead of the
 | 
			
		||||
      // adjoint in the last multiplication
 | 
			
		||||
      PokeIndex<LorentzIndex>(der,  -2.0 * factor * der_mu, mu);
 | 
			
		||||
    } 
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  void Minv(const GaugeField& in, GaugeField& inverted){
 | 
			
		||||
@@ -134,14 +149,12 @@ class LaplacianAdjointField: public Metric<typename Impl::Field> {
 | 
			
		||||
    Solver(HermOp, in, inverted);
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  void MInvSquareRoot(GaugeField& P){
 | 
			
		||||
    // Takes a gaussian gauge field and multiplies by the metric
 | 
			
		||||
    // need the rational approximation for the square root 
 | 
			
		||||
  void MSquareRoot(GaugeField& P){
 | 
			
		||||
    GaugeField Gp(P._grid);
 | 
			
		||||
    HermitianLinearOperator<LaplacianAdjointField<Impl>,GaugeField> HermOp(*this);
 | 
			
		||||
    ConjugateGradientMultiShift<GaugeField> msCG(param.MaxIter,PowerNegHalf);
 | 
			
		||||
    ConjugateGradientMultiShift<GaugeField> msCG(param.MaxIter,PowerHalf);
 | 
			
		||||
    msCG(HermOp,P,Gp);
 | 
			
		||||
    P = Gp; // now P has the correct distribution
 | 
			
		||||
    P = Gp; 
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
 
 | 
			
		||||
@@ -39,7 +39,8 @@ public:
 | 
			
		||||
  virtual void ImportGauge(const Field&)   = 0;
 | 
			
		||||
  virtual void M(const Field&, Field&)     = 0;
 | 
			
		||||
  virtual void Minv(const Field&, Field&)  = 0;
 | 
			
		||||
  virtual void MInvSquareRoot(Field&) = 0;
 | 
			
		||||
  virtual void MSquareRoot(Field&) = 0;
 | 
			
		||||
  virtual void MDeriv(const Field&, Field&) = 0;
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
@@ -54,9 +55,12 @@ public:
 | 
			
		||||
  virtual void Minv(const Field& in, Field& out){
 | 
			
		||||
    out = in;
 | 
			
		||||
  }
 | 
			
		||||
  virtual void MInvSquareRoot(Field& P){
 | 
			
		||||
  virtual void MSquareRoot(Field& P){
 | 
			
		||||
    // do nothing
 | 
			
		||||
  }
 | 
			
		||||
  virtual void MDeriv(const Field& in, Field& out){
 | 
			
		||||
    out = zero;
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
@@ -64,17 +68,17 @@ public:
 | 
			
		||||
// Generalised momenta
 | 
			
		||||
///////////////////////////////
 | 
			
		||||
 | 
			
		||||
template <typename Implementation, typename Metric>
 | 
			
		||||
template <typename Implementation>
 | 
			
		||||
class GeneralisedMomenta{
 | 
			
		||||
public:
 | 
			
		||||
  typedef typename Implementation::Field MomentaField;  //for readability
 | 
			
		||||
 | 
			
		||||
  Metric M;
 | 
			
		||||
  typedef typename Implementation::GaugeLinkField MomentaLinkField;  //for readability
 | 
			
		||||
  Metric<MomentaField>& M;
 | 
			
		||||
  MomentaField Mom;
 | 
			
		||||
 | 
			
		||||
  GeneralisedMomenta(GridBase* grid): Mom(grid){}
 | 
			
		||||
 | 
			
		||||
  GeneralisedMomenta(GridBase* grid, Metric<MomentaField>& M): M(M), Mom(grid){}
 | 
			
		||||
 | 
			
		||||
  // Correct
 | 
			
		||||
  void MomentaDistribution(GridParallelRNG& pRNG){
 | 
			
		||||
    // Generate a distribution for
 | 
			
		||||
    // 1/2 P^dag G P
 | 
			
		||||
@@ -83,26 +87,45 @@ public:
 | 
			
		||||
    // Generate gaussian momenta
 | 
			
		||||
    Implementation::generate_momenta(Mom, pRNG);
 | 
			
		||||
    // Modify the distribution with the metric
 | 
			
		||||
    M.MInvSquareRoot(Mom);
 | 
			
		||||
    M.MSquareRoot(Mom);
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  void Derivative(MomentaField& in, MomentaField& der){
 | 
			
		||||
  // Correct
 | 
			
		||||
  RealD MomentaAction(){
 | 
			
		||||
    MomentaField inv(Mom._grid);
 | 
			
		||||
    inv = zero;
 | 
			
		||||
    M.Minv(Mom, inv);
 | 
			
		||||
    LatticeComplex Hloc(Mom._grid);
 | 
			
		||||
    Hloc = zero;
 | 
			
		||||
    for (int mu = 0; mu < Nd; mu++) {
 | 
			
		||||
      // This is not very general
 | 
			
		||||
      // hide in the operators
 | 
			
		||||
      auto Mom_mu = PeekIndex<LorentzIndex>(Mom, mu);
 | 
			
		||||
      auto inv_mu = PeekIndex<LorentzIndex>(inv, mu);
 | 
			
		||||
      Hloc += trace(Mom_mu * inv_mu);
 | 
			
		||||
    }
 | 
			
		||||
    Complex Hsum = sum(Hloc);
 | 
			
		||||
    return Hsum.real();
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  // Correct
 | 
			
		||||
  void DerivativeU(MomentaField& in, MomentaField& der){
 | 
			
		||||
    // Compute the derivative of the kinetic term
 | 
			
		||||
    // with respect to the gauge field
 | 
			
		||||
    MomentaField MomDer(in._grid);
 | 
			
		||||
    MomentaField MDer(in._grid);
 | 
			
		||||
    MomentaField X(in._grid);
 | 
			
		||||
 | 
			
		||||
    M.Minv(in, X); // X = G in
 | 
			
		||||
    M.MDeriv(X, MomDer, DaggerNo); // MomDer = dM/dU X
 | 
			
		||||
    // MomDer is just the derivative
 | 
			
		||||
    MomDer = adj(X)* MomDer;
 | 
			
		||||
    // Traceless Antihermitian
 | 
			
		||||
    // assuming we are in the algebra
 | 
			
		||||
    der = Implementation::projectForce(MomDer);
 | 
			
		||||
    X = zero;
 | 
			
		||||
    M.Minv(in, X);  // X = G in
 | 
			
		||||
    M.MDeriv(X, MDer);  // MDer = U * dS/dU
 | 
			
		||||
    der = Implementation::projectForce(MDer);  // Ta if gauge fields
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
  //   
 | 
			
		||||
  void DerivativeP(MomentaField& der){
 | 
			
		||||
    der = zero;
 | 
			
		||||
    M.Minv(Mom, der);
 | 
			
		||||
    der = Implementation::projectForce(der); 
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
};
 | 
			
		||||
 
 | 
			
		||||
							
								
								
									
										178
									
								
								tests/forces/Test_laplacian_force.cc
									
									
									
									
									
										Normal file
									
								
							
							
						
						
									
										178
									
								
								tests/forces/Test_laplacian_force.cc
									
									
									
									
									
										Normal file
									
								
							@@ -0,0 +1,178 @@
 | 
			
		||||
    /*************************************************************************************
 | 
			
		||||
 | 
			
		||||
    Grid physics library, www.github.com/paboyle/Grid 
 | 
			
		||||
 | 
			
		||||
    Source file: ./tests/Test_rect_force.cc
 | 
			
		||||
 | 
			
		||||
    Copyright (C) 2015
 | 
			
		||||
 | 
			
		||||
Author: Azusa Yamaguchi <ayamaguc@staffmail.ed.ac.uk>
 | 
			
		||||
 | 
			
		||||
    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 std;
 | 
			
		||||
using namespace Grid;
 | 
			
		||||
using namespace Grid::QCD;
 | 
			
		||||
 | 
			
		||||
#define parallel_for PARALLEL_FOR_LOOP for
 | 
			
		||||
 | 
			
		||||
int main (int argc, char ** argv)
 | 
			
		||||
{
 | 
			
		||||
  Grid_init(&argc,&argv);
 | 
			
		||||
 | 
			
		||||
  std::vector<int> latt_size   = GridDefaultLatt();
 | 
			
		||||
  std::vector<int> simd_layout = GridDefaultSimd(Nd,vComplex::Nsimd());
 | 
			
		||||
  std::vector<int> mpi_layout  = GridDefaultMpi();
 | 
			
		||||
 | 
			
		||||
  GridCartesian               Grid(latt_size,simd_layout,mpi_layout);
 | 
			
		||||
  GridRedBlackCartesian     RBGrid(latt_size,simd_layout,mpi_layout);
 | 
			
		||||
 | 
			
		||||
  int threads = GridThread::GetThreads();
 | 
			
		||||
  std::cout<<GridLogMessage << "Grid is setup to use "<<threads<<" threads"<<std::endl;
 | 
			
		||||
 | 
			
		||||
  std::vector<int> seeds({1,2,3,4});
 | 
			
		||||
 | 
			
		||||
  GridParallelRNG          pRNG(&Grid);
 | 
			
		||||
  pRNG.SeedRandomDevice();
 | 
			
		||||
 | 
			
		||||
  LatticeGaugeField U(&Grid);
 | 
			
		||||
  LatticeGaugeField P(&Grid);
 | 
			
		||||
  LatticeColourMatrix P_mu(&Grid);
 | 
			
		||||
  // Matrix in the algebra
 | 
			
		||||
  for (int mu = 0; mu < Nd; mu++) {
 | 
			
		||||
    SU<Nc>::GaussianFundamentalLieAlgebraMatrix(pRNG, P_mu);
 | 
			
		||||
    PokeIndex<LorentzIndex>(P, P_mu, mu);
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  SU3::HotConfiguration(pRNG,U);
 | 
			
		||||
  
 | 
			
		||||
 | 
			
		||||
  ConjugateGradient<LatticeGaugeField> CG(1.0e-8, 10000);
 | 
			
		||||
  LaplacianParams LapPar(0.001, 1.0, 1000, 1e-8, 10, 64);
 | 
			
		||||
  RealD Kappa = 0.99;
 | 
			
		||||
  LaplacianAdjointField<PeriodicGimplR> Laplacian(&Grid, CG, LapPar, Kappa);
 | 
			
		||||
  GeneralisedMomenta<PeriodicGimplR> LaplacianMomenta(&Grid, Laplacian);
 | 
			
		||||
  LaplacianMomenta.M.ImportGauge(U);
 | 
			
		||||
  LaplacianMomenta.MomentaDistribution(pRNG);// fills the Momenta with the correct distr
 | 
			
		||||
  
 | 
			
		||||
 | 
			
		||||
  std::cout << std::setprecision(15);
 | 
			
		||||
  std::cout << GridLogMessage << "MomentaAction" << std::endl;
 | 
			
		||||
  ComplexD S    = LaplacianMomenta.MomentaAction();
 | 
			
		||||
 | 
			
		||||
  // get the deriv with respect to "U"
 | 
			
		||||
  LatticeGaugeField UdSdU(&Grid);
 | 
			
		||||
  std::cout << GridLogMessage<< "DerivativeU" << std::endl;
 | 
			
		||||
  LaplacianMomenta.DerivativeU(LaplacianMomenta.Mom, UdSdU);
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
  ////////////////////////////////////
 | 
			
		||||
  // Modify the gauge field a little 
 | 
			
		||||
  ////////////////////////////////////
 | 
			
		||||
  RealD dt = 0.001;
 | 
			
		||||
 | 
			
		||||
  LatticeColourMatrix mommu(&Grid); 
 | 
			
		||||
  LatticeColourMatrix forcemu(&Grid); 
 | 
			
		||||
  LatticeGaugeField mom(&Grid); 
 | 
			
		||||
  LatticeGaugeField Uprime(&Grid); 
 | 
			
		||||
 | 
			
		||||
  std::cout << GridLogMessage << "Update the U " << std::endl;
 | 
			
		||||
  for(int mu=0;mu<Nd;mu++){
 | 
			
		||||
  // Traceless antihermitian momentum; gaussian in lie algebra
 | 
			
		||||
    SU3::GaussianFundamentalLieAlgebraMatrix(pRNG, mommu); 
 | 
			
		||||
    auto Umu = PeekIndex<LorentzIndex>(U, mu);
 | 
			
		||||
    PokeIndex<LorentzIndex>(mom,mommu,mu);
 | 
			
		||||
    Umu = expMat(mommu, dt, 12) * Umu;
 | 
			
		||||
    PokeIndex<LorentzIndex>(Uprime, ProjectOnGroup(Umu), mu);
 | 
			
		||||
  
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  std::cout << GridLogMessage << "New action " << std::endl;
 | 
			
		||||
  LaplacianMomenta.M.ImportGauge(Uprime);
 | 
			
		||||
  ComplexD Sprime    = LaplacianMomenta.MomentaAction();
 | 
			
		||||
 | 
			
		||||
  //////////////////////////////////////////////
 | 
			
		||||
  // Use derivative to estimate dS
 | 
			
		||||
  //////////////////////////////////////////////
 | 
			
		||||
 | 
			
		||||
  LatticeComplex dS(&Grid); dS = zero;
 | 
			
		||||
 | 
			
		||||
  for(int mu=0;mu<Nd;mu++){
 | 
			
		||||
    auto UdSdUmu = PeekIndex<LorentzIndex>(UdSdU,mu);
 | 
			
		||||
         mommu   = PeekIndex<LorentzIndex>(mom,mu);
 | 
			
		||||
 | 
			
		||||
    // Update gauge action density
 | 
			
		||||
    // U = exp(p dt) U
 | 
			
		||||
    // dU/dt = p U
 | 
			
		||||
    // so dSdt = trace( dUdt dSdU) = trace( p UdSdUmu ) 
 | 
			
		||||
 | 
			
		||||
    dS = dS + trace(mommu*UdSdUmu)*dt*2.0;
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  Complex dSpred    = sum(dS);
 | 
			
		||||
 | 
			
		||||
  std::cout << GridLogMessage << " S      "<<S<<std::endl;
 | 
			
		||||
  std::cout << GridLogMessage << " Sprime "<<Sprime<<std::endl;
 | 
			
		||||
  std::cout << GridLogMessage << "dS      "<<Sprime-S<<std::endl;
 | 
			
		||||
  std::cout << GridLogMessage << "pred dS "<< dSpred <<std::endl;
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
  // P derivative
 | 
			
		||||
  // Increment p 
 | 
			
		||||
  dt = 0.0001;
 | 
			
		||||
  LaplacianMomenta.M.ImportGauge(U);
 | 
			
		||||
  LatticeGaugeField UdSdP(&Grid);
 | 
			
		||||
  LaplacianMomenta.DerivativeP(UdSdP);
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
  LaplacianMomenta.Mom += dt*P;
 | 
			
		||||
  
 | 
			
		||||
  Sprime    = LaplacianMomenta.MomentaAction();
 | 
			
		||||
 | 
			
		||||
  // Prediciton
 | 
			
		||||
 | 
			
		||||
  dS = zero;
 | 
			
		||||
   for(int mu=0;mu<Nd;mu++){
 | 
			
		||||
    auto dSdPmu = PeekIndex<LorentzIndex>(UdSdP,mu);
 | 
			
		||||
    auto Pmu = PeekIndex<LorentzIndex>(P,mu);
 | 
			
		||||
    // Update gauge action density
 | 
			
		||||
    // 
 | 
			
		||||
    // dMom/dt = P
 | 
			
		||||
    // so dSdt = trace( dPdt dSdP) = trace( P dSdP ) 
 | 
			
		||||
 | 
			
		||||
    dS = dS + trace(Pmu*dSdPmu)*dt*2.0;
 | 
			
		||||
  } 
 | 
			
		||||
  dSpred    = sum(dS);
 | 
			
		||||
 | 
			
		||||
  std::cout << GridLogMessage << " S      "<<S<<std::endl;
 | 
			
		||||
  std::cout << GridLogMessage << " Sprime "<<Sprime<<std::endl;
 | 
			
		||||
  std::cout << GridLogMessage << "dS      "<<Sprime-S<<std::endl;
 | 
			
		||||
  std::cout << GridLogMessage << "pred dS "<< dSpred <<std::endl; 
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
  std::cout<< GridLogMessage << "Done" <<std::endl;
 | 
			
		||||
  Grid_finalize();
 | 
			
		||||
}
 | 
			
		||||
							
								
								
									
										97
									
								
								tests/hmc/Test_hmc_WilsonGauge_Implicit.cc
									
									
									
									
									
										Normal file
									
								
							
							
						
						
									
										97
									
								
								tests/hmc/Test_hmc_WilsonGauge_Implicit.cc
									
									
									
									
									
										Normal file
									
								
							@@ -0,0 +1,97 @@
 | 
			
		||||
/*************************************************************************************
 | 
			
		||||
 | 
			
		||||
Grid physics library, www.github.com/paboyle/Grid
 | 
			
		||||
 | 
			
		||||
Source file: ./tests/Test_hmc_WilsonFermionGauge.cc
 | 
			
		||||
 | 
			
		||||
Copyright (C) 2015
 | 
			
		||||
 | 
			
		||||
Author: Peter Boyle <pabobyle@ph.ed.ac.uk>
 | 
			
		||||
Author: neo <cossu@post.kek.jp>
 | 
			
		||||
Author: Guido Cossu <guido.cossu@ed.ac.uk>
 | 
			
		||||
 | 
			
		||||
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>
 | 
			
		||||
 | 
			
		||||
int main(int argc, char **argv) {
 | 
			
		||||
  using namespace Grid;
 | 
			
		||||
  using namespace Grid::QCD;
 | 
			
		||||
 | 
			
		||||
  Grid_init(&argc, &argv);
 | 
			
		||||
  int threads = GridThread::GetThreads();
 | 
			
		||||
  // here make a routine to print all the relevant information on the run
 | 
			
		||||
  std::cout << GridLogMessage << "Grid is setup to use " << threads << " threads" << std::endl;
 | 
			
		||||
 | 
			
		||||
   // Typedefs to simplify notation
 | 
			
		||||
  typedef GenericHMCRunner<ImplicitLeapFrog> HMCWrapper;  // Uses the default minimum norm
 | 
			
		||||
  HMCWrapper TheHMC;
 | 
			
		||||
 | 
			
		||||
  // Grid from the command line
 | 
			
		||||
  TheHMC.Resources.AddFourDimGrid("gauge");
 | 
			
		||||
  // Possibile to create the module by hand 
 | 
			
		||||
  // hardcoding parameters or using a Reader
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
  // Checkpointer definition
 | 
			
		||||
  CheckpointerParameters CPparams;  
 | 
			
		||||
  CPparams.config_prefix = "ckpoint_lat";
 | 
			
		||||
  CPparams.rng_prefix = "ckpoint_rng";
 | 
			
		||||
  CPparams.saveInterval = 20;
 | 
			
		||||
  CPparams.format = "IEEE64BIG";
 | 
			
		||||
  
 | 
			
		||||
  TheHMC.Resources.LoadBinaryCheckpointer(CPparams);
 | 
			
		||||
 | 
			
		||||
  RNGModuleParameters RNGpar;
 | 
			
		||||
  RNGpar.serial_seeds = "1 2 3 4 5";
 | 
			
		||||
  RNGpar.parallel_seeds = "6 7 8 9 10 12";
 | 
			
		||||
  TheHMC.Resources.SetRNGSeeds(RNGpar);
 | 
			
		||||
 | 
			
		||||
  // Construct observables
 | 
			
		||||
  // here there is too much indirection 
 | 
			
		||||
  PlaquetteObsParameters PlPar;
 | 
			
		||||
  PlPar.output_prefix = "Plaquette";
 | 
			
		||||
  PlaquetteMod<HMCWrapper::ImplPolicy> PlaqModule(PlPar);
 | 
			
		||||
  TheHMC.Resources.AddObservable(&PlaqModule);
 | 
			
		||||
  //////////////////////////////////////////////
 | 
			
		||||
 | 
			
		||||
  /////////////////////////////////////////////////////////////
 | 
			
		||||
  // Collect actions, here use more encapsulation
 | 
			
		||||
  // need wrappers of the fermionic classes 
 | 
			
		||||
  // that have a complex construction
 | 
			
		||||
  // standard
 | 
			
		||||
  RealD beta = 5.6;
 | 
			
		||||
  WilsonGaugeActionR Waction(beta);
 | 
			
		||||
  
 | 
			
		||||
  ActionLevel<HMCWrapper::Field> Level1(1);
 | 
			
		||||
  Level1.push_back(&Waction);
 | 
			
		||||
  //Level1.push_back(WGMod.getPtr());
 | 
			
		||||
  TheHMC.TheAction.push_back(Level1);
 | 
			
		||||
  /////////////////////////////////////////////////////////////
 | 
			
		||||
 | 
			
		||||
  // HMC parameters are serialisable 
 | 
			
		||||
  TheHMC.Parameters.MD.MDsteps = 60;
 | 
			
		||||
  TheHMC.Parameters.MD.trajL   = 1.0;
 | 
			
		||||
 | 
			
		||||
  TheHMC.ReadCommandLine(argc, argv); // these can be parameters from file
 | 
			
		||||
  TheHMC.Run();  // no smearing
 | 
			
		||||
 | 
			
		||||
  Grid_finalize();
 | 
			
		||||
 | 
			
		||||
} // main
 | 
			
		||||
@@ -70,13 +70,13 @@ int main (int argc, char ** argv)
 | 
			
		||||
  LatticeColourMatrix src_mu(&Grid);
 | 
			
		||||
  for (int mu = 0; mu < Nd; mu++) {
 | 
			
		||||
    SU<Nc>::GaussianFundamentalLieAlgebraMatrix(pRNG, src_mu);
 | 
			
		||||
    PokeIndex<LorentzIndex>(src_f, src_mu, mu);
 | 
			
		||||
    PokeIndex<LorentzIndex>(src_f, timesI(src_mu), mu);
 | 
			
		||||
  }
 | 
			
		||||
  LatticeGaugeField result_f(&Grid);
 | 
			
		||||
 | 
			
		||||
  // Definition of the Laplacian operator
 | 
			
		||||
  ConjugateGradient<LatticeGaugeField> CG(1.0e-8,10000);
 | 
			
		||||
  LaplacianParams LapPar(0.001, 1.0, 1000, 1e-8, 10, 64);
 | 
			
		||||
  LaplacianParams LapPar(0.00001, 1.0, 1000, 1e-8, 10, 64);
 | 
			
		||||
  LaplacianAdjointField<PeriodicGimplR> Laplacian(&Grid, CG, LapPar, Kappa);
 | 
			
		||||
  Laplacian.ImportGauge(Umu);
 | 
			
		||||
  std::cout << GridLogMessage << "Testing the Laplacian using the full matrix" <<std::endl;
 | 
			
		||||
@@ -84,7 +84,7 @@ int main (int argc, char ** argv)
 | 
			
		||||
  
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
  Laplacian.MomentaDistribution(src_f);
 | 
			
		||||
  Laplacian.MSquareRoot(src_f);
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
  // Tests also the version using the algebra decomposition
 | 
			
		||||
 
 | 
			
		||||
		Reference in New Issue
	
	Block a user