AmpGen 2.1
Loading...
Searching...
No Matches
EventType.h
Go to the documentation of this file.
1#ifndef AMPGEN_EVENTTYPE_H
2#define AMPGEN_EVENTTYPE_H
3
4#include <functional>
5#include <iostream>
6#include <map>
7#include <vector>
8#include <initializer_list>
9#include "AmpGen/Event.h"
10
11namespace AmpGen
12{
13 class Projection;
18 class EventType;
19 std::ostream& operator<<( std::ostream& os, const EventType& type );
20
22 {
23 public:
25 EventType() = default;
26
29 EventType( const std::vector<std::string>&, const bool& isTD = false );
30
33 EventType conj( const bool& headOnly = 0, const bool& dontConjHead = 0 ) const;
34
36 std::map<std::string, unsigned> getEventFormat( const bool& outputNames = false ) const;
37
40 std::pair<unsigned,unsigned> count(const unsigned& index) const;
41 std::pair<double, double> minmax( const std::vector<unsigned>& indices) const;
42 std::vector<double> masses() const;
43 std::string mother() const;
44 double mass( const unsigned& index ) const;
45 double motherMass() const;
46 std::vector<std::string> finalStates() const;
47 bool isTimeDependent() const;
48 unsigned eventSize() const;
49 unsigned size() const;
50 unsigned dof() const;
51 std::string operator[]( const unsigned& index ) const;
52 std::string decayDescriptor() const;
53 std::string label( const unsigned& index, bool isRoot = true ) const;
54 std::string label( const std::vector<unsigned>& index, bool isRoot = true ) const;
55 std::vector<Projection> defaultProjections(const unsigned& nBins=100) const;
56 Projection projection(const unsigned& nBins, const std::vector<unsigned>& indices, const std::string& observable = "mass2") const;
57
58 bool operator==( const EventType& other ) const;
59 bool has( const std::string& name ) const;
60
61 void extendEventType( const std::string& branch );
62
64 std::function<void( Event& )> symmetriser() const;
65 std::function<bool( Event&, const std::vector<int>& ids)> automaticOrdering() const;
66
68 std::pair<unsigned, unsigned> dim() const;
69
70 friend std::ostream& AmpGen::operator<<( std::ostream& os, const EventType& type );
71
72 Event makeEvent() const;
73
74 std::string constructor_string() const;
75 private:
76 std::string m_mother;
77 double m_motherMass;
78 std::vector<std::string> m_particleNames;
79 std::vector<std::string> m_particleNamesPickled;
80 std::vector<double> m_particleMasses;
81 std::vector<bool> m_ignore;
82 bool m_timeDependent;
83 std::vector<std::string> m_eventTypeExtensions;
84 std::pair<unsigned, unsigned> m_dim;
85 bool m_alt_part_names;
86 };
87
88
89 extern "C" unsigned python__EventType__dim( const char* eventType );
90} // namespace AmpGen
91
92#endif
Encapsulates the final state particles of a single event.
Definition Event.h:18
Deals with final state configuration of events, specifically dealing with the ordering of particles i...
Definition EventType.h:22
bool has(const std::string &name) const
std::pair< double, double > minmax(const std::vector< unsigned > &indices) const
EventType()=default
Default constructor.
unsigned eventSize() const
std::string label(const unsigned &index, bool isRoot=true) const
void extendEventType(const std::string &branch)
std::string mother() const
unsigned dof() const
std::vector< std::string > finalStates() const
bool operator==(const EventType &other) const
std::pair< unsigned, unsigned > dim() const
Calculates the number of spin indices associated with the initial and final state,...
std::function< void(Event &)> symmetriser() const
Functor to randomly symmetrise data of this event type, using the Fisher-Yates shuffle.
std::string decayDescriptor() const
std::string operator[](const unsigned &index) const
std::vector< Projection > defaultProjections(const unsigned &nBins=100) const
std::pair< unsigned, unsigned > count(const unsigned &index) const
Counts the number of particles in this event type with the same name as the index'th name.
std::string label(const std::vector< unsigned > &index, bool isRoot=true) const
std::string constructor_string() const
Event makeEvent() const
std::map< std::string, unsigned > getEventFormat(const bool &outputNames=false) const
Returns the event format, that matches between expressions and the names used in Particle.
Projection projection(const unsigned &nBins, const std::vector< unsigned > &indices, const std::string &observable="mass2") const
std::function< bool(Event &, const std::vector< int > &ids)> automaticOrdering() const
unsigned size() const
double mass(const unsigned &index) const
EventType(const std::vector< std::string > &, const bool &isTD=false)
Takes a list of particles, beginning with the head of the decay and then then final state particles,...
EventType conj(const bool &headOnly=0, const bool &dontConjHead=0) const
Returns the CP-conjugated event type.
std::vector< double > masses() const
bool isTimeDependent() const
double motherMass() const
std::ostream & operator<<(std::ostream &os, const CompiledExpressionBase &expression)
unsigned python__EventType__dim(const char *eventType)