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