CMS 3D CMS Logo

Public Member Functions | Private Types | Private Member Functions | Private Attributes

MatcherByPullsAlgorithm Class Reference

#include <MuonAnalysis/MuonAssociators/interface/MatcherByPullsAlgorithm.cc>

List of all members.

Public Member Functions

void fillInvCov (const reco::Track &tk, AlgebraicSymMatrix55 &invCov) const
 Fill the inverse covariance matrix for the match(track, candidate, invCov) method.
std::pair< bool, float > match (const reco::Track &tk, const reco::Candidate &mc, const AlgebraicSymMatrix55 &invertedCovariance) const
std::pair< int, float > match (const reco::RecoCandidate &src, const std::vector< reco::GenParticle > &cands, const std::vector< uint8_t > &good) const
std::pair< int, float > match (const reco::Track &src, const std::vector< reco::GenParticle > &cands, const std::vector< uint8_t > &good) const
 MatcherByPullsAlgorithm (const edm::ParameterSet &)
void matchMany (const reco::Track &src, const std::vector< reco::GenParticle > &cands, const std::vector< uint8_t > &good, std::vector< std::pair< double, int > > &matchesToFill) const
void matchMany (const reco::RecoCandidate &src, const std::vector< reco::GenParticle > &cands, const std::vector< uint8_t > &good, std::vector< std::pair< double, int > > &matchesToFill) const
 ~MatcherByPullsAlgorithm ()

Private Types

enum  TrackChoice { StaTrack, TrkTrack, GlbTrack }
 

Enum to define which track to use.

More...

Private Member Functions

const reco::Tracktrack (const reco::RecoCandidate &src) const
 Get track out of Candidate, NULL if missing.

Private Attributes

double cut_
 Cut on the pull.
bool diagOnly_
 Use only the diagonal terms of the covariance matrix.
double dr2_
 DeltaR of the matching cone.
TrackChoice track_
 Track to be used in matching.
bool useVertex_
 Use also dxy / dsz in the matching.

Detailed Description

Description: Matches a RecoCandidate to a GenParticle (or any other Candidate) using the pulls of the helix parameters

Implementation: <Notes on="" implementation>="">

Definition at line 35 of file MatcherByPullsAlgorithm.h.


Member Enumeration Documentation

Enum to define which track to use.

Enumerator:
StaTrack 
TrkTrack 
GlbTrack 

Definition at line 84 of file MatcherByPullsAlgorithm.h.


Constructor & Destructor Documentation

MatcherByPullsAlgorithm::MatcherByPullsAlgorithm ( const edm::ParameterSet iConfig) [explicit]

Definition at line 13 of file MatcherByPullsAlgorithm.cc.

References Exception, edm::ParameterSet::getParameter(), GlbTrack, StaTrack, track(), track_, and TrkTrack.

                                                                               :
    dr2_(std::pow(iConfig.getParameter<double>("maxDeltaR"),2)),
    cut_(iConfig.getParameter<double>("maxPull")),
    diagOnly_(iConfig.getParameter<bool>("diagonalElementsOnly")),
    useVertex_(iConfig.getParameter<bool>("useVertexVariables"))
{
    std::string track = iConfig.getParameter<std::string>("track");
    if      (track == "standAloneMuon") track_ = StaTrack;
    else if (track == "combinedMuon")   track_ = GlbTrack;
    else if (track == "track")          track_ = TrkTrack;
    else throw cms::Exception("Configuration") << "MatcherByPullsAlgorithm: track '"<<track<<"' is not known\n" << 
               "Allowed values are: 'track', 'combinedMuon', 'standAloneMuon' (as per reco::RecoCandidate object)\n";
}
MatcherByPullsAlgorithm::~MatcherByPullsAlgorithm ( )

Definition at line 27 of file MatcherByPullsAlgorithm.cc.

{
}

Member Function Documentation

void MatcherByPullsAlgorithm::fillInvCov ( const reco::Track tk,
AlgebraicSymMatrix55 invCov 
) const

Fill the inverse covariance matrix for the match(track, candidate, invCov) method.

Definition at line 147 of file MatcherByPullsAlgorithm.cc.

References reco::TrackBase::covariance(), diagOnly_, i, j, and useVertex_.

Referenced by match(), and matchMany().

