CMS 3D CMS Logo

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

QualityCutsAnalyzer Class Reference

Inheritance diagram for QualityCutsAnalyzer:
edm::EDAnalyzer

List of all members.

Classes

struct  histogram_element_t
class  histogram_t

Public Member Functions

 QualityCutsAnalyzer (const edm::ParameterSet &)
 ~QualityCutsAnalyzer ()

Private Types

typedef std::vector
< std::vector
< histogram_element_t > > 
histogram_data_t
typedef std::vector< int > vint
typedef std::vector< std::string > vstring

Private Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void beginJob ()
virtual void endJob ()
void LoopOverJetTracksAssociation (const edm::ESHandle< TransientTrackBuilder > &, const edm::Handle< reco::VertexCollection > &, const edm::Handle< reco::JetTracksAssociationCollection > &)

Private Attributes

TrackClassifier classifier_
histogram_data_t histogram_data_
edm::InputTag jetTracksAssociation_
double maximumChiSquared_
int minimumNumberOfHits_
int minimumNumberOfPixelHits_
double minimumTransverseMomentum_
edm::InputTag primaryVertexProducer_
std::string rootFile_
edm::InputTag trackProducer_
reco::TrackBase::TrackQuality trackQuality_
bool useAllQualities_

Detailed Description

Definition at line 36 of file QualityCutsAnalyzer.cc.


Member Typedef Documentation

typedef std::vector<std::vector<histogram_element_t> > QualityCutsAnalyzer::histogram_data_t [private]

Definition at line 115 of file QualityCutsAnalyzer.cc.

typedef std::vector<int> QualityCutsAnalyzer::vint [private]

Definition at line 52 of file QualityCutsAnalyzer.cc.

typedef std::vector<std::string> QualityCutsAnalyzer::vstring [private]

Definition at line 53 of file QualityCutsAnalyzer.cc.


Constructor & Destructor Documentation

QualityCutsAnalyzer::QualityCutsAnalyzer ( const edm::ParameterSet config) [explicit]

Definition at line 222 of file QualityCutsAnalyzer.cc.

References gather_cfg::cout, edm::ParameterSet::getUntrackedParameter(), jetTracksAssociation_, maximumChiSquared_, minimumNumberOfHits_, minimumNumberOfPixelHits_, minimumTransverseMomentum_, primaryVertexProducer_, reco::TrackBase::qualityByName(), rootFile_, trackProducer_, trackQuality_, and useAllQualities_.

                                                                      : classifier_(config)
{
    trackProducer_         = config.getUntrackedParameter<edm::InputTag> ( "trackProducer" );
    primaryVertexProducer_ = config.getUntrackedParameter<edm::InputTag> ( "primaryVertexProducer" );
    jetTracksAssociation_  = config.getUntrackedParameter<edm::InputTag> ( "jetTracksAssociation" );

    rootFile_ = config.getUntrackedParameter<std::string> ( "rootFile" );

    minimumNumberOfHits_       = config.getUntrackedParameter<int> ( "minimumNumberOfHits" );
    minimumNumberOfPixelHits_  = config.getUntrackedParameter<int> ( "minimumNumberOfPixelHits" );
    minimumTransverseMomentum_ = config.getUntrackedParameter<double> ( "minimumTransverseMomentum" );
    maximumChiSquared_         = config.getUntrackedParameter<double> ( "maximumChiSquared" );

    std::string trackQualityType = config.getUntrackedParameter<std::string>("trackQualityClass"); //used
    trackQuality_ =  reco::TrackBase::qualityByName(trackQualityType);
    useAllQualities_ = false;

    std::transform(trackQualityType.begin(), trackQualityType.end(), trackQualityType.begin(), (int(*)(int)) std::tolower);
    if (trackQualityType == "any")
    {
        std::cout << "Using any" << std::endl;
        useAllQualities_ = true;
    }
}
QualityCutsAnalyzer::~QualityCutsAnalyzer ( )

Definition at line 247 of file QualityCutsAnalyzer.cc.

{}

Member Function Documentation

void QualityCutsAnalyzer::analyze ( const edm::Event event,
const edm::EventSetup setup 
) [private, virtual]

Implements edm::EDAnalyzer.

Definition at line 255 of file QualityCutsAnalyzer.cc.

References classifier_, edm::EventSetup::get(), jetTracksAssociation_, LoopOverJetTracksAssociation(), TrackClassifier::newEvent(), primaryVertexProducer_, and trackProducer_.

{
    // Track collection
    edm::Handle<edm::View<reco::Track> > trackCollection;
    event.getByLabel(trackProducer_,trackCollection);
    // Primary vertex
    edm::Handle<reco::VertexCollection> primaryVertexCollection;
    event.getByLabel(primaryVertexProducer_, primaryVertexCollection);
    // Jet to tracks associator
    edm::Handle<reco::JetTracksAssociationCollection> jetTracks;
    event.getByLabel(jetTracksAssociation_, jetTracks);
    // Trasient track builder
    edm::ESHandle<TransientTrackBuilder> TTbuilder;
    setup.get<TransientTrackRecord>().get("TransientTrackBuilder", TTbuilder);

    // Setting up event information for the track categories.
    classifier_.newEvent(event, setup);

    LoopOverJetTracksAssociation(
        TTbuilder,
        primaryVertexCollection,
        jetTracks
    );
}
void QualityCutsAnalyzer::beginJob ( void  ) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 283 of file QualityCutsAnalyzer.cc.

References histogram_data_.

{
    histogram_data_.resize(6);
}
void QualityCutsAnalyzer::endJob ( void  ) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 291 of file QualityCutsAnalyzer.cc.

References mergeVDriftHistosByStation::file, QualityCutsAnalyzer::histogram_t::Fill(), histogram_data_, i, j, rootFile_, and QualityCutsAnalyzer::histogram_t::Write().

{
    TFile file(rootFile_.c_str(), "RECREATE");
    file.cd();

    // saving the histograms
    for (std::size_t i=0; i<6; i++)
    {
        std::string particle;
        if (i == 0)
            particle = std::string("B_tracks");
        else if (i == 1)
            particle = std::string("C_tracks");
        else if (i == 2)
            particle = std::string("nonB_tracks");
        else if (i == 3)
            particle = std::string("displaced_tracks");
        else if (i == 4)
            particle = std::string("bad_tracks");
        else
            particle = std::string("fake_tracks");

        histogram_t histogram(particle);

        for (std::size_t j=0; j<histogram_data_[i].size(); j++)
            histogram.Fill(histogram_data_[i][j]);

        histogram.Write();
    }

    file.Flush();
}
void QualityCutsAnalyzer::LoopOverJetTracksAssociation ( const edm::ESHandle< TransientTrackBuilder > &  TTbuilder,
const edm::Handle< reco::VertexCollection > &  primaryVertexProducer_,
const edm::Handle< reco::JetTracksAssociationCollection > &  jetTracksAssociation 
) [private]

Definition at line 326 of file QualityCutsAnalyzer.cc.

References TrackCategories::Bad, TransientTrackBuilder::build(), TrackCategories::BWeakDecay, classifier_, TrackCategories::CWeakDecay, debug_cff::d0, TrackClassifier::evaluate(), TrackCategories::Fake, histogram_data_, i, getHLTprescales::index, TrackCategories::is(), IPTools::jetTrackDistance(), maximumChiSquared_, minimumNumberOfHits_, minimumNumberOfPixelHits_, minimumTransverseMomentum_, AlCaHLTBitMon_ParallelJobs::p, TrackCategories::PrimaryVertex, edm::ESHandle< T >::product(), edm::Handle< T >::product(), IPTools::signedDecayLength3D(), IPTools::signedImpactParameter3D(), IPTools::signedTransverseImpactParameter(), PrimaryVertexSorter::sortedList(), trackQuality_, testEve_cfg::tracks, useAllQualities_, and reco::Vertex::z().

Referenced by analyze().

