CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
reco::V0Filter Class Reference

#include <V0Filter.h>

Public Member Functions

bool operator() (const reco::TrackRef *tracks, unsigned int n) const
 
bool operator() (const reco::Track *tracks, unsigned int n) const
 
bool operator() (const std::vector< reco::CandidatePtr > &tracks) const
 
bool operator() (const std::vector< const Track * > &tracks) const
 
bool operator() (const std::vector< reco::TrackRef > &tracks) const
 
bool operator() (const std::vector< reco::Track > &tracks) const
 
bool operator() (const reco::Track *const *tracks, unsigned int n) const
 
 V0Filter (const edm::ParameterSet &params)
 
 ~V0Filter ()
 

Private Attributes

double k0sMassWindow
 

Detailed Description

Definition at line 14 of file V0Filter.h.

Constructor & Destructor Documentation

V0Filter::V0Filter ( const edm::ParameterSet params)

Definition at line 10 of file V0Filter.cc.

10  :
11  k0sMassWindow(params.getParameter<double>("k0sMassWindow"))
12 {
13 }
T getParameter(std::string const &) const
double k0sMassWindow
Definition: V0Filter.h:37
reco::V0Filter::~V0Filter ( )
inline

Definition at line 17 of file V0Filter.h.

References gen::n, operator()(), and l1t::tracks.

17 {}

Member Function Documentation

bool V0Filter::operator() ( const reco::TrackRef tracks,
unsigned int  n 
) const

Definition at line 47 of file V0Filter.cc.

References mps_fire::i, and gen::n.

Referenced by operator()(), and ~V0Filter().

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 }
bool V0Filter::operator() ( const reco::Track tracks,
unsigned int  n 
) const

Definition at line 80 of file V0Filter.cc.

References mps_fire::i, and gen::n.

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 }
bool V0Filter::operator() ( const std::vector< reco::CandidatePtr > &  tracks) const

Definition at line 57 of file V0Filter.cc.

References funct::abs(), ALCARECOTkAlJpsiMuMu_cff::charge, reco::ParticleMasses::k0, and k0sMassWindow.

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 }
double k0sMassWindow
Definition: V0Filter.h:37
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool V0Filter::operator() ( const std::vector< const Track * > &  tracks) const

Definition at line 73 of file V0Filter.cc.

74 {
75  return (*this)(&tracks[0], tracks.size());
76 }
bool reco::V0Filter::operator() ( const std::vector< reco::TrackRef > &  tracks) const
inline

Definition at line 26 of file V0Filter.h.

27  { return (*this)(&tracks[0], tracks.size()); }
bool reco::V0Filter::operator() ( const std::vector< reco::Track > &  tracks) const
inline

Definition at line 30 of file V0Filter.h.

References operator()().

31  { return (*this)(&tracks[0], tracks.size()); }
bool V0Filter::operator() ( const reco::Track *const *  tracks,
unsigned int  n 
) const

Definition at line 16 of file V0Filter.cc.

References funct::abs(), ALCARECOTkAlJpsiMuMu_cff::charge, reco::ParticleMasses::k0, k0sMassWindow, and reco::ParticleMasses::piPlus.

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 }
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
std::vector< vec1 > vec2
Definition: HCALResponse.h:16

Member Data Documentation

double reco::V0Filter::k0sMassWindow
private

Definition at line 37 of file V0Filter.h.

Referenced by operator()().