{
    if (useVertex_) {
        invCov = tk.covariance();
        if (diagOnly_) {
            for (size_t i = 0; i < 5; ++i) { for (size_t j = i+1; j < 5; ++j) { 
                invCov(i,j) = 0;
            } }
        }
        invCov.Invert();
    } else {
        AlgebraicSymMatrix33 momCov = tk.covariance().Sub<AlgebraicSymMatrix33>(0,0); // get 3x3 matrix
        if (diagOnly_) { momCov(0,1) = 0; momCov(0,2) = 0; momCov(1,2) = 0; }
        momCov.Invert();
        invCov.Place_at(momCov,0,0);
    }
}
std::pair< int, float > MatcherByPullsAlgorithm::match ( const reco::RecoCandidate src,
const std::vector< reco::GenParticle > &  cands,
const std::vector< uint8_t > &  good 
) const

Match Reco Candidate to MC Candidates, skipping the ones which are not good Return index of matchin and pull, or (-1,9e9)

Definition at line 73 of file MatcherByPullsAlgorithm.cc.

References match(), and track().

{
    const reco::Track * tk = track(src);
    return (tk == 0 ? 
            std::pair<int,float>(-1,9e9) : 
            match(*tk, cands, good));
}
std::pair< int, float > MatcherByPullsAlgorithm::match ( const reco::Track src,
const std::vector< reco::GenParticle > &  cands,
const std::vector< uint8_t > &  good 
) const

Match Reco Track to MC Tracks, skipping the ones which are not good Return index of matchin and pull, or (-1,9e9)

Definition at line 84 of file MatcherByPullsAlgorithm.cc.

References fillInvCov(), i, m, match(), and n.

{
    std::pair<int,float> best(-1,9e9);

    AlgebraicSymMatrix55 invCov; fillInvCov(tk, invCov);
    for (int i = 0, n = cands.size(); i < n; ++i) {
        if (!good[i]) continue;
        std::pair<bool,float> m = match(tk, cands[i], invCov);
        if (m.first && (m.second < best.second)) {
            best.first  = i;
            best.second = m.second;
        } 
    }
    return best;
}
std::pair< bool, float > MatcherByPullsAlgorithm::match ( const reco::Track tk,
const reco::Candidate mc,
const AlgebraicSymMatrix55 invertedCovariance 
) const

Match Track to MC Candidate, using already inverted covariance matrix Return status of match and pull, or (-1,9e9)

Definition at line 39 of file MatcherByPullsAlgorithm.cc.

References reco::Candidate::charge(), reco::TrackBase::charge(), gather_cfg::cout, cut_, SiPixelRawToDigiRegional_cfi::deltaPhi, Geom::deltaR2(), diffTreeTool::diff, dr2_, reco::TrackBase::dsz(), reco::TrackBase::dxy(), reco::TrackBase::eta(), reco::Candidate::eta(), i, reco::Candidate::p(), reco::Candidate::phi(), reco::TrackBase::phi(), reco::TrackBase::pt(), reco::Candidate::pt(), reco::TrackBase::qoverp(), mathSSE::sqrt(), reco::Candidate::theta(), reco::TrackBase::theta(), reco::Candidate::vertex(), reco::Candidate::vx(), reco::TrackBase::vx(), reco::TrackBase::vy(), reco::Candidate::vy(), reco::TrackBase::vz(), and reco::Candidate::vz().

