CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

VertexHistoryAnalyzer Class Reference

Inheritance diagram for VertexHistoryAnalyzer:
edm::EDAnalyzer

List of all members.

Public Member Functions

 VertexHistoryAnalyzer (const edm::ParameterSet &)

Private Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void beginJob ()
virtual void beginRun (const edm::Run &, const edm::EventSetup &)
std::string particleString (int) const
std::string vertexString (HepMC::GenVertex::particles_in_const_iterator, HepMC::GenVertex::particles_in_const_iterator, HepMC::GenVertex::particles_out_const_iterator, HepMC::GenVertex::particles_out_const_iterator) const
std::string vertexString (TrackingParticleRefVector, TrackingParticleRefVector) const

Private Attributes

VertexClassifier classifier_
edm::ESHandle< ParticleDataTablepdt_
edm::InputTag vertexProducer_

Detailed Description

Definition at line 34 of file VertexHistoryAnalyzer.cc.


Constructor & Destructor Documentation

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

Definition at line 70 of file VertexHistoryAnalyzer.cc.

References edm::ParameterSet::getUntrackedParameter(), and vertexProducer_.

                                                                          : classifier_(config)
{
    vertexProducer_ = config.getUntrackedParameter<edm::InputTag> ("vertexProducer");
}

Member Function Documentation

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

Implements edm::EDAnalyzer.

Definition at line 76 of file VertexHistoryAnalyzer.cc.

References classifier_, gather_cfg::cout, VertexClassifier::evaluate(), VertexCategories::Fake, genParticleCandidates2GenParticles_cfi::genParticles, HistoryBase::genParticleTrail(), HistoryBase::genVertexTrail(), VertexClassifier::history(), getHLTprescales::index, VertexCategories::is(), VertexClassifier::newEvent(), particleString(), benchmark_cfg::pdgId, HistoryBase::simParticleTrail(), HistoryBase::simVertexTrail(), ExpressReco_HICollisions_FallBack::vertexCollection, vertexProducer_, and vertexString().

{
    // Set the classifier for a new event
    classifier_.newEvent(event, setup);

    // Vertex collection
    edm::Handle<edm::View<reco::Vertex> > vertexCollection;
    event.getByLabel(vertexProducer_, vertexCollection);

    // Get a constant reference to the track history associated to the classifier
    VertexHistory const & tracer = classifier_.history();

    // Loop over the track collection.
    for (std::size_t index = 0; index < vertexCollection->size(); index++)
    {
        std::cout << std::endl << "History for vertex #" << index << " : " << std::endl;

        // Classify the track and detect for fakes
        if ( !classifier_.evaluate( reco::VertexBaseRef(vertexCollection, index) ).is(VertexClassifier::Fake) )
        {
            // Get the list of TrackingParticles associated to
            VertexHistory::SimParticleTrail simParticles(tracer.simParticleTrail());

            // Loop over all simParticles
            for (std::size_t hindex=0; hindex<simParticles.size(); hindex++)
            {
                std::cout << "  simParticles [" << hindex << "] : "
                          << particleString(simParticles[hindex]->pdgId())
                          << std::endl;
            }

            // Get the list of TrackingVertexes associated to
            VertexHistory::SimVertexTrail simVertexes(tracer.simVertexTrail());

            // Loop over all simVertexes
            if ( !simVertexes.empty() )
            {
                for (std::size_t hindex=0; hindex<simVertexes.size(); hindex++)
                {
                    std::cout << "  simVertex    [" << hindex << "] : "
                              << vertexString(
                                  simVertexes[hindex]->sourceTracks(),
                                  simVertexes[hindex]->daughterTracks()
                              )
                              << std::endl;
                }
            }
            else
                std::cout << "  simVertex no found" << std::endl;

            // Get the list of GenParticles associated to
            VertexHistory::GenParticleTrail genParticles(tracer.genParticleTrail());

            // Loop over all genParticles
            for (std::size_t hindex=0; hindex<genParticles.size(); hindex++)
            {
                std::cout << "  genParticles [" << hindex << "] : "
                          << particleString(genParticles[hindex]->pdg_id())
                          << std::endl;
            }

            // Get the list of TrackingVertexes associated to
            VertexHistory::GenVertexTrail genVertexes(tracer.genVertexTrail());

            // Loop over all simVertexes
            if ( !genVertexes.empty() )
            {
                for (std::size_t hindex=0; hindex<genVertexes.size(); hindex++)
                {
                    std::cout << "  genVertex    [" << hindex << "] : "
                              << vertexString(
                                  genVertexes[hindex]->particles_in_const_begin(),
                                  genVertexes[hindex]->particles_in_const_end(),
                                  genVertexes[hindex]->particles_out_const_begin(),
                                  genVertexes[hindex]->particles_out_const_end()
                              )
                              << std::endl;
                }
            }
            else
                std::cout << "  genVertex no found" << std::endl;
        }
        else
            std::cout << "  fake vertex" << std::endl;

        std::cout << "  vertex categories : " << classifier_;
        std::cout << std::endl;
    }
}
void VertexHistoryAnalyzer::beginJob ( void  ) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 175 of file VertexHistoryAnalyzer.cc.

