CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Discriminants.cc
Go to the documentation of this file.
2 
3 namespace PFTauDiscriminants {
4 using namespace std;
5 
7 
8 void
10 {
11  // convert to int for TTree
12  result.push_back(static_cast<int>(input->getDecayMode()->getDecayMode()));
13 }
14 
15 void
17 {
18  size_t output = 0;
19  const std::vector<const reco::Candidate*>& theOutlierObjects = input->outlierObjectsSortedByPt();
20  for(std::vector<const reco::Candidate*>::const_iterator iObject = theOutlierObjects.begin();
21  iObject != theOutlierObjects.end();
22  ++iObject)
23  {
24  const reco::Candidate* currentObject = *iObject;
25  if (currentObject && currentObject->charge() != 0 )
26  output++;
27  }
28  // convert to int for TTree
29  result.push_back(static_cast<int>(output));
30 }
31 
32 void
34 {
35  const std::vector<const reco::Candidate*>& theOutlierObjects = input->outlierObjectsSortedByPt();
36  size_t output = theOutlierObjects.size();
37  // convert to int for TTree
38  result.push_back(static_cast<int>(output));
39 }
40 
41 void
43 {
44  result.push_back(input->getDecayMode()->pt());
45 }
46 
47 void
49 {
50  result.push_back(std::abs(input->getDecayMode()->eta()));
51 }
52 
53 void
55 {
56  const reco::Candidate* theMainTrack = input->mainTrack();
57  if (theMainTrack)
58  result.push_back(theMainTrack->pt());
59 }
60 
61 void
63 {
64  math::XYZVector signalObjectsAxis = input->getDecayMode()->momentum();
65 
66  const reco::Candidate* theMainTrack = input->mainTrack();
67 
68  DeltaR<math::XYZVector> myDRComputer;
69 
70  if (theMainTrack)
71  result.push_back(myDRComputer(theMainTrack->momentum(), signalObjectsAxis));
72 }
73 
74 void
76 {
77  const std::vector<const reco::Candidate*>& theSignalObjects = input->signalObjectsSortedByPt();
78 
79  const reco::Candidate* theMainTrack = input->mainTrack();
80 
81  for(std::vector<const reco::Candidate*>::const_iterator iObject = theSignalObjects.begin();
82  iObject != theSignalObjects.end();
83  ++iObject)
84  {
85  const reco::Candidate* currentObject = *iObject;
86  if (currentObject->charge() && currentObject != theMainTrack)
87  result.push_back(currentObject->pt());
88  }
89 }
90 
91 void
93 {
94  const std::vector<const reco::Candidate*>& theSignalObjects = input->signalObjectsSortedByPt();
95 
96  for(std::vector<const reco::Candidate*>::const_iterator iObject = theSignalObjects.begin();
97  iObject != theSignalObjects.end();
98  ++iObject)
99  {
100  const reco::Candidate* currentObject = *iObject;
101  if (!currentObject->charge())
102  result.push_back(currentObject->pt());
103  }
104 }
105 
106 void
108 {
109  const reco::PFTauDecayMode* theDecayMode = input->getDecayMode();
110  if (!theDecayMode)
111  return;
112 
113  const std::vector<const reco::Candidate*> theFilteredObjects = theDecayMode->filteredObjectCandidates();
114 
115  for(std::vector<const reco::Candidate*>::const_iterator iObject = theFilteredObjects.begin();
116  iObject != theFilteredObjects.end();
117  ++iObject)
118  {
119  const reco::Candidate* currentObject = *iObject;
120  result.push_back(currentObject->pt());
121  }
122 }
123 
124 void
126 {
127  const std::vector<const reco::Candidate*>& theSignalObjects = input->signalObjectsSortedByPt();
128 
129  for(std::vector<const reco::Candidate*>::const_iterator iObject = theSignalObjects.begin();
130  iObject != theSignalObjects.end();
131  ++iObject)
132  {
133  const reco::Candidate* currentObject = *iObject;
134  if (!currentObject->charge())
135  result.push_back(input->getLeafDaughters(currentObject).size());
136  }
137 }
138 
139 void
141 {
142  const std::vector<const reco::Candidate*>& theSignalObjects = input->signalObjectsSortedByPt();
143 
144  for(std::vector<const reco::Candidate*>::const_iterator iObject = theSignalObjects.begin();
145  iObject != theSignalObjects.end();
146  ++iObject)
147  {
148  const reco::Candidate* currentObject = *iObject;
149  if (!currentObject->charge())
150  {
151  std::vector<const reco::Candidate*> daughters = input->getLeafDaughters(currentObject);
152  for(std::vector<const reco::Candidate*>::const_iterator iDaughter = daughters.begin();
153  iDaughter != daughters.end();
154  ++iDaughter)
155  {
156  result.push_back((*iDaughter)->pt());
157  }
158  }
159  }
160 }
161 
162 
163 void
165 {
166  const std::vector<const reco::Candidate*>& theSignalObjects = input->signalObjectsSortedByPt();
167 
168  math::XYZVector signalObjectsAxis = input->getDecayMode()->momentum();
169 
170  const reco::Candidate* theMainTrack = input->mainTrack();
171 
172  DeltaR<math::XYZVector> myDRComputer;
173 
174  for(std::vector<const reco::Candidate*>::const_iterator iObject = theSignalObjects.begin();
175  iObject != theSignalObjects.end();
176  ++iObject)
177  {
178  const reco::Candidate* currentObject = *iObject;
179  if (currentObject->charge() && currentObject != theMainTrack)
180  result.push_back(myDRComputer(currentObject->momentum(), signalObjectsAxis));
181  }
182 }
183 
184 void
186 {
187  const std::vector<const reco::Candidate*>& theSignalObjects = input->signalObjectsSortedByPt();
188 
189  math::XYZVector signalObjectsAxis = input->getDecayMode()->momentum();
190 
191  DeltaR<math::XYZVector> myDRComputer;
192 
193  for(std::vector<const reco::Candidate*>::const_iterator iObject = theSignalObjects.begin();
194  iObject != theSignalObjects.end();
195  ++iObject)
196  {
197  const reco::Candidate* currentObject = *iObject;
198  if (!currentObject->charge())
199  result.push_back(myDRComputer(currentObject->momentum(), signalObjectsAxis));
200  }
201 }
202 
203 void
205 {
206  const std::vector<const reco::Candidate*>& theSignalObjects = input->signalObjectsSortedByPt();
207  const reco::Candidate* theMainTrack = input->mainTrack();
208  if (!theMainTrack)
209  return;
210  LorentzVector mainTrackFourVector = theMainTrack->p4();
211 
212  for(std::vector<const reco::Candidate*>::const_iterator iObject = theSignalObjects.begin();
213  iObject != theSignalObjects.end();
214  ++iObject)
215  {
216  const reco::Candidate* currentObject = *iObject;
217  if (currentObject != theMainTrack)
218  {
219  LorentzVector combinedFourVector = mainTrackFourVector + currentObject->p4();
220  result.push_back(combinedFourVector.mass());
221  }
222  }
223 }
224 
225 // takes invariant mass of all objects in signal cone
226 void
228 {
229  result.push_back(input->getDecayMode()->mass());
230 }
231 
232 // takes invariant mass of all objects in signal cone + Filtered objects
233 void
235 {
236  LorentzVector signalObjects = input->getDecayMode()->p4();
237  signalObjects += input->getDecayMode()->filteredObjects().p4();
238  result.push_back(signalObjects.M());
239 }
240 
241 // returns vector of invariant masses of larger and larger subsets of all signal objects e.g. result[2] is
242 // the invariant mass of the lead track with the next highest Pt object
243 
244 void
246 {
247  const std::vector<const reco::Candidate*>& theSignalObjects = input->signalObjectsSortedByPt();
248  const reco::Candidate* theMainTrack = input->mainTrack();
249  if (!theMainTrack)
250  return;
251  LorentzVector fourVectorSoFar = theMainTrack->p4();
252 
253  for(std::vector<const reco::Candidate*>::const_iterator iObject = theSignalObjects.begin();
254  iObject != theSignalObjects.end();
255  ++iObject)
256  {
257  const reco::Candidate* currentObject = *iObject;
258  if (currentObject != theMainTrack)
259  {
260  fourVectorSoFar += currentObject->p4();
261  result.push_back(fourVectorSoFar.mass());
262  }
263  }
264 }
265 
266 void
268 {
269  const std::vector<const reco::Candidate*>& theOutlierObjects = input->outlierObjectsSortedByPt();
270  for(std::vector<const reco::Candidate*>::const_iterator iObject = theOutlierObjects.begin();
271  iObject != theOutlierObjects.end();
272  ++iObject)
273  {
274  const reco::Candidate* currentObject = *iObject;
275  if (currentObject)
276  result.push_back(currentObject->pt());
277  }
278 }
279 
280 void
282 {
283  LorentzVector totalFourVector;
284  const std::vector<const reco::Candidate*>& theOutlierObjects = input->outlierObjectsSortedByPt();
285  for(std::vector<const reco::Candidate*>::const_iterator iObject = theOutlierObjects.begin();
286  iObject != theOutlierObjects.end();
287  ++iObject)
288  {
289  const reco::Candidate* currentObject = *iObject;
290  if (currentObject)
291  totalFourVector += currentObject->p4();
292  }
293  result.push_back(totalFourVector.pt());
294 }
295 
296 void
298 {
299  LorentzVector totalFourVector;
300  const std::vector<const reco::Candidate*>& theOutlierObjects = input->outlierObjectsSortedByPt();
301  for(std::vector<const reco::Candidate*>::const_iterator iObject = theOutlierObjects.begin();
302  iObject != theOutlierObjects.end();
303  ++iObject)
304  {
305  const reco::Candidate* currentObject = *iObject;
306  if (currentObject)
307  totalFourVector += currentObject->p4();
308  }
309  result.push_back(totalFourVector.M());
310 }
311 
312 void
314 {
315  const std::vector<const reco::Candidate*>& theoutlierObjects = input->outlierObjectsSortedByPt();
316  math::XYZVector signalObjectsAxis = input->getDecayMode()->momentum();
317  DeltaR<math::XYZVector> myDRComputer;
318  for(std::vector<const reco::Candidate*>::const_iterator iObject = theoutlierObjects.begin();
319  iObject != theoutlierObjects.end();
320  ++iObject)
321  {
322  const reco::Candidate* currentObject = *iObject;
323  if (currentObject)
324  result.push_back(myDRComputer(currentObject->momentum(), signalObjectsAxis));
325  }
326 }
327 
328 void
330 {
331  const std::vector<const reco::Candidate*>& theOutlierObjects = input->outlierObjectsSortedByPt();
332  for(std::vector<const reco::Candidate*>::const_iterator iObject = theOutlierObjects.begin();
333  iObject != theOutlierObjects.end();
334  ++iObject)
335  {
336  const reco::Candidate* currentObject = *iObject;
337  if (currentObject && currentObject->charge())
338  result.push_back(currentObject->pt());
339  }
340 }
341 
342 void
344 {
345  const std::vector<const reco::Candidate*>& theOutlierObjects = input->outlierObjectsSortedByPt();
346  double output = 0.0;
347  for(std::vector<const reco::Candidate*>::const_iterator iObject = theOutlierObjects.begin();
348  iObject != theOutlierObjects.end();
349  ++iObject)
350  {
351  const reco::Candidate* currentObject = *iObject;
352  if (currentObject && currentObject->charge())
353  output += currentObject->pt();
354  }
355  result.push_back(output);
356 }
357 
358 void
360 {
361  const std::vector<const reco::Candidate*>& theoutlierObjects = input->outlierObjectsSortedByPt();
362  math::XYZVector signalObjectsAxis = input->getDecayMode()->momentum();
363  DeltaR<math::XYZVector> myDRComputer;
364  for(std::vector<const reco::Candidate*>::const_iterator iObject = theoutlierObjects.begin();
365  iObject != theoutlierObjects.end();
366  ++iObject)
367  {
368  const reco::Candidate* currentObject = *iObject;
369  if (currentObject && currentObject->charge())
370  result.push_back(myDRComputer(currentObject->momentum(), signalObjectsAxis));
371  }
372 }
373 
374 void
376 {
377  const std::vector<const reco::Candidate*>& theOutlierObjects = input->outlierObjectsSortedByPt();
378  for(std::vector<const reco::Candidate*>::const_iterator iObject = theOutlierObjects.begin();
379  iObject != theOutlierObjects.end();
380  ++iObject)
381  {
382  const reco::Candidate* currentObject = *iObject;
383  if (currentObject && !currentObject->charge())
384  result.push_back(currentObject->pt());
385  }
386 }
387 
388 void
390 {
391  const std::vector<const reco::Candidate*>& theOutlierObjects = input->outlierObjectsSortedByPt();
392  double output = 0.0;
393  for(std::vector<const reco::Candidate*>::const_iterator iObject = theOutlierObjects.begin();
394  iObject != theOutlierObjects.end();
395  ++iObject)
396  {
397  const reco::Candidate* currentObject = *iObject;
398  if (currentObject && !currentObject->charge())
399  output += currentObject->pt();
400  }
401  result.push_back(output);
402 }
403 
404 void
406 {
407  const std::vector<const reco::Candidate*>& theoutlierObjects = input->outlierObjectsSortedByPt();
408  math::XYZVector signalObjectsAxis = input->getDecayMode()->momentum();
409  DeltaR<math::XYZVector> myDRComputer;
410  for(std::vector<const reco::Candidate*>::const_iterator iObject = theoutlierObjects.begin();
411  iObject != theoutlierObjects.end();
412  ++iObject)
413  {
414  const reco::Candidate* currentObject = *iObject;
415  if (currentObject && !currentObject->charge())
416  result.push_back(myDRComputer(currentObject->momentum(), signalObjectsAxis));
417  }
418 }
419 
420 
421 }
422 
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
const CompositeCandidate & filteredObjects() const
returns references to PF objects that were filtered
reco::Particle::LorentzVector LorentzVector
Definition: Discriminants.h:32
virtual double pt() const =0
transverse momentum
void doComputation(PFTauDiscriminantManager *input, std::vector< int > &result)
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
virtual Vector momentum() const
spatial momentum vector
#define abs(x)
Definition: mlp_lapack.h:159
static std::vector< const reco::Candidate * > getLeafDaughters(const reco::Candidate *input)
return the lowest level constituent candidates of a composite candidate
const std::vector< const reco::Candidate * > & outlierObjectsSortedByPt()
virtual double eta() const
momentum pseudorapidity
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
std::vector< const Candidate * > filteredObjectCandidates(int absCharge=-2) const
returns pointers to filtered objects (i.e. those not included in signal objects)
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
virtual double mass() const
mass
virtual Vector momentum() const =0
spatial momentum vector
hadronicTauDecayModes getDecayMode() const
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
tuple result
Definition: query.py:137
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
void doComputation(PFTauDiscriminantManager *input, std::vector< int > &result)
const std::vector< const reco::Candidate * > & signalObjectsSortedByPt()
accessed by Discriminant classes (caches to prevent multiple sorts)
virtual int charge() const =0
electric charge
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
const reco::Candidate * mainTrack()
get the &#39;main&#39; track (track computed for relevancy to tau decay resonances) (ie pi- in pi+pi+pi-) ...
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
const reco::PFTauDecayMode * getDecayMode() const
returns associated PFTauDecayMode
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
virtual double pt() const
transverse momentum
void doComputation(PFTauDiscriminantManager *input, std::vector< int > &result)
Definition: Discriminants.cc:9
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:25
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
void doComputation(PFTauDiscriminantManager *input, std::vector< double > &result)
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector