CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/SimTracker/TrackHistory/src/TrackClassifier.cc

Go to the documentation of this file.
00001 
00002 #include <math.h>
00003 #include <cstdlib>
00004 #include <iostream>
00005 
00006 #include "HepPDT/ParticleID.hh"
00007 
00008 #include "SimTracker/TrackHistory/interface/TrackClassifier.h"
00009 
00010 #define update(a, b) do { (a) = (a) | (b); } while(0)
00011 
00012 TrackClassifier::TrackClassifier(edm::ParameterSet const & config) : TrackCategories(),
00013         hepMCLabel_( config.getUntrackedParameter<edm::InputTag>("hepMC") ),
00014         beamSpotLabel_( config.getUntrackedParameter<edm::InputTag>("beamSpot") ),
00015         tracer_(config),
00016         quality_(config)
00017 {
00018     // Set the history depth after hadronization
00019     tracer_.depth(-2);
00020 
00021     // Set the maximum d0pull for the bad category
00022     badPull_ = config.getUntrackedParameter<double>("badPull");
00023 
00024     // Set the minimum decay length for detecting long decays
00025     longLivedDecayLength_ = config.getUntrackedParameter<double>("longLivedDecayLength");
00026 
00027     // Set the distance for clustering vertices
00028     float vertexClusteringDistance = config.getUntrackedParameter<double>("vertexClusteringDistance");
00029     vertexClusteringSqDistance_ = vertexClusteringDistance * vertexClusteringDistance;
00030 
00031     // Set the number of innermost layers to check for bad hits
00032     numberOfInnerLayers_ = config.getUntrackedParameter<unsigned int>("numberOfInnerLayers");
00033 
00034     // Set the minimum number of simhits in the tracker
00035     minTrackerSimHits_ = config.getUntrackedParameter<unsigned int>("minTrackerSimHits");
00036 }
00037 
00038 
00039 void TrackClassifier::newEvent ( edm::Event const & event, edm::EventSetup const & setup )
00040 {
00041     // Get the new event information for the tracer
00042     tracer_.newEvent(event, setup);
00043 
00044     // Get the new event information for the track quality analyser
00045     quality_.newEvent(event, setup);
00046 
00047     // Get hepmc of the event
00048     event.getByLabel(hepMCLabel_, mcInformation_);
00049 
00050     // Magnetic field
00051     setup.get<IdealMagneticFieldRecord>().get(magneticField_);
00052 
00053     // Get the partivle data table
00054     setup.getData(particleDataTable_);
00055 
00056     // get the beam spot
00057     event.getByLabel(beamSpotLabel_, beamSpot_);
00058 
00059     // Transient track builder
00060     setup.get<TransientTrackRecord>().get("TransientTrackBuilder", transientTrackBuilder_);
00061 
00062     // Create the list of primary vertices associated to the event
00063     genPrimaryVertices();
00064 }
00065 
00066 
00067 TrackClassifier const & TrackClassifier::evaluate (reco::TrackBaseRef const & track)
00068 {
00069     // Initializing the category vector
00070     reset();
00071 
00072     // Associate and evaluate the track history (check for fakes)
00073     if ( tracer_.evaluate(track) )
00074     {
00075         // Classify all the tracks by their association and reconstruction information
00076         reconstructionInformation(track);
00077 
00078         // Get all the information related to the simulation details
00079         simulationInformation();
00080 
00081         // Analyse the track reconstruction quality
00082         qualityInformation(track);
00083 
00084         // Get hadron flavor of the initial hadron
00085         hadronFlavor();
00086 
00087         // Get all the information related to decay process
00088         processesAtGenerator();
00089 
00090         // Get information about conversion and other interactions
00091         processesAtSimulation();
00092 
00093         // Get geometrical information about the vertices
00094         vertexInformation();
00095 
00096         // Check for unkown classification
00097         unknownTrack();
00098     }
00099     else
00100         flags_[Fake] = true;
00101 
00102     return *this;
00103 }
00104 
00105 
00106 TrackClassifier const & TrackClassifier::evaluate (TrackingParticleRef const & track)
00107 {
00108     // Initializing the category vector
00109     reset();
00110 
00111     // Trace the history for the given TP
00112     tracer_.evaluate(track);
00113 
00114     // Collect the associated reco track
00115     const reco::TrackBaseRef & recotrack = tracer_.recoTrack();
00116 
00117     // If there is a reco truck then evaluate the simulated history
00118     if ( recotrack.isNonnull() )
00119     {
00120         flags_[Reconstructed] = true;
00121         // Classify all the tracks by their association and reconstruction information
00122         reconstructionInformation(recotrack);
00123         // Analyse the track reconstruction quality
00124         qualityInformation(recotrack);
00125     }
00126     else
00127         flags_[Reconstructed] = false;
00128 
00129     // Get all the information related to the simulation details
00130     simulationInformation();
00131 
00132     // Get hadron flavor of the initial hadron
00133     hadronFlavor();
00134 
00135     // Get all the information related to decay process
00136     processesAtGenerator();
00137 
00138     // Get information about conversion and other interactions
00139     processesAtSimulation();
00140 
00141     // Get geometrical information about the vertices
00142     vertexInformation();
00143 
00144     // Check for unkown classification
00145     unknownTrack();
00146 
00147     return *this;
00148 }
00149 
00150 
00151 void TrackClassifier::reconstructionInformation(reco::TrackBaseRef const & track)
00152 {
00153     TrackingParticleRef tpr = tracer_.simParticle();
00154 
00155     // Compute tracking particle parameters at point of closest approach to the beamline
00156 
00157     const SimTrack * assocTrack = &(*tpr->g4Track_begin());
00158 
00159     FreeTrajectoryState ftsAtProduction(
00160         GlobalPoint(
00161             tpr->vertex().x(),
00162             tpr->vertex().y(),
00163             tpr->vertex().z()
00164         ),
00165         GlobalVector(
00166             assocTrack->momentum().x(),
00167             assocTrack->momentum().y(),
00168             assocTrack->momentum().z()
00169         ),
00170         TrackCharge(track->charge()),
00171         magneticField_.product()
00172     );
00173 
00174     try
00175     {
00176         TSCPBuilderNoMaterial tscpBuilder;
00177         TrajectoryStateClosestToPoint tsAtClosestApproach = tscpBuilder(
00178                     ftsAtProduction,
00179                     GlobalPoint(beamSpot_->x0(), beamSpot_->y0(), beamSpot_->z0())
00180                 );
00181 
00182         GlobalVector v = tsAtClosestApproach.theState().position()
00183                          - GlobalPoint(beamSpot_->x0(), beamSpot_->y0(), beamSpot_->z0());
00184         GlobalVector p = tsAtClosestApproach.theState().momentum();
00185 
00186         // Simulated dxy
00187         double dxySim = -v.x()*sin(p.phi()) + v.y()*cos(p.phi());
00188 
00189         // Simulated dz
00190         double dzSim = v.z() - (v.x()*p.x() + v.y()*p.y())*p.z()/p.perp2();
00191 
00192         // Calculate the dxy pull
00193         double dxyPull = std::abs(
00194                              track->dxy( reco::TrackBase::Point(beamSpot_->x0(), beamSpot_->y0(), beamSpot_->z0()) ) - dxySim
00195                          ) / track->dxyError();
00196 
00197         // Calculate the dx pull
00198         double dzPull = std::abs(
00199                             track->dz( reco::TrackBase::Point(beamSpot_->x0(), beamSpot_->y0(), beamSpot_->z0()) ) - dzSim
00200                         ) / track->dzError();
00201 
00202         // Return true if d0Pull > badD0Pull sigmas
00203         flags_[Bad] = (dxyPull > badPull_ || dzPull > badPull_);
00204 
00205     }
00206     catch (cms::Exception exception)
00207     {
00208         flags_[Bad] = true;
00209     }
00210 }
00211 
00212 
00213 void TrackClassifier::simulationInformation()
00214 {
00215     // Get the event id for the initial TP.
00216     EncodedEventId eventId = tracer_.simParticle()->eventId();
00217     // Check for signal events
00218     flags_[SignalEvent] = !eventId.bunchCrossing() && !eventId.event();
00219     // Check for muons
00220     flags_[Muon] = (abs(tracer_.simParticle()->pdgId()) == 13);
00221     // Check for the number of psimhit in tracker
00222     flags_[TrackerSimHits] = tracer_.simParticle()->matchedHit() >= (int)minTrackerSimHits_;
00223 }
00224 
00225 
00226 void TrackClassifier::qualityInformation(reco::TrackBaseRef const & track)
00227 {
00228     // run the hit-by-hit reconstruction quality analysis
00229     quality_.evaluate(tracer_.simParticleTrail(), track);
00230 
00231     unsigned int maxLayers = std::min(numberOfInnerLayers_, quality_.numberOfLayers());
00232 
00233     // check the innermost layers for bad hits
00234     for (unsigned int i = 0; i < maxLayers; i++)
00235     {
00236         const TrackQuality::Layer &layer = quality_.layer(i);
00237 
00238         // check all hits in that layer
00239         for (unsigned int j = 0; j < layer.hits.size(); j++)
00240         {
00241             const TrackQuality::Layer::Hit &hit = layer.hits[j];
00242 
00243             // In those cases the bad hit was used by track reconstruction
00244             if (hit.state == TrackQuality::Layer::Noise ||
00245                     hit.state == TrackQuality::Layer::Misassoc)
00246                 flags_[BadInnerHits] = true;
00247             else if (hit.state == TrackQuality::Layer::Shared)
00248                 flags_[SharedInnerHits] = true;
00249         }
00250     }
00251 }
00252 
00253 
00254 void TrackClassifier::hadronFlavor()
00255 {
00256     // Get the initial hadron
00257     const HepMC::GenParticle * particle = tracer_.genParticle();
00258 
00259     // Check for the initial hadron
00260     if (particle)
00261     {
00262         HepPDT::ParticleID pid(particle->pdg_id());
00263         flags_[Bottom] = pid.hasBottom();
00264         flags_[Charm] =  pid.hasCharm();
00265         flags_[Light] = !pid.hasCharm() && !pid.hasBottom();
00266     }
00267 }
00268 
00269 
00270 void TrackClassifier::processesAtGenerator()
00271 {
00272     // pdgid of the "in" particle to the production vertex
00273     int pdgid = 0;
00274 
00275     // Get the generated particles from track history
00276     TrackHistory::GenParticleTrail const & genParticleTrail = tracer_.genParticleTrail();
00277 
00278     // Loop over the generated particles
00279     for (TrackHistory::GenParticleTrail::const_iterator iparticle = genParticleTrail.begin(); iparticle != genParticleTrail.end(); ++iparticle)
00280     {
00281         // Get the source vertex for the particle
00282         HepMC::GenVertex * productionVertex = (*iparticle)->production_vertex();
00283 
00284         // Get the pointer to the vertex by removing the const-ness (no const methos in HepMC::GenVertex)
00285         // HepMC::GenVertex * vertex = const_cast<HepMC::GenVertex *>(*ivertex);
00286 
00287         // Check for a non-null pointer to the production vertex
00288         if (productionVertex)
00289         {
00290             // Only case track history will navegate (one in or source particle per vertex)
00291             if ( productionVertex->particles_in_size() == 1 )
00292             {
00293                 // Look at the pdgid of the first "in" particle to the vertex
00294                 pdgid = std::abs((*productionVertex->particles_in_const_begin())->pdg_id());
00295                 // Get particle type
00296                 HepPDT::ParticleID particleID(pdgid);
00297 
00298                 // Check if the particle type is valid one
00299                 if (particleID.isValid())
00300                 {
00301                     // Get particle data
00302                     ParticleData const * particleData = particleDataTable_->particle(particleID);
00303                     // Check if the particle exist in the table
00304                     if (particleData)
00305                     {
00306                         // Check if their life time is bigger than longLivedDecayLength_
00307                         if ( particleData->lifetime() > longLivedDecayLength_ )
00308                             update(flags_[LongLivedDecay], true);
00309                         // Check for B and C weak decays
00310                         update(flags_[BWeakDecay], particleID.hasBottom());
00311                         update(flags_[CWeakDecay], particleID.hasCharm());
00312                         // Check for B and C pure leptonic decay
00313                         int daughterId = abs((*iparticle)->pdg_id());
00314                         update(flags_[FromBWeakDecayMuon], particleID.hasBottom() && daughterId == 13);
00315                         update(flags_[FromCWeakDecayMuon], particleID.hasCharm() && daughterId == 13);
00316                     }
00317                     // Check Tau, Ks and Lambda decay
00318                     update(flags_[ChargePionDecay], pdgid == 211);
00319                     update(flags_[ChargeKaonDecay], pdgid == 321);
00320                     update(flags_[TauDecay], pdgid == 15);
00321                     update(flags_[KsDecay], pdgid == 310);
00322                     update(flags_[LambdaDecay], pdgid == 3122);
00323                     update(flags_[JpsiDecay], pdgid == 443);
00324                     update(flags_[XiDecay], pdgid == 3312);
00325                     update(flags_[SigmaPlusDecay], pdgid == 3222);
00326                     update(flags_[SigmaMinusDecay], pdgid == 3112);
00327                 }
00328             }
00329         }
00330     }
00331     // Decays in flight
00332     update(flags_[FromChargePionMuon], flags_[Muon] && flags_[ChargePionDecay]);
00333     update(flags_[FromChargeKaonMuon], flags_[Muon] && flags_[ChargeKaonDecay]);
00334     update(flags_[DecayOnFlightMuon], (flags_[FromChargePionMuon] || flags_[FromChargeKaonMuon]));
00335 }
00336 
00337 
00338 void TrackClassifier::processesAtSimulation()
00339 {
00340     TrackHistory::SimParticleTrail const & simParticleTrail = tracer_.simParticleTrail();
00341 
00342     // Loop over the simulated particles
00343     for (
00344         TrackHistory::SimParticleTrail::const_iterator iparticle = simParticleTrail.begin();
00345         iparticle != simParticleTrail.end();
00346         ++iparticle
00347     )
00348     {
00349         // pdgid of the real source parent vertex
00350         int pdgid = 0;
00351 
00352         // Get a reference to the TP's parent vertex
00353         TrackingVertexRef const & parentVertex = (*iparticle)->parentVertex();
00354 
00355         // Look for the original source track
00356         if ( parentVertex.isNonnull() )
00357         {
00358             // select the original source in case of combined vertices
00359             bool flag = false;
00360             TrackingVertex::tp_iterator itd, its;
00361 
00362             for (its = parentVertex->sourceTracks_begin(); its != parentVertex->sourceTracks_end(); ++its)
00363             {
00364                 for (itd = parentVertex->daughterTracks_begin(); itd != parentVertex->daughterTracks_end(); ++itd)
00365                     if (itd != its)
00366                     {
00367                         flag = true;
00368                         break;
00369                     }
00370                 if (flag)
00371                     break;
00372             }
00373 
00374             // Collect the pdgid of the original source track
00375             if ( its != parentVertex->sourceTracks_end() )
00376                 pdgid = std::abs((*its)->pdgId());
00377             else
00378                 pdgid = 0;
00379         }
00380 
00381         // Check for the number of psimhit if different from zero
00382         if ((*iparticle)->trackPSimHit().empty()) continue;
00383 
00384         // Collect the G4 process of the first psimhit (it should be the same for all of them)
00385         unsigned short process = (*iparticle)->pSimHit_begin()->processType();
00386 
00387         // Flagging all the different processes
00388 
00389         update(
00390             flags_[KnownProcess],
00391             process != G4::Undefined &&
00392             process != G4::Unknown &&
00393             process != G4::Primary
00394         );
00395 
00396         update(flags_[UndefinedProcess], process == G4::Undefined);
00397         update(flags_[UnknownProcess], process == G4::Unknown);
00398         update(flags_[PrimaryProcess], process == G4::Primary);
00399         update(flags_[HadronicProcess], process == G4::Hadronic);
00400         update(flags_[DecayProcess], process == G4::Decay);
00401         update(flags_[ComptonProcess], process == G4::Compton);
00402         update(flags_[AnnihilationProcess], process == G4::Annihilation);
00403         update(flags_[EIoniProcess], process == G4::EIoni);
00404         update(flags_[HIoniProcess], process == G4::HIoni);
00405         update(flags_[MuIoniProcess], process == G4::MuIoni);
00406         update(flags_[PhotonProcess], process == G4::Photon);
00407         update(flags_[MuPairProdProcess], process == G4::MuPairProd);
00408         update(flags_[ConversionsProcess], process == G4::Conversions);
00409         update(flags_[EBremProcess], process == G4::EBrem);
00410         update(flags_[SynchrotronRadiationProcess], process == G4::SynchrotronRadiation);
00411         update(flags_[MuBremProcess], process == G4::MuBrem);
00412         update(flags_[MuNuclProcess], process == G4::MuNucl);
00413 
00414         // Get particle type
00415         HepPDT::ParticleID particleID(pdgid);
00416 
00417         // Check if the particle type is valid one
00418         if (particleID.isValid())
00419         {
00420             // Get particle data
00421             ParticleData const * particleData = particleDataTable_->particle(particleID);
00422             // Special treatment for decays
00423             if (process == G4::Decay)
00424             {
00425                 // Check if the particle exist in the table
00426                 if (particleData)
00427                 {
00428                     // Check if their life time is bigger than 1e-14
00429                     if ( particleDataTable_->particle(particleID)->lifetime() > longLivedDecayLength_ )
00430                         update(flags_[LongLivedDecay], true);
00431 
00432                     // Check for B and C weak decays
00433                     update(flags_[BWeakDecay], particleID.hasBottom());
00434                     update(flags_[CWeakDecay], particleID.hasCharm());
00435 
00436                     // Check for B or C pure leptonic decays
00437                     int daughtId = abs((*iparticle)->pdgId());
00438                     update(flags_[FromBWeakDecayMuon], particleID.hasBottom() && daughtId == 13);
00439                     update(flags_[FromCWeakDecayMuon], particleID.hasCharm() && daughtId == 13);
00440                 }
00441                 // Check decays
00442                 update(flags_[ChargePionDecay], pdgid == 211);
00443                 update(flags_[ChargeKaonDecay], pdgid == 321);
00444                 update(flags_[TauDecay], pdgid == 15);
00445                 update(flags_[KsDecay], pdgid == 310);
00446                 update(flags_[LambdaDecay], pdgid == 3122);
00447                 update(flags_[JpsiDecay], pdgid == 443);
00448                 update(flags_[XiDecay], pdgid == 3312);
00449                 update(flags_[OmegaDecay], pdgid == 3334);
00450                 update(flags_[SigmaPlusDecay], pdgid == 3222);
00451                 update(flags_[SigmaMinusDecay], pdgid == 3112);
00452             }
00453         }
00454     }
00455     // Decays in flight
00456     update(flags_[FromChargePionMuon], flags_[Muon] && flags_[ChargePionDecay]);
00457     update(flags_[FromChargeKaonMuon], flags_[Muon] && flags_[ChargeKaonDecay]);
00458     update(flags_[DecayOnFlightMuon], flags_[FromChargePionMuon] || flags_[FromChargeKaonMuon]);
00459 }
00460 
00461 
00462 void TrackClassifier::vertexInformation()
00463 {
00464     // Get the main primary vertex from the list
00465     GeneratedPrimaryVertex const & genpv = genpvs_.back();
00466 
00467     // Get the generated history of the tracks
00468     TrackHistory::GenParticleTrail & genParticleTrail = const_cast<TrackHistory::GenParticleTrail &> (tracer_.genParticleTrail());
00469 
00470     // Vertex counter
00471     int counter = 0;
00472 
00473     // Unit transformation from mm to cm
00474     double const mm = 0.1;
00475 
00476     double oldX = genpv.x;
00477     double oldY = genpv.y;
00478     double oldZ = genpv.z;
00479 
00480     // Loop over the generated particles
00481     for (
00482         TrackHistory::GenParticleTrail::reverse_iterator iparticle = genParticleTrail.rbegin();
00483         iparticle != genParticleTrail.rend();
00484         ++iparticle
00485     )
00486     {
00487         // Look for those with production vertex
00488         HepMC::GenVertex * parent = (*iparticle)->production_vertex();
00489         if (parent)
00490         {
00491             HepMC::ThreeVector p = parent->point3d();
00492 
00493             double distance2   = pow(p.x() * mm - genpv.x, 2) + pow(p.y() * mm - genpv.y, 2) + pow(p.z() * mm - genpv.z, 2);
00494             double difference2 = pow(p.x() * mm - oldX, 2)    + pow(p.y() * mm - oldY, 2)    + pow(p.z() * mm - oldZ, 2);
00495 
00496             // std::cout << "Distance2 : " << distance2 << " (" << p.x() * mm << "," << p.y() * mm << "," << p.z() * mm << ")" << std::endl;
00497             // std::cout << "Difference2 : " << difference2 << std::endl;
00498 
00499             if ( difference2 > vertexClusteringSqDistance_ )
00500             {
00501                 if ( distance2 > vertexClusteringSqDistance_ ) counter++;
00502                 oldX = p.x() * mm;
00503                 oldY = p.y() * mm;
00504                 oldZ = p.z() * mm;
00505             }
00506         }
00507     }
00508 
00509     TrackHistory::SimParticleTrail & simParticleTrail = const_cast<TrackHistory::SimParticleTrail &> (tracer_.simParticleTrail());
00510 
00511     // Loop over the generated particles
00512     for (
00513         TrackHistory::SimParticleTrail::reverse_iterator iparticle = simParticleTrail.rbegin();
00514         iparticle != simParticleTrail.rend();
00515         ++iparticle
00516     )
00517     {
00518         // Look for those with production vertex
00519         TrackingParticle::Point p = (*iparticle)->vertex();
00520 
00521         double distance2   = pow(p.x() - genpv.x, 2) + pow(p.y() - genpv.y, 2) + pow(p.z() - genpv.z, 2);
00522         double difference2 = pow(p.x() - oldX, 2)    + pow(p.y() - oldY, 2)    + pow(p.z() - oldZ, 2);
00523 
00524         // std::cout << "Distance2 : " << distance2 << " (" << p.x() << "," << p.y() << "," << p.z() << ")" << std::endl;
00525         // std::cout << "Difference2 : " << difference2 << std::endl;
00526 
00527         if ( difference2 > vertexClusteringSqDistance_ )
00528         {
00529             if ( distance2 > vertexClusteringSqDistance_ ) counter++;
00530             oldX = p.x();
00531             oldY = p.y();
00532             oldZ = p.z();
00533         }
00534     }
00535 
00536     if ( !counter )
00537         flags_[PrimaryVertex] = true;
00538     else if ( counter == 1 )
00539         flags_[SecondaryVertex] = true;
00540     else
00541         flags_[TertiaryVertex] = true;
00542 }
00543 
00544 
00545 bool TrackClassifier::isFinalstateParticle(const HepMC::GenParticle * p)
00546 {
00547     return !p->end_vertex() && p->status() == 1;
00548 }
00549 
00550 
00551 bool TrackClassifier::isCharged(const HepMC::GenParticle * p)
00552 {
00553     const ParticleData * part = particleDataTable_->particle( p->pdg_id() );
00554     if (part)
00555         return part->charge()!=0;
00556     else
00557     {
00558         // the new/improved particle table doesn't know anti-particles
00559         return  particleDataTable_->particle( -p->pdg_id() ) != 0;
00560     }
00561 }
00562 
00563 
00564 void TrackClassifier::genPrimaryVertices()
00565 {
00566     genpvs_.clear();
00567 
00568     const HepMC::GenEvent * event = mcInformation_->GetEvent();
00569 
00570     if (event)
00571     {
00572         int idx = 0;
00573 
00574         // Loop over the different GenVertex
00575         for ( HepMC::GenEvent::vertex_const_iterator ivertex = event->vertices_begin(); ivertex != event->vertices_end(); ++ivertex )
00576         {
00577             bool hasParentVertex = false;
00578 
00579             // Loop over the parents looking to see if they are coming from a production vertex
00580             for (
00581                 HepMC::GenVertex::particle_iterator iparent = (*ivertex)->particles_begin(HepMC::parents);
00582                 iparent != (*ivertex)->particles_end(HepMC::parents);
00583                 ++iparent
00584             )
00585                 if ( (*iparent)->production_vertex() )
00586                 {
00587                     hasParentVertex = true;
00588                     break;
00589                 }
00590 
00591             // Reject those vertices with parent vertices
00592             if (hasParentVertex) continue;
00593 
00594             // Get the position of the vertex
00595             HepMC::FourVector pos = (*ivertex)->position();
00596 
00597             double const mm = 0.1;
00598 
00599             GeneratedPrimaryVertex pv(pos.x()*mm, pos.y()*mm, pos.z()*mm);
00600 
00601             std::vector<GeneratedPrimaryVertex>::iterator ientry = genpvs_.begin();
00602 
00603             // Search for a VERY close vertex in the list
00604             for (; ientry != genpvs_.end(); ++ientry)
00605             {
00606                 double distance2 = pow(pv.x - ientry->x, 2) + pow(pv.y - ientry->y, 2) + pow(pv.z - ientry->z, 2);
00607                 if ( distance2 < vertexClusteringSqDistance_ )
00608                     break;
00609             }
00610 
00611             // Check if there is not a VERY close vertex and added to the list
00612             if (ientry == genpvs_.end())
00613                 ientry = genpvs_.insert(ientry,pv);
00614 
00615             // Add the vertex barcodes to the new or existent vertices
00616             ientry->genVertex.push_back((*ivertex)->barcode());
00617 
00618             // Collect final state descendants
00619             for (
00620                 HepMC::GenVertex::particle_iterator idecendants  = (*ivertex)->particles_begin(HepMC::descendants);
00621                 idecendants != (*ivertex)->particles_end(HepMC::descendants);
00622                 ++idecendants
00623             )
00624             {
00625                 if (isFinalstateParticle(*idecendants))
00626                     if ( find(ientry->finalstateParticles.begin(), ientry->finalstateParticles.end(), (*idecendants)->barcode()) == ientry->finalstateParticles.end() )
00627                     {
00628                         ientry->finalstateParticles.push_back((*idecendants)->barcode());
00629                         HepMC::FourVector m = (*idecendants)->momentum();
00630 
00631                         ientry->ptot.setPx(ientry->ptot.px() + m.px());
00632                         ientry->ptot.setPy(ientry->ptot.py() + m.py());
00633                         ientry->ptot.setPz(ientry->ptot.pz() + m.pz());
00634                         ientry->ptot.setE(ientry->ptot.e() + m.e());
00635                         ientry->ptsq += m.perp() * m.perp();
00636 
00637                         if ( m.perp() > 0.8 && std::abs(m.pseudoRapidity()) < 2.5 && isCharged(*idecendants) ) ientry->nGenTrk++;
00638                     }
00639             }
00640             idx++;
00641         }
00642     }
00643 
00644     std::sort(genpvs_.begin(), genpvs_.end());
00645 }