{

}
void VertexHistoryAnalyzer::beginRun ( const edm::Run run,
const edm::EventSetup setup 
) [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 168 of file VertexHistoryAnalyzer.cc.

References edm::EventSetup::getData(), and pdt_.

{
    // Get the particles table.
    setup.getData( pdt_ );
}
std::string VertexHistoryAnalyzer::particleString ( int  pdgId) const [private]

Definition at line 181 of file VertexHistoryAnalyzer.cc.

References RecoTau_DiTaus_pt_20-420_cfg::ParticleID, ExpressReco_HICollisions_FallBack::particleType, benchmark_cfg::pdgId, pdt_, and evf::utils::pid.

Referenced by analyze().

{
    ParticleData const * pid;

    std::ostringstream vDescription;

    HepPDT::ParticleID particleType(pdgId);

    if (particleType.isValid())
    {
        pid = pdt_->particle(particleType);
        if (pid)
            vDescription << pid->name();
        else
            vDescription << pdgId;
    }
    else
        vDescription << pdgId;

    return vDescription.str();
}
std::string VertexHistoryAnalyzer::vertexString ( HepMC::GenVertex::particles_in_const_iterator  in_begin,
HepMC::GenVertex::particles_in_const_iterator  in_end,
HepMC::GenVertex::particles_out_const_iterator  out_begin,
HepMC::GenVertex::particles_out_const_iterator  out_end 
) const [private]

Definition at line 261 of file VertexHistoryAnalyzer.cc.

References recoMuon::in, j, dbtoconf::out, RecoTau_DiTaus_pt_20-420_cfg::ParticleID, ExpressReco_HICollisions_FallBack::particleType, pdt_, and evf::utils::pid.

{
    ParticleData const * pid;

    std::ostringstream vDescription;

    std::size_t j = 0;

    HepMC::GenVertex::particles_in_const_iterator in, itmp;

    for (in = in_begin; in != in_end; in++, j++)
    {
        if (!j) vDescription << "(";

        HepPDT::ParticleID particleType((*in)->pdg_id());

        if (particleType.isValid())
        {
            pid = pdt_->particle(particleType);
            if (pid)
                vDescription << pid->name();
            else
                vDescription << (*in)->pdg_id();
        }
        else
            vDescription << (*in)->pdg_id();

        itmp = in;

        if (++itmp == in_end) vDescription << ")";
        else vDescription << ",";
    }

    vDescription << "->";
    j = 0;

    HepMC::GenVertex::particles_out_const_iterator out, otmp;

    for (out = out_begin; out != out_end; out++, j++)
    {
        if (!j) vDescription << "(";

        HepPDT::ParticleID particleType((*out)->pdg_id());

        if (particleType.isValid())
        {
            pid = pdt_->particle(particleType);
            if (pid)
                vDescription << pid->name();
            else
                vDescription << (*out)->pdg_id();
        }
        else
            vDescription << (*out)->pdg_id();

        otmp = out;

        if (++otmp == out_end) vDescription << ")";
        else vDescription << ",";
    }

    return vDescription.str();
}
std::string VertexHistoryAnalyzer::vertexString ( TrackingParticleRefVector  in,
TrackingParticleRefVector  out 
) const [private]

Definition at line 204 of file VertexHistoryAnalyzer.cc.

References j, RecoTau_DiTaus_pt_20-420_cfg::ParticleID, ExpressReco_HICollisions_FallBack::particleType, benchmark_cfg::pdgId, pdt_, evf::utils::pid, and edm::RefVector< C, T, F >::size().

Referenced by analyze().

{
    ParticleData const * pid;

    std::ostringstream vDescription;

    for (std::size_t j = 0; j < in.size(); j++)
    {
        if (!j) vDescription << "(";

        HepPDT::ParticleID particleType(in[j]->pdgId());

        if (particleType.isValid())
        {
            pid = pdt_->particle(particleType);
            if (pid)
                vDescription << pid->name();
            else
                vDescription << in[j]->pdgId();
        }
        else
            vDescription << in[j]->pdgId();

        if (j == in.size() - 1) vDescription << ")";
        else vDescription << ",";
    }

    vDescription << "->";

    for (std::size_t j = 0; j < out.size(); j++)
    {
        if (!j) vDescription << "(";

        HepPDT::ParticleID particleType(out[j]->pdgId());

        if (particleType.isValid())
        {
            pid = pdt_->particle(particleType);
            if (pid)
                vDescription << pid->name();
            else
                vDescription << out[j]->pdgId();
        }
        else
            vDescription << out[j]->pdgId();

        if (j == out.size() - 1) vDescription << ")";
        else vDescription << ",";
    }

    return vDescription.str();
}

Member Data Documentation

Definition at line 48 of file VertexHistoryAnalyzer.cc.

Referenced by analyze().

Definition at line 52 of file VertexHistoryAnalyzer.cc.

Referenced by beginRun(), particleString(), and vertexString().

Definition at line 50 of file VertexHistoryAnalyzer.cc.

Referenced by analyze(), and VertexHistoryAnalyzer().