AmpGen 2.1
Loading...
Searching...
No Matches
EventListSIMD.h
Go to the documentation of this file.
1#ifndef AMPGEN_EVENTLIST2_H
2#define AMPGEN_EVENTLIST2_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/EventList.h"
12#include <chrono>
13#include <functional>
14#include <numeric>
15#include <cstddef>
16#include <algorithm>
17
18#include <TH1D.h>
19#include <TH2D.h>
20#include <TTree.h>
21
22#ifdef _OPENMP
23 #include <omp.h>
24#endif
25
27#include "AmpGen/simd/utils.h"
28#include "AmpGen/Store.h"
29
30namespace AmpGen
31{
34 {
35 private:
37 std::vector<real_v> m_weights {};
38 std::vector<real_v> m_genPDF {};
39 EventType m_eventType {};
40 public:
42 EventListSIMD() = default;
43 EventListSIMD( const EventType& type );
44 template < class ... ARGS > EventListSIMD( const std::string& fname, const EventType& evtType, const ARGS&... args ) : EventListSIMD(evtType)
45 {
46 loadFromFile( fname, ArgumentPack(args...) );
47 }
48 template < class ... ARGS > EventListSIMD( const std::string& fname, const ARGS&... args ) : EventListSIMD()
49 {
50 loadFromFile( fname, ArgumentPack(args...) );
51 }
52 template < class ... ARGS > EventListSIMD( const std::vector<std::string>& fname, const EventType& evtType, const ARGS&... args ) : EventListSIMD(evtType)
53 {
54 for( auto& f : fname ) loadFromFile( f, ArgumentPack(args...) );
55 }
56 template < class ... ARGS > EventListSIMD( TTree* tree, const EventType& evtType, const ARGS&... args ) : EventListSIMD(evtType)
57 {
58 loadFromTree( tree, ArgumentPack(args...) );
59 }
60 EventListSIMD( const EventList& other );
61 const real_v* data() const { return m_data.data(); }
62 operator Store<real_v, Alignment::AoS> () const { return m_data ; }
63 const auto& store() const { return m_data; }
64 const Event at(const unsigned& p) const { return EventListSIMD::operator[](p) ; }
65 const real_v* block(const unsigned& p) const { return m_data.data() + p * m_data.nFields(); }
66 real_v* block(const unsigned& p) { return m_data.data() + p * m_data.nFields(); }
67 real_v weight(const unsigned& p) const { return m_weights[p]; }
68 real_v genPDF(const unsigned& p) const { return m_genPDF [p]; }
69 const auto nFields() const { return m_data.nFields(); }
70 void setWeight( const unsigned& block, const real_v& w, const real_v& g=1.f)
71 {
72 m_weights[block] = w;
73 m_genPDF[block] = g;
74 }
75 void setGenPDF( const unsigned& block, const real_v& g)
76 {
77 m_genPDF[block] = g;
78 }
79 void resize( const unsigned nEvents )
80 {
81 m_data = Store<real_v, Alignment::AoS>(nEvents, m_eventType.eventSize());
82 m_weights.resize( aligned_size(), 1.f);
83 m_genPDF.resize( aligned_size(), 1.f );
84 }
85 const Event operator[]( const size_t&) const;
86 std::array<Event, utils::size<real_v>::value> scatter(unsigned) const;
87 void gather(const std::array<Event, utils::size<real_v>::value>&, unsigned);
89 auto end() const { return make_scatter_iterator<utils::size<real_v>::value>(size(), (const EventListSIMD*)(nullptr) ); }
92 EventType eventType() const { return m_eventType; }
93 size_t aligned_size() const { return m_data.aligned_size(); }
94 double integral() const;
95 size_t eventSize() const { return m_data.nFields(); }
96 size_t size() const { return m_data.size(); }
97 size_t nBlocks() const { return m_data.nBlocks(); }
98 void setEventType( const EventType& type ) { m_eventType = type; }
99 void add( const EventListSIMD& evts );
100 void loadFromTree( TTree* tree, const ArgumentPack& args );
101 void loadFromFile( const std::string& fname, const ArgumentPack& args );
102 void clear();
103
104 TTree* tree( const std::string& name, const std::vector<std::string>& extraBranches = {} ) const;
105
106 TH1D* makeProjection(const Projection& projection , const ArgumentPack& args = ArgumentPack()) const;
107 TH2D* makeProjection(const Projection2D& projection, const ArgumentPack& args = ArgumentPack()) const;
108 std::vector<TH1D*> makeProjections( const std::vector<Projection>& projections, const ArgumentPack& args );
109
110 template <class... ARGS> std::vector<TH1D*> makeDefaultProjections( const ARGS&... args )
111 {
112 auto argPack = ArgumentPack( args... );
113 size_t nBins = argPack.getArg<PlotOptions::Bins>(100);
114 auto proj = eventType().defaultProjections(nBins);
115 return makeProjections( proj , argPack );
116 }
117
118 template <typename... ARGS> std::vector<TH1D*> makeProjections( const std::vector<Projection>& projections, const ARGS&... args )
119 {
120 return makeProjections( projections, ArgumentPack( args... ) );
121 }
122
123 template <typename... ARGS,
124 typename = std::enable_if_t< ! std::is_same<zeroType<ARGS...>, ArgumentPack>::value > >
125 TH1D* makeProjection( const Projection& projection, const ARGS&... args ) const
126 {
127 return makeProjection( projection, ArgumentPack(args...) );
128 }
129
130 template <typename... ARGS,
131 typename = std::enable_if_t< ! std::is_same<zeroType<ARGS...>, ArgumentPack>::value > >
132 TH2D* makeProjection( const Projection2D& projection, const ARGS&... args )
133 {
134 return makeProjection( projection, ArgumentPack(args...) );
135 }
136
137 template <typename functor> EventListSIMD& transform( functor&& fcn )
138 {
139 for ( auto& event : *this ) fcn( event );
140 return *this;
141 }
142 static std::vector<real_v> makeEvent( const Event& event )
143 {
144 std::vector<real_v> rt( event.size() );
145 for( unsigned i = 0 ; i != event.size(); ++i ) rt[i] = event[i];
146 return rt;
147 }
148 };
149
150} // namespace AmpGen
151#endif
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
unsigned size() const
Definition Event.h:30
void resize(const unsigned nEvents)
EventListSIMD & transform(functor &&fcn)
std::vector< TH1D * > makeProjections(const std::vector< Projection > &projections, const ARGS &... args)
const Event at(const unsigned &p) const
EventListSIMD(const EventList &other)
EventListSIMD(const std::string &fname, const EventType &evtType, const ARGS &... args)
const Event operator[](const size_t &) const
EventListSIMD(TTree *tree, const EventType &evtType, const ARGS &... args)
const auto & store() const
real_v genPDF(const unsigned &p) const
void gather(const std::array< Event, utils::size< real_v >::value > &, unsigned)
TH2D * makeProjection(const Projection2D &projection, const ARGS &... args)
TH2D * makeProjection(const Projection2D &projection, const ArgumentPack &args=ArgumentPack()) const
EventType eventType() const
double integral() const
TH1D * makeProjection(const Projection &projection, const ARGS &... args) const
EventListSIMD(const std::string &fname, const ARGS &... args)
TH1D * makeProjection(const Projection &projection, const ArgumentPack &args=ArgumentPack()) const
size_t eventSize() const
void loadFromFile(const std::string &fname, const ArgumentPack &args)
size_t size() const
std::array< Event, utils::size< real_v >::value > scatter(unsigned) const
EventListSIMD(const std::vector< std::string > &fname, const EventType &evtType, const ARGS &... args)
size_t nBlocks() const
std::vector< TH1D * > makeProjections(const std::vector< Projection > &projections, const ArgumentPack &args)
void add(const EventListSIMD &evts)
real_v weight(const unsigned &p) const
EventListSIMD(const EventType &type)
size_t aligned_size() const
void loadFromTree(TTree *tree, const ArgumentPack &args)
void setEventType(const EventType &type)
const real_v * data() const
static std::vector< real_v > makeEvent(const Event &event)
TTree * tree(const std::string &name, const std::vector< std::string > &extraBranches={}) const
std::vector< TH1D * > makeDefaultProjections(const ARGS &... args)
const auto nFields() const
void setGenPDF(const unsigned &block, const real_v &g)
void setWeight(const unsigned &block, const real_v &w, const real_v &g=1.f)
const real_v * block(const unsigned &p) const
real_v * block(const unsigned &p)
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
AVX::real_v real_v
Definition utils.h:46
typename detail::zeroType< args... >::type zeroType
Definition MetaUtils.h:35
auto make_scatter_iterator(const unsigned &pos, store_type *store)
Definition iterator.h:43
static constexpr unsigned value
Definition utils.h:51