{
    const TransientTrackBuilder * bproduct = TTbuilder.product();

    // getting the primary vertex
    // use first pv of the collection
    reco::Vertex pv;

    if (primaryVertexProducer_->size() != 0)
    {
        PrimaryVertexSorter pvs;
        std::vector<reco::Vertex> sortedList = pvs.sortedList(*(primaryVertexProducer_.product()));
        pv = (sortedList.front());
    }
    else
    { // create a dummy PV
        // cout << "NO PV FOUND" << endl;
        reco::Vertex::Error e;
        e(0,0)=0.0015*0.0015;
        e(1,1)=0.0015*0.0015;
        e(2,2)=15.*15.;
        reco::Vertex::Point p(0,0,0);
        pv = reco::Vertex(p,e,1,1,1);
    }

    reco::JetTracksAssociationCollection::const_iterator it = jetTracksAssociation->begin();

    int i=0;

    for (; it != jetTracksAssociation->end(); it++, i++)
    {
        // get jetTracks object
        reco::JetTracksAssociationRef jetTracks(jetTracksAssociation, i);

        double pvZ = pv.z();
        GlobalVector direction(jetTracks->first->px(), jetTracks->first->py(), jetTracks->first->pz());

        // get the tracks associated to the jet
        reco::TrackRefVector tracks = jetTracks->second;
        for (std::size_t index = 0; index < tracks.size(); index++)
        {
            edm::RefToBase<reco::Track> track(tracks[index]);

            double pt = tracks[index]->pt();
            double chi2 = tracks[index]->normalizedChi2();
            int hits = tracks[index]->hitPattern().numberOfValidHits();
            int pixelHits = tracks[index]->hitPattern().numberOfValidPixelHits();

            if (
                hits < minimumNumberOfHits_ ||
                pixelHits < minimumNumberOfPixelHits_ ||
                pt < minimumTransverseMomentum_ ||
                chi2 >  maximumChiSquared_ ||
                (!useAllQualities_ && !tracks[index]->quality(trackQuality_))
            ) continue;

            const reco::TransientTrack transientTrack = bproduct->build(&(*tracks[index]));
            double dta = - IPTools::jetTrackDistance(transientTrack, direction, pv).second.value();
            double sdl = IPTools::signedDecayLength3D(transientTrack, direction, pv).second.value();
            double ips = IPTools::signedImpactParameter3D(transientTrack, direction, pv).second.value();
            double d0 = IPTools::signedTransverseImpactParameter(transientTrack, direction, pv).second.value();
            double dz = tracks[index]->dz() - pvZ;

            // Classify the reco track;
            classifier_.evaluate( edm::RefToBase<reco::Track>(tracks[index]) );

            // Check for the different categories
            if ( classifier_.is(TrackClassifier::Fake) )
                histogram_data_[5].push_back(histogram_element_t(sdl, dta, d0, dz, ips, pt, chi2, hits, pixelHits));
            else if ( classifier_.is(TrackClassifier::BWeakDecay) )
                histogram_data_[0].push_back(histogram_element_t(sdl, dta, d0, dz, ips, pt, chi2, hits, pixelHits));
            else if ( classifier_.is(TrackClassifier::Bad) )
                histogram_data_[4].push_back(histogram_element_t(sdl, dta, d0, dz, ips, pt, chi2, hits, pixelHits));
            else if (
                !classifier_.is(TrackClassifier::CWeakDecay) &&
                !classifier_.is(TrackClassifier::PrimaryVertex)
            )
                histogram_data_[3].push_back(histogram_element_t(sdl, dta, d0, dz, ips, pt, chi2, hits, pixelHits));
            else if ( classifier_.is(TrackClassifier::CWeakDecay) )
                histogram_data_[1].push_back(histogram_element_t(sdl, dta, d0, dz, ips, pt, chi2, hits, pixelHits));
            else
                histogram_data_[2].push_back(histogram_element_t(sdl, dta, d0, dz, ips, pt, chi2, hits, pixelHits));

        }
    }
}

Member Data Documentation

Definition at line 214 of file QualityCutsAnalyzer.cc.

Referenced by analyze(), and LoopOverJetTracksAssociation().

Definition at line 116 of file QualityCutsAnalyzer.cc.

Referenced by beginJob(), endJob(), and LoopOverJetTracksAssociation().

Definition at line 57 of file QualityCutsAnalyzer.cc.

Referenced by analyze(), and QualityCutsAnalyzer().

Definition at line 62 of file QualityCutsAnalyzer.cc.

Referenced by LoopOverJetTracksAssociation(), and QualityCutsAnalyzer().

Definition at line 61 of file QualityCutsAnalyzer.cc.

Referenced by LoopOverJetTracksAssociation(), and QualityCutsAnalyzer().

Definition at line 61 of file QualityCutsAnalyzer.cc.

Referenced by LoopOverJetTracksAssociation(), and QualityCutsAnalyzer().

Definition at line 62 of file QualityCutsAnalyzer.cc.

Referenced by LoopOverJetTracksAssociation(), and QualityCutsAnalyzer().

Definition at line 56 of file QualityCutsAnalyzer.cc.

Referenced by analyze(), and QualityCutsAnalyzer().

std::string QualityCutsAnalyzer::rootFile_ [private]

Definition at line 59 of file QualityCutsAnalyzer.cc.

Referenced by endJob(), and QualityCutsAnalyzer().

Definition at line 55 of file QualityCutsAnalyzer.cc.

Referenced by analyze(), and QualityCutsAnalyzer().

Definition at line 65 of file QualityCutsAnalyzer.cc.

Referenced by LoopOverJetTracksAssociation(), and QualityCutsAnalyzer().

Definition at line 64 of file QualityCutsAnalyzer.cc.

Referenced by LoopOverJetTracksAssociation(), and QualityCutsAnalyzer().