FairRoot/PandaRoot
FairEvtFilterOnSingleParticleCounts.h
Go to the documentation of this file.
1 // -------------------------------------------------------------------------
2 // ----- FairEvtFilterOnSingleParticleCounts header file -----
3 // -------------------------------------------------------------------------
4 
5 
35 #ifndef FairEvtFilterOnSingleParticleCounts_H_
36 #define FairEvtFilterOnSingleParticleCounts_H_
37 
38 #include "FairEvtFilter.h"
39 #include "TParticle.h"
40 #include "TMath.h"
41 #include<map>
42 #include<vector>
43 #include<iostream>
44 
45 
46 typedef std::map<Int_t, std::vector<Int_t> > PdgGroupId;
47 typedef std::pair<Int_t,std::vector<Int_t> > PdgGroupIdPair;
48 typedef std::map<Int_t,std::vector<Int_t> >::iterator PdgGroupIdIterator;
49 
50 std::ostream& operator <<(std::ostream& os, const std::vector<Int_t>& v);
51 std::ostream& operator <<(std::ostream& os, const std::vector< std::pair<Double_t, Double_t> >& vpair);
52 std::ostream& operator <<(std::ostream& os, const std::vector< std::pair<Int_t, Int_t> >& vpair);
53 std::ostream& operator <<(std::ostream& os, const std::map<Int_t, std::vector<Int_t> >& PdgGroupId);
54 
55 
57 {
58 
59 public:
60 
61  //Default constructor
63 
64  //Constructor with name and title
65  FairEvtFilterOnSingleParticleCounts(const char* name, const char* title="FairEvtFilterOnSingleParticleCounts");
66 
67  // Destructor
69 
71  // User interfaces -- Pdg Code Min and Max
73  // Use this for grouping up to 8 pdgCodes into 1 groupId
74  // all particles belonging to the groupId are regarded as being indistinguishable
75  // min defines how many particles you want in your events AT LEAST
76  // max defines how many particles you want in your events AT MOST
77  // the min and max numbers are used for all particles with one of the above pdgCodes
78  // returns kTRUE if the filter was added, otherwise returns kFALSE
79  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 );
80 
81 
82  // use if you want to group more than 8 pdgCodes into 1 groupId
83  // all particles belonging to the groupId are regarded as being indistinguishable
84  // min defines how many particles you want in your events AT LEAST
85  // max defines how many particles you want in your events AT MOST
86  // the min and max numbers are used for all particles with pdgCodes in the vector pdgCodes
87  // returns kTRUE if the filter was added, otherwise returns kFALSE
88  // all entries containing invalidPdgCode in the pdgCodes vector are skipped (as there is no pdgCode invalidPdgCode)
89  Bool_t AndMinMaxPdgCodes( Int_t min, Int_t max, std::vector<Int_t> &pdgCodes );
90 
91 
93  // User interfaces -- Pdg Code Min only
95  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){
96  return AndMinMaxPdgCodes( min, 9999, pdgCode1, pdgCode2, pdgCode3, pdgCode4, pdgCode5, pdgCode6, pdgCode7, pdgCode8 );
97  };
98  Bool_t AndMinPdgCodes(Int_t min, std::vector<Int_t> &pdgCodes ){
99  return AndMinMaxPdgCodes( min, 9999, pdgCodes );
100  };
101 
102 
104  // User interfaces -- Pdg Code Max only
106  Bool_t AndMaxPdgCodes(Int_t max, std::vector<Int_t> &pdgCodes ){
107  return AndMinMaxPdgCodes( 0, max, pdgCodes );
108  };
109  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){
110  return AndMinMaxPdgCodes( 0, max, pdgCode1, pdgCode2, pdgCode3, pdgCode4, pdgCode5, pdgCode6, pdgCode7, pdgCode8 );
111  };
112 
113 
115  // User interfaces -- Charge
117  // min defines how many particles you want in your events AT LEAST
118  // max defines how many particles you want in your events AT MOST
119  // Use for setting multiplicity constraints based on particle charges
120  // ChargeState
121  // kNeutral - neutral particles
122  // kPlus - positively charged particles
123  // kMinus - negatively charged particles
124  // kCharged - positively or negatively charged particles
125  // kAll - all particles (neutral or charged)
126  Bool_t AndMinMaxCharge( Int_t min, Int_t max, ChargeState charge );
127 
128 
129  //use the methods above with fixed minimum or maximum number of particles
130  Bool_t AndMinCharge(Int_t min, ChargeState charge){return AndMinMaxCharge( min, 9999, charge);};
131  Bool_t AndMaxCharge(Int_t max, ChargeState charge){return AndMinMaxCharge( 0, max, charge);};
132 
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 
142  //use the methods above with fixed minimum or maximum number of particles
143  Bool_t AndMinAllParticles(Int_t min){ return AndMinMaxCharge( min, 9999, kAll); };
144  Bool_t AndMaxAllParticles(Int_t max){ return AndMinMaxCharge( 0, max, kAll); };
145 
146 
148  // User interfaces -- Momentum constraints
150  //use if you want to consider particles which fulfill certain momentum constraints
151  //could be used in combination with any other filter, if you want to set momentum limits for several particles
152  Bool_t AndPRange(Double_t min, Double_t max){return AndMinMaxMom( min, max, kMomTotal);};
153  Bool_t AndPtRange(Double_t min, Double_t max){return AndMinMaxMom( min, max, kMomTrans);};
154  Bool_t AndPzRange(Double_t min, Double_t max){return AndMinMaxMom( min, max, kMomZ);};
155 
156 
158  // User interfaces -- Angular constraints
160  //use if you want to consider particles in certain angle or vertex ranges
161  //could be used in combination with any other filter, if you want to set geometric limits for several particles
162  Bool_t AndThetaRange(Double_t min, Double_t max){return AndMinMaxGeom( min, max, kTheta);};
163  Bool_t AndPhiRange(Double_t min, Double_t max){return AndMinMaxGeom( min, max, kPhi);};
164 
165 
167  // User interfaces -- Vertex constraints
169  //use if you want to consider particles in certain angle or vertex ranges
170  //could be used in combination with any other filter, if you want to set geometric limits for several particles
171  Bool_t AndVzRange(Double_t min, Double_t max){return AndMinMaxGeom( min, max, kVertexZ);};
172  Bool_t AndRhoRange(Double_t min, Double_t max){return AndMinMaxGeom( min, max, kVertexRho);};
174 
175 
176  // checks if the particles in the event (fParticleList) suit the filter settings
177  // kTRUE if event matches, kFALSE otherwise
178  Bool_t EventMatches(Int_t evtNr);
179 
180 
181  // returns kTRUE if the filter on pdg codes or charge states (or all particles) is turned on
182  // Momentum and angular constraints do not play any role here
184  return ( fFilterPdg || fFilterCharge );
185  };
186 
187 
188 private:
189  // in these methods the functions of the user interfaces from above are realized
190  // all user interfaces are only calling this two methods with different arguments
193 
194 
195  // initialize the counters for pdg codes or charge states
196  void InitCounters();
197 
198 
199  // Set default boundaries for momenta, angles and vertices
200  // that are always fulfilled
201  // this implements the default case of not filtering on momenta, angles or vertices
202  void SetDefaultBoundaries();
203 
204 
205  // count up the appropriate element of fCountCharge
206  void CountCharge(TParticle* particle);
207  // count up the appropriate element of fCountGroupId
208  void CountPdg(TParticle* particle);
209 
210  // check whether the particle fulfills the momentum constraints
211  Bool_t AcceptMomentum(TParticle* particle);
212  // check whether the particle fulfills the geometry constraints
213  Bool_t AcceptGeometry(TParticle* particle);
214 
215  // check whether the event satisfies the constraints on pdg multiplicities (fCountGroupId vector has to be filled already)
217  //check whether the event satisfies the constraints on charge state multiplicities (fCountCharge vector has to be filled already)
219 
220 
221  // first input = pdg code (as used by event generators)
222  // second input = "vector of group IDs that the pdg code is assigned to"
223  // one pdg code can be assigned to mutliple group IDs
224  // all pdg codes within the same group ID are treated as indistinguishable
225  // 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
227 
228  // for counting the multiplicities based on group IDs
229  std::vector<Int_t> fCountGroupId;
230  // for counts the multiplicities based on electrical charge
231  std::vector<Int_t> fCountCharge;
232 
233 
234  // vector defining the minimum and maximum acceptable multiplicities for each group ID
235  std::vector< std::pair<Int_t,Int_t> > fGroupIdCountsMinMax;
236  // vector containing the minimum and maximum acceptable multiplicities for each electrical charge state
237  // check ChargeState in FairEvtFilter.h for translating the index to an electrical charge
238  std::vector< std::pair<Int_t,Int_t> > fChargeCountsMinMax;
239 
240 
241  // only particles with a momentum within the min and max value will be counted
242  // check MomState in FairEvtFilter.h for translating the index to an electrical charge
243  std::vector< std::pair<Double_t,Double_t> > fMomMinMax;
244  // only particles with angles and vertices within the min and max value will be counted
245  // check GeomState in FairEvtFilter.h for translating the index to an electrical charge
246  std::vector< std::pair<Double_t,Double_t> > fGeomMinMax;
247 
248 
253 
254 
257  static const Int_t kInvalidPdgCode = 0;
258 
259 
261 
262 };
263 
264 
265 #endif /* FairEvtFilterOnSingleParticleCounts_H_ */
static const Int_t kInvalidPdgCode
constant holding an integer number which is not used as a pdg code this serves as a place holder when...
Bool_t AndMaxPdgCodes(Int_t max, std::vector< Int_t > &pdgCodes)
Bool_t fFilterMom
is kTRUE if any momentum filter is set
std::vector< std::pair< Int_t, Int_t > > fChargeCountsMinMax
Bool_t fFilterGeom
is kTRUE if any geometry filter is set
Bool_t fFilterCharge
is kTRUE if any filter on electrical charge is set (also kTRUE if filter on the total number of parti...
std::map< Int_t, std::vector< Int_t > >::iterator PdgGroupIdIterator
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)
Double_t mom
Definition: plot_dirc.C:14
friend F32vec4 max(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:26
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:4
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)
const int particle
std::vector< std::pair< Double_t, Double_t > > fGeomMinMax
Bool_t AndMinMaxMom(Double_t min, Double_t max, MomState mom)
Bool_t AndMinCharge(Int_t min, ChargeState charge)
Bool_t AndMinMaxCharge(Int_t min, Int_t max, ChargeState charge)
Double_t
std::vector< std::pair< Int_t, Int_t > > fGroupIdCountsMinMax
Bool_t AndMaxCharge(Int_t max, ChargeState charge)
Bool_t AndMinPdgCodes(Int_t min, std::vector< Int_t > &pdgCodes)
Bool_t AndMinMaxGeom(Double_t min, Double_t max, GeomState geom)
Bool_t fFilterPdg
is kTRUE if any pdg / group ID filter is set
friend F32vec4 min(const F32vec4 &a, const F32vec4 &b)
Definition: P4_F32vec4.h:25
std::map< Int_t, std::vector< Int_t > > PdgGroupId
TString name
ClassDef(FairEvtFilterOnSingleParticleCounts, 1)
std::vector< std::pair< Double_t, Double_t > > fMomMinMax
std::pair< Int_t, std::vector< Int_t > > PdgGroupIdPair