CMS 3D CMS Logo

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

HistoryBase Class Reference

Base class to all the history types. More...

#include <HistoryBase.h>

Inheritance diagram for HistoryBase:
TrackHistory VertexHistory

List of all members.

Public Types

typedef std::vector< const
HepMC::GenParticle * > 
GenParticleTrail
 GenParticle trail type.
typedef std::vector< const
HepMC::GenVertex * > 
GenVertexTrail
 GenVertex trail type.
typedef std::set< const
HepMC::GenVertex * > 
GenVertexTrailHelper
 GenVertex trail helper type.
typedef std::vector
< TrackingParticleRef
SimParticleTrail
 SimParticle trail type.
typedef std::vector
< TrackingVertexRef
SimVertexTrail
 SimVertex trail type.

Public Member Functions

void depth (int d)
 Set the depth of the history.
const HepMC::GenParticle * genParticle () const
 Returns a pointer to most primitive status 1 or 2 particle.
GenParticleTrail const & genParticleTrail () const
 Return all generated particle in the history.
GenVertexTrail const & genVertexTrail () const
 Return all generated vertex in the history.
 HistoryBase ()
const TrackingParticleRefsimParticle () const
 Return the initial tracking particle from the history.
SimParticleTrail const & simParticleTrail () const
 Return all the simulated particle in the history.
const TrackingVertexRefsimVertex () const
 Return the initial tracking vertex from the history.
SimVertexTrail const & simVertexTrail () const
 Return all the simulated vertices in the history.

Protected Member Functions

bool evaluate (TrackingParticleRef tpr)
 Evaluate track history using a TrackingParticleRef.
bool evaluate (TrackingVertexRef tvr)
 Evaluate track history using a TrackingParticleRef.

Protected Attributes

GenParticleTrail genParticleTrail_
GenVertexTrail genVertexTrail_
GenVertexTrailHelper genVertexTrailHelper_
SimParticleTrail simParticleTrail_
SimVertexTrail simVertexTrail_

Private Member Functions

void resetTrails ()
 Reset trail functions.
void resetTrails (TrackingParticleRef tpr)
void traceGenHistory (HepMC::GenVertex const *)
 Trace all the simulated information for a given pointer to a GenVertex.
void traceGenHistory (HepMC::GenParticle const *)
 Trace all the simulated information for a given pointer to a GenParticle.
bool traceSimHistory (TrackingParticleRef const &, int)
 Trace all the simulated information for a given reference to a TrackingParticle.
bool traceSimHistory (TrackingVertexRef const &, int)
 Trace all the simulated information for a given reference to a TrackingVertex.

Private Attributes

int depth_

Detailed Description

Base class to all the history types.

Definition at line 12 of file HistoryBase.h.


Member Typedef Documentation

typedef std::vector<const HepMC::GenParticle *> HistoryBase::GenParticleTrail

GenParticle trail type.

Definition at line 18 of file HistoryBase.h.

typedef std::vector<const HepMC::GenVertex *> HistoryBase::GenVertexTrail

GenVertex trail type.

Definition at line 21 of file HistoryBase.h.

typedef std::set<const HepMC::GenVertex *> HistoryBase::GenVertexTrailHelper

GenVertex trail helper type.

Definition at line 24 of file HistoryBase.h.

SimParticle trail type.

Definition at line 27 of file HistoryBase.h.

SimVertex trail type.

Definition at line 30 of file HistoryBase.h.


Constructor & Destructor Documentation

HistoryBase::HistoryBase ( ) [inline]

Definition at line 33 of file HistoryBase.h.

References depth_.

    {
        // Default depth
        depth_ = -1;
    }

Member Function Documentation

void HistoryBase::depth ( int  d) [inline]

Set the depth of the history.

Definition at line 47 of file HistoryBase.h.

References depth_.

Referenced by traceSimHistory(), TrackClassifier::TrackClassifier(), and VertexClassifier::VertexClassifier().

    {
        depth_ = d;
    }
bool HistoryBase::evaluate ( TrackingParticleRef  tpr) [inline, protected]

Evaluate track history using a TrackingParticleRef.

Reimplemented in TrackHistory.

Definition at line 114 of file HistoryBase.h.

References depth_, resetTrails(), and traceSimHistory().

    {
        resetTrails(tpr);
        return traceSimHistory(tpr, depth_);
    }
bool HistoryBase::evaluate ( TrackingVertexRef  tvr) [inline, protected]

Evaluate track history using a TrackingParticleRef.

Reimplemented in VertexHistory.

Definition at line 128 of file HistoryBase.h.

References depth_, resetTrails(), and traceSimHistory().

    {
        resetTrails();
        return traceSimHistory(tvr, depth_);
    }
const HepMC::GenParticle* HistoryBase::genParticle ( ) const [inline]

Returns a pointer to most primitive status 1 or 2 particle.

Definition at line 89 of file HistoryBase.h.

References genParticleTrail_.

Referenced by TrackClassifier::hadronFlavor(), and GenTrackMatcher::produce().

    {
        if ( genParticleTrail_.empty() ) return 0;
        return genParticleTrail_[genParticleTrail_.size()-1];
    }
GenParticleTrail const& HistoryBase::genParticleTrail ( ) const [inline]

Return all generated particle in the history.

Definition at line 71 of file HistoryBase.h.

References genParticleTrail_.

Referenced by VertexHistoryAnalyzer::analyze(), TrackHistoryAnalyzer::analyze(), TrackClassifier::processesAtGenerator(), and TrackClassifier::vertexInformation().

    {
        return genParticleTrail_;
    }
GenVertexTrail const& HistoryBase::genVertexTrail ( ) const [inline]

Return all generated vertex in the history.

Definition at line 65 of file HistoryBase.h.

References genVertexTrail_.

Referenced by VertexHistoryAnalyzer::analyze(), TrackHistoryAnalyzer::analyze(), VertexClassifier::processesAtGenerator(), and VertexClassifier::vertexInformation().

    {
        return genVertexTrail_;
    }
void HistoryBase::resetTrails ( TrackingParticleRef  tpr) [inline, private]

Definition at line 160 of file HistoryBase.h.

References resetTrails(), and simParticleTrail_.

    {
        resetTrails();
        simParticleTrail_.push_back(tpr);
    }
void HistoryBase::resetTrails ( ) [inline, private]

Reset trail functions.

Definition at line 151 of file HistoryBase.h.

References genParticleTrail_, genVertexTrail_, genVertexTrailHelper_, simParticleTrail_, and simVertexTrail_.

Referenced by evaluate(), and resetTrails().

const TrackingParticleRef& HistoryBase::simParticle ( ) const [inline]

Return the initial tracking particle from the history.

Definition at line 77 of file HistoryBase.h.

References simParticleTrail_.

Referenced by TrackClassifier::reconstructionInformation(), and TrackClassifier::simulationInformation().

    {
        return simParticleTrail_[0];
    }
SimParticleTrail const& HistoryBase::simParticleTrail ( ) const [inline]
const TrackingVertexRef& HistoryBase::simVertex ( ) const [inline]

Return the initial tracking vertex from the history.

Definition at line 83 of file HistoryBase.h.

References simVertexTrail_.

Referenced by recoBSVTagInfoValidationAnalyzer::analyze(), SVTagInfoValidationAnalyzer::analyze(), and VertexClassifier::simulationInformation().

    {
        return simVertexTrail_[0];
    }
SimVertexTrail const& HistoryBase::simVertexTrail ( ) const [inline]

Return all the simulated vertices in the history.

Definition at line 53 of file HistoryBase.h.

References simVertexTrail_.

Referenced by VertexHistoryAnalyzer::analyze(), TrackHistoryAnalyzer::analyze(), VertexClassifier::processesAtSimulation(), and VertexClassifier::vertexInformation().

    {
        return simVertexTrail_;
    }
void HistoryBase::traceGenHistory ( HepMC::GenVertex const *  genVertex) [private]

Trace all the simulated information for a given pointer to a GenVertex.

Definition at line 21 of file HistoryBase.cc.

References genVertexTrail_, genVertexTrailHelper_, and traceGenHistory().

{
    // Verify if has a vertex associated
    if (genVertex)
    {
        // Skip if already exist in the collection
        if ( genVertexTrailHelper_.find(genVertex) != genVertexTrailHelper_.end() )
            return;
        // Add vertex to the history
        genVertexTrail_.push_back(genVertex);
        genVertexTrailHelper_.insert(genVertex);
        // Verify if the vertex has incoming particles
        if ( genVertex->particles_in_size() )
            traceGenHistory( *(genVertex->particles_in_const_begin()) );
    }
}
void HistoryBase::traceGenHistory ( HepMC::GenParticle const *  genParticle) [private]

Trace all the simulated information for a given pointer to a GenParticle.

Definition at line 8 of file HistoryBase.cc.

References abs, depth_, and genParticleTrail_.

Referenced by traceGenHistory(), and traceSimHistory().

{
    // Third stop criteria: status abs(depth_) particles after the hadronization.
    // The after hadronization is done by detecting the pdg_id pythia code from 88 to 99
    if ( genParticle->status() <= abs(depth_) && (genParticle->pdg_id() < 88 || genParticle->pdg_id() > 99) )
    {
        genParticleTrail_.push_back(genParticle);
        // Get the producer vertex and trace it history
        traceGenHistory( genParticle->production_vertex() );
    }
}
bool HistoryBase::traceSimHistory ( TrackingVertexRef const &  trackingVertex,
int  depth 
) [private]

Trace all the simulated information for a given reference to a TrackingVertex.

Definition at line 59 of file HistoryBase.cc.

References spr::find(), edm::Ref< C, T, F >::isNonnull(), LogDebug, simParticleTrail_, simVertexTrail_, traceGenHistory(), and traceSimHistory().

{
    // verify if the parent vertex exists
    if ( trackingVertex.isNonnull() )
    {
        // save the vertex in the trail
        simVertexTrail_.push_back(trackingVertex);

        if ( !trackingVertex->sourceTracks().empty() )
        {
            LogDebug("TrackHistory") << "Moving on to the parent particle." << std::endl;

            // select the original source in case of combined vertices
            bool flag = false;
            TrackingVertex::tp_iterator itd, its;

            for (its = trackingVertex->sourceTracks_begin(); its != trackingVertex->sourceTracks_end(); its++)
            {
                for (itd = trackingVertex->daughterTracks_begin(); itd != trackingVertex->daughterTracks_end(); itd++)
                    if (itd != its)
                    {
                        flag = true;
                        break;
                    }
                if (flag)
                    break;
            }

            // verify if the new particle is not in the trail (looping partiles)
            if (
                std::find(
                    simParticleTrail_.begin(),
                    simParticleTrail_.end(),
                    *its
                ) != simParticleTrail_.end()
            )
            {
                LogDebug("TrackHistory") <<  "WARNING: Looping track found." << std::endl;
                return false;
            }

            // save particle in the trail
            simParticleTrail_.push_back(*its);
            return traceSimHistory (*its, --depth);
        }
        else if ( !trackingVertex->genVertices().empty() )
        {
            // navigate over all the associated generated vertexes
            LogDebug("TrackHistory") << "Vertex has " << trackingVertex->genVertices().size() << "GenVertex image." << std::endl;
            for (
                TrackingVertex::genv_iterator ivertex = trackingVertex->genVertices_begin();
                ivertex != trackingVertex->genVertices_end();
                ++ivertex
            )
                traceGenHistory(&(**(ivertex)));
            return true;
        }
        else
        {
            LogDebug("TrackHistory") <<  "WARNING: Source track for tracking vertex cannot be found." << std::endl;
        }
    }
    else
    {
        LogDebug("TrackHistory") << " WARNING: Vertex cannot be found.";
    }

    return false;
}
bool HistoryBase::traceSimHistory ( TrackingParticleRef const &  trackingParticle,
int  depth 
) [private]

Trace all the simulated information for a given reference to a TrackingParticle.

Definition at line 39 of file HistoryBase.cc.

References depth(), depth_, LogDebug, and traceGenHistory().

Referenced by evaluate(), and traceSimHistory().

{
    // first stop condition: if the required depth is reached
    if ( depth == depth_ && depth_ >= 0 ) return true;

    // sencond stop condition: if a gen particle is associated to the TP
    if ( !trackingParticle->genParticle().empty() )
    {
        LogDebug("TrackHistory") << "Particle " << trackingParticle->pdgId() << " has a GenParicle image." << std::endl;
        traceGenHistory(&(**(trackingParticle->genParticle_begin())));
        return true;
    }

    LogDebug("TrackHistory") << "No GenParticle image for " << trackingParticle->pdgId() << std::endl;

    // get a reference to the TP's parent vertex and trace it history
    return traceSimHistory( trackingParticle->parentVertex(), depth );
}

Member Data Documentation

int HistoryBase::depth_ [private]

Definition at line 136 of file HistoryBase.h.

Referenced by depth(), evaluate(), HistoryBase(), traceGenHistory(), and traceSimHistory().

Definition at line 99 of file HistoryBase.h.

Referenced by genParticle(), genParticleTrail(), resetTrails(), and traceGenHistory().

Definition at line 98 of file HistoryBase.h.

Referenced by genVertexTrail(), resetTrails(), and traceGenHistory().

Definition at line 104 of file HistoryBase.h.

Referenced by resetTrails(), and traceGenHistory().

Definition at line 101 of file HistoryBase.h.

Referenced by resetTrails(), simParticle(), simParticleTrail(), and traceSimHistory().

Definition at line 100 of file HistoryBase.h.

Referenced by resetTrails(), simVertex(), simVertexTrail(), and traceSimHistory().