CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

reco::V0Filter Class Reference

#include <V0Filter.h>

List of all members.

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::TrackRef > &tracks) const
bool operator() (const std::vector< reco::Track > &tracks) const
 V0Filter (const edm::ParameterSet &params)
 ~V0Filter ()

Private Member Functions

bool operator() (const reco::Track **tracks, unsigned int n) const

Private Attributes

double k0sMassWindow

Detailed Description

Definition at line 12 of file V0Filter.h.


Constructor & Destructor Documentation

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

Definition at line 16 of file V0Filter.cc.

                                                :
        k0sMassWindow(params.getParameter<double>("k0sMassWindow"))
{
}
reco::V0Filter::~V0Filter ( ) [inline]

Definition at line 15 of file V0Filter.h.

{}

Member Function Documentation

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

Definition at line 53 of file V0Filter.cc.

References i, and n.

{
        std::vector<const reco::Track*> trackPtrs(n);
        for(unsigned int i = 0; i < n; i++)
                trackPtrs[i] = &*tracks[i];

        return (*this)(&trackPtrs[0], n);
}
bool V0Filter::operator() ( const reco::Track **  tracks,
unsigned int  n 
) const [inline, private]

Definition at line 22 of file V0Filter.cc.

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

{
        // only check for K0s for now

        if (n != 2)
                return true;

        if (tracks[0]->charge() * tracks[1]->charge() > 0)
                return true;

        ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > vec1;
        ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > vec2;

        vec1.SetPx(tracks[0]->px());
        vec1.SetPy(tracks[0]->py());
        vec1.SetPz(tracks[0]->pz());
        vec1.SetM(ParticleMasses::piPlus);

        vec2.SetPx(tracks[1]->px());
        vec2.SetPy(tracks[1]->py());
        vec2.SetPz(tracks[1]->pz());
        vec2.SetM(ParticleMasses::piPlus);

        double invariantMass = (vec1 + vec2).M();
        if (std::abs(invariantMass - ParticleMasses::k0) < k0sMassWindow)
                return false;

        return true;
}
bool reco::V0Filter::operator() ( const std::vector< reco::Track > &  tracks) const [inline]

Definition at line 25 of file V0Filter.h.

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

Definition at line 21 of file V0Filter.h.

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

Definition at line 63 of file V0Filter.cc.

References i, and n.

{
        std::vector<const reco::Track*> trackPtrs(n);
        for(unsigned int i = 0; i < n; i++)
                trackPtrs[i] = &tracks[i];

        return (*this)(&trackPtrs[0], n);
}

Member Data Documentation

Definition at line 32 of file V0Filter.h.

Referenced by operator()().