mirror of
				https://github.com/paboyle/Grid.git
				synced 2025-11-03 21:44:33 +00:00 
			
		
		
		
	Auxiliary fields
This commit is contained in:
		@@ -144,7 +144,7 @@ class HMCWrapperTemplate: public HMCRunnerBase<ReaderClass> {
 | 
			
		||||
    LaplacianParams LapPar(0.0001, 1.0, 1000, 1e-8, 12, 64);
 | 
			
		||||
    RealD Kappa = 0.9;
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
    // Better to pass the generalised momenta to the integrator
 | 
			
		||||
    LaplacianAdjointField<PeriodicGimplR> Laplacian(UGrid, CG, LapPar, Kappa);
 | 
			
		||||
    TheIntegrator MDynamics(UGrid, Parameters.MD, TheAction, Smearing, Laplacian);
 | 
			
		||||
 | 
			
		||||
 
 | 
			
		||||
@@ -128,17 +128,23 @@ class Integrator {
 | 
			
		||||
      Mom -= force * ep; 
 | 
			
		||||
    }
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
    // Generalised momenta
 | 
			
		||||
    MomentaField MomDer(P.Mom._grid);
 | 
			
		||||
    P.M.ImportGauge(U);
 | 
			
		||||
    P.DerivativeU(P.Mom, MomDer);
 | 
			
		||||
    Mom -= MomDer * ep;
 | 
			
		||||
 | 
			
		||||
    // Auxiliary fields
 | 
			
		||||
    //P.update_auxiliary_momenta(ep*0.5);
 | 
			
		||||
    //P.AuxiliaryFieldsDerivative(MomDer);
 | 
			
		||||
    //Mom -= MomDer * ep;
 | 
			
		||||
    //P.update_auxiliary_momenta(ep*0.5);
 | 
			
		||||
 | 
			
		||||
    // Force from the other representations
 | 
			
		||||
    as[level].apply(update_P_hireps, Representations, Mom, U, ep);
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  void implicit_update_P(Field& U, int level, double ep) {
 | 
			
		||||
  void implicit_update_P(Field& U, int level, double ep, bool intermediate = false) {
 | 
			
		||||
    t_P[level] += ep;
 | 
			
		||||
 | 
			
		||||
    std::cout << GridLogIntegrator << "[" << level << "] P "
 | 
			
		||||
@@ -170,16 +176,18 @@ class Integrator {
 | 
			
		||||
    P.M.ImportGauge(U);
 | 
			
		||||
    MomentaField MomDer(P.Mom._grid);
 | 
			
		||||
    MomentaField MomDer1(P.Mom._grid);
 | 
			
		||||
    MomentaField AuxDer(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)
 | 
			
		||||
    if (intermediate)
 | 
			
		||||
      P.DerivativeU(P.Mom, MomDer1);
 | 
			
		||||
 | 
			
		||||
    call = 1;
 | 
			
		||||
    // Auxiliary fields
 | 
			
		||||
    //P.update_auxiliary_momenta(ep*0.5);
 | 
			
		||||
    //P.AuxiliaryFieldsDerivative(AuxDer);
 | 
			
		||||
    //Msum += AuxDer;
 | 
			
		||||
    
 | 
			
		||||
 | 
			
		||||
    // Here run recursively
 | 
			
		||||
    int counter = 1;
 | 
			
		||||
@@ -194,7 +202,7 @@ class Integrator {
 | 
			
		||||
      std::cout << GridLogIntegrator << "|Force| laplacian site average: " << force_abs
 | 
			
		||||
                << std::endl;
 | 
			
		||||
 | 
			
		||||
      NewMom = P.Mom - ep* 0.5 * (2.0*Msum + MomDer + MomDer1);
 | 
			
		||||
      NewMom = P.Mom - ep* 0.5 * (2.0*Msum + MomDer + MomDer1);// simplify
 | 
			
		||||
      diff = NewMom - OldMom;
 | 
			
		||||
      counter++;
 | 
			
		||||
      RelativeError = std::sqrt(norm2(diff))/std::sqrt(norm2(NewMom));
 | 
			
		||||
@@ -205,7 +213,7 @@ class Integrator {
 | 
			
		||||
    P.Mom = NewMom;
 | 
			
		||||
 | 
			
		||||
    // update the auxiliary fields momenta    
 | 
			
		||||
    // todo
 | 
			
		||||
    //P.update_auxiliary_momenta(ep*0.5);
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
@@ -239,7 +247,7 @@ class Integrator {
 | 
			
		||||
    Field diff(U._grid);
 | 
			
		||||
    Real threshold = 1e-6;
 | 
			
		||||
    int counter = 1;
 | 
			
		||||
    int MaxCounter = 1000;
 | 
			
		||||
    int MaxCounter = 100;
 | 
			
		||||
 | 
			
		||||
    Field OldU = U;
 | 
			
		||||
    Field NewU = U;
 | 
			
		||||
@@ -247,6 +255,10 @@ class Integrator {
 | 
			
		||||
    P.M.ImportGauge(U);
 | 
			
		||||
    P.DerivativeP(Mom1); // first term in the derivative 
 | 
			
		||||
 | 
			
		||||
 
 | 
			
		||||
    //P.update_auxiliary_fields(ep*0.5);
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
    do {
 | 
			
		||||
      std::cout << GridLogIntegrator << "UpdateU implicit step "<< counter << std::endl;
 | 
			
		||||
      
 | 
			
		||||
@@ -271,6 +283,9 @@ class Integrator {
 | 
			
		||||
    } while (RelativeError > threshold && counter < MaxCounter);
 | 
			
		||||
 | 
			
		||||
    U = NewU;
 | 
			
		||||
 | 
			
		||||
    //P.update_auxiliary_fields(ep*0.5);
 | 
			
		||||
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
 
 | 
			
		||||
@@ -344,7 +344,7 @@ class ImplicitLeapFrog : public Integrator<FieldImplementation, SmearingPolicy,
 | 
			
		||||
      if (last_step){
 | 
			
		||||
        this->update_P(U, level, eps / 2.0);
 | 
			
		||||
      } else {
 | 
			
		||||
      this->implicit_update_P(U, level, eps);
 | 
			
		||||
      this->implicit_update_P(U, level, eps, true);// intermediate step
 | 
			
		||||
      }
 | 
			
		||||
    }
 | 
			
		||||
  }
 | 
			
		||||
 
 | 
			
		||||
@@ -144,6 +144,25 @@ class LaplacianAdjointField: public Metric<typename Impl::Field> {
 | 
			
		||||
    } 
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  // separating this temporarily
 | 
			
		||||
  void MDeriv(const GaugeField& left, const GaugeField& right,
 | 
			
		||||
              GaugeField& der) {
 | 
			
		||||
    // in is anti-hermitian
 | 
			
		||||
    RealD factor = -kappa / (double(4 * Nd));
 | 
			
		||||
 | 
			
		||||
    for (int mu = 0; mu < Nd; mu++) {
 | 
			
		||||
      GaugeLinkField der_mu(der._grid);
 | 
			
		||||
      der_mu = zero;
 | 
			
		||||
      for (int nu = 0; nu < Nd; nu++) {
 | 
			
		||||
        GaugeLinkField left_nu = PeekIndex<LorentzIndex>(left, nu);
 | 
			
		||||
        GaugeLinkField right_nu = PeekIndex<LorentzIndex>(right, nu);
 | 
			
		||||
        der_mu += U[mu] * Cshift(left_nu, mu, 1) * adj(U[mu]) * right_nu;
 | 
			
		||||
        der_mu += U[mu] * Cshift(right_nu, mu, 1) * adj(U[mu]) * left_nu;
 | 
			
		||||
      }
 | 
			
		||||
      PokeIndex<LorentzIndex>(der, -factor * der_mu, mu);
 | 
			
		||||
    }
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  void Minv(const GaugeField& in, GaugeField& inverted){
 | 
			
		||||
    HermitianLinearOperator<LaplacianAdjointField<Impl>,GaugeField> HermOp(*this);
 | 
			
		||||
    Solver(HermOp, in, inverted);
 | 
			
		||||
@@ -157,6 +176,14 @@ class LaplacianAdjointField: public Metric<typename Impl::Field> {
 | 
			
		||||
    P = Gp; 
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  void MInvSquareRoot(GaugeField& P){
 | 
			
		||||
    GaugeField Gp(P._grid);
 | 
			
		||||
    HermitianLinearOperator<LaplacianAdjointField<Impl>,GaugeField> HermOp(*this);
 | 
			
		||||
    ConjugateGradientMultiShift<GaugeField> msCG(param.MaxIter,PowerInvHalf);
 | 
			
		||||
    msCG(HermOp,P,Gp);
 | 
			
		||||
    P = Gp; 
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
 private:
 | 
			
		||||
 
 | 
			
		||||
@@ -40,7 +40,9 @@ public:
 | 
			
		||||
  virtual void M(const Field&, Field&)     = 0;
 | 
			
		||||
  virtual void Minv(const Field&, Field&)  = 0;
 | 
			
		||||
  virtual void MSquareRoot(Field&) = 0;
 | 
			
		||||
  virtual void MInvSquareRoot(Field&) = 0;
 | 
			
		||||
  virtual void MDeriv(const Field&, Field&) = 0;
 | 
			
		||||
  virtual void MDeriv(const Field&, const Field&, Field&) = 0;
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
@@ -58,9 +60,15 @@ public:
 | 
			
		||||
  virtual void MSquareRoot(Field& P){
 | 
			
		||||
    // do nothing
 | 
			
		||||
  }
 | 
			
		||||
  virtual void MInvSquareRoot(Field& P){
 | 
			
		||||
    // do nothing
 | 
			
		||||
  }
 | 
			
		||||
  virtual void MDeriv(const Field& in, Field& out){
 | 
			
		||||
    out = zero;
 | 
			
		||||
  }
 | 
			
		||||
  virtual void MDeriv(const Field& left, const Field& right, Field& out){
 | 
			
		||||
    out = zero;
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
@@ -76,18 +84,37 @@ public:
 | 
			
		||||
  Metric<MomentaField>& M;
 | 
			
		||||
  MomentaField Mom;
 | 
			
		||||
 | 
			
		||||
  GeneralisedMomenta(GridBase* grid, Metric<MomentaField>& M): M(M), Mom(grid){}
 | 
			
		||||
  // Auxiliary fields
 | 
			
		||||
  // not hard coded but inherit the type from the metric
 | 
			
		||||
  // created Nd new fields
 | 
			
		||||
  // hide these in the metric?
 | 
			
		||||
  //typedef Lattice<iVector<iScalar<iMatrix<vComplex, Nc> >, Nd/2 > > AuxiliaryMomentaType;
 | 
			
		||||
  MomentaField AuxMom;
 | 
			
		||||
  MomentaField AuxField;
 | 
			
		||||
 | 
			
		||||
  GeneralisedMomenta(GridBase* grid, Metric<MomentaField>& M): M(M), Mom(grid), AuxMom(grid), AuxField(grid){}
 | 
			
		||||
 | 
			
		||||
  // Correct
 | 
			
		||||
  void MomentaDistribution(GridParallelRNG& pRNG){
 | 
			
		||||
    // Generate a distribution for
 | 
			
		||||
    // 1/2 P^dag G P
 | 
			
		||||
    // P^dag G P
 | 
			
		||||
    // where G = M^-1
 | 
			
		||||
 | 
			
		||||
    // Generate gaussian momenta
 | 
			
		||||
    Implementation::generate_momenta(Mom, pRNG);
 | 
			
		||||
    // Modify the distribution with the metric
 | 
			
		||||
    M.MSquareRoot(Mom);
 | 
			
		||||
 | 
			
		||||
    if (0) {
 | 
			
		||||
      // Auxiliary momenta
 | 
			
		||||
      // do nothing if trivial, so hide in the metric
 | 
			
		||||
      MomentaField AuxMomTemp(Mom._grid);
 | 
			
		||||
      Implementation::generate_momenta(AuxMom, pRNG);
 | 
			
		||||
      Implementation::generate_momenta(AuxField, pRNG);
 | 
			
		||||
      // Modify the distribution with the metric
 | 
			
		||||
      // Aux^dag M Aux
 | 
			
		||||
      M.MInvSquareRoot(AuxMom);  // AuxMom = M^{-1/2} AuxMomTemp
 | 
			
		||||
    }
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  // Correct
 | 
			
		||||
@@ -99,17 +126,34 @@ public:
 | 
			
		||||
    Hloc = zero;
 | 
			
		||||
    for (int mu = 0; mu < Nd; mu++) {
 | 
			
		||||
      // This is not very general
 | 
			
		||||
      // hide in the operators
 | 
			
		||||
      // hide in the metric
 | 
			
		||||
      auto Mom_mu = PeekIndex<LorentzIndex>(Mom, mu);
 | 
			
		||||
      auto inv_mu = PeekIndex<LorentzIndex>(inv, mu);
 | 
			
		||||
      Hloc += trace(Mom_mu * inv_mu);
 | 
			
		||||
    }
 | 
			
		||||
 | 
			
		||||
    if (0) {
 | 
			
		||||
      // Auxiliary Fields
 | 
			
		||||
      // hide in the metric
 | 
			
		||||
      M.M(AuxMom, inv);
 | 
			
		||||
      for (int mu = 0; mu < Nd; mu++) {
 | 
			
		||||
        // This is not very general
 | 
			
		||||
        // hide in the operators
 | 
			
		||||
        auto inv_mu = PeekIndex<LorentzIndex>(inv, mu);
 | 
			
		||||
        auto am_mu = PeekIndex<LorentzIndex>(AuxMom, mu);
 | 
			
		||||
        auto af_mu = PeekIndex<LorentzIndex>(AuxField, mu);
 | 
			
		||||
        Hloc += trace(am_mu * inv_mu);// p M p
 | 
			
		||||
        Hloc += trace(af_mu * af_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 MDer(in._grid);
 | 
			
		||||
@@ -118,16 +162,54 @@ public:
 | 
			
		||||
    M.Minv(in, X);  // X = G in
 | 
			
		||||
    M.MDeriv(X, MDer);  // MDer = U * dS/dU
 | 
			
		||||
    der = Implementation::projectForce(MDer);  // Ta if gauge fields
 | 
			
		||||
    
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  void AuxiliaryFieldsDerivative(MomentaField& der){
 | 
			
		||||
    der = zero;
 | 
			
		||||
    if (0){
 | 
			
		||||
    // Auxiliary fields
 | 
			
		||||
    MomentaField der_temp(der._grid);
 | 
			
		||||
    MomentaField X(der._grid);
 | 
			
		||||
    X=zero;
 | 
			
		||||
    //M.M(AuxMom, X); // X = M Aux
 | 
			
		||||
    // Two derivative terms
 | 
			
		||||
    // the Mderiv need separation of left and right terms
 | 
			
		||||
    M.MDeriv(AuxMom, der); 
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
    // this one should not be necessary (identical to the previous one)
 | 
			
		||||
    //M.MDeriv(X, AuxMom, der_temp); der += der_temp;
 | 
			
		||||
 | 
			
		||||
    der = -1.0*Implementation::projectForce(der);
 | 
			
		||||
    }
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  //   
 | 
			
		||||
  void DerivativeP(MomentaField& der){
 | 
			
		||||
    der = zero;
 | 
			
		||||
    M.Minv(Mom, der);
 | 
			
		||||
    // is the projection necessary here?
 | 
			
		||||
    // no for fields in the algebra
 | 
			
		||||
    der = Implementation::projectForce(der); 
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  void update_auxiliary_momenta(RealD ep){
 | 
			
		||||
    if(0){
 | 
			
		||||
      AuxMom -= ep * AuxField;
 | 
			
		||||
    }
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
  void update_auxiliary_fields(RealD ep){
 | 
			
		||||
    if (0) {
 | 
			
		||||
      MomentaField tmp(AuxMom._grid);
 | 
			
		||||
      MomentaField tmp2(AuxMom._grid);
 | 
			
		||||
      M.M(AuxMom, tmp);
 | 
			
		||||
      // M.M(tmp, tmp2);
 | 
			
		||||
      AuxField += ep * tmp;  // M^2 AuxMom
 | 
			
		||||
      // factor of 2?
 | 
			
		||||
    }
 | 
			
		||||
  }
 | 
			
		||||
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
 | 
			
		||||
 
 | 
			
		||||
@@ -79,14 +79,16 @@ int main (int argc, char ** argv)
 | 
			
		||||
 | 
			
		||||
  // get the deriv with respect to "U"
 | 
			
		||||
  LatticeGaugeField UdSdU(&Grid);
 | 
			
		||||
  LatticeGaugeField AuxDer(&Grid);
 | 
			
		||||
  std::cout << GridLogMessage<< "DerivativeU" << std::endl;
 | 
			
		||||
  LaplacianMomenta.DerivativeU(LaplacianMomenta.Mom, UdSdU);
 | 
			
		||||
 | 
			
		||||
  LaplacianMomenta.AuxiliaryFieldsDerivative(AuxDer);
 | 
			
		||||
  UdSdU += AuxDer;
 | 
			
		||||
 | 
			
		||||
  ////////////////////////////////////
 | 
			
		||||
  // Modify the gauge field a little 
 | 
			
		||||
  ////////////////////////////////////
 | 
			
		||||
  RealD dt = 0.001;
 | 
			
		||||
  RealD dt = 0.0001;
 | 
			
		||||
 | 
			
		||||
  LatticeColourMatrix mommu(&Grid); 
 | 
			
		||||
  LatticeColourMatrix forcemu(&Grid); 
 | 
			
		||||
 
 | 
			
		||||
@@ -76,7 +76,7 @@ int main(int argc, char **argv) {
 | 
			
		||||
  // need wrappers of the fermionic classes 
 | 
			
		||||
  // that have a complex construction
 | 
			
		||||
  // standard
 | 
			
		||||
  RealD beta = 5.6;
 | 
			
		||||
  RealD beta = 0.0;
 | 
			
		||||
  WilsonGaugeActionR Waction(beta);
 | 
			
		||||
  
 | 
			
		||||
  ActionLevel<HMCWrapper::Field> Level1(1);
 | 
			
		||||
 
 | 
			
		||||
		Reference in New Issue
	
	Block a user