PandaRoot
FairEvtFilterOnSingleParticleCounts.h
Go to the documentation of this file.
1 // -------------------------------------------------------------------------
2 // ----- FairEvtFilterOnSingleParticleCounts header file -----
3 // -------------------------------------------------------------------------
4 
32 #ifndef FairEvtFilterOnSingleParticleCounts_H_
33 #define FairEvtFilterOnSingleParticleCounts_H_
34 
35 #include "FairEvtFilter.h"
36 #include "TParticle.h"
37 #include "TMath.h"
38 #include <map>
39 #include <vector>
40 #include <iostream>
41 
42 typedef std::map<Int_t, std::vector<Int_t>> PdgGroupId;
43 typedef std::pair<Int_t, std::vector<Int_t>> PdgGroupIdPair;
44 typedef std::map<Int_t, std::vector<Int_t>>::iterator PdgGroupIdIterator;
45 
46 std::ostream &operator<<(std::ostream &os, const std::vector<Int_t> &v);
47 std::ostream &operator<<(std::ostream &os, const std::vector<std::pair<Double_t, Double_t>> &vpair);
48 std::ostream &operator<<(std::ostream &os, const std::vector<std::pair<Int_t, Int_t>> &vpair);
49 std::ostream &operator<<(std::ostream &os, const std::map<Int_t, std::vector<Int_t>> &PdgGroupId);
50 
52 
53  public:
54  // Default constructor
56 
57  // Constructor with name and title
58  FairEvtFilterOnSingleParticleCounts(const char *name, const char *title = "FairEvtFilterOnSingleParticleCounts");
59 
60  // Destructor
62 
64  // User interfaces -- Pdg Code Min and Max
66  // Use this for grouping up to 8 pdgCodes into 1 groupId
67  // all particles belonging to the groupId are regarded as being indistinguishable
68  // min defines how many particles you want in your events AT LEAST
69  // max defines how many particles you want in your events AT MOST
70  // the min and max numbers are used for all particles with one of the above pdgCodes
71  // returns kTRUE if the filter was added, otherwise returns kFALSE
72  Bool_t AndMinMaxPdgCodes(Int_t min, Int_t max, Int_t pdgCode1, Int_t pdgCode2 = kInvalidPdgCode, Int_t pdgCode3 = kInvalidPdgCode, Int_t pdgCode4 = kInvalidPdgCode,
73  Int_t pdgCode5 = kInvalidPdgCode, Int_t pdgCode6 = kInvalidPdgCode, Int_t pdgCode7 = kInvalidPdgCode, Int_t pdgCode8 = kInvalidPdgCode);
74 
75  // use if you want to group more than 8 pdgCodes into 1 groupId
76  // all particles belonging to the groupId are regarded as being indistinguishable
77  // min defines how many particles you want in your events AT LEAST
78  // max defines how many particles you want in your events AT MOST
79  // the min and max numbers are used for all particles with pdgCodes in the vector pdgCodes
80  // returns kTRUE if the filter was added, otherwise returns kFALSE
81  // all entries containing invalidPdgCode in the pdgCodes vector are skipped (as there is no pdgCode invalidPdgCode)
82  Bool_t AndMinMaxPdgCodes(Int_t min, Int_t max, std::vector<Int_t> &pdgCodes);
83 
85  // User interfaces -- Pdg Code Min only
87  Bool_t AndMinPdgCodes(Int_t min, Int_t pdgCode1, Int_t pdgCode2 = kInvalidPdgCode, Int_t pdgCode3 = kInvalidPdgCode, Int_t pdgCode4 = kInvalidPdgCode,
88  Int_t pdgCode5 = kInvalidPdgCode, Int_t pdgCode6 = kInvalidPdgCode, Int_t pdgCode7 = kInvalidPdgCode, Int_t pdgCode8 = kInvalidPdgCode)
89  {
90  return AndMinMaxPdgCodes(min, 9999, pdgCode1, pdgCode2, pdgCode3, pdgCode4, pdgCode5, pdgCode6, pdgCode7, pdgCode8);
91  };
92  Bool_t AndMinPdgCodes(Int_t min, std::vector<Int_t> &pdgCodes) { return AndMinMaxPdgCodes(min, 9999, pdgCodes); };
93 
95  // User interfaces -- Pdg Code Max only
97  Bool_t AndMaxPdgCodes(Int_t max, std::vector<Int_t> &pdgCodes) { return AndMinMaxPdgCodes(0, max, pdgCodes); };
98  Bool_t AndMaxPdgCodes(Int_t max, Int_t pdgCode1, Int_t pdgCode2 = kInvalidPdgCode, Int_t pdgCode3 = kInvalidPdgCode, Int_t pdgCode4 = kInvalidPdgCode,
99  Int_t pdgCode5 = kInvalidPdgCode, Int_t pdgCode6 = kInvalidPdgCode, Int_t pdgCode7 = kInvalidPdgCode, Int_t pdgCode8 = kInvalidPdgCode)
100  {
101  return AndMinMaxPdgCodes(0, max, pdgCode1, pdgCode2, pdgCode3, pdgCode4, pdgCode5, pdgCode6, pdgCode7, pdgCode8);
102  };
103 
105  // User interfaces -- Charge
107  // min defines how many particles you want in your events AT LEAST
108  // max defines how many particles you want in your events AT MOST
109  // Use for setting multiplicity constraints based on particle charges
110  // ChargeState
111  // kNeutral - neutral particles
112  // kPlus - positively charged particles
113  // kMinus - negatively charged particles
114  // kCharged - positively or negatively charged particles
115  // kAll - all particles (neutral or charged)
116  Bool_t AndMinMaxCharge(Int_t min, Int_t max, ChargeState charge);
117 
118  // use the methods above with fixed minimum or maximum number of particles
119  Bool_t AndMinCharge(Int_t min, ChargeState charge) { return AndMinMaxCharge(min, 9999, charge); };
120  Bool_t AndMaxCharge(Int_t max, ChargeState charge) { return AndMinMaxCharge(0, max, charge); };
121 
123  // User interfaces -- All particles in event
125  // Use if you want to put constraints on the number of particles in the event
126  // Functionality is implemented in AndMinMaxCharge
127  Bool_t AndMinMaxAllParticles(Int_t min, Int_t max) { return AndMinMaxCharge(min, max, kAll); };
128 
129  // use the methods above with fixed minimum or maximum number of particles
130  Bool_t AndMinAllParticles(Int_t min) { return AndMinMaxCharge(min, 9999, kAll); };
131  Bool_t AndMaxAllParticles(Int_t max) { return AndMinMaxCharge(0, max, kAll); };
132 
134  // User interfaces -- Momentum constraints
136  // use if you want to consider particles which fulfill certain momentum constraints
137  // could be used in combination with any other filter, if you want to set momentum limits for several particles
138  Bool_t AndPRange(Double_t min, Double_t max) { return AndMinMaxMom(min, max, kMomTotal); };
139  Bool_t AndPtRange(Double_t min, Double_t max) { return AndMinMaxMom(min, max, kMomTrans); };
140  Bool_t AndPzRange(Double_t min, Double_t max) { return AndMinMaxMom(min, max, kMomZ); };
141 
143  // User interfaces -- Angular constraints
145  // use if you want to consider particles in certain angle or vertex ranges
146  // could be used in combination with any other filter, if you want to set geometric limits for several particles
147  Bool_t AndThetaRange(Double_t min, Double_t max) { return AndMinMaxGeom(min, max, kTheta); };
148  Bool_t AndPhiRange(Double_t min, Double_t max) { return AndMinMaxGeom(min, max, kPhi); };
149 
151  // User interfaces -- Vertex constraints
153  // use if you want to consider particles in certain angle or vertex ranges
154  // could be used in combination with any other filter, if you want to set geometric limits for several particles
155  Bool_t AndVzRange(Double_t min, Double_t max) { return AndMinMaxGeom(min, max, kVertexZ); };
156  Bool_t AndRhoRange(Double_t min, Double_t max) { return AndMinMaxGeom(min, max, kVertexRho); };
157  Bool_t AndRadiusRange(Double_t min, Double_t max) { return AndMinMaxGeom(min, max, kVertexRadius); };
158 
159  // checks if the particles in the event (fParticleList) suit the filter settings
160  // kTRUE if event matches, kFALSE otherwise
161  Bool_t EventMatches(Int_t evtNr);
162 
163  // returns kTRUE if the filter on pdg codes or charge states (or all particles) is turned on
164  // Momentum and angular constraints do not play any role here
165  Bool_t FilterActive() { return (fFilterPdg || fFilterCharge); };
166 
167  private:
168  // in these methods the functions of the user interfaces from above are realized
169  // all user interfaces are only calling this two methods with different arguments
170  Bool_t AndMinMaxMom(Double_t min, Double_t max, MomState mom);
171  Bool_t AndMinMaxGeom(Double_t min, Double_t max, GeomState geom);
172 
173  // initialize the counters for pdg codes or charge states
174  void InitCounters();
175 
176  // Set default boundaries for momenta, angles and vertices
177  // that are always fulfilled
178  // this implements the default case of not filtering on momenta, angles or vertices
179  void SetDefaultBoundaries();
180 
181  // count up the appropriate element of fCountCharge
182  void CountCharge(TParticle *particle);
183  // count up the appropriate element of fCountGroupId
184  void CountPdg(TParticle *particle);
185 
186  // check whether the particle fulfills the momentum constraints
187  Bool_t AcceptMomentum(TParticle *particle);
188  // check whether the particle fulfills the geometry constraints
189  Bool_t AcceptGeometry(TParticle *particle);
190 
191  // check whether the event satisfies the constraints on pdg multiplicities (fCountGroupId vector has to be filled already)
192  Bool_t AcceptPdgCounter();
193  // check whether the event satisfies the constraints on charge state multiplicities (fCountCharge vector has to be filled already)
194  Bool_t AcceptChargeCounter();
195 
196  // first input = pdg code (as used by event generators)
197  // second input = "vector of group IDs that the pdg code is assigned to"
198  // one pdg code can be assigned to mutliple group IDs
199  // all pdg codes within the same group ID are treated as indistinguishable
200  // group IDs are automatically assigned; they are not known to the user and have to start with 0 and be consecutive as they are used as the index for fPartCounts
201  PdgGroupId fPdgGroupId;
202 
203  // for counting the multiplicities based on group IDs
204  std::vector<Int_t> fCountGroupId;
205  // for counts the multiplicities based on electrical charge
206  std::vector<Int_t> fCountCharge;
207 
208  // vector defining the minimum and maximum acceptable multiplicities for each group ID
209  std::vector<std::pair<Int_t, Int_t>> fGroupIdCountsMinMax;
210  // vector containing the minimum and maximum acceptable multiplicities for each electrical charge state
211  // check ChargeState in FairEvtFilter.h for translating the index to an electrical charge
212  std::vector<std::pair<Int_t, Int_t>> fChargeCountsMinMax;
213 
214  // only particles with a momentum within the min and max value will be counted
215  // check MomState in FairEvtFilter.h for translating the index to an electrical charge
216  std::vector<std::pair<Double_t, Double_t>> fMomMinMax;
217  // only particles with angles and vertices within the min and max value will be counted
218  // check GeomState in FairEvtFilter.h for translating the index to an electrical charge
219  std::vector<std::pair<Double_t, Double_t>> fGeomMinMax;
220 
221  Bool_t fFilterPdg;
222  Bool_t fFilterCharge;
223  Bool_t fFilterMom;
224  Bool_t fFilterGeom;
225 
228  static const Int_t kInvalidPdgCode = 0;
229 
231 };
232 
233 #endif /* FairEvtFilterOnSingleParticleCounts_H_ */
Bool_t AndMaxPdgCodes(Int_t max, std::vector< Int_t > &pdgCodes)
Bool_t AndMinMaxPdgCodes(Int_t min, Int_t max, Int_t pdgCode1, Int_t pdgCode2=kInvalidPdgCode, Int_t pdgCode3=kInvalidPdgCode, Int_t pdgCode4=kInvalidPdgCode, Int_t pdgCode5=kInvalidPdgCode, Int_t pdgCode6=kInvalidPdgCode, Int_t pdgCode7=kInvalidPdgCode, Int_t pdgCode8=kInvalidPdgCode)
std::pair< Int_t, std::vector< Int_t > > PdgGroupIdPair
friend F32vec4 max(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:25
Bool_t AndMinPdgCodes(Int_t min, Int_t pdgCode1, Int_t pdgCode2=kInvalidPdgCode, Int_t pdgCode3=kInvalidPdgCode, Int_t pdgCode4=kInvalidPdgCode, Int_t pdgCode5=kInvalidPdgCode, Int_t pdgCode6=kInvalidPdgCode, Int_t pdgCode7=kInvalidPdgCode, Int_t pdgCode8=kInvalidPdgCode)
__m128 v
Definition: P4_F32vec4.h:3
Bool_t AndMaxPdgCodes(Int_t max, Int_t pdgCode1, Int_t pdgCode2=kInvalidPdgCode, Int_t pdgCode3=kInvalidPdgCode, Int_t pdgCode4=kInvalidPdgCode, Int_t pdgCode5=kInvalidPdgCode, Int_t pdgCode6=kInvalidPdgCode, Int_t pdgCode7=kInvalidPdgCode, Int_t pdgCode8=kInvalidPdgCode)
std::map< Int_t, std::vector< Int_t > >::iterator PdgGroupIdIterator
Bool_t AndMinCharge(Int_t min, ChargeState charge)
Bool_t AndMinMaxCharge(Int_t min, Int_t max, ChargeState charge)
Bool_t AndMaxCharge(Int_t max, ChargeState charge)
Bool_t AndMinPdgCodes(Int_t min, std::vector< Int_t > &pdgCodes)
friend F32vec4 min(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:24
std::map< Int_t, std::vector< Int_t > > PdgGroupId