CMS 3D CMS Logo

V0Filter.cc
Go to the documentation of this file.
1 #include <cmath>
2 
3 #include <Math/GenVector/PxPyPzM4D.h>
7 
8 using namespace reco;
9 
11  k0sMassWindow(params.getParameter<double>("k0sMassWindow"))
12 {
13 }
14 
15 bool
16 V0Filter::operator () (const reco::Track *const *tracks, unsigned int n) const
17 {
18  // only check for K0s for now
19 
20  if (n != 2)
21  return true;
22 
23  if (tracks[0]->charge() * tracks[1]->charge() > 0)
24  return true;
25 
26  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > vec1;
27  ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > vec2;
28 
29  vec1.SetPx(tracks[0]->px());
30  vec1.SetPy(tracks[0]->py());
31  vec1.SetPz(tracks[0]->pz());
32  vec1.SetM(ParticleMasses::piPlus);
33 
34  vec2.SetPx(tracks[1]->px());
35  vec2.SetPy(tracks[1]->py());
36  vec2.SetPz(tracks[1]->pz());
37  vec2.SetM(ParticleMasses::piPlus);
38 
39  double invariantMass = (vec1 + vec2).M();
40  if (std::abs(invariantMass - ParticleMasses::k0) < k0sMassWindow)
41  return false;
42 
43  return true;
44 }
45 
46 bool
47 V0Filter::operator () (const reco::TrackRef *tracks, unsigned int n) const
48 {
49  std::vector<const reco::Track*> trackPtrs(n);
50  for(unsigned int i = 0; i < n; i++)
51  trackPtrs[i] = &*tracks[i];
52 
53  return (*this)(&trackPtrs[0], n);
54 }
55 
56 bool
57 V0Filter::operator () (const std::vector<reco::CandidatePtr> & tracks) const
58 {
59  if (tracks.size() != 2)
60  return true;
61 
62  if (tracks[0]->charge() * tracks[1]->charge() > 0)
63  return true;
64 
65  double invariantMass = (tracks[0]->p4()+tracks[1]->p4()).M();
66  if (std::abs(invariantMass - ParticleMasses::k0) < k0sMassWindow)
67  return false;
68 
69  return true;
70 }
71 
72 bool
73 V0Filter::operator () (const std::vector<const reco::Track *> & tracks) const
74 {
75  return (*this)(&tracks[0], tracks.size());
76 }
77 
78 
79 bool
80 V0Filter::operator () (const reco::Track *tracks, unsigned int n) const
81 {
82  std::vector<const reco::Track*> trackPtrs(n);
83  for(unsigned int i = 0; i < n; i++)
84  trackPtrs[i] = &tracks[i];
85 
86  return (*this)(&trackPtrs[0], n);
87 }
88 
const double piPlus
Definition: ParticleMasses.h:9
double k0sMassWindow
Definition: V0Filter.h:37
std::vector< double > vec1
Definition: HCALResponse.h:15
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
V0Filter(const edm::ParameterSet &params)
Definition: V0Filter.cc:10
fixed size matrix
bool operator()(const reco::TrackRef *tracks, unsigned int n) const
Definition: V0Filter.cc:47
std::vector< vec1 > vec2
Definition: HCALResponse.h:16