Get track history and classify it in function of their . More...
#include <TrackClassifier.h>
Classes | |
struct | G4 |
struct | GeneratedPrimaryVertex |
Auxiliary class holding simulated primary vertices. More... | |
Public Types | |
typedef TrackCategories | Categories |
Type to the associate category. | |
Public Member Functions | |
TrackClassifier const & | evaluate (reco::TrackBaseRef const &) |
Classify the RecoTrack in categories. | |
TrackClassifier const & | evaluate (TrackingParticleRef const &) |
Classify the TrackingParticle in categories. | |
TrackClassifier const & | evaluate (reco::TrackRef const &track) |
Classify the RecoTrack in categories. | |
TrackHistory const & | history () const |
Returns a reference to the track history used in the classification. | |
void | newEvent (edm::Event const &, edm::EventSetup const &) |
Pre-process event information (for accessing reconstraction information) | |
TrackQuality const & | quality () const |
Returns a reference to the track quality used in the classification. | |
TrackClassifier (edm::ParameterSet const &) | |
Constructor by ParameterSet. | |
Private Member Functions | |
void | genPrimaryVertices () |
void | hadronFlavor () |
Get hadron flavor of the initial hadron. | |
bool | isCharged (const HepMC::GenParticle *) |
bool | isFinalstateParticle (const HepMC::GenParticle *) |
void | processesAtGenerator () |
Get all the information related to decay process. | |
void | processesAtSimulation () |
Get information about conversion and other interactions. | |
void | qualityInformation (reco::TrackBaseRef const &) |
Classify all the tracks by their reconstruction quality. | |
void | reconstructionInformation (reco::TrackBaseRef const &) |
Classify all the tracks by their association and reconstruction information. | |
void | simulationInformation () |
Get all the information related to the simulation details. | |
void | vertexInformation () |
Get geometrical information about the vertices. | |
Private Attributes | |
double | badPull_ |
edm::Handle< reco::BeamSpot > | beamSpot_ |
const edm::InputTag | beamSpotLabel_ |
std::vector < GeneratedPrimaryVertex > | genpvs_ |
const edm::InputTag | hepMCLabel_ |
double | longLivedDecayLength_ |
edm::ESHandle< MagneticField > | magneticField_ |
edm::Handle< edm::HepMCProduct > | mcInformation_ |
unsigned int | minTrackerSimHits_ |
unsigned int | numberOfInnerLayers_ |
edm::ESHandle< ParticleDataTable > | particleDataTable_ |
TrackQuality | quality_ |
TrackHistory | tracer_ |
edm::ESHandle < TransientTrackBuilder > | transientTrackBuilder_ |
double | vertexClusteringSqDistance_ |
Get track history and classify it in function of their .
Definition at line 26 of file TrackClassifier.h.
Type to the associate category.
Definition at line 32 of file TrackClassifier.h.
TrackClassifier::TrackClassifier | ( | edm::ParameterSet const & | config | ) |
Constructor by ParameterSet.
Definition at line 12 of file TrackClassifier.cc.
References badPull_, HistoryBase::depth(), edm::ParameterSet::getUntrackedParameter(), longLivedDecayLength_, minTrackerSimHits_, numberOfInnerLayers_, tracer_, and vertexClusteringSqDistance_.
: TrackCategories(), hepMCLabel_( config.getUntrackedParameter<edm::InputTag>("hepMC") ), beamSpotLabel_( config.getUntrackedParameter<edm::InputTag>("beamSpot") ), tracer_(config), quality_(config) { // Set the history depth after hadronization tracer_.depth(-2); // Set the maximum d0pull for the bad category badPull_ = config.getUntrackedParameter<double>("badPull"); // Set the minimum decay length for detecting long decays longLivedDecayLength_ = config.getUntrackedParameter<double>("longLivedDecayLength"); // Set the distance for clustering vertices float vertexClusteringDistance = config.getUntrackedParameter<double>("vertexClusteringDistance"); vertexClusteringSqDistance_ = vertexClusteringDistance * vertexClusteringDistance; // Set the number of innermost layers to check for bad hits numberOfInnerLayers_ = config.getUntrackedParameter<unsigned int>("numberOfInnerLayers"); // Set the minimum number of simhits in the tracker minTrackerSimHits_ = config.getUntrackedParameter<unsigned int>("minTrackerSimHits"); }
TrackClassifier const & TrackClassifier::evaluate | ( | reco::TrackBaseRef const & | track | ) |
Classify the RecoTrack in categories.
Definition at line 67 of file TrackClassifier.cc.
References TrackHistory::evaluate(), TrackCategories::Fake, TrackCategories::flags_, hadronFlavor(), processesAtGenerator(), processesAtSimulation(), qualityInformation(), reconstructionInformation(), TrackCategories::reset(), simulationInformation(), tracer_, TrackCategories::unknownTrack(), and vertexInformation().
Referenced by TrackCategoriesAnalyzer::analyze(), TrackHistoryAnalyzer::analyze(), TrackingParticleCategoriesAnalyzer::analyze(), evaluate(), QualityCutsAnalyzer::LoopOverJetTracksAssociation(), and JetVetoedTracksAssociationDRVertex::produce().
{ // Initializing the category vector reset(); // Associate and evaluate the track history (check for fakes) if ( tracer_.evaluate(track) ) { // Classify all the tracks by their association and reconstruction information reconstructionInformation(track); // Get all the information related to the simulation details simulationInformation(); // Analyse the track reconstruction quality qualityInformation(track); // Get hadron flavor of the initial hadron hadronFlavor(); // Get all the information related to decay process processesAtGenerator(); // Get information about conversion and other interactions processesAtSimulation(); // Get geometrical information about the vertices vertexInformation(); // Check for unkown classification unknownTrack(); } else flags_[Fake] = true; return *this; }
TrackClassifier const& TrackClassifier::evaluate | ( | reco::TrackRef const & | track | ) | [inline] |
Classify the RecoTrack in categories.
Definition at line 47 of file TrackClassifier.h.
References evaluate().
{ return evaluate( reco::TrackBaseRef(track) ); }
TrackClassifier const & TrackClassifier::evaluate | ( | TrackingParticleRef const & | track | ) |
Classify the TrackingParticle in categories.
Reimplemented in TrackClassifierByProxy< Collection >.
Definition at line 106 of file TrackClassifier.cc.
References TrackHistory::evaluate(), TrackCategories::flags_, hadronFlavor(), edm::RefToBase< T >::isNonnull(), processesAtGenerator(), processesAtSimulation(), qualityInformation(), TrackCategories::Reconstructed, reconstructionInformation(), TrackHistory::recoTrack(), Reconstruction_cff::recotrack, TrackCategories::reset(), simulationInformation(), tracer_, TrackCategories::unknownTrack(), and vertexInformation().
{ // Initializing the category vector reset(); // Trace the history for the given TP tracer_.evaluate(track); // Collect the associated reco track const reco::TrackBaseRef & recotrack = tracer_.recoTrack(); // If there is a reco truck then evaluate the simulated history if ( recotrack.isNonnull() ) { flags_[Reconstructed] = true; // Classify all the tracks by their association and reconstruction information reconstructionInformation(recotrack); // Analyse the track reconstruction quality qualityInformation(recotrack); } else flags_[Reconstructed] = false; // Get all the information related to the simulation details simulationInformation(); // Get hadron flavor of the initial hadron hadronFlavor(); // Get all the information related to decay process processesAtGenerator(); // Get information about conversion and other interactions processesAtSimulation(); // Get geometrical information about the vertices vertexInformation(); // Check for unkown classification unknownTrack(); return *this; }
void TrackClassifier::genPrimaryVertices | ( | ) | [private] |
Definition at line 564 of file TrackClassifier.cc.
References abs, event(), spr::find(), genpvs_, UserOptions_cff::idx, isCharged(), isFinalstateParticle(), m, mcInformation_, parents, pos, funct::pow(), python::multivaluedict::sort(), and vertexClusteringSqDistance_.
Referenced by newEvent().
{ genpvs_.clear(); const HepMC::GenEvent * event = mcInformation_->GetEvent(); if (event) { int idx = 0; // Loop over the different GenVertex for ( HepMC::GenEvent::vertex_const_iterator ivertex = event->vertices_begin(); ivertex != event->vertices_end(); ++ivertex ) { bool hasParentVertex = false; // Loop over the parents looking to see if they are coming from a production vertex for ( HepMC::GenVertex::particle_iterator iparent = (*ivertex)->particles_begin(HepMC::parents); iparent != (*ivertex)->particles_end(HepMC::parents); ++iparent ) if ( (*iparent)->production_vertex() ) { hasParentVertex = true; break; } // Reject those vertices with parent vertices if (hasParentVertex) continue; // Get the position of the vertex HepMC::FourVector pos = (*ivertex)->position(); double const mm = 0.1; GeneratedPrimaryVertex pv(pos.x()*mm, pos.y()*mm, pos.z()*mm); std::vector<GeneratedPrimaryVertex>::iterator ientry = genpvs_.begin(); // Search for a VERY close vertex in the list for (; ientry != genpvs_.end(); ++ientry) { double distance2 = pow(pv.x - ientry->x, 2) + pow(pv.y - ientry->y, 2) + pow(pv.z - ientry->z, 2); if ( distance2 < vertexClusteringSqDistance_ ) break; } // Check if there is not a VERY close vertex and added to the list if (ientry == genpvs_.end()) ientry = genpvs_.insert(ientry,pv); // Add the vertex barcodes to the new or existent vertices ientry->genVertex.push_back((*ivertex)->barcode()); // Collect final state descendants for ( HepMC::GenVertex::particle_iterator idecendants = (*ivertex)->particles_begin(HepMC::descendants); idecendants != (*ivertex)->particles_end(HepMC::descendants); ++idecendants ) { if (isFinalstateParticle(*idecendants)) if ( find(ientry->finalstateParticles.begin(), ientry->finalstateParticles.end(), (*idecendants)->barcode()) == ientry->finalstateParticles.end() ) { ientry->finalstateParticles.push_back((*idecendants)->barcode()); HepMC::FourVector m = (*idecendants)->momentum(); ientry->ptot.setPx(ientry->ptot.px() + m.px()); ientry->ptot.setPy(ientry->ptot.py() + m.py()); ientry->ptot.setPz(ientry->ptot.pz() + m.pz()); ientry->ptot.setE(ientry->ptot.e() + m.e()); ientry->ptsq += m.perp() * m.perp(); if ( m.perp() > 0.8 && std::abs(m.pseudoRapidity()) < 2.5 && isCharged(*idecendants) ) ientry->nGenTrk++; } } idx++; } } std::sort(genpvs_.begin(), genpvs_.end()); }
void TrackClassifier::hadronFlavor | ( | ) | [private] |
Get hadron flavor of the initial hadron.
Definition at line 254 of file TrackClassifier.cc.
References TrackCategories::Bottom, TrackCategories::Charm, TrackCategories::flags_, HistoryBase::genParticle(), configurableAnalysis::GenParticle, TrackCategories::Light, RecoTau_DiTaus_pt_20-420_cfg::ParticleID, evf::utils::pid, and tracer_.
Referenced by evaluate().
{ // Get the initial hadron const HepMC::GenParticle * particle = tracer_.genParticle(); // Check for the initial hadron if (particle) { HepPDT::ParticleID pid(particle->pdg_id()); flags_[Bottom] = pid.hasBottom(); flags_[Charm] = pid.hasCharm(); flags_[Light] = !pid.hasCharm() && !pid.hasBottom(); } }
TrackHistory const& TrackClassifier::history | ( | ) | const [inline] |
Returns a reference to the track history used in the classification.
Definition at line 53 of file TrackClassifier.h.
References tracer_.
Referenced by TrackHistoryAnalyzer::analyze().
{ return tracer_; }
bool TrackClassifier::isCharged | ( | const HepMC::GenParticle * | p | ) | [private] |
Definition at line 551 of file TrackClassifier.cc.
References particleDataTable_.
Referenced by genPrimaryVertices().
{ const ParticleData * part = particleDataTable_->particle( p->pdg_id() ); if (part) return part->charge()!=0; else { // the new/improved particle table doesn't know anti-particles return particleDataTable_->particle( -p->pdg_id() ) != 0; } }
bool TrackClassifier::isFinalstateParticle | ( | const HepMC::GenParticle * | p | ) | [private] |
Definition at line 545 of file TrackClassifier.cc.
Referenced by genPrimaryVertices().
void TrackClassifier::newEvent | ( | edm::Event const & | event, |
edm::EventSetup const & | setup | ||
) |
Pre-process event information (for accessing reconstraction information)
Reimplemented in TrackClassifierByProxy< Collection >.
Definition at line 39 of file TrackClassifier.cc.
References beamSpot_, beamSpotLabel_, genPrimaryVertices(), edm::EventSetup::get(), edm::EventSetup::getData(), hepMCLabel_, magneticField_, mcInformation_, TrackQuality::newEvent(), TrackHistory::newEvent(), particleDataTable_, quality_, tracer_, and transientTrackBuilder_.
Referenced by QualityCutsAnalyzer::analyze(), TrackCategoriesAnalyzer::analyze(), TrackHistoryAnalyzer::analyze(), TrackingParticleCategoriesAnalyzer::analyze(), and JetVetoedTracksAssociatorAtVertex::produce().
{ // Get the new event information for the tracer tracer_.newEvent(event, setup); // Get the new event information for the track quality analyser quality_.newEvent(event, setup); // Get hepmc of the event event.getByLabel(hepMCLabel_, mcInformation_); // Magnetic field setup.get<IdealMagneticFieldRecord>().get(magneticField_); // Get the partivle data table setup.getData(particleDataTable_); // get the beam spot event.getByLabel(beamSpotLabel_, beamSpot_); // Transient track builder setup.get<TransientTrackRecord>().get("TransientTrackBuilder", transientTrackBuilder_); // Create the list of primary vertices associated to the event genPrimaryVertices(); }
void TrackClassifier::processesAtGenerator | ( | ) | [private] |
Get all the information related to decay process.
Definition at line 270 of file TrackClassifier.cc.
References abs, TrackCategories::BWeakDecay, TrackCategories::ChargeKaonDecay, TrackCategories::ChargePionDecay, TrackCategories::CWeakDecay, TrackCategories::DecayOnFlightMuon, TrackCategories::flags_, TrackCategories::FromBWeakDecayMuon, TrackCategories::FromChargeKaonMuon, TrackCategories::FromChargePionMuon, TrackCategories::FromCWeakDecayMuon, HistoryBase::genParticleTrail(), TrackCategories::JpsiDecay, TrackCategories::KsDecay, TrackCategories::LambdaDecay, TrackCategories::LongLivedDecay, longLivedDecayLength_, TrackCategories::Muon, particleDataTable_, RecoTau_DiTaus_pt_20-420_cfg::ParticleID, TrackCategories::SigmaMinusDecay, TrackCategories::SigmaPlusDecay, TrackCategories::TauDecay, tracer_, update, and TrackCategories::XiDecay.
Referenced by evaluate().
{ // pdgid of the "in" particle to the production vertex int pdgid = 0; // Get the generated particles from track history TrackHistory::GenParticleTrail const & genParticleTrail = tracer_.genParticleTrail(); // Loop over the generated particles for (TrackHistory::GenParticleTrail::const_iterator iparticle = genParticleTrail.begin(); iparticle != genParticleTrail.end(); ++iparticle) { // Get the source vertex for the particle HepMC::GenVertex * productionVertex = (*iparticle)->production_vertex(); // Get the pointer to the vertex by removing the const-ness (no const methos in HepMC::GenVertex) // HepMC::GenVertex * vertex = const_cast<HepMC::GenVertex *>(*ivertex); // Check for a non-null pointer to the production vertex if (productionVertex) { // Only case track history will navegate (one in or source particle per vertex) if ( productionVertex->particles_in_size() == 1 ) { // Look at the pdgid of the first "in" particle to the vertex pdgid = std::abs((*productionVertex->particles_in_const_begin())->pdg_id()); // Get particle type HepPDT::ParticleID particleID(pdgid); // Check if the particle type is valid one if (particleID.isValid()) { // Get particle data ParticleData const * particleData = particleDataTable_->particle(particleID); // Check if the particle exist in the table if (particleData) { // Check if their life time is bigger than longLivedDecayLength_ if ( particleData->lifetime() > longLivedDecayLength_ ) update(flags_[LongLivedDecay], true); // Check for B and C weak decays update(flags_[BWeakDecay], particleID.hasBottom()); update(flags_[CWeakDecay], particleID.hasCharm()); // Check for B and C pure leptonic decay int daughterId = abs((*iparticle)->pdg_id()); update(flags_[FromBWeakDecayMuon], particleID.hasBottom() && daughterId == 13); update(flags_[FromCWeakDecayMuon], particleID.hasCharm() && daughterId == 13); } // Check Tau, Ks and Lambda decay update(flags_[ChargePionDecay], pdgid == 211); update(flags_[ChargeKaonDecay], pdgid == 321); update(flags_[TauDecay], pdgid == 15); update(flags_[KsDecay], pdgid == 310); update(flags_[LambdaDecay], pdgid == 3122); update(flags_[JpsiDecay], pdgid == 443); update(flags_[XiDecay], pdgid == 3312); update(flags_[SigmaPlusDecay], pdgid == 3222); update(flags_[SigmaMinusDecay], pdgid == 3112); } } } } // Decays in flight update(flags_[FromChargePionMuon], flags_[Muon] && flags_[ChargePionDecay]); update(flags_[FromChargeKaonMuon], flags_[Muon] && flags_[ChargeKaonDecay]); update(flags_[DecayOnFlightMuon], (flags_[FromChargePionMuon] || flags_[FromChargeKaonMuon])); }
void TrackClassifier::processesAtSimulation | ( | ) | [private] |
Get information about conversion and other interactions.
Definition at line 338 of file TrackClassifier.cc.
References abs, TrackClassifier::G4::Annihilation, TrackCategories::AnnihilationProcess, TrackCategories::BWeakDecay, TrackCategories::ChargeKaonDecay, TrackCategories::ChargePionDecay, TrackClassifier::G4::Compton, TrackCategories::ComptonProcess, TrackClassifier::G4::Conversions, TrackCategories::ConversionsProcess, TrackCategories::CWeakDecay, TrackClassifier::G4::Decay, TrackCategories::DecayOnFlightMuon, TrackCategories::DecayProcess, TrackClassifier::G4::EBrem, TrackCategories::EBremProcess, TrackClassifier::G4::EIoni, TrackCategories::EIoniProcess, TrackCategories::flags_, TrackCategories::FromBWeakDecayMuon, TrackCategories::FromChargeKaonMuon, TrackCategories::FromChargePionMuon, TrackCategories::FromCWeakDecayMuon, TrackClassifier::G4::Hadronic, TrackCategories::HadronicProcess, TrackClassifier::G4::HIoni, TrackCategories::HIoniProcess, edm::Ref< C, T, F >::isNonnull(), TrackCategories::JpsiDecay, TrackCategories::KnownProcess, TrackCategories::KsDecay, TrackCategories::LambdaDecay, TrackCategories::LongLivedDecay, longLivedDecayLength_, TrackClassifier::G4::MuBrem, TrackCategories::MuBremProcess, TrackClassifier::G4::MuIoni, TrackCategories::MuIoniProcess, TrackClassifier::G4::MuNucl, TrackCategories::MuNuclProcess, TrackCategories::Muon, TrackClassifier::G4::MuPairProd, TrackCategories::MuPairProdProcess, TrackCategories::OmegaDecay, particleDataTable_, RecoTau_DiTaus_pt_20-420_cfg::ParticleID, TrackClassifier::G4::Photon, TrackCategories::PhotonProcess, TrackClassifier::G4::Primary, TrackCategories::PrimaryProcess, LaserDQM_cfg::process, TrackCategories::SigmaMinusDecay, TrackCategories::SigmaPlusDecay, HistoryBase::simParticleTrail(), TrackClassifier::G4::SynchrotronRadiation, TrackCategories::SynchrotronRadiationProcess, TrackCategories::TauDecay, tracer_, TrackClassifier::G4::Undefined, TrackCategories::UndefinedProcess, TrackClassifier::G4::Unknown, TrackCategories::UnknownProcess, update, and TrackCategories::XiDecay.
Referenced by evaluate().
{ TrackHistory::SimParticleTrail const & simParticleTrail = tracer_.simParticleTrail(); // Loop over the simulated particles for ( TrackHistory::SimParticleTrail::const_iterator iparticle = simParticleTrail.begin(); iparticle != simParticleTrail.end(); ++iparticle ) { // pdgid of the real source parent vertex int pdgid = 0; // Get a reference to the TP's parent vertex TrackingVertexRef const & parentVertex = (*iparticle)->parentVertex(); // Look for the original source track if ( parentVertex.isNonnull() ) { // select the original source in case of combined vertices bool flag = false; TrackingVertex::tp_iterator itd, its; for (its = parentVertex->sourceTracks_begin(); its != parentVertex->sourceTracks_end(); ++its) { for (itd = parentVertex->daughterTracks_begin(); itd != parentVertex->daughterTracks_end(); ++itd) if (itd != its) { flag = true; break; } if (flag) break; } // Collect the pdgid of the original source track if ( its != parentVertex->sourceTracks_end() ) pdgid = std::abs((*its)->pdgId()); else pdgid = 0; } // Check for the number of psimhit if different from zero if ((*iparticle)->trackPSimHit().empty()) continue; // Collect the G4 process of the first psimhit (it should be the same for all of them) unsigned short process = (*iparticle)->pSimHit_begin()->processType(); // Flagging all the different processes update( flags_[KnownProcess], process != G4::Undefined && process != G4::Unknown && process != G4::Primary ); update(flags_[UndefinedProcess], process == G4::Undefined); update(flags_[UnknownProcess], process == G4::Unknown); update(flags_[PrimaryProcess], process == G4::Primary); update(flags_[HadronicProcess], process == G4::Hadronic); update(flags_[DecayProcess], process == G4::Decay); update(flags_[ComptonProcess], process == G4::Compton); update(flags_[AnnihilationProcess], process == G4::Annihilation); update(flags_[EIoniProcess], process == G4::EIoni); update(flags_[HIoniProcess], process == G4::HIoni); update(flags_[MuIoniProcess], process == G4::MuIoni); update(flags_[PhotonProcess], process == G4::Photon); update(flags_[MuPairProdProcess], process == G4::MuPairProd); update(flags_[ConversionsProcess], process == G4::Conversions); update(flags_[EBremProcess], process == G4::EBrem); update(flags_[SynchrotronRadiationProcess], process == G4::SynchrotronRadiation); update(flags_[MuBremProcess], process == G4::MuBrem); update(flags_[MuNuclProcess], process == G4::MuNucl); // Get particle type HepPDT::ParticleID particleID(pdgid); // Check if the particle type is valid one if (particleID.isValid()) { // Get particle data ParticleData const * particleData = particleDataTable_->particle(particleID); // Special treatment for decays if (process == G4::Decay) { // Check if the particle exist in the table if (particleData) { // Check if their life time is bigger than 1e-14 if ( particleDataTable_->particle(particleID)->lifetime() > longLivedDecayLength_ ) update(flags_[LongLivedDecay], true); // Check for B and C weak decays update(flags_[BWeakDecay], particleID.hasBottom()); update(flags_[CWeakDecay], particleID.hasCharm()); // Check for B or C pure leptonic decays int daughtId = abs((*iparticle)->pdgId()); update(flags_[FromBWeakDecayMuon], particleID.hasBottom() && daughtId == 13); update(flags_[FromCWeakDecayMuon], particleID.hasCharm() && daughtId == 13); } // Check decays update(flags_[ChargePionDecay], pdgid == 211); update(flags_[ChargeKaonDecay], pdgid == 321); update(flags_[TauDecay], pdgid == 15); update(flags_[KsDecay], pdgid == 310); update(flags_[LambdaDecay], pdgid == 3122); update(flags_[JpsiDecay], pdgid == 443); update(flags_[XiDecay], pdgid == 3312); update(flags_[OmegaDecay], pdgid == 3334); update(flags_[SigmaPlusDecay], pdgid == 3222); update(flags_[SigmaMinusDecay], pdgid == 3112); } } } // Decays in flight update(flags_[FromChargePionMuon], flags_[Muon] && flags_[ChargePionDecay]); update(flags_[FromChargeKaonMuon], flags_[Muon] && flags_[ChargeKaonDecay]); update(flags_[DecayOnFlightMuon], flags_[FromChargePionMuon] || flags_[FromChargeKaonMuon]); }
TrackQuality const& TrackClassifier::quality | ( | void | ) | const [inline] |
Returns a reference to the track quality used in the classification.
Definition at line 59 of file TrackClassifier.h.
References quality_.
{ return quality_; }
void TrackClassifier::qualityInformation | ( | reco::TrackBaseRef const & | track | ) | [private] |
Classify all the tracks by their reconstruction quality.
Definition at line 226 of file TrackClassifier.cc.
References TrackCategories::BadInnerHits, TrackQuality::evaluate(), TrackCategories::flags_, TrackQuality::Layer::hits, i, j, TrackQuality::layer(), min, TrackQuality::Layer::Misassoc, TrackQuality::Layer::Noise, numberOfInnerLayers_, TrackQuality::numberOfLayers(), quality_, TrackQuality::Layer::Shared, TrackCategories::SharedInnerHits, HistoryBase::simParticleTrail(), TrackQuality::Layer::Hit::state, and tracer_.
Referenced by evaluate().
{ // run the hit-by-hit reconstruction quality analysis quality_.evaluate(tracer_.simParticleTrail(), track); unsigned int maxLayers = std::min(numberOfInnerLayers_, quality_.numberOfLayers()); // check the innermost layers for bad hits for (unsigned int i = 0; i < maxLayers; i++) { const TrackQuality::Layer &layer = quality_.layer(i); // check all hits in that layer for (unsigned int j = 0; j < layer.hits.size(); j++) { const TrackQuality::Layer::Hit &hit = layer.hits[j]; // In those cases the bad hit was used by track reconstruction if (hit.state == TrackQuality::Layer::Noise || hit.state == TrackQuality::Layer::Misassoc) flags_[BadInnerHits] = true; else if (hit.state == TrackQuality::Layer::Shared) flags_[SharedInnerHits] = true; } } }
void TrackClassifier::reconstructionInformation | ( | reco::TrackBaseRef const & | track | ) | [private] |
Classify all the tracks by their association and reconstruction information.
Definition at line 151 of file TrackClassifier.cc.
References abs, TrackCategories::Bad, badPull_, beamSpot_, funct::cos(), exception, TrackCategories::flags_, magneticField_, CoreSimTrack::momentum(), FreeTrajectoryState::momentum(), AlCaHLTBitMon_ParallelJobs::p, FreeTrajectoryState::position(), edm::ESHandle< T >::product(), HistoryBase::simParticle(), funct::sin(), TrajectoryStateClosestToPoint::theState(), tracer_, v, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().
Referenced by evaluate().
{ TrackingParticleRef tpr = tracer_.simParticle(); // Compute tracking particle parameters at point of closest approach to the beamline const SimTrack * assocTrack = &(*tpr->g4Track_begin()); FreeTrajectoryState ftsAtProduction( GlobalPoint( tpr->vertex().x(), tpr->vertex().y(), tpr->vertex().z() ), GlobalVector( assocTrack->momentum().x(), assocTrack->momentum().y(), assocTrack->momentum().z() ), TrackCharge(track->charge()), magneticField_.product() ); try { TSCPBuilderNoMaterial tscpBuilder; TrajectoryStateClosestToPoint tsAtClosestApproach = tscpBuilder( ftsAtProduction, GlobalPoint(beamSpot_->x0(), beamSpot_->y0(), beamSpot_->z0()) ); GlobalVector v = tsAtClosestApproach.theState().position() - GlobalPoint(beamSpot_->x0(), beamSpot_->y0(), beamSpot_->z0()); GlobalVector p = tsAtClosestApproach.theState().momentum(); // Simulated dxy double dxySim = -v.x()*sin(p.phi()) + v.y()*cos(p.phi()); // Simulated dz double dzSim = v.z() - (v.x()*p.x() + v.y()*p.y())*p.z()/p.perp2(); // Calculate the dxy pull double dxyPull = std::abs( track->dxy( reco::TrackBase::Point(beamSpot_->x0(), beamSpot_->y0(), beamSpot_->z0()) ) - dxySim ) / track->dxyError(); // Calculate the dx pull double dzPull = std::abs( track->dz( reco::TrackBase::Point(beamSpot_->x0(), beamSpot_->y0(), beamSpot_->z0()) ) - dzSim ) / track->dzError(); // Return true if d0Pull > badD0Pull sigmas flags_[Bad] = (dxyPull > badPull_ || dzPull > badPull_); } catch (cms::Exception exception) { flags_[Bad] = true; } }
void TrackClassifier::simulationInformation | ( | ) | [private] |
Get all the information related to the simulation details.
Definition at line 213 of file TrackClassifier.cc.
References abs, EncodedEventId::bunchCrossing(), EncodedEventId::event(), TrackCategories::flags_, minTrackerSimHits_, TrackCategories::Muon, TrackCategories::SignalEvent, HistoryBase::simParticle(), tracer_, and TrackCategories::TrackerSimHits.
Referenced by evaluate().
{ // Get the event id for the initial TP. EncodedEventId eventId = tracer_.simParticle()->eventId(); // Check for signal events flags_[SignalEvent] = !eventId.bunchCrossing() && !eventId.event(); // Check for muons flags_[Muon] = (abs(tracer_.simParticle()->pdgId()) == 13); // Check for the number of psimhit in tracker flags_[TrackerSimHits] = tracer_.simParticle()->matchedHit() >= (int)minTrackerSimHits_; }
void TrackClassifier::vertexInformation | ( | ) | [private] |
Get geometrical information about the vertices.
Definition at line 462 of file TrackClassifier.cc.
References TrackCategories::flags_, HistoryBase::genParticleTrail(), genpvs_, AlCaHLTBitMon_ParallelJobs::p, dbtoconf::parent, funct::pow(), TrackCategories::PrimaryVertex, TrackCategories::SecondaryVertex, HistoryBase::simParticleTrail(), TrackCategories::TertiaryVertex, tracer_, vertexClusteringSqDistance_, TrackClassifier::GeneratedPrimaryVertex::x, TrackClassifier::GeneratedPrimaryVertex::y, and TrackClassifier::GeneratedPrimaryVertex::z.
Referenced by evaluate().
{ // Get the main primary vertex from the list GeneratedPrimaryVertex const & genpv = genpvs_.back(); // Get the generated history of the tracks TrackHistory::GenParticleTrail & genParticleTrail = const_cast<TrackHistory::GenParticleTrail &> (tracer_.genParticleTrail()); // Vertex counter int counter = 0; // Unit transformation from mm to cm double const mm = 0.1; double oldX = genpv.x; double oldY = genpv.y; double oldZ = genpv.z; // Loop over the generated particles for ( TrackHistory::GenParticleTrail::reverse_iterator iparticle = genParticleTrail.rbegin(); iparticle != genParticleTrail.rend(); ++iparticle ) { // Look for those with production vertex HepMC::GenVertex * parent = (*iparticle)->production_vertex(); if (parent) { HepMC::ThreeVector p = parent->point3d(); double distance2 = pow(p.x() * mm - genpv.x, 2) + pow(p.y() * mm - genpv.y, 2) + pow(p.z() * mm - genpv.z, 2); double difference2 = pow(p.x() * mm - oldX, 2) + pow(p.y() * mm - oldY, 2) + pow(p.z() * mm - oldZ, 2); // std::cout << "Distance2 : " << distance2 << " (" << p.x() * mm << "," << p.y() * mm << "," << p.z() * mm << ")" << std::endl; // std::cout << "Difference2 : " << difference2 << std::endl; if ( difference2 > vertexClusteringSqDistance_ ) { if ( distance2 > vertexClusteringSqDistance_ ) counter++; oldX = p.x() * mm; oldY = p.y() * mm; oldZ = p.z() * mm; } } } TrackHistory::SimParticleTrail & simParticleTrail = const_cast<TrackHistory::SimParticleTrail &> (tracer_.simParticleTrail()); // Loop over the generated particles for ( TrackHistory::SimParticleTrail::reverse_iterator iparticle = simParticleTrail.rbegin(); iparticle != simParticleTrail.rend(); ++iparticle ) { // Look for those with production vertex TrackingParticle::Point p = (*iparticle)->vertex(); double distance2 = pow(p.x() - genpv.x, 2) + pow(p.y() - genpv.y, 2) + pow(p.z() - genpv.z, 2); double difference2 = pow(p.x() - oldX, 2) + pow(p.y() - oldY, 2) + pow(p.z() - oldZ, 2); // std::cout << "Distance2 : " << distance2 << " (" << p.x() << "," << p.y() << "," << p.z() << ")" << std::endl; // std::cout << "Difference2 : " << difference2 << std::endl; if ( difference2 > vertexClusteringSqDistance_ ) { if ( distance2 > vertexClusteringSqDistance_ ) counter++; oldX = p.x(); oldY = p.y(); oldZ = p.z(); } } if ( !counter ) flags_[PrimaryVertex] = true; else if ( counter == 1 ) flags_[SecondaryVertex] = true; else flags_[TertiaryVertex] = true; }
double TrackClassifier::badPull_ [private] |
Definition at line 69 of file TrackClassifier.h.
Referenced by reconstructionInformation(), and TrackClassifier().
edm::Handle<reco::BeamSpot> TrackClassifier::beamSpot_ [private] |
Definition at line 111 of file TrackClassifier.h.
Referenced by newEvent(), and reconstructionInformation().
const edm::InputTag TrackClassifier::beamSpotLabel_ [private] |
Definition at line 67 of file TrackClassifier.h.
Referenced by newEvent().
std::vector<GeneratedPrimaryVertex> TrackClassifier::genpvs_ [private] |
Definition at line 156 of file TrackClassifier.h.
Referenced by genPrimaryVertices(), and vertexInformation().
const edm::InputTag TrackClassifier::hepMCLabel_ [private] |
Definition at line 66 of file TrackClassifier.h.
Referenced by newEvent().
double TrackClassifier::longLivedDecayLength_ [private] |
Definition at line 70 of file TrackClassifier.h.
Referenced by processesAtGenerator(), processesAtSimulation(), and TrackClassifier().
Definition at line 103 of file TrackClassifier.h.
Referenced by newEvent(), and reconstructionInformation().
Definition at line 105 of file TrackClassifier.h.
Referenced by genPrimaryVertices(), and newEvent().
unsigned int TrackClassifier::minTrackerSimHits_ [private] |
Definition at line 73 of file TrackClassifier.h.
Referenced by simulationInformation(), and TrackClassifier().
unsigned int TrackClassifier::numberOfInnerLayers_ [private] |
Definition at line 72 of file TrackClassifier.h.
Referenced by qualityInformation(), and TrackClassifier().
Definition at line 107 of file TrackClassifier.h.
Referenced by isCharged(), newEvent(), processesAtGenerator(), and processesAtSimulation().
TrackQuality TrackClassifier::quality_ [private] |
Definition at line 101 of file TrackClassifier.h.
Referenced by newEvent(), quality(), and qualityInformation().
TrackHistory TrackClassifier::tracer_ [private] |
Definition at line 99 of file TrackClassifier.h.
Referenced by evaluate(), hadronFlavor(), history(), newEvent(), processesAtGenerator(), processesAtSimulation(), qualityInformation(), reconstructionInformation(), simulationInformation(), TrackClassifier(), and vertexInformation().
Definition at line 109 of file TrackClassifier.h.
Referenced by newEvent().
double TrackClassifier::vertexClusteringSqDistance_ [private] |
Definition at line 71 of file TrackClassifier.h.
Referenced by genPrimaryVertices(), TrackClassifier(), and vertexInformation().