Referenced by match(), and matchMany().

                                                                                                                    {
    if (::deltaR2(tk,c) <= dr2_) {
        AlgebraicVector5 diff(tk.qoverp() - c.charge()/c.p(), 
                              tk.theta() - c.theta(), 
                              ::deltaPhi(tk.phi(), c.phi()),
                              tk.dxy(c.vertex()),
                              tk.dsz(c.vertex())); 
        double pull = ROOT::Math::Similarity(diff, invCov);
#if 0
        std::cout << "Tk charge/pt/eta/phi/vx/vy/vz " << tk.charge() << "\t" << tk.pt() << "\t" << tk.eta() << "\t" << tk.phi() << "\t" << tk.vx() << "\t" << tk.vy() << "\t" << tk.vz() << std::endl;
        std::cout << "MC charge/pt/eta/phi/vx/vy/vz " << c.charge() << "\t" << c.pt() << "\t" << c.eta() << "\t" << c.phi() << "\t" << c.vx() << "\t" << c.vy() << "\t" << c.vz() << std::endl;
        std::cout << "Delta: " << diff << std::endl;
        std::cout << "Sigmas: ";
        for (size_t i = 0; i < 5; ++i) {
            if (invCov(i,i) == 0) std::cout << "---\t";
            else std::cout << std::sqrt(1.0/invCov(i,i)) << "\t";
        }
        std::cout << std::endl;
        std::cout << "Items: ";
        for (size_t i = 0; i < 5; ++i) {
            if (invCov(i,i) == 0) std::cout << "---\t";
            else std::cout << diff(i)*std::sqrt(invCov(i,i)) << "\t";
        }
        std::cout << std::endl;
        std::cout << "Pull: "  << pull << std::endl;
#endif
        return std::pair<bool,float>(pull < cut_, pull);
    }
    return std::pair<bool,float>(false,9e9);
}
void MatcherByPullsAlgorithm::matchMany ( const reco::RecoCandidate src,
const std::vector< reco::GenParticle > &  cands,
const std::vector< uint8_t > &  good,
std::vector< std::pair< double, int > > &  matchesToFill 
) const

Match Reco Candidate to MC Candidates, allowing multiple matches and skipping the ones which are not good It will fill in the vector of <double,int> with pull and index for all matching candidates, already sorted by pulls. This method assumes that matchesToFill is empty when the method is called

Definition at line 103 of file MatcherByPullsAlgorithm.cc.

References track().

{
    const reco::Track * tk = track(src);
    if (tk != 0) matchMany(*tk, cands, good, matchesToFill);
}
void MatcherByPullsAlgorithm::matchMany ( const reco::Track src,
const std::vector< reco::GenParticle > &  cands,
const std::vector< uint8_t > &  good,
std::vector< std::pair< double, int > > &  matchesToFill 
) const

Match Reco Track to MC Tracks, allowing multiple matches and skipping the ones which are not good It will fill in the vector of <double,int> with pull and index for all matching candidates, already sorted by pulls. This method assumes that matchesToFill is empty when the method is called

Definition at line 113 of file MatcherByPullsAlgorithm.cc.

References fillInvCov(), i, m, match(), n, and python::multivaluedict::sort().

{

    AlgebraicSymMatrix55 invCov; fillInvCov(tk, invCov);
    for (int i = 0, n = cands.size(); i < n; ++i) {
        if (!good[i]) continue;
        std::pair<bool,double> m = match(tk, cands[i], invCov);
        if (m.first) matchesToFill.push_back(std::make_pair(m.second,i));
    }
    std::sort(matchesToFill.begin(),matchesToFill.end());
}
const reco::Track * MatcherByPullsAlgorithm::track ( const reco::RecoCandidate src) const [private]

Get track out of Candidate, NULL if missing.

Definition at line 136 of file MatcherByPullsAlgorithm.cc.

References reco::RecoCandidate::combinedMuon(), edm::Ref< C, T, F >::get(), GlbTrack, edm::Ref< C, T, F >::isNonnull(), reco::RecoCandidate::standAloneMuon(), StaTrack, reco::RecoCandidate::track(), track_, and TrkTrack.

Referenced by match(), MatcherByPullsAlgorithm(), and matchMany().

                                                                  {
    switch (track_) {
        case StaTrack : return muon.standAloneMuon().isNonnull() ? muon.standAloneMuon().get() : 0;
        case GlbTrack : return muon.combinedMuon().isNonnull()   ? muon.combinedMuon().get()   : 0;
        case TrkTrack : return muon.track().isNonnull()          ? muon.track().get()          : 0;
    }
    assert(false);
}

Member Data Documentation

Cut on the pull.

Definition at line 93 of file MatcherByPullsAlgorithm.h.

Referenced by match().

Use only the diagonal terms of the covariance matrix.

Definition at line 96 of file MatcherByPullsAlgorithm.h.

Referenced by fillInvCov().

DeltaR of the matching cone.

Definition at line 90 of file MatcherByPullsAlgorithm.h.

Referenced by match().

Track to be used in matching.

Definition at line 87 of file MatcherByPullsAlgorithm.h.

Referenced by MatcherByPullsAlgorithm(), and track().

Use also dxy / dsz in the matching.

Definition at line 99 of file MatcherByPullsAlgorithm.h.

Referenced by fillInvCov().