1
0
mirror of https://github.com/paboyle/Grid.git synced 2024-11-10 07:55:35 +00:00

Hadrons: progress on the interface, genetic algorithm freezing

This commit is contained in:
Antonin Portelli 2017-12-01 19:38:23 +00:00
parent a3fe874a5b
commit 514993ed17
14 changed files with 164 additions and 93 deletions

View File

@ -98,6 +98,8 @@ void Application::run(void)
{
parseParameterFile(parameterFileName_);
}
env().checkGraph();
env().printContent();
if (!scheduled_)
{
schedule();
@ -124,8 +126,14 @@ void Application::parseParameterFile(const std::string parameterFileName)
LOG(Message) << "Building application from '" << parameterFileName << "'..." << std::endl;
read(reader, "parameters", par);
setPar(par);
push(reader, "modules");
push(reader, "module");
if (!push(reader, "modules"))
{
HADRON_ERROR("Cannot open node 'modules' in parameter file '" + parameterFileName + "'");
}
if (!push(reader, "module"))
{
HADRON_ERROR("Cannot open node 'modules/module' in parameter file '" + parameterFileName + "'");
}
do
{
read(reader, "id", id);
@ -186,6 +194,8 @@ void Application::schedule(void)
// build module dependency graph
LOG(Message) << "Building module graph..." << std::endl;
auto graph = env().makeModuleGraph();
LOG(Debug) << "Module graph:" << std::endl;
LOG(Debug) << graph << std::endl;
auto con = graph.getConnectedComponents();
// constrained topological sort using a genetic algorithm
@ -329,4 +339,3 @@ void Application::memoryProfile(void)
HadronsLogMessage.Active(msg);
}

View File

@ -333,6 +333,17 @@ Graph<unsigned int> Environment::makeModuleGraph(void) const
return moduleGraph;
}
void Environment::checkGraph(void) const
{
for (auto &o: object_)
{
if (o.module < 0)
{
HADRON_ERROR("object '" + o.name + "' does not have a creator");
}
}
}
#define BIG_SEP "==============="
#define SEP "---------------"
#define MEM_MSG(size)\
@ -346,6 +357,7 @@ Environment::executeProgram(const std::vector<unsigned int> &p)
bool continueCollect, nothingFreed;
// build garbage collection schedule
LOG(Debug) << "Building garbage collection schedule..." << std::endl;
freeProg.resize(p.size());
for (unsigned int i = 0; i < object_.size(); ++i)
{
@ -359,11 +371,12 @@ Environment::executeProgram(const std::vector<unsigned int> &p)
auto it = std::find_if(p.rbegin(), p.rend(), pred);
if (it != p.rend())
{
freeProg[std::distance(p.rend(), it) - 1].insert(i);
freeProg[std::distance(it, p.rend()) - 1].insert(i);
}
}
// program execution
LOG(Debug) << "Executing program..." << std::endl;
for (unsigned int i = 0; i < p.size(); ++i)
{
// execute module
@ -712,16 +725,16 @@ void Environment::freeAll(void)
void Environment::printContent(void)
{
LOG(Message) << "Modules: " << std::endl;
LOG(Debug) << "Modules: " << std::endl;
for (unsigned int i = 0; i < module_.size(); ++i)
{
LOG(Message) << std::setw(4) << i << ": "
<< getModuleName(i) << std::endl;
LOG(Debug) << std::setw(4) << i << ": "
<< getModuleName(i) << std::endl;
}
LOG(Message) << "Objects: " << std::endl;
LOG(Debug) << "Objects: " << std::endl;
for (unsigned int i = 0; i < object_.size(); ++i)
{
LOG(Message) << std::setw(4) << i << ": "
<< getObjectName(i) << std::endl;
LOG(Debug) << std::setw(4) << i << ": "
<< getObjectName(i) << std::endl;
}
}

View File

@ -76,6 +76,7 @@ public:
typedef std::unique_ptr<GridRedBlackCartesian> GridRbPt;
typedef std::unique_ptr<GridParallelRNG> RngPt;
typedef std::unique_ptr<LatticeBase> LatticePt;
enum class Storage {object, cache, temporary};
private:
struct ModuleInfo
{
@ -88,6 +89,7 @@ private:
struct ObjInfo
{
Size size{0};
Storage storage{Storage::object};
unsigned int Ls{0};
const std::type_info *type{nullptr};
std::string name;
@ -140,18 +142,17 @@ public:
bool hasModule(const unsigned int address) const;
bool hasModule(const std::string name) const;
Graph<unsigned int> makeModuleGraph(void) const;
void checkGraph(void) const;
Size executeProgram(const std::vector<unsigned int> &p);
Size executeProgram(const std::vector<std::string> &p);
// general memory management
void addObject(const std::string name,
const int moduleAddress = -1);
template <typename T, typename ... Ts>
template <typename T, typename P>
void createObject(const std::string name,
const Storage storage,
const unsigned int Ls,
Ts ... args);
template <typename T>
void createLattice(const std::string name,
const unsigned int Ls = 1);
P &&pt);
template <typename T>
T * getObject(const unsigned int address) const;
template <typename T>
@ -203,6 +204,7 @@ private:
// module and related maps
std::vector<ModuleInfo> module_;
std::map<std::string, unsigned int> moduleAddress_;
std::string currentModule_{""};
// lattice store
std::map<unsigned int, LatticePt> lattice_;
// object store
@ -281,9 +283,11 @@ M * Environment::getModule(const std::string name) const
return getModule<M>(getModuleAddress(name));
}
template <typename T, typename ... Ts>
void Environment::createObject(const std::string name, const unsigned int Ls,
Ts ... args)
template <typename T, typename P>
void Environment::createObject(const std::string name,
const Environment::Storage storage,
const unsigned int Ls,
P &&pt)
{
if (!hasObject(name))
{
@ -296,11 +300,13 @@ void Environment::createObject(const std::string name, const unsigned int Ls,
{
MemoryStats memStats;
MemoryProfiler::stats = &memStats;
object_[address].Ls = Ls;
object_[address].data.reset(new Holder<T>(new T(args...)));
object_[address].size = memStats.totalAllocated;
object_[address].type = &typeid(T);
MemoryProfiler::stats = &memStats;
object_[address].storage = storage;
object_[address].Ls = Ls;
object_[address].data.reset(new Holder<T>(pt));
object_[address].size = memStats.totalAllocated;
object_[address].type = &typeid(T);
MemoryProfiler::stats = nullptr;
}
else
{
@ -308,14 +314,6 @@ void Environment::createObject(const std::string name, const unsigned int Ls,
}
}
template <typename T>
void Environment::createLattice(const std::string name, const unsigned int Ls)
{
GridCartesian *g = getGrid(Ls);
createObject<T>(name, Ls, g);
}
template <typename T>
T * Environment::getObject(const unsigned int address) const
{

View File

@ -87,13 +87,54 @@ public:\
static ns##mod##ModuleRegistrar ns##mod##ModuleRegistrarInstance;
#define ARG(...) __VA_ARGS__
#define MACRO_REDIRECT(arg1, arg2, arg3, macro, ...) macro
#define mCreateObj(type, name, Ls, ...)\
env().template createObject<type>(name, Ls, __VA_ARGS__)
#define mGetObj(type, name)\
#define envGet(type, name)\
*env().template getObject<type>(name)
#define envGetTmp(type, name)\
*env().template getObject<type>(getName() + "_tmp_" + name)
#define envIsType(type, name)\
env().template getObject<type>(name)
#define envCreate(type, name, Ls, pt)\
env().template createObject<type>(name, Environment::Storage::object, Ls, pt)
#define envCreateLat4(type, name)\
envCreate(type, name, 1, new type(env().getGrid()))
#define envCreateLat5(type, name, Ls)\
envCreate(type, name, Ls, new type(env().getGrid(Ls)))
#define envCreateLat(...)\
MACRO_REDIRECT(__VA_ARGS__, envCreateLat5, envCreateLat4)(__VA_ARGS__)
#define envCache(type, name, Ls, pt)\
env().template createObject<type>(name, Environment::Storage::cache, Ls, pt)
#define envCacheLat4(type, name)\
envCache(type, name, 1, new type(env().getGrid()))
#define envCacheLat5(type, name, Ls)\
envCache(type, name, Ls, new type(env().getGrid(Ls)))
#define envCacheLat(...)\
MACRO_REDIRECT(__VA_ARGS__, envCacheLat5, envCacheLat4)(__VA_ARGS__)
#define envTmp(type, name, Ls, pt)\
env().template createObject<type>(getName() + "_tmp_" + name, \
Environment::Storage::temporary, Ls, pt)
#define envTmpLat4(type, name)\
envTmp(type, name, 1, new type(env().getGrid()))
#define envTmpLat5(type, name, Ls)\
envTmp(type, name, Ls, new type(env().getGrid(Ls)))
#define envTmpLat(...)\
MACRO_REDIRECT(__VA_ARGS__, envTmpLat5, envTmpLat4)(__VA_ARGS__)
/******************************************************************************
* Module class *
******************************************************************************/

View File

@ -33,13 +33,13 @@ See the full license in the file "LICENSE" in the top level distribution directo
// #include <Grid/Hadrons/Modules/MContraction/Baryon.hpp>
// #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/Meson.hpp>
// #include <Grid/Hadrons/Modules/MContraction/WardIdentity.hpp>
// #include <Grid/Hadrons/Modules/MContraction/WeakHamiltonian.hpp>
// #include <Grid/Hadrons/Modules/MContraction/WeakHamiltonianEye.hpp>
// #include <Grid/Hadrons/Modules/MContraction/WeakHamiltonianNonEye.hpp>
// #include <Grid/Hadrons/Modules/MContraction/WeakNeutral4ptDisc.hpp>
// #include <Grid/Hadrons/Modules/MFermion/GaugeProp.hpp>
#include <Grid/Hadrons/Modules/MFermion/GaugeProp.hpp>
// #include <Grid/Hadrons/Modules/MGauge/Load.hpp>
// #include <Grid/Hadrons/Modules/MGauge/Random.hpp>
// #include <Grid/Hadrons/Modules/MGauge/StochEm.hpp>
@ -48,10 +48,10 @@ See the full license in the file "LICENSE" in the top level distribution directo
// #include <Grid/Hadrons/Modules/MScalar/ChargedProp.hpp>
// #include <Grid/Hadrons/Modules/MScalar/FreeProp.hpp>
// #include <Grid/Hadrons/Modules/MScalar/Scalar.hpp>
// #include <Grid/Hadrons/Modules/MSink/Point.hpp>
#include <Grid/Hadrons/Modules/MSink/Point.hpp>
// #include <Grid/Hadrons/Modules/MSink/Smear.hpp>
#include <Grid/Hadrons/Modules/MSolver/RBPrecCG.hpp>
// #include <Grid/Hadrons/Modules/MSource/Point.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>

View File

@ -110,15 +110,16 @@ void TDWF<FImpl>::setup(void)
LOG(Message) << "Fermion boundary conditions: " << par().boundary
<< std::endl;
env().createGrid(par().Ls);
auto &U = mGetObj(LatticeGaugeField, par().gauge);
auto &U = envGet(LatticeGaugeField, par().gauge);
auto &g4 = *env().getGrid();
auto &grb4 = *env().getRbGrid();
auto &g5 = *env().getGrid(par().Ls);
auto &grb5 = *env().getRbGrid(par().Ls);
std::vector<Complex> boundary = strToVec<Complex>(par().boundary);
typename DomainWallFermion<FImpl>::ImplParams implParams(boundary);
mCreateObj(DomainWallFermion<FImpl>, getName(), par().Ls,
U, g5, grb5, g4, grb4, par().mass, par().M5, implParams);
envCreate(FMat, getName(), par().Ls,
new DomainWallFermion<FImpl>(U, g5, grb5, g4, grb4, par().mass,
par().M5, implParams));
}
// execution ///////////////////////////////////////////////////////////////////

View File

@ -105,13 +105,14 @@ void TWilson<FImpl>::setup(void)
<< " using gauge field '" << par().gauge << "'" << std::endl;
LOG(Message) << "Fermion boundary conditions: " << par().boundary
<< std::endl;
auto &U = mGetObj(LatticeGaugeField, par().gauge);
auto &U = envGet(LatticeGaugeField, par().gauge);
auto &grid = *env().getGrid();
auto &gridRb = *env().getRbGrid();
std::vector<Complex> boundary = strToVec<Complex>(par().boundary);
typename WilsonFermion<FImpl>::ImplParams implParams(boundary);
mCreateObj(WilsonFermion<FImpl>, getName(), 1, U, grid, gridRb, par().mass,
implParams);
envCreate(FMat, getName(), 1, new WilsonFermion<FImpl>(U, grid, gridRb,
par().mass,
implParams));
}
// execution ///////////////////////////////////////////////////////////////////

View File

@ -153,7 +153,6 @@ void TMeson<FImpl1, FImpl2>::parseGammaString(std::vector<GammaPair> &gammaList)
}
}
// execution ///////////////////////////////////////////////////////////////////
#define mesonConnected(q1, q2, gSnk, gSrc) \
(g5*(gSnk))*(q1)*(adj(gSrc)*g5)*adj(q2)
@ -180,11 +179,11 @@ void TMeson<FImpl1, FImpl2>::execute(void)
result[i].gamma_src = gammaList[i].second;
result[i].corr.resize(nt);
}
if (env().template isObjectOfType<SlicedPropagator1>(par().q1) and
env().template isObjectOfType<SlicedPropagator2>(par().q2))
if (envIsType(SlicedPropagator1, par().q1) and
envIsType(SlicedPropagator2, par().q2))
{
SlicedPropagator1 &q1 = *env().template getObject<SlicedPropagator1>(par().q1);
SlicedPropagator2 &q2 = *env().template getObject<SlicedPropagator2>(par().q2);
SlicedPropagator1 &q1 = envGet(SlicedPropagator1, par().q1);
SlicedPropagator2 &q2 = envGet(SlicedPropagator2, par().q2);
LOG(Message) << "(propagator already sinked)" << std::endl;
for (unsigned int i = 0; i < result.size(); ++i)
@ -200,8 +199,8 @@ void TMeson<FImpl1, FImpl2>::execute(void)
}
else
{
PropagatorField1 &q1 = *env().template getObject<PropagatorField1>(par().q1);
PropagatorField2 &q2 = *env().template getObject<PropagatorField2>(par().q2);
PropagatorField1 &q1 = envGet(PropagatorField1, par().q1);
PropagatorField2 &q2 = envGet(PropagatorField2, par().q2);
LatticeComplex c(env().getGrid());
LOG(Message) << "(using sink '" << par().sink << "')" << std::endl;
@ -214,15 +213,14 @@ void TMeson<FImpl1, FImpl2>::execute(void)
ns = env().getModuleNamespace(env().getObjectModule(par().sink));
if (ns == "MSource")
{
PropagatorField1 &sink =
*env().template getObject<PropagatorField1>(par().sink);
PropagatorField1 &sink = envGet(PropagatorField1, par().sink);
c = trace(mesonConnected(q1, q2, gSnk, gSrc)*sink);
sliceSum(c, buf, Tp);
}
else if (ns == "MSink")
{
SinkFnScalar &sink = *env().template getObject<SinkFnScalar>(par().sink);
SinkFnScalar &sink = envGet(SinkFnScalar, par().sink);
c = trace(mesonConnected(q1, q2, gSnk, gSrc));
buf = sink(c);

View File

@ -127,10 +127,13 @@ template <typename FImpl>
void TGaugeProp<FImpl>::setup(void)
{
Ls_ = env().getObjectLs(par().solver);
env().template registerLattice<PropagatorField>(getName());
envCreateLat(PropagatorField, getName());
envTmpLat(FermionField, "source", Ls_);
envTmpLat(FermionField, "sol", Ls_);
envTmpLat(FermionField, "tmp");
if (Ls_ > 1)
{
env().template registerLattice<PropagatorField>(getName() + "_5d", Ls_);
envCreateLat(PropagatorField, getName() + "_5d", Ls_);
}
}
@ -139,21 +142,18 @@ template <typename FImpl>
void TGaugeProp<FImpl>::execute(void)
{
LOG(Message) << "Computing quark propagator '" << getName() << "'"
<< std::endl;
<< std::endl;
FermionField source(env().getGrid(Ls_)), sol(env().getGrid(Ls_)),
tmp(env().getGrid());
FermionField &source = envGetTmp(FermionField, "source");
FermionField &sol = envGetTmp(FermionField, "sol");
FermionField &tmp = envGetTmp(FermionField, "tmp");
std::string propName = (Ls_ == 1) ? getName() : (getName() + "_5d");
PropagatorField &prop = *env().template createLattice<PropagatorField>(propName);
PropagatorField &fullSrc = *env().template getObject<PropagatorField>(par().source);
SolverFn &solver = *env().template getObject<SolverFn>(par().solver);
if (Ls_ > 1)
{
env().template createLattice<PropagatorField>(getName());
}
PropagatorField &prop = envGet(PropagatorField, propName);
PropagatorField &fullSrc = envGet(PropagatorField, par().source);
SolverFn &solver = envGet(SolverFn, par().solver);
LOG(Message) << "Inverting using solver '" << par().solver
<< "' on source '" << par().source << "'" << std::endl;
<< "' on source '" << par().source << "'" << std::endl;
for (unsigned int s = 0; s < Ns; ++s)
for (unsigned int c = 0; c < Nc; ++c)
{
@ -190,8 +190,7 @@ void TGaugeProp<FImpl>::execute(void)
// create 4D propagators from 5D one if necessary
if (Ls_ > 1)
{
PropagatorField &p4d =
*env().template getObject<PropagatorField>(getName());
PropagatorField &p4d = envGet(PropagatorField, getName());
make_4D(sol, tmp, Ls_);
FermToProp(p4d, tmp, s, c);
}

View File

@ -57,13 +57,13 @@ std::vector<std::string> TUnit::getOutput(void)
// setup ///////////////////////////////////////////////////////////////////////
void TUnit::setup(void)
{
mCreateObj(LatticeGaugeField, getName(), 1, env().getGrid());
envCreateLat(LatticeGaugeField, getName());
}
// execution ///////////////////////////////////////////////////////////////////
void TUnit::execute(void)
{
LOG(Message) << "Creating unit gauge configuration" << std::endl;
auto &U = mGetObj(LatticeGaugeField, getName());
auto &U = envGet(LatticeGaugeField, getName());
SU3::ColdConfiguration(*env().get4dRng(), U);
}

View File

@ -65,6 +65,9 @@ public:
virtual void setup(void);
// execution
virtual void execute(void);
private:
bool hasPhase_{false};
std::string momphName_;
};
MODULE_REGISTER_NS(Point, TPoint<FIMPL>, MSink);
@ -77,6 +80,7 @@ MODULE_REGISTER_NS(ScalarPoint, TPoint<ScalarImplCR>, MSink);
template <typename FImpl>
TPoint<FImpl>::TPoint(const std::string name)
: Module<PointPar>(name)
, momphName_ (name + "_momph")
{}
// dependencies/products ///////////////////////////////////////////////////////
@ -100,30 +104,36 @@ std::vector<std::string> TPoint<FImpl>::getOutput(void)
template <typename FImpl>
void TPoint<FImpl>::setup(void)
{
unsigned int size;
size = env().template lattice4dSize<LatticeComplex>();
env().registerObject(getName(), size);
envTmpLat(LatticeComplex, "coor");
envCacheLat(LatticeComplex, momphName_);
envCreate(SinkFn, getName(), 1, nullptr);
}
// execution ///////////////////////////////////////////////////////////////////
template <typename FImpl>
void TPoint<FImpl>::execute(void)
{
std::vector<Real> p = strToVec<Real>(par().mom);
LatticeComplex ph(env().getGrid()), coor(env().getGrid());
std::vector<Real> p = strToVec<Real>(par().mom);
LatticeComplex &ph = envGet(LatticeComplex, momphName_);
Complex i(0.0,1.0);
LOG(Message) << "Setting up point sink function for momentum ["
<< par().mom << "]" << std::endl;
ph = zero;
for(unsigned int mu = 0; mu < env().getNd(); mu++)
if (!hasPhase_)
{
LatticeCoordinate(coor, mu);
ph = ph + (p[mu]/env().getGrid()->_fdimensions[mu])*coor;
LatticeComplex &coor = envGetTmp(LatticeComplex, "coor");
ph = zero;
for(unsigned int mu = 0; mu < env().getNd(); mu++)
{
LatticeCoordinate(coor, mu);
ph = ph + (p[mu]/env().getGrid()->_fdimensions[mu])*coor;
}
ph = exp((Real)(2*M_PI)*i*ph);
hasPhase_ = true;
}
ph = exp((Real)(2*M_PI)*i*ph);
auto sink = [ph](const PropagatorField &field)
auto sink = [&ph](const PropagatorField &field)
{
SlicedPropagator res;
PropagatorField tmp = ph*field;
@ -132,7 +142,7 @@ void TPoint<FImpl>::execute(void)
return res;
};
env().setObject(getName(), new SinkFn(sink));
envGet(SinkFn, getName()) = sink;
}
END_MODULE_NAMESPACE

View File

@ -105,7 +105,7 @@ void TRBPrecCG<FImpl>::setup(void)
<< par().residual << std::endl;
auto Ls = env().getObjectLs(par().action);
auto &mat = mGetObj(FMat, par().action);
auto &mat = envGet(FMat, par().action);
auto solver = [&mat, this](FermionField &sol, const FermionField &source)
{
ConjugateGradient<FermionField> cg(par().residual, 10000);
@ -113,7 +113,7 @@ void TRBPrecCG<FImpl>::setup(void)
schurSolver(mat, source, sol);
};
mCreateObj(SolverFn, getName(), Ls, solver);
envCreate(SolverFn, getName(), Ls, new SolverFn(solver));
env().addOwnership(getName(), par().action);
}

View File

@ -111,7 +111,7 @@ std::vector<std::string> TPoint<FImpl>::getOutput(void)
template <typename FImpl>
void TPoint<FImpl>::setup(void)
{
env().template registerLattice<PropagatorField>(getName());
envCreateLat(PropagatorField, getName());
}
// execution ///////////////////////////////////////////////////////////////////
@ -123,7 +123,7 @@ void TPoint<FImpl>::execute(void)
LOG(Message) << "Creating point source at position [" << par().position
<< "]" << std::endl;
PropagatorField &src = *env().template createLattice<PropagatorField>(getName());
PropagatorField &src = envGet(PropagatorField, getName());
id = 1.;
src = zero;
pokeSite(id, src, position);

View File

@ -12,13 +12,16 @@ modules_cc =\
modules_hpp =\
Modules/MAction/DWF.hpp \
Modules/MAction/Wilson.hpp \
Modules/MSink/Point.hpp \
Modules/MSource/Point.hpp \
Modules/MGauge/Unit.hpp \
Modules/MSolver/RBPrecCG.hpp
Modules/MSolver/RBPrecCG.hpp \
Modules/MFermion/GaugeProp.hpp \
Modules/MContraction/Meson.hpp
# Modules/MContraction/Baryon.hpp \
# Modules/MContraction/DiscLoop.hpp \
# Modules/MContraction/Gamma3pt.hpp \
# Modules/MContraction/Meson.hpp \
# Modules/MContraction/WardIdentity.hpp \
# Modules/MContraction/WeakHamiltonian.hpp \
# Modules/MContraction/WeakHamiltonianEye.hpp \
@ -32,10 +35,8 @@ modules_hpp =\
# Modules/MScalar/ChargedProp.hpp \
# Modules/MScalar/FreeProp.hpp \
# Modules/MScalar/Scalar.hpp \
# Modules/MSink/Point.hpp \
# Modules/MSink/Smear.hpp \
# Modules/MSolver/RBPrecCG.hpp \
# Modules/MSource/Point.hpp \
# Modules/MSource/SeqConserved.hpp \
# Modules/MSource/SeqGamma.hpp \
# Modules/MSource/Wall.hpp \