AmpGen 2.1
Loading...
Searching...
No Matches
EventList.h
Go to the documentation of this file.
1#ifndef AMPGEN_EVENTLIST_H
2#define AMPGEN_EVENTLIST_H
3
5#include "AmpGen/EventType.h"
6#include "AmpGen/MsgService.h"
7#include "AmpGen/Event.h"
8#include "AmpGen/Projection.h"
9#include "AmpGen/Utilities.h"
10#include "AmpGen/MetaUtils.h"
11#include "AmpGen/Units.h"
12
13#include <chrono>
14#include <functional>
15#include <numeric>
16#include <cstddef>
17#include <algorithm>
18
19#include <TH1D.h>
20#include <TH2D.h>
21#include <TTree.h>
22
23#ifdef _OPENMP
24 #include <omp.h>
25#endif
26
27namespace AmpGen
28{
29 namespace PlotOptions { DECLARE_ARGUMENT(Bins, size_t); }
32 {
33 private:
34 std::vector<Event> m_data = {};
35 EventType m_eventType = {};
36 std::map<std::string, unsigned int> m_extensions = {};
37
38 public:
40 EventList() = default;
41 EventList( const EventType& type );
42 template < class ... ARGS > EventList( const std::string& fname, const EventType& evtType, const ARGS&... args ) : EventList(evtType)
43 {
44 loadFromFile( fname, ArgumentPack(args...) );
45 }
46 template < class ... ARGS > EventList( const std::string& fname, const ARGS&... args ) : EventList()
47 {
48 loadFromFile( fname, ArgumentPack(args...) );
49 }
50 template < class ... ARGS > EventList( const std::vector<std::string>& fname, const EventType& evtType, const ARGS&... args ) : EventList(evtType)
51 {
52 for( auto& f : fname ) loadFromFile( f, ArgumentPack(args...) );
53 }
54 template < class ... ARGS > EventList( TTree* tree, const EventType& evtType, const ARGS&... args ) : EventList(evtType)
55 {
56 loadFromTree( tree, ArgumentPack(args...) );
57 }
58 const EventList& store() const { return *this;}
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(); }
63 Event& operator[]( const size_t& pos ) { return m_data[pos]; }
64 real_t* getEvent( const size_t& index ) { return (( *this )[index] ); }
65 const real_t* getEvent( const size_t& index ) const { return (const real_t*)(( *this )[index] ); }
66 std::vector<Event>::const_iterator begin() const { return m_data.cbegin(); }
67 std::vector<Event>::const_iterator end() const { return m_data.cend(); }
68 const Event& operator[]( const size_t& pos ) const { return m_data[pos]; }
69 EventType eventType() const { return m_eventType; }
70 const Event& at( const size_t& pos ) const { return m_data[pos]; }
71 size_t size() const { return m_data.size(); }
72 size_t aligned_size() const { return m_data.size() ; }
73 size_t nBlocks() const { return m_data.size() ; }
74 double integral() const;
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
79 {
80 auto it = m_extensions.find(key);
81 if( it == m_extensions.end() ) return m_data[0].size() - 1;
82 return it->second;
83 }
84 void reserve( const size_t& size );
85 void resize ( const size_t& size );
86 void push_back( const Event& evt );
87 void emplace_back( const Event& evt);
88 void setEventType( const EventType& type ) { m_eventType = type; }
89 void add( const EventList& evts );
90 void loadFromTree( TTree* tree, const ArgumentPack& args );
91 void loadFromFile( const std::string& fname, const ArgumentPack& args );
92 void clear();
93 void extend(const std::string& key, unsigned pos)
94 {
95 m_extensions[key] = pos;
96 }
97
98 void setWeight( const unsigned int& pos, const double& w, const double&g=+1)
99 {
100 m_data[pos].setWeight(w);
101 m_data[pos].setGenPdf(g);
102 }
103 void setGenPDF( const unsigned int& pos, const double& g)
104 {
105 m_data[pos].setGenPdf(g);
106 }
107 void erase( const std::vector<Event>::iterator& begin, const std::vector<Event>::iterator& end );
108
109 TTree* tree( const std::string& name, const std::vector<std::string>& extraBranches = {} ) const;
110
111 TH1D* makeProjection(const Projection& projection , const ArgumentPack& args = ArgumentPack()) const;
112 TH2D* makeProjection(const Projection2D& projection, const ArgumentPack& args = ArgumentPack()) const;
113 std::vector<TH1D*> makeProjections( const std::vector<Projection>& projections, const ArgumentPack& args );
114
115 template <class... ARGS> std::vector<TH1D*> makeDefaultProjections( const ARGS&... args )
116 {
117 auto argPack = ArgumentPack( args... );
118 size_t nBins = argPack.getArg<PlotOptions::Bins>(100);
119 auto proj = eventType().defaultProjections(nBins);
120 return makeProjections( proj , argPack );
121 }
122
123 template <typename... ARGS> std::vector<TH1D*> makeProjections( const std::vector<Projection>& projections, const ARGS&... args )
124 {
125 return makeProjections( projections, ArgumentPack( args... ) );
126 }
127
128 template <typename... ARGS,
129 typename = std::enable_if_t< ! std::is_same<zeroType<ARGS...>, ArgumentPack>::value > >
130 TH1D* makeProjection( const Projection& projection, const ARGS&... args ) const
131 {
132 return makeProjection( projection, ArgumentPack(args...) );
133 }
134
135 template <typename... ARGS,
136 typename = std::enable_if_t< ! std::is_same<zeroType<ARGS...>, ArgumentPack>::value > >
137 TH2D* makeProjection( const Projection2D& projection, const ARGS&... args )
138 {
139 return makeProjection( projection, ArgumentPack(args...) );
140 }
141
142 template <typename functor> EventList& transform( functor&& fcn )
143 {
144 for ( auto& event : m_data ) fcn( event );
145 return *this;
146 }
147
148 template <typename functor> void filter( functor&& fcn )
149 {
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");
153 }
154
155 template <typename functor> unsigned count( functor&& fcn ) const
156 {
157 return std::count_if( std::begin(*this), std::end(*this), fcn );
158 }
159 };
160 DECLARE_ARGUMENT(Branches, std::vector<std::string>);
161 DECLARE_ARGUMENT(ExtraBranches, std::vector<std::string>);
162 DECLARE_ARGUMENT(IdBranches, std::vector<std::string>);
163 DECLARE_ARGUMENT(EntryList, std::vector<size_t>);
164 DECLARE_ARGUMENT(GetGenPdf, bool);
165 DECLARE_ARGUMENT(Filter, std::string);
166 DECLARE_ARGUMENT(WeightBranch, std::string);
167 DECLARE_ARGUMENT(ApplySym, bool);
168 DECLARE_ARGUMENT(WeightFunction, std::function<double( const Event& )>);
169 DECLARE_ARGUMENT(InputUnits, AmpGen::Units);
170} // namespace AmpGen
171
172#endif
#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.
Definition Event.h:18
double * block(const unsigned pos)
Definition EventList.h:75
EventList(const EventType &type)
std::vector< Event >::iterator begin()
Definition EventList.h:61
EventList & transform(functor &&fcn)
Definition EventList.h:142
const Event & at(const size_t &pos) const
Definition EventList.h:70
void setEventType(const EventType &type)
Definition EventList.h:88
EventList(const std::string &fname, const ARGS &... args)
Definition EventList.h:46
real_t * getEvent(const size_t &index)
Definition EventList.h:64
std::vector< TH1D * > makeProjections(const std::vector< Projection > &projections, const ARGS &... args)
Definition EventList.h:123
TTree * tree(const std::string &name, const std::vector< std::string > &extraBranches={}) const
unsigned key(const std::string &key) const
Definition EventList.h:78
std::vector< TH1D * > makeDefaultProjections(const ARGS &... args)
Definition EventList.h:115
std::vector< TH1D * > makeProjections(const std::vector< Projection > &projections, const ArgumentPack &args)
void setGenPDF(const unsigned int &pos, const double &g)
Definition EventList.h:103
void reserve(const size_t &size)
real_t weight(const size_t &pos) const
Definition EventList.h:76
std::vector< Event >::reverse_iterator rbegin()
Definition EventList.h:59
Event & operator[](const size_t &pos)
Definition EventList.h:63
EventType eventType() const
Definition EventList.h:69
TH1D * makeProjection(const Projection &projection, const ArgumentPack &args=ArgumentPack()) const
EventList(const std::vector< std::string > &fname, const EventType &evtType, const ARGS &... args)
Definition EventList.h:50
std::vector< Event >::const_iterator end() const
Definition EventList.h:67
double integral() const
void loadFromFile(const std::string &fname, const ArgumentPack &args)
TH2D * makeProjection(const Projection2D &projection, const ARGS &... args)
Definition EventList.h:137
TH1D * makeProjection(const Projection &projection, const ARGS &... args) const
Definition EventList.h:130
void emplace_back(const Event &evt)
void filter(functor &&fcn)
Definition EventList.h:148
std::vector< Event >::const_iterator begin() const
Definition EventList.h:66
size_t aligned_size() const
Definition EventList.h:72
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()
Definition EventList.h:62
std::vector< Event >::reverse_iterator rend()
Definition EventList.h:60
void loadFromTree(TTree *tree, const ArgumentPack &args)
size_t nBlocks() const
Definition EventList.h:73
EventList(const std::string &fname, const EventType &evtType, const ARGS &... args)
Definition EventList.h:42
EventList(TTree *tree, const EventType &evtType, const ARGS &... args)
Definition EventList.h:54
void add(const EventList &evts)
void extend(const std::string &key, unsigned pos)
Definition EventList.h:93
unsigned count(functor &&fcn) const
Definition EventList.h:155
void resize(const size_t &size)
size_t size() const
Definition EventList.h:71
EventList()=default
const EventList & store() const
Definition EventList.h:58
void setWeight(const unsigned int &pos, const double &w, const double &g=+1)
Definition EventList.h:98
const Event & operator[](const size_t &pos) const
Definition EventList.h:68
void push_back(const Event &evt)
const real_t * getEvent(const size_t &index) const
Definition EventList.h:65
real_t genPDF(const size_t &pos) const
Definition EventList.h:77
Deals with final state configuration of events, specifically dealing with the ordering of particles i...
Definition EventType.h:22
std::vector< Projection > defaultProjections(const unsigned &nBins=100) const
#define INFO(X)
Used for printing information messages, and will always be printed.
Definition MsgService.h:75
double real_t
Definition Types.h:6
typename detail::zeroType< args... >::type zeroType
Definition MetaUtils.h:35