1#ifndef AMPGEN_EVENTLIST_H
2#define AMPGEN_EVENTLIST_H
34 std::vector<Event> m_data = {};
36 std::map<std::string, unsigned int> m_extensions = {};
46 template <
class ... ARGS >
EventList(
const std::string& fname,
const ARGS&... args ) :
EventList()
50 template <
class ... ARGS >
EventList(
const std::vector<std::string>& fname,
const EventType& evtType,
const ARGS&... args ) :
EventList(evtType)
59 std::vector<Event>::reverse_iterator
rbegin() {
return m_data.rbegin(); }
60 std::vector<Event>::reverse_iterator
rend() {
return m_data.rend(); }
61 std::vector<Event>::iterator
begin() {
return m_data.begin(); }
62 std::vector<Event>::iterator
end() {
return m_data.end(); }
66 std::vector<Event>::const_iterator
begin()
const {
return m_data.cbegin(); }
67 std::vector<Event>::const_iterator
end()
const {
return m_data.cend(); }
70 const Event&
at(
const size_t& pos )
const {
return m_data[pos]; }
71 size_t size()
const {
return m_data.size(); }
73 size_t nBlocks()
const {
return m_data.size() ; }
75 double*
block(
const unsigned pos) {
return m_data[pos].address(); }
76 real_t weight(
const size_t& pos)
const {
return m_data[pos].weight(); }
77 real_t genPDF(
const size_t& pos)
const {
return m_data[pos].genPdf(); }
78 unsigned key(
const std::string&
key)
const
80 auto it = m_extensions.find(
key);
81 if( it == m_extensions.end() )
return m_data[0].size() - 1;
95 m_extensions[
key] = pos;
98 void setWeight(
const unsigned int& pos,
const double& w,
const double&g=+1)
100 m_data[pos].setWeight(w);
101 m_data[pos].setGenPdf(g);
103 void setGenPDF(
const unsigned int& pos,
const double& g)
105 m_data[pos].setGenPdf(g);
107 void erase(
const std::vector<Event>::iterator&
begin,
const std::vector<Event>::iterator&
end );
109 TTree*
tree(
const std::string& name,
const std::vector<std::string>& extraBranches = {} )
const;
118 size_t nBins = argPack.getArg<PlotOptions::Bins>(100);
123 template <
typename... ARGS> std::vector<TH1D*>
makeProjections(
const std::vector<Projection>& projections,
const ARGS&... args )
128 template <
typename... ARGS,
135 template <
typename... ARGS,
144 for (
auto& event : m_data )
fcn( event );
148 template <
typename functor>
void filter( functor&&
fcn )
150 unsigned currentSize =
size();
151 m_data.erase( std::remove_if( m_data.begin(), m_data.end(),
fcn ) , m_data.end() );
152 INFO(
"Filter retains " <<
size() <<
" / " << currentSize <<
" events");
155 template <
typename functor>
unsigned count( functor&&
fcn )
const
157 return std::count_if( std::begin(*
this), std::end(*
this),
fcn );
#define DECLARE_ARGUMENT(X, Y)
Container for a set of arguments Contains a set of arguments packed from a variadic constructor,...
Base class for compiled expressions, i.e.
Encapsulates the final state particles of a single event.
double * block(const unsigned pos)
EventList(const EventType &type)
std::vector< Event >::iterator begin()
EventList & transform(functor &&fcn)
const Event & at(const size_t &pos) const
void setEventType(const EventType &type)
EventList(const std::string &fname, const ARGS &... args)
real_t * getEvent(const size_t &index)
std::vector< TH1D * > makeProjections(const std::vector< Projection > &projections, const ARGS &... args)
TTree * tree(const std::string &name, const std::vector< std::string > &extraBranches={}) const
unsigned key(const std::string &key) const
std::vector< TH1D * > makeDefaultProjections(const ARGS &... args)
std::vector< TH1D * > makeProjections(const std::vector< Projection > &projections, const ArgumentPack &args)
void setGenPDF(const unsigned int &pos, const double &g)
void reserve(const size_t &size)
real_t weight(const size_t &pos) const
std::vector< Event >::reverse_iterator rbegin()
Event & operator[](const size_t &pos)
EventType eventType() const
TH1D * makeProjection(const Projection &projection, const ArgumentPack &args=ArgumentPack()) const
EventList(const std::vector< std::string > &fname, const EventType &evtType, const ARGS &... args)
std::vector< Event >::const_iterator end() const
void loadFromFile(const std::string &fname, const ArgumentPack &args)
TH2D * makeProjection(const Projection2D &projection, const ARGS &... args)
TH1D * makeProjection(const Projection &projection, const ARGS &... args) const
void emplace_back(const Event &evt)
void filter(functor &&fcn)
std::vector< Event >::const_iterator begin() const
size_t aligned_size() const
TH2D * makeProjection(const Projection2D &projection, const ArgumentPack &args=ArgumentPack()) const
void erase(const std::vector< Event >::iterator &begin, const std::vector< Event >::iterator &end)
std::vector< Event >::iterator end()
std::vector< Event >::reverse_iterator rend()
void loadFromTree(TTree *tree, const ArgumentPack &args)
EventList(const std::string &fname, const EventType &evtType, const ARGS &... args)
EventList(TTree *tree, const EventType &evtType, const ARGS &... args)
void add(const EventList &evts)
void extend(const std::string &key, unsigned pos)
unsigned count(functor &&fcn) const
void resize(const size_t &size)
const EventList & store() const
void setWeight(const unsigned int &pos, const double &w, const double &g=+1)
const Event & operator[](const size_t &pos) const
void push_back(const Event &evt)
const real_t * getEvent(const size_t &index) const
real_t genPDF(const size_t &pos) const
Deals with final state configuration of events, specifically dealing with the ordering of particles i...
std::vector< Projection > defaultProjections(const unsigned &nBins=100) const
#define INFO(X)
Used for printing information messages, and will always be printed.
typename detail::zeroType< args... >::type zeroType