mirror of
				https://github.com/paboyle/Grid.git
				synced 2025-11-03 21:44:33 +00:00 
			
		
		
		
	Hadrons: modules for testing conserved current contractions and sequential insertion.
This commit is contained in:
		@@ -32,6 +32,8 @@ See the full license in the file "LICENSE" in the top level distribution directo
 | 
			
		||||
#include <Grid/Hadrons/Modules/MContraction/DiscLoop.hpp>
 | 
			
		||||
#include <Grid/Hadrons/Modules/MContraction/Gamma3pt.hpp>
 | 
			
		||||
#include <Grid/Hadrons/Modules/MContraction/Meson.hpp>
 | 
			
		||||
#include <Grid/Hadrons/Modules/MContraction/WardIdentity.hpp>
 | 
			
		||||
#include <Grid/Hadrons/Modules/MContraction/WardIdentitySeq.hpp>
 | 
			
		||||
#include <Grid/Hadrons/Modules/MContraction/WeakHamiltonian.hpp>
 | 
			
		||||
#include <Grid/Hadrons/Modules/MContraction/WeakHamiltonianEye.hpp>
 | 
			
		||||
#include <Grid/Hadrons/Modules/MContraction/WeakHamiltonianNonEye.hpp>
 | 
			
		||||
@@ -42,6 +44,7 @@ See the full license in the file "LICENSE" in the top level distribution directo
 | 
			
		||||
#include <Grid/Hadrons/Modules/MLoop/NoiseLoop.hpp>
 | 
			
		||||
#include <Grid/Hadrons/Modules/MSolver/RBPrecCG.hpp>
 | 
			
		||||
#include <Grid/Hadrons/Modules/MSource/Point.hpp>
 | 
			
		||||
#include <Grid/Hadrons/Modules/MSource/SeqConserved.hpp>
 | 
			
		||||
#include <Grid/Hadrons/Modules/MSource/SeqGamma.hpp>
 | 
			
		||||
#include <Grid/Hadrons/Modules/MSource/Wall.hpp>
 | 
			
		||||
#include <Grid/Hadrons/Modules/MSource/Z2.hpp>
 | 
			
		||||
 
 | 
			
		||||
							
								
								
									
										151
									
								
								extras/Hadrons/Modules/MContraction/WardIdentity.hpp
									
									
									
									
									
										Normal file
									
								
							
							
						
						
									
										151
									
								
								extras/Hadrons/Modules/MContraction/WardIdentity.hpp
									
									
									
									
									
										Normal file
									
								
							@@ -0,0 +1,151 @@
 | 
			
		||||
#ifndef Hadrons_WardIdentity_hpp_
 | 
			
		||||
#define Hadrons_WardIdentity_hpp_
 | 
			
		||||
 | 
			
		||||
#include <Grid/Hadrons/Global.hpp>
 | 
			
		||||
#include <Grid/Hadrons/Module.hpp>
 | 
			
		||||
#include <Grid/Hadrons/ModuleFactory.hpp>
 | 
			
		||||
 | 
			
		||||
BEGIN_HADRONS_NAMESPACE
 | 
			
		||||
 | 
			
		||||
/*
 | 
			
		||||
  Ward Identity contractions
 | 
			
		||||
 -----------------------------
 | 
			
		||||
 
 | 
			
		||||
 * options:
 | 
			
		||||
 - q:      propagator, 5D if available (string)
 | 
			
		||||
 - q4d:    4D propagator, duplicate of q if q is not 5D (string)
 | 
			
		||||
 - action: action module used for propagator solution (string)
 | 
			
		||||
 - mass:   mass of quark (double)
 | 
			
		||||
*/
 | 
			
		||||
 | 
			
		||||
/******************************************************************************
 | 
			
		||||
 *                              WardIdentity                                  *
 | 
			
		||||
 ******************************************************************************/
 | 
			
		||||
BEGIN_MODULE_NAMESPACE(MContraction)
 | 
			
		||||
 | 
			
		||||
class WardIdentityPar: Serializable
 | 
			
		||||
{
 | 
			
		||||
public:
 | 
			
		||||
    GRID_SERIALIZABLE_CLASS_MEMBERS(WardIdentityPar,
 | 
			
		||||
                                    std::string, q,
 | 
			
		||||
                                    std::string, q4d,
 | 
			
		||||
                                    std::string, action,
 | 
			
		||||
                                    double,      mass);
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
template <typename FImpl>
 | 
			
		||||
class TWardIdentity: public Module<WardIdentityPar>
 | 
			
		||||
{
 | 
			
		||||
public:
 | 
			
		||||
    TYPE_ALIASES(FImpl,);
 | 
			
		||||
public:
 | 
			
		||||
    // constructor
 | 
			
		||||
    TWardIdentity(const std::string name);
 | 
			
		||||
    // destructor
 | 
			
		||||
    virtual ~TWardIdentity(void) = default;
 | 
			
		||||
    // dependency relation
 | 
			
		||||
    virtual std::vector<std::string> getInput(void);
 | 
			
		||||
    virtual std::vector<std::string> getOutput(void);
 | 
			
		||||
    // setup
 | 
			
		||||
    virtual void setup(void);
 | 
			
		||||
    // execution
 | 
			
		||||
    virtual void execute(void);
 | 
			
		||||
private:
 | 
			
		||||
    unsigned int Ls_;
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
MODULE_REGISTER_NS(WardIdentity, TWardIdentity<FIMPL>, MContraction);
 | 
			
		||||
 | 
			
		||||
/******************************************************************************
 | 
			
		||||
 *                     TWardIdentity implementation                           *
 | 
			
		||||
 ******************************************************************************/
 | 
			
		||||
// constructor /////////////////////////////////////////////////////////////////
 | 
			
		||||
template <typename FImpl>
 | 
			
		||||
TWardIdentity<FImpl>::TWardIdentity(const std::string name)
 | 
			
		||||
: Module<WardIdentityPar>(name)
 | 
			
		||||
{}
 | 
			
		||||
 | 
			
		||||
// dependencies/products ///////////////////////////////////////////////////////
 | 
			
		||||
template <typename FImpl>
 | 
			
		||||
std::vector<std::string> TWardIdentity<FImpl>::getInput(void)
 | 
			
		||||
{
 | 
			
		||||
    std::vector<std::string> in = {par().q, par().q4d, par().action};
 | 
			
		||||
    
 | 
			
		||||
    return in;
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
template <typename FImpl>
 | 
			
		||||
std::vector<std::string> TWardIdentity<FImpl>::getOutput(void)
 | 
			
		||||
{
 | 
			
		||||
    std::vector<std::string> out = {getName()};
 | 
			
		||||
    
 | 
			
		||||
    return out;
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
// setup ///////////////////////////////////////////////////////////////////////
 | 
			
		||||
template <typename FImpl>
 | 
			
		||||
void TWardIdentity<FImpl>::setup(void)
 | 
			
		||||
{
 | 
			
		||||
    Ls_ = env().getObjectLs(par().q);
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
// execution ///////////////////////////////////////////////////////////////////
 | 
			
		||||
template <typename FImpl>
 | 
			
		||||
void TWardIdentity<FImpl>::execute(void)
 | 
			
		||||
{
 | 
			
		||||
    LOG(Message) << "Performing Ward Identity checks for quark '" << par().q
 | 
			
		||||
                 << "'." << std::endl;
 | 
			
		||||
 | 
			
		||||
    PropagatorField psi(env().getGrid()), tmp(env().getGrid());
 | 
			
		||||
    PropagatorField q    = *env().template getObject<PropagatorField>(par().q);
 | 
			
		||||
    PropagatorField q4d  = *env().template getObject<PropagatorField>(par().q4d);
 | 
			
		||||
    FMat            &act = *(env().template getObject<FMat>(par().action));
 | 
			
		||||
    Gamma           g5(Gamma::Algebra::Gamma5);
 | 
			
		||||
    LatticeComplex  PP(env().getGrid()), PA(env().getGrid()),
 | 
			
		||||
                    c(env().getGrid()), PJ5q(env().getGrid()),
 | 
			
		||||
                    vector_WI(env().getGrid()), defect(env().getGrid());
 | 
			
		||||
    c = zero; PJ5q = zero; vector_WI = zero; defect = zero;
 | 
			
		||||
    std::vector<LatticeComplex> Vmu(Nd, c);
 | 
			
		||||
    std::vector<LatticeComplex> Amu(Nd, c);
 | 
			
		||||
 | 
			
		||||
    // Get PP, PA, V_mu, A_mu for 4D.
 | 
			
		||||
    PP = trace(adj(q4d)*q4d);
 | 
			
		||||
    PA = trace(adj(q4d)*g5*q4d);
 | 
			
		||||
    for (unsigned int mu = 0; mu < Nd; ++mu)
 | 
			
		||||
    {
 | 
			
		||||
        act.ContractConservedCurrent(q, q, tmp, Current::Vector, mu);
 | 
			
		||||
        Vmu[mu] = trace(g5*tmp);
 | 
			
		||||
        act.ContractConservedCurrent(q, q, tmp, Current::Axial, mu);
 | 
			
		||||
        Amu[mu] = trace(g5*tmp);
 | 
			
		||||
    }
 | 
			
		||||
 | 
			
		||||
    // Get PJ5q for 5D (zero for 4D).
 | 
			
		||||
    if (Ls_ > 1)
 | 
			
		||||
    {
 | 
			
		||||
        ExtractSlice(psi, q, 0, Ls_/2 - 1);
 | 
			
		||||
        psi  = 0.5 * (psi + g5*psi);
 | 
			
		||||
        ExtractSlice(tmp, q, 0, Ls_/2);
 | 
			
		||||
        psi += 0.5 * (tmp - g5*tmp);
 | 
			
		||||
        PJ5q = trace(adj(psi)*psi);
 | 
			
		||||
    }
 | 
			
		||||
 | 
			
		||||
    // Test ward identities, D_mu V_mu = 0; D_mu A_mu = 2m<PP> + 2 PJ5q
 | 
			
		||||
    for (unsigned int mu = 0; mu < Nd; ++mu)
 | 
			
		||||
    {
 | 
			
		||||
        vector_WI += Vmu[mu] - Cshift(Vmu[mu], mu, -1); 
 | 
			
		||||
        defect    += Amu[mu] - Cshift(Amu[mu], mu, -1);
 | 
			
		||||
    }
 | 
			
		||||
    defect -= 2.*PJ5q;
 | 
			
		||||
    defect -= 2.*(par().mass)*PP;
 | 
			
		||||
 | 
			
		||||
    LOG(Message) << "Vector Ward Identity check Delta_mu V_mu = " 
 | 
			
		||||
                 << norm2(vector_WI) << std::endl;
 | 
			
		||||
    LOG(Message) << "Axial Ward Identity defect Delta_mu A_mu = "
 | 
			
		||||
                 << norm2(defect) << std::endl;
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
END_MODULE_NAMESPACE
 | 
			
		||||
 | 
			
		||||
END_HADRONS_NAMESPACE
 | 
			
		||||
 | 
			
		||||
#endif // Hadrons_WardIdentity_hpp_
 | 
			
		||||
							
								
								
									
										117
									
								
								extras/Hadrons/Modules/MContraction/WardIdentitySeq.hpp
									
									
									
									
									
										Normal file
									
								
							
							
						
						
									
										117
									
								
								extras/Hadrons/Modules/MContraction/WardIdentitySeq.hpp
									
									
									
									
									
										Normal file
									
								
							@@ -0,0 +1,117 @@
 | 
			
		||||
#ifndef Hadrons_WardIdentitySeq_hpp_
 | 
			
		||||
#define Hadrons_WardIdentitySeq_hpp_
 | 
			
		||||
 | 
			
		||||
#include <Grid/Hadrons/Global.hpp>
 | 
			
		||||
#include <Grid/Hadrons/Module.hpp>
 | 
			
		||||
#include <Grid/Hadrons/ModuleFactory.hpp>
 | 
			
		||||
 | 
			
		||||
BEGIN_HADRONS_NAMESPACE
 | 
			
		||||
 | 
			
		||||
/*
 | 
			
		||||
  Ward Identity contractions using sequential propagators.
 | 
			
		||||
 -----------------------------
 | 
			
		||||
 
 | 
			
		||||
 * options:
 | 
			
		||||
 - q_x: propagator, mu = x current insertion (string).
 | 
			
		||||
 - q_y: propagator, mu = y current insertion (string). 
 | 
			
		||||
 - q_z: propagator, mu = z current insertion (string).
 | 
			
		||||
 - q_t: propagator, mu = t current insertion (string).
 | 
			
		||||
*/
 | 
			
		||||
 | 
			
		||||
/******************************************************************************
 | 
			
		||||
 *                            WardIdentitySeq                                 *
 | 
			
		||||
 ******************************************************************************/
 | 
			
		||||
BEGIN_MODULE_NAMESPACE(MContraction)
 | 
			
		||||
 | 
			
		||||
class WardIdentitySeqPar: Serializable
 | 
			
		||||
{
 | 
			
		||||
public:
 | 
			
		||||
    GRID_SERIALIZABLE_CLASS_MEMBERS(WardIdentitySeqPar,
 | 
			
		||||
                                    std::string, q_x,
 | 
			
		||||
                                    std::string, q_y,
 | 
			
		||||
                                    std::string, q_z,
 | 
			
		||||
                                    std::string, q_t);
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
template <typename FImpl>
 | 
			
		||||
class TWardIdentitySeq: public Module<WardIdentitySeqPar>
 | 
			
		||||
{
 | 
			
		||||
public:
 | 
			
		||||
    TYPE_ALIASES(FImpl,);
 | 
			
		||||
public:
 | 
			
		||||
    // constructor
 | 
			
		||||
    TWardIdentitySeq(const std::string name);
 | 
			
		||||
    // destructor
 | 
			
		||||
    virtual ~TWardIdentitySeq(void) = default;
 | 
			
		||||
    // dependency relation
 | 
			
		||||
    virtual std::vector<std::string> getInput(void);
 | 
			
		||||
    virtual std::vector<std::string> getOutput(void);
 | 
			
		||||
    // setup
 | 
			
		||||
    virtual void setup(void);
 | 
			
		||||
    // execution
 | 
			
		||||
    virtual void execute(void);
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
MODULE_REGISTER_NS(WardIdentitySeq, TWardIdentitySeq<FIMPL>, MContraction);
 | 
			
		||||
 | 
			
		||||
/******************************************************************************
 | 
			
		||||
 *                 TWardIdentitySeq implementation                             *
 | 
			
		||||
 ******************************************************************************/
 | 
			
		||||
// constructor /////////////////////////////////////////////////////////////////
 | 
			
		||||
template <typename FImpl>
 | 
			
		||||
TWardIdentitySeq<FImpl>::TWardIdentitySeq(const std::string name)
 | 
			
		||||
: Module<WardIdentitySeqPar>(name)
 | 
			
		||||
{}
 | 
			
		||||
 | 
			
		||||
// dependencies/products ///////////////////////////////////////////////////////
 | 
			
		||||
template <typename FImpl>
 | 
			
		||||
std::vector<std::string> TWardIdentitySeq<FImpl>::getInput(void)
 | 
			
		||||
{
 | 
			
		||||
    std::vector<std::string> in = {par().q_x, par().q_y, par().q_z, par().q_t};
 | 
			
		||||
    
 | 
			
		||||
    return in;
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
template <typename FImpl>
 | 
			
		||||
std::vector<std::string> TWardIdentitySeq<FImpl>::getOutput(void)
 | 
			
		||||
{
 | 
			
		||||
    std::vector<std::string> out = {getName()};
 | 
			
		||||
    
 | 
			
		||||
    return out;
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
// setup ///////////////////////////////////////////////////////////////////////
 | 
			
		||||
template <typename FImpl>
 | 
			
		||||
void TWardIdentitySeq<FImpl>::setup(void)
 | 
			
		||||
{
 | 
			
		||||
    
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
// execution ///////////////////////////////////////////////////////////////////
 | 
			
		||||
template <typename FImpl>
 | 
			
		||||
void TWardIdentitySeq<FImpl>::execute(void)
 | 
			
		||||
{
 | 
			
		||||
    LatticeComplex  vector_WI(env().getGrid()), c(env().getGrid());
 | 
			
		||||
    PropagatorField q_x = *env().template getObject<PropagatorField>(par().q_x);
 | 
			
		||||
    PropagatorField q_y = *env().template getObject<PropagatorField>(par().q_y);
 | 
			
		||||
    PropagatorField q_z = *env().template getObject<PropagatorField>(par().q_z);
 | 
			
		||||
    PropagatorField q_t = *env().template getObject<PropagatorField>(par().q_t);
 | 
			
		||||
    PropagatorField *q[Nd] = {&q_x, &q_y, &q_z, &q_t};
 | 
			
		||||
    Gamma           g5(Gamma::Algebra::Gamma5);
 | 
			
		||||
 | 
			
		||||
    // Check D_mu V_mu = 0
 | 
			
		||||
    for (unsigned int mu = 0; mu < Nd; ++mu)
 | 
			
		||||
    {
 | 
			
		||||
        c = trace(g5*(*q[mu]));
 | 
			
		||||
        vector_WI += c - Cshift(c, mu, -1);
 | 
			
		||||
    }
 | 
			
		||||
 | 
			
		||||
    LOG(Message) << "Ward Identity checks for sequential vector current "
 | 
			
		||||
                 << "insertion = " << norm2(vector_WI) << std::endl;
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
END_MODULE_NAMESPACE
 | 
			
		||||
 | 
			
		||||
END_HADRONS_NAMESPACE
 | 
			
		||||
 | 
			
		||||
#endif // Hadrons_WardIdentitySeq_hpp_
 | 
			
		||||
							
								
								
									
										129
									
								
								extras/Hadrons/Modules/MSource/SeqConserved.hpp
									
									
									
									
									
										Normal file
									
								
							
							
						
						
									
										129
									
								
								extras/Hadrons/Modules/MSource/SeqConserved.hpp
									
									
									
									
									
										Normal file
									
								
							@@ -0,0 +1,129 @@
 | 
			
		||||
#ifndef Hadrons_SeqConserved_hpp_
 | 
			
		||||
#define Hadrons_SeqConserved_hpp_
 | 
			
		||||
 | 
			
		||||
#include <Grid/Hadrons/Global.hpp>
 | 
			
		||||
#include <Grid/Hadrons/Module.hpp>
 | 
			
		||||
#include <Grid/Hadrons/ModuleFactory.hpp>
 | 
			
		||||
 | 
			
		||||
BEGIN_HADRONS_NAMESPACE
 | 
			
		||||
 | 
			
		||||
/*
 | 
			
		||||
 
 | 
			
		||||
 Sequential source
 | 
			
		||||
 -----------------------------
 | 
			
		||||
 * src_x = q_x * theta(x_3 - tA) * theta(tB - x_3) * J_mu * exp(i x.mom)
 | 
			
		||||
 
 | 
			
		||||
 * options:
 | 
			
		||||
 - q: input propagator (string)
 | 
			
		||||
 - action: fermion action used for propagator q (string)
 | 
			
		||||
 - tA: begin timeslice (integer)
 | 
			
		||||
 - tB: end timesilce (integer)
 | 
			
		||||
 - curr_type: type of conserved current to insert (Current)
 | 
			
		||||
 - mu: Lorentz index of current to insert (integer)
 | 
			
		||||
 - mom: momentum insertion, space-separated float sequence (e.g ".1 .2 1. 0.")
 | 
			
		||||
 
 | 
			
		||||
 */
 | 
			
		||||
 | 
			
		||||
/******************************************************************************
 | 
			
		||||
 *                              SeqConserved                                  *
 | 
			
		||||
 ******************************************************************************/
 | 
			
		||||
BEGIN_MODULE_NAMESPACE(MSource)
 | 
			
		||||
 | 
			
		||||
class SeqConservedPar: Serializable
 | 
			
		||||
{
 | 
			
		||||
public:
 | 
			
		||||
    GRID_SERIALIZABLE_CLASS_MEMBERS(SeqConservedPar,
 | 
			
		||||
                                    std::string,  q,
 | 
			
		||||
                                    std::string,  action,
 | 
			
		||||
                                    unsigned int, tA,
 | 
			
		||||
                                    unsigned int, tB,
 | 
			
		||||
                                    Current,      curr_type,
 | 
			
		||||
                                    unsigned int, mu,
 | 
			
		||||
                                    std::string,  mom);
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
template <typename FImpl>
 | 
			
		||||
class TSeqConserved: public Module<SeqConservedPar>
 | 
			
		||||
{
 | 
			
		||||
public:
 | 
			
		||||
    TYPE_ALIASES(FImpl,);
 | 
			
		||||
public:
 | 
			
		||||
    // constructor
 | 
			
		||||
    TSeqConserved(const std::string name);
 | 
			
		||||
    // destructor
 | 
			
		||||
    virtual ~TSeqConserved(void) = default;
 | 
			
		||||
    // dependency relation
 | 
			
		||||
    virtual std::vector<std::string> getInput(void);
 | 
			
		||||
    virtual std::vector<std::string> getOutput(void);
 | 
			
		||||
    // setup
 | 
			
		||||
    virtual void setup(void);
 | 
			
		||||
    // execution
 | 
			
		||||
    virtual void execute(void);
 | 
			
		||||
};
 | 
			
		||||
 | 
			
		||||
MODULE_REGISTER_NS(SeqConserved, TSeqConserved<FIMPL>, MSource);
 | 
			
		||||
 | 
			
		||||
/******************************************************************************
 | 
			
		||||
 *                      TSeqConserved implementation                          *
 | 
			
		||||
 ******************************************************************************/
 | 
			
		||||
// constructor /////////////////////////////////////////////////////////////////
 | 
			
		||||
template <typename FImpl>
 | 
			
		||||
TSeqConserved<FImpl>::TSeqConserved(const std::string name)
 | 
			
		||||
: Module<SeqConservedPar>(name)
 | 
			
		||||
{}
 | 
			
		||||
 | 
			
		||||
// dependencies/products ///////////////////////////////////////////////////////
 | 
			
		||||
template <typename FImpl>
 | 
			
		||||
std::vector<std::string> TSeqConserved<FImpl>::getInput(void)
 | 
			
		||||
{
 | 
			
		||||
    std::vector<std::string> in;
 | 
			
		||||
    
 | 
			
		||||
    return in;
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
template <typename FImpl>
 | 
			
		||||
std::vector<std::string> TSeqConserved<FImpl>::getOutput(void)
 | 
			
		||||
{
 | 
			
		||||
    std::vector<std::string> out = {getName()};
 | 
			
		||||
    
 | 
			
		||||
    return out;
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
// setup ///////////////////////////////////////////////////////////////////////
 | 
			
		||||
template <typename FImpl>
 | 
			
		||||
void TSeqConserved<FImpl>::setup(void)
 | 
			
		||||
{
 | 
			
		||||
    
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
// execution ///////////////////////////////////////////////////////////////////
 | 
			
		||||
template <typename FImpl>
 | 
			
		||||
void TSeqConserved<FImpl>::execute(void)
 | 
			
		||||
{
 | 
			
		||||
    if (par().tA == par().tB)
 | 
			
		||||
    {
 | 
			
		||||
        LOG(Message) << "Generating sequential source with conserved "
 | 
			
		||||
                     << par().curr_type << " current insertion (mu = " 
 | 
			
		||||
                     << par().mu << ") at " << "t = " << par().tA << std::endl;
 | 
			
		||||
    }
 | 
			
		||||
    else
 | 
			
		||||
    {
 | 
			
		||||
        LOG(Message) << "Generating sequential source with conserved "
 | 
			
		||||
                     << par().curr_type << " current insertion (mu = " 
 | 
			
		||||
                     << par().mu << ") for " << par().tA << " <= t <= " 
 | 
			
		||||
                     << par().tB << std::endl;
 | 
			
		||||
    }
 | 
			
		||||
    PropagatorField &src = *env().template createLattice<PropagatorField>(getName());
 | 
			
		||||
    PropagatorField &q   = *env().template getObject<PropagatorField>(par().q);
 | 
			
		||||
    FMat            &mat = *(env().template getObject<FMat>(par().action));
 | 
			
		||||
 | 
			
		||||
    std::vector<Real> mom = strToVec<Real>(par().mom);
 | 
			
		||||
    mat.SeqConservedCurrent(q, src, par().curr_type, par().mu, 
 | 
			
		||||
                            mom, par().tA, par().tB);
 | 
			
		||||
}
 | 
			
		||||
 | 
			
		||||
END_MODULE_NAMESPACE
 | 
			
		||||
 | 
			
		||||
END_HADRONS_NAMESPACE
 | 
			
		||||
 | 
			
		||||
#endif // Hadrons_SeqConserved_hpp_
 | 
			
		||||
@@ -13,6 +13,8 @@ modules_hpp =\
 | 
			
		||||
  Modules/MContraction/DiscLoop.hpp \
 | 
			
		||||
  Modules/MContraction/Gamma3pt.hpp \
 | 
			
		||||
  Modules/MContraction/Meson.hpp \
 | 
			
		||||
  Modules/MContraction/WardIdentity.hpp \
 | 
			
		||||
  Modules/MContraction/WardIdentitySeq.hpp \
 | 
			
		||||
  Modules/MContraction/WeakHamiltonian.hpp \
 | 
			
		||||
  Modules/MContraction/WeakHamiltonianEye.hpp \
 | 
			
		||||
  Modules/MContraction/WeakHamiltonianNonEye.hpp \
 | 
			
		||||
@@ -23,6 +25,7 @@ modules_hpp =\
 | 
			
		||||
  Modules/MLoop/NoiseLoop.hpp \
 | 
			
		||||
  Modules/MSolver/RBPrecCG.hpp \
 | 
			
		||||
  Modules/MSource/Point.hpp \
 | 
			
		||||
  Modules/MSource/SeqConserved.hpp \
 | 
			
		||||
  Modules/MSource/SeqGamma.hpp \
 | 
			
		||||
  Modules/MSource/Wall.hpp \
 | 
			
		||||
  Modules/MSource/Z2.hpp \
 | 
			
		||||
 
 | 
			
		||||
		Reference in New Issue
	
	Block a user