CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/SimTracker/TrackHistory/src/VertexClassifier.cc

Go to the documentation of this file.
00001 /*
00002  *  VertexClassifier.C
00003  */
00004 
00005 #include <math.h>
00006 #include <cstdlib>
00007 #include <iostream>
00008 
00009 #include "HepPDT/ParticleID.hh"
00010 
00011 #include "SimTracker/TrackHistory/interface/VertexClassifier.h"
00012 
00013 
00014 #define update(a, b) do { (a) = (a) | (b); } while(0)
00015 
00016 
00017 VertexClassifier::VertexClassifier(edm::ParameterSet const & config) : VertexCategories(), tracer_(config),
00018         hepMCLabel_( config.getUntrackedParameter<edm::InputTag>("hepMC") )
00019 {
00020     // Set the history depth after hadronization
00021     tracer_.depth(-2);
00022 
00023     // Set the minimum decay length for detecting long decays
00024     longLivedDecayLength_ = config.getUntrackedParameter<double>("longLivedDecayLength");
00025 
00026     // Set the distance for clustering vertices
00027     vertexClusteringDistance_ = config.getUntrackedParameter<double>("vertexClusteringDistance");
00028 }
00029 
00030 
00031 void VertexClassifier::newEvent(edm::Event const & event, edm::EventSetup const & setup)
00032 {
00033     // Get the new event information for the tracer
00034     tracer_.newEvent(event, setup);
00035 
00036     // Get hepmc of the event
00037     event.getByLabel(hepMCLabel_, mcInformation_);
00038 
00039     // Get the partivle data table
00040     setup.getData(particleDataTable_);
00041 
00042     // Create the list of primary vertices associated to the event
00043     genPrimaryVertices();
00044 }
00045 
00046 
00047 VertexClassifier const & VertexClassifier::evaluate (reco::VertexBaseRef const & vertex)
00048 {
00049     // Initializing the category vector
00050     reset();
00051 
00052     // Associate and evaluate the vertex history (check for fakes)
00053     if ( tracer_.evaluate(vertex) )
00054     {
00055         // Get all the information related to the simulation details
00056         simulationInformation();
00057 
00058         // Get all the information related to decay process
00059         processesAtGenerator();
00060 
00061         // Get information about conversion and other interactions
00062         processesAtSimulation();
00063 
00064         // Get geometrical information about the vertices
00065         vertexInformation();
00066 
00067         // Check for unkown classification
00068         unknownVertex();
00069     }
00070     else
00071         flags_[Fake] = true;
00072 
00073     return *this;
00074 }
00075 
00076 
00077 VertexClassifier const & VertexClassifier::evaluate (TrackingVertexRef const & vertex)
00078 {
00079     // Initializing the category vector
00080     reset();
00081 
00082     // Trace the history for the given TP
00083     tracer_.evaluate(vertex);
00084 
00085     // Check for a reconstructed track
00086     if ( tracer_.recoVertex().isNonnull() )
00087         flags_[Reconstructed] = true;
00088     else
00089         flags_[Reconstructed] = false;
00090 
00091     // Get all the information related to the simulation details
00092     simulationInformation();
00093 
00094     // Get all the information related to decay process
00095     processesAtGenerator();
00096 
00097     // Get information about conversion and other interactions
00098     processesAtSimulation();
00099 
00100     // Get geometrical information about the vertices
00101     vertexInformation();
00102 
00103     // Check for unkown classification
00104     unknownVertex();
00105 
00106     return *this;
00107 }
00108 
00109 
00110 void VertexClassifier::simulationInformation()
00111 {
00112     // Get the event id for the initial TP.
00113     EncodedEventId eventId = tracer_.simVertex()->eventId();
00114     // Check for signal events
00115     flags_[SignalEvent] = !eventId.bunchCrossing() && !eventId.event();
00116 }
00117 
00118 
00119 void VertexClassifier::processesAtGenerator()
00120 {
00121     // Get the generated vetices from track history
00122     VertexHistory::GenVertexTrail const & genVertexTrail = tracer_.genVertexTrail();
00123 
00124     // Loop over the generated vertices
00125     for (
00126         VertexHistory::GenVertexTrail::const_iterator ivertex = genVertexTrail.begin();
00127         ivertex != genVertexTrail.end();
00128         ++ivertex
00129     )
00130     {
00131         // Get the pointer to the vertex by removing the const-ness (no const methos in HepMC::GenVertex)
00132         HepMC::GenVertex * vertex = const_cast<HepMC::GenVertex *>(*ivertex);
00133 
00134         // Loop over the sources looking for specific decays
00135         for (
00136             HepMC::GenVertex::particle_iterator iparent = vertex->particles_begin(HepMC::parents);
00137             iparent != vertex->particles_end(HepMC::parents);
00138             ++iparent
00139         )
00140         {
00141             // Collect the pdgid of the parent
00142             int pdgid = std::abs((*iparent)->pdg_id());
00143             // Get particle type
00144             HepPDT::ParticleID particleID(pdgid);
00145 
00146             // Check if the particle type is valid one
00147             if (particleID.isValid())
00148             {
00149                 // Get particle data
00150                 ParticleData const * particleData = particleDataTable_->particle(particleID);
00151                 // Check if the particle exist in the table
00152                 if (particleData)
00153                 {
00154                     // Check if their life time is bigger than longLivedDecayLength_
00155                     if ( particleData->lifetime() > longLivedDecayLength_ )
00156                     {
00157                         // Check for B, C weak decays and long lived decays
00158                         update(flags_[BWeakDecay], particleID.hasBottom());
00159                         update(flags_[CWeakDecay], particleID.hasCharm());
00160                         update(flags_[LongLivedDecay], true);
00161                     }
00162                     // Check Tau, Ks and Lambda decay
00163                     update(flags_[TauDecay], pdgid == 15);
00164                     update(flags_[KsDecay], pdgid == 310);
00165                     update(flags_[LambdaDecay], pdgid == 3122);
00166                     update(flags_[JpsiDecay], pdgid == 443);
00167                     update(flags_[XiDecay], pdgid == 3312);
00168                     update(flags_[OmegaDecay], pdgid == 3334);
00169                     update(flags_[SigmaPlusDecay], pdgid == 3222);
00170                     update(flags_[SigmaMinusDecay], pdgid == 3112);
00171                 }
00172             }
00173         }
00174     }
00175 }
00176 
00177 
00178 void VertexClassifier::processesAtSimulation()
00179 {
00180     VertexHistory::SimVertexTrail const & simVertexTrail = tracer_.simVertexTrail();
00181 
00182     for (
00183         VertexHistory::SimVertexTrail::const_iterator ivertex = simVertexTrail.begin();
00184         ivertex != simVertexTrail.end();
00185         ++ivertex
00186     )
00187     {
00188         // pdgid of the real source parent vertex
00189         int pdgid = 0;
00190 
00191         // select the original source in case of combined vertices
00192         bool flag = false;
00193         TrackingVertex::tp_iterator itd, its;
00194 
00195         for (its = (*ivertex)->sourceTracks_begin(); its != (*ivertex)->sourceTracks_end(); ++its)
00196         {
00197             for (itd = (*ivertex)->daughterTracks_begin(); itd != (*ivertex)->daughterTracks_end(); ++itd)
00198                 if (itd != its)
00199                 {
00200                     flag = true;
00201                     break;
00202                 }
00203             if (flag)
00204                 break;
00205         }
00206         // Collect the pdgid of the original source track
00207         if ( its != (*ivertex)->sourceTracks_end() )
00208             pdgid = std::abs((*its)->pdgId());
00209         else
00210             pdgid = 0;
00211 
00212         // Loop over the simulated particles
00213         for (
00214             TrackingVertex::tp_iterator iparticle = (*ivertex)->daughterTracks_begin();
00215             iparticle != (*ivertex)->daughterTracks_end();
00216             ++iparticle
00217         )
00218         {
00219 
00220             if ( (*iparticle)->matchedHit() )
00221             {
00222                 // Collect the G4 process of the first psimhit (it should be the same for all of them)
00223                 unsigned short process = (*iparticle)->pSimHit_begin()->processType();
00224 
00225                 // Flagging all the different processes
00226 
00227                 update(
00228                     flags_[KnownProcess],
00229                     process != G4::Undefined &&
00230                     process != G4::Unknown &&
00231                     process != G4::Primary
00232                 );
00233 
00234                 update(flags_[UndefinedProcess], process == G4::Undefined);
00235                 update(flags_[UnknownProcess], process == G4::Unknown);
00236                 update(flags_[PrimaryProcess], process == G4::Primary);
00237                 update(flags_[HadronicProcess], process == G4::Hadronic);
00238                 update(flags_[DecayProcess], process == G4::Decay);
00239                 update(flags_[ComptonProcess], process == G4::Compton);
00240                 update(flags_[AnnihilationProcess], process == G4::Annihilation);
00241                 update(flags_[EIoniProcess], process == G4::EIoni);
00242                 update(flags_[HIoniProcess], process == G4::HIoni);
00243                 update(flags_[MuIoniProcess], process == G4::MuIoni);
00244                 update(flags_[PhotonProcess], process == G4::Photon);
00245                 update(flags_[MuPairProdProcess], process == G4::MuPairProd);
00246                 update(flags_[ConversionsProcess], process == G4::Conversions);
00247                 update(flags_[EBremProcess], process == G4::EBrem);
00248                 update(flags_[SynchrotronRadiationProcess], process == G4::SynchrotronRadiation);
00249                 update(flags_[MuBremProcess], process == G4::MuBrem);
00250                 update(flags_[MuNuclProcess], process == G4::MuNucl);
00251 
00252                 // Special treatment for decays
00253                 if (process == G4::Decay)
00254                 {
00255                     // Get particle type
00256                     HepPDT::ParticleID particleID(pdgid);
00257                     // Check if the particle type is valid one
00258                     if (particleID.isValid())
00259                     {
00260                         // Get particle data
00261                         ParticleData const * particleData = particleDataTable_->particle(particleID);
00262                         // Check if the particle exist in the table
00263                         if (particleData)
00264                         {
00265                             // Check if their life time is bigger than 1e-14
00266                             if ( particleDataTable_->particle(particleID)->lifetime() > longLivedDecayLength_ )
00267                             {
00268                                 // Check for B, C weak decays and long lived decays
00269                                 update(flags_[BWeakDecay], particleID.hasBottom());
00270                                 update(flags_[CWeakDecay], particleID.hasCharm());
00271                                 update(flags_[LongLivedDecay], true);
00272                             }
00273                             // Check Tau, Ks and Lambda decay
00274                             update(flags_[TauDecay], pdgid == 15);
00275                             update(flags_[KsDecay], pdgid == 310);
00276                             update(flags_[LambdaDecay], pdgid == 3122);
00277                             update(flags_[JpsiDecay], pdgid == 443);
00278                             update(flags_[XiDecay], pdgid == 3312);
00279                             update(flags_[OmegaDecay], pdgid == 3334);
00280                             update(flags_[SigmaPlusDecay], pdgid == 3222);
00281                             update(flags_[SigmaMinusDecay], pdgid == 3112);
00282                         }
00283                     }
00284                 }
00285             }
00286         }
00287     }
00288 }
00289 
00290 
00291 void VertexClassifier::vertexInformation()
00292 {
00293     // Helper class for clusterization
00294     typedef std::multimap<double, HepMC::ThreeVector> Clusters;
00295     typedef std::pair<double, HepMC::ThreeVector> ClusterPair;
00296 
00297     Clusters clusters;
00298 
00299     // Get the main primary vertex from the list
00300     GeneratedPrimaryVertex const & genpv = genpvs_.back();
00301 
00302     // Get the generated history of the tracks
00303     VertexHistory::GenVertexTrail & genVertexTrail = const_cast<VertexHistory::GenVertexTrail &> (tracer_.genVertexTrail());
00304 
00305     // Unit transformation from mm to cm
00306     double const mm = 0.1;
00307 
00308     // Loop over the generated vertexes
00309     for (
00310         VertexHistory::GenVertexTrail::const_iterator ivertex = genVertexTrail.begin();
00311         ivertex != genVertexTrail.end();
00312         ++ivertex
00313     )
00314     {
00315         // Check vertex exist
00316         if (*ivertex)
00317         {
00318             // Measure the distance2 respecto the primary vertex
00319             HepMC::ThreeVector p = (*ivertex)->point3d();
00320             double distance = sqrt( pow(p.x() * mm - genpv.x, 2) + pow(p.y() * mm - genpv.y, 2) + pow(p.z() * mm - genpv.z, 2) );
00321 
00322             // If there is not any clusters add the first vertex.
00323             if ( clusters.empty() )
00324             {
00325                 clusters.insert( ClusterPair(distance, HepMC::ThreeVector(p.x() * mm, p.y() * mm, p.z() * mm)) );
00326                 continue;
00327             }
00328 
00329             // Check if there is already a cluster in the given distance from primary vertex
00330             Clusters::const_iterator icluster = clusters.lower_bound(distance - vertexClusteringDistance_);
00331 
00332             if ( icluster == clusters.upper_bound(distance + vertexClusteringDistance_) )
00333             {
00334                 clusters.insert ( ClusterPair(distance, HepMC::ThreeVector(p.x() * mm, p.y() * mm, p.z() * mm)) );
00335                 continue;
00336             }
00337 
00338             bool cluster = false;
00339 
00340             // Looping over the vertex clusters of a given distance from primary vertex
00341             for (;
00342                     icluster != clusters.upper_bound(distance + vertexClusteringDistance_);
00343                     ++icluster
00344                 )
00345             {
00346                 double difference = sqrt (
00347                                         pow(p.x() * mm - icluster->second.x(), 2) +
00348                                         pow(p.y() * mm - icluster->second.y(), 2) +
00349                                         pow(p.z() * mm - icluster->second.z(), 2)
00350                                     );
00351 
00352                 if ( difference < vertexClusteringDistance_ )
00353                 {
00354                     cluster = true;
00355                     break;
00356                 }
00357             }
00358 
00359             if (!cluster) clusters.insert ( ClusterPair(distance, HepMC::ThreeVector(p.x() * mm, p.y() * mm, p.z() * mm)) );
00360         }
00361     }
00362 
00363     VertexHistory::SimVertexTrail & simVertexTrail = const_cast<VertexHistory::SimVertexTrail &> (tracer_.simVertexTrail());
00364 
00365     // Loop over the generated particles
00366     for (
00367         VertexHistory::SimVertexTrail::reverse_iterator ivertex = simVertexTrail.rbegin();
00368         ivertex != simVertexTrail.rend();
00369         ++ivertex
00370     )
00371     {
00372         // Look for those with production vertex
00373         TrackingVertex::LorentzVector p = (*ivertex)->position();
00374 
00375         double distance = sqrt( pow(p.x() - genpv.x, 2) + pow(p.y() - genpv.y, 2) + pow(p.z() - genpv.z, 2) );
00376 
00377         // If there is not any clusters add the first vertex.
00378         if ( clusters.empty() )
00379         {
00380             clusters.insert( ClusterPair(distance, HepMC::ThreeVector(p.x(), p.y(), p.z())) );
00381             continue;
00382         }
00383 
00384         // Check if there is already a cluster in the given distance from primary vertex
00385         Clusters::const_iterator icluster = clusters.lower_bound(distance - vertexClusteringDistance_);
00386 
00387         if ( icluster == clusters.upper_bound(distance + vertexClusteringDistance_) )
00388         {
00389             clusters.insert ( ClusterPair(distance, HepMC::ThreeVector(p.x(), p.y(), p.z())) );
00390             continue;
00391         }
00392 
00393         bool cluster = false;
00394 
00395         // Looping over the vertex clusters of a given distance from primary vertex
00396         for (;
00397                 icluster != clusters.upper_bound(distance + vertexClusteringDistance_);
00398                 ++icluster
00399             )
00400         {
00401             double difference = sqrt (
00402                                     pow(p.x() - icluster->second.x(), 2) +
00403                                     pow(p.y() - icluster->second.y(), 2) +
00404                                     pow(p.z() - icluster->second.z(), 2)
00405                                 );
00406 
00407             if ( difference < vertexClusteringDistance_ )
00408             {
00409                 cluster = true;
00410                 break;
00411             }
00412         }
00413 
00414         if (!cluster) clusters.insert ( ClusterPair(distance, HepMC::ThreeVector(p.x(), p.y(), p.z())) );
00415     }
00416 
00417     if ( clusters.size() == 1 )
00418         flags_[PrimaryVertex] = true;
00419     else if ( clusters.size() == 2 )
00420         flags_[SecondaryVertex] = true;
00421     else
00422         flags_[TertiaryVertex] = true;
00423 }
00424 
00425 
00426 bool VertexClassifier::isFinalstateParticle(const HepMC::GenParticle * p)
00427 {
00428     return !p->end_vertex() && p->status() == 1;
00429 }
00430 
00431 
00432 bool VertexClassifier::isCharged(const HepMC::GenParticle * p)
00433 {
00434     const ParticleData * part = particleDataTable_->particle( p->pdg_id() );
00435     if (part)
00436         return part->charge()!=0;
00437     else
00438     {
00439         // the new/improved particle table doesn't know anti-particles
00440         return  particleDataTable_->particle( -p->pdg_id() ) != 0;
00441     }
00442 }
00443 
00444 
00445 void VertexClassifier::genPrimaryVertices()
00446 {
00447     genpvs_.clear();
00448 
00449     const HepMC::GenEvent * event = mcInformation_->GetEvent();
00450 
00451     if (event)
00452     {
00453         int idx = 0;
00454 
00455         // Loop over the different GenVertex
00456         for ( HepMC::GenEvent::vertex_const_iterator ivertex = event->vertices_begin(); ivertex != event->vertices_end(); ++ivertex )
00457         {
00458             bool hasParentVertex = false;
00459 
00460             // Loop over the parents looking to see if they are coming from a production vertex
00461             for (
00462                 HepMC::GenVertex::particle_iterator iparent = (*ivertex)->particles_begin(HepMC::parents);
00463                 iparent != (*ivertex)->particles_end(HepMC::parents);
00464                 ++iparent
00465             )
00466                 if ( (*iparent)->production_vertex() )
00467                 {
00468                     hasParentVertex = true;
00469                     break;
00470                 }
00471 
00472             // Reject those vertices with parent vertices
00473             if (hasParentVertex) continue;
00474 
00475             // Get the position of the vertex
00476             HepMC::FourVector pos = (*ivertex)->position();
00477 
00478             double const mm = 0.1;
00479 
00480             GeneratedPrimaryVertex pv(pos.x()*mm, pos.y()*mm, pos.z()*mm);
00481 
00482             std::vector<GeneratedPrimaryVertex>::iterator ientry = genpvs_.begin();
00483 
00484             // Search for a VERY close vertex in the list
00485             for (; ientry != genpvs_.end(); ++ientry)
00486             {
00487                 double distance = sqrt( pow(pv.x - ientry->x, 2) + pow(pv.y - ientry->y, 2) + pow(pv.z - ientry->z, 2) );
00488                 if ( distance < vertexClusteringDistance_ )
00489                     break;
00490             }
00491 
00492             // Check if there is not a VERY close vertex and added to the list
00493             if (ientry == genpvs_.end())
00494                 ientry = genpvs_.insert(ientry,pv);
00495 
00496             // Add the vertex barcodes to the new or existent vertices
00497             ientry->genVertex.push_back((*ivertex)->barcode());
00498 
00499             // Collect final state descendants
00500             for (
00501                 HepMC::GenVertex::particle_iterator idecendants  = (*ivertex)->particles_begin(HepMC::descendants);
00502                 idecendants != (*ivertex)->particles_end(HepMC::descendants);
00503                 ++idecendants
00504             )
00505             {
00506                 if (isFinalstateParticle(*idecendants))
00507                     if ( find(ientry->finalstateParticles.begin(), ientry->finalstateParticles.end(), (*idecendants)->barcode()) == ientry->finalstateParticles.end() )
00508                     {
00509                         ientry->finalstateParticles.push_back((*idecendants)->barcode());
00510                         HepMC::FourVector m = (*idecendants)->momentum();
00511 
00512                         ientry->ptot.setPx(ientry->ptot.px() + m.px());
00513                         ientry->ptot.setPy(ientry->ptot.py() + m.py());
00514                         ientry->ptot.setPz(ientry->ptot.pz() + m.pz());
00515                         ientry->ptot.setE(ientry->ptot.e() + m.e());
00516                         ientry->ptsq += m.perp() * m.perp();
00517 
00518                         if ( m.perp() > 0.8 && std::abs(m.pseudoRapidity()) < 2.5 && isCharged(*idecendants) ) ientry->nGenTrk++;
00519                     }
00520             }
00521             idx++;
00522         }
00523     }
00524 
00525     std::sort(genpvs_.begin(), genpvs_.end());
00526 }