CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10/src/SimGeneral/TrackingAnalysis/plugins/TrackingTruthProducer.cc

Go to the documentation of this file.
00001 
00002 #include "DataFormats/DetId/interface/DetId.h"
00003 #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
00004 #include "DataFormats/SiPixelDetId/interface/PXFDetId.h"
00005 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00006 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
00007 #include "DataFormats/SiStripDetId/interface/TECDetId.h"
00008 #include "DataFormats/SiStripDetId/interface/TIBDetId.h"
00009 #include "DataFormats/SiStripDetId/interface/TIDDetId.h"
00010 #include "DataFormats/SiStripDetId/interface/TOBDetId.h"
00011 
00012 #include "FWCore/Framework/interface/MakerMacros.h"
00013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00014 
00015 #include "SimDataFormats/EncodedEventId/interface/EncodedEventId.h"
00016 #include "SimDataFormats/TrackingAnalysis/interface/TrackingParticle.h"
00017 #include "SimDataFormats/TrackingAnalysis/interface/TrackingVertex.h"
00018 
00019 #include "SimGeneral/TrackingAnalysis/interface/TrackingTruthProducer.h"
00020 
00021 typedef edm::Ref<edm::HepMCProduct, HepMC::GenParticle > GenParticleRef;
00022 typedef edm::Ref<edm::HepMCProduct, HepMC::GenVertex >   GenVertexRef;
00023 typedef math::XYZTLorentzVectorD    LorentzVector;
00024 typedef math::XYZPoint Vector;
00025 
00026 TrackingTruthProducer::TrackingTruthProducer(const edm::ParameterSet & config) :
00027         pSimHitSelector_(config), pixelPSimHitSelector_(config), trackerPSimHitSelector_(config), muonPSimHitSelector_(config)
00028 {
00029     // Initialize global parameters
00030     dataLabels_             = config.getParameter<std::vector<std::string> >("HepMCDataLabels");
00031     useMultipleHepMCLabels_ = config.getParameter<bool>("useMultipleHepMCLabels");
00032 
00033     distanceCut_            = config.getParameter<double>("vertexDistanceCut");
00034     volumeRadius_           = config.getParameter<double>("volumeRadius");
00035     volumeZ_                = config.getParameter<double>("volumeZ");
00036     mergedBremsstrahlung_   = config.getParameter<bool>("mergedBremsstrahlung");
00037     removeDeadModules_      = config.getParameter<bool>("removeDeadModules");
00038     simHitLabel_            = config.getParameter<std::string>("simHitLabel");
00039 
00040     // Initialize selection for building TrackingParticles
00041     if ( config.exists("select") )
00042     {
00043         edm::ParameterSet param = config.getParameter<edm::ParameterSet>("select");
00044         selector_ = TrackingParticleSelector(
00045                         param.getParameter<double>("ptMinTP"),
00046                         param.getParameter<double>("minRapidityTP"),
00047                         param.getParameter<double>("maxRapidityTP"),
00048                         param.getParameter<double>("tipTP"),
00049                         param.getParameter<double>("lipTP"),
00050                         param.getParameter<int>("minHitTP"),
00051                         param.getParameter<bool>("signalOnlyTP"),
00052                         param.getParameter<bool>("chargedOnlyTP"),
00053                         param.getParameter<bool>("stableOnlyTP"),
00054                         param.getParameter<std::vector<int> >("pdgIdTP")
00055                     );
00056         selectorFlag_ = true;
00057     }
00058     else
00059         selectorFlag_ = false;
00060 
00061     MessageCategory_       = "TrackingTruthProducer";
00062 
00063     edm::LogInfo (MessageCategory_) << "Setting up TrackingTruthProducer";
00064     edm::LogInfo (MessageCategory_) << "Vertex distance cut set to " << distanceCut_  << " mm";
00065     edm::LogInfo (MessageCategory_) << "Volume radius set to "       << volumeRadius_ << " mm";
00066     edm::LogInfo (MessageCategory_) << "Volume Z      set to "       << volumeZ_      << " mm";
00067 
00068     if (useMultipleHepMCLabels_) edm::LogInfo (MessageCategory_) << "Collecting generator information from pileup.";
00069     if (mergedBremsstrahlung_) edm::LogInfo (MessageCategory_) << "Merging electrom bremsstralung";
00070     if (removeDeadModules_) edm::LogInfo (MessageCategory_) << "Removing psimhit from dead modules";
00071 
00072     if (!mergedBremsstrahlung_)
00073     {
00074         produces<TrackingVertexCollection>();
00075         produces<TrackingParticleCollection>();
00076     }
00077     else
00078     {
00079         produces<TrackingVertexCollection>();
00080         produces<TrackingParticleCollection>();
00081         produces<TrackingVertexCollection>("MergedTrackTruth");
00082         produces<TrackingParticleCollection>("MergedTrackTruth");
00083     }
00084 }
00085 
00086 
00087 void TrackingTruthProducer::produce(edm::Event &event, const edm::EventSetup & setup)
00088 {
00089     // Clean the list of hepmc products
00090     hepMCProducts_.clear();
00091 
00092     // Collect all the HepMCProducts
00093     edm::Handle<edm::HepMCProduct> hepMCHandle;
00094 
00095     for (std::vector<std::string>::const_iterator source = dataLabels_.begin(); source != dataLabels_.end(); ++source)
00096     {
00097         if ( event.getByLabel(*source, hepMCHandle) )
00098         {
00099             hepMCProducts_.push_back(hepMCHandle);
00100             edm::LogInfo (MessageCategory_) << "Using HepMC source " << *source;
00101             if (!useMultipleHepMCLabels_) break;
00102         }
00103     }
00104 
00105     if ( hepMCProducts_.empty() )
00106     {
00107         edm::LogWarning (MessageCategory_) << "No HepMC source found";
00108         return;
00109     }
00110     else if (hepMCProducts_.size() > 1 || useMultipleHepMCLabels_)
00111     {
00112         edm::LogInfo (MessageCategory_) << "You are using more than one HepMC source.";
00113         edm::LogInfo (MessageCategory_) << "If the labels are not in the same order as the events in the crossing frame (i.e. signal, pileup(s) ) ";
00114         edm::LogInfo (MessageCategory_) << "or there are fewer labels than events in the crossing frame";
00115         edm::LogInfo (MessageCategory_) << MessageCategory_ << " may try to access data in the wrong HepMCProduct and crash.";
00116     }
00117 
00118     // Collect all the simtracks from the crossing frame
00119     edm::Handle<CrossingFrame<SimTrack> > cfSimTracks;
00120     event.getByLabel("mix", simHitLabel_, cfSimTracks);
00121 
00122     // Create a mix collection from one simtrack collection
00123     simTracks_ = std::auto_ptr<MixCollection<SimTrack> >( new MixCollection<SimTrack>(cfSimTracks.product()) );
00124 
00125     // Collect all the simvertex from the crossing frame
00126     edm::Handle<CrossingFrame<SimVertex> > cfSimVertexes;
00127     event.getByLabel("mix", simHitLabel_, cfSimVertexes);
00128 
00129     // Create a mix collection from one simvertex collection
00130     simVertexes_ = std::auto_ptr<MixCollection<SimVertex> >( new MixCollection<SimVertex>(cfSimVertexes.product()) );
00131 
00132     // Create collections of things we will put in event
00133     trackingParticles_ = std::auto_ptr<TrackingParticleCollection>( new TrackingParticleCollection );
00134     trackingVertexes_ = std::auto_ptr<TrackingVertexCollection>( new TrackingVertexCollection );
00135 
00136     // Get references before put so we can cross reference
00137     refTrackingParticles_ = event.getRefBeforePut<TrackingParticleCollection>();
00138     refTrackingVertexes_ = event.getRefBeforePut<TrackingVertexCollection>();
00139 
00140     // Create a list of psimhits
00141     if (removeDeadModules_)
00142     {
00143         pSimHits_.clear();
00144         pixelPSimHitSelector_.select(pSimHits_, event, setup);
00145         trackerPSimHitSelector_.select(pSimHits_, event, setup);
00146         muonPSimHitSelector_.select(pSimHits_, event, setup);
00147     }
00148     else
00149     {
00150         pSimHits_.clear();
00151         pSimHitSelector_.select(pSimHits_, event, setup);
00152     }
00153 
00154     // Create a multimap between trackId and hit indices
00155     associator(pSimHits_, trackIdToHits_);
00156 
00157     // Create a map between trackId and track index
00158     associator(simTracks_, trackIdToIndex_);
00159 
00160     // Create a map between vertexId and vertex index
00161     associator(simVertexes_, vertexIdToIndex_);
00162 
00163     createTrackingTruth();
00164 
00165     if (mergedBremsstrahlung_)
00166     {
00167         // Create collections of things we will put in event,
00168         mergedTrackingParticles_ = std::auto_ptr<TrackingParticleCollection>( new TrackingParticleCollection );
00169         mergedTrackingVertexes_ = std::auto_ptr<TrackingVertexCollection>( new TrackingVertexCollection );
00170 
00171         // Get references before put so we can cross reference
00172         refMergedTrackingParticles_ = event.getRefBeforePut<TrackingParticleCollection>("MergedTrackTruth");
00173         refMergedTrackingVertexes_ = event.getRefBeforePut<TrackingVertexCollection>("MergedTrackTruth");
00174 
00175         // Merged brem electrons
00176         mergeBremsstrahlung();
00177 
00178         // Put TrackingParticles and TrackingVertices in event
00179         event.put(mergedTrackingParticles_, "MergedTrackTruth");
00180         event.put(mergedTrackingVertexes_, "MergedTrackTruth");
00181         event.put(trackingParticles_);
00182         event.put(trackingVertexes_);
00183     }
00184     else
00185     {
00186         // Put TrackingParticles and TrackingVertices in event
00187         event.put(trackingParticles_);
00188         event.put(trackingVertexes_);
00189     }
00190 }
00191 
00192 
00193 void TrackingTruthProducer::associator(
00194     std::vector<PSimHit> const & pSimHits,
00195     EncodedTruthIdToIndexes & association
00196 )
00197 {
00198     // Clear the association map
00199     association.clear();
00200 
00201     // Create a association from simtracks to overall index in the mix collection
00202     for (std::size_t i = 0; i < pSimHits.size(); ++i)
00203     {
00204         EncodedTruthIdToIndexes::key_type objectId = EncodedTruthIdToIndexes::key_type(pSimHits[i].eventId(), pSimHits[i].trackId());
00205         association.insert( std::make_pair(objectId, i) );
00206     }
00207 }
00208 
00209 
00210 void TrackingTruthProducer::associator(
00211     std::auto_ptr<MixCollection<SimTrack> > const & mixCollection,
00212     EncodedTruthIdToIndex & association
00213 )
00214 {
00215     int index = 0;
00216     // Clear the association map
00217     association.clear();
00218     // Create a association from simtracks to overall index in the mix collection
00219     for (MixCollection<SimTrack>::MixItr iterator = mixCollection->begin(); iterator != mixCollection->end(); ++iterator, ++index)
00220     {
00221         EncodedTruthId objectId = EncodedTruthId(iterator->eventId(), iterator->trackId());
00222         association.insert( std::make_pair(objectId, index) );
00223     }
00224 }
00225 
00226 
00227 void TrackingTruthProducer::associator(
00228     std::auto_ptr<MixCollection<SimVertex> > const & mixCollection,
00229     EncodedTruthIdToIndex & association
00230 )
00231 {
00232     int index = 0;
00233 
00234     // Solution to the problem of not having vertexId
00235     bool useVertexId = true;
00236     EncodedEventIdToIndex vertexId;
00237     EncodedEventId oldEventId;
00238     unsigned int oldVertexId = 0;
00239 
00240     // Loop for finding repeated vertexId (vertexId problem hack)
00241     for (MixCollection<SimVertex>::MixItr iterator = mixCollection->begin(); iterator != mixCollection->end(); ++iterator, ++index)
00242     {
00243         if (!index || iterator->eventId() != oldEventId)
00244         {
00245             oldEventId = iterator->eventId();
00246             oldVertexId = iterator->vertexId();
00247             continue;
00248         }
00249 
00250         if ( iterator->vertexId() == oldVertexId )
00251         {
00252             edm::LogWarning(MessageCategory_) << "Multiple vertexId found, no using vertexId.";
00253             useVertexId = false;
00254             break;
00255         }
00256     }
00257 
00258     // Reset the index
00259     index = 0;
00260 
00261     // Clear the association map
00262     association.clear();
00263 
00264     // Create a association from simvertexes to overall index in the mix collection
00265     for (MixCollection<SimVertex>::MixItr iterator = mixCollection->begin(); iterator != mixCollection->end(); ++iterator, ++index)
00266     {
00267         EncodedTruthId objectId;
00268         if (useVertexId)
00269             objectId = EncodedTruthId(iterator->eventId(), iterator->vertexId());
00270         else
00271             objectId = EncodedTruthId(iterator->eventId(), vertexId[iterator->eventId()]++);
00272         association.insert( std::make_pair(objectId, index) );
00273     }
00274 }
00275 
00276 
00277 void TrackingTruthProducer::mergeBremsstrahlung()
00278 {
00279     unsigned int index = 0;
00280 
00281     std::set<unsigned int> excludedTV, excludedTP;
00282 
00283     // Merge Bremsstrahlung vertexes
00284     for (TrackingVertexCollection::iterator iVC = trackingVertexes_->begin(); iVC != trackingVertexes_->end(); ++iVC, ++index)
00285     {
00286         // Check Bremsstrahlung vertex
00287         if ( isBremsstrahlungVertex(*iVC, trackingParticles_) )
00288         {
00289             // Get a pointer to the source track (A Ref<> cannot be use with a product!)
00290             TrackingParticle * track = &trackingParticles_->at(iVC->sourceTracks_begin()->key());
00291             // Get a Ref<> to the source track
00292             TrackingParticleRef trackRef = *iVC->sourceTracks_begin();
00293             // Pointer to electron daughter
00294             TrackingParticle * daughter = 0;
00295             // Ref<> to electron daughter
00296             TrackingParticleRef daughterRef;
00297 
00298             // Select the electron daughter and redirect the photon
00299             for (TrackingVertex::tp_iterator idaughter = iVC->daughterTracks_begin(); idaughter != iVC->daughterTracks_end(); ++idaughter)
00300             {
00301                 TrackingParticle * pointer = &trackingParticles_->at(idaughter->key());
00302                 if ( std::abs( pointer->pdgId() ) == 11 )
00303                 {
00304                     // Set pointer to the electron daughter
00305                     daughter = pointer;
00306                     // Set Ref<> to the electron daughter
00307                     daughterRef = *idaughter;
00308                 }
00309                 else if ( pointer->pdgId() == 22 )
00310                 {
00311                     // Delete the photon original parent vertex
00312                     pointer->clearParentVertex();
00313                     // Set the new parent vertex to the vertex of the source track
00314                     pointer->setParentVertex( track->parentVertex() );
00315                     // Get a non-const pointer to the parent vertex
00316                     TrackingVertex * vertex = &trackingVertexes_->at( track->parentVertex().key() );
00317                     // Add the photon to the doughter list of the parent vertex
00318                     vertex->addDaughterTrack( *idaughter );
00319                 }
00320             }
00321 
00322             // Add the electron segments from the electron daughter
00323             // track must not be the same particle as daughter
00324             // if (track != daughter)
00325             for (TrackingParticle::g4t_iterator isegment = daughter->g4Track_begin(); isegment != daughter->g4Track_end(); ++isegment)
00326                 track->addG4Track(*isegment);
00327 
00328             // Copy all the simhits to the new track
00329             for (std::vector<PSimHit>::const_iterator ihit = daughter->pSimHit_begin(); ihit != daughter->pSimHit_end(); ++ihit)
00330                 track->addPSimHit(*ihit);
00331 
00332             // Make a copy of the decay vertexes of the track
00333             TrackingVertexRefVector decayVertices( track->decayVertices() );
00334 
00335             // Clear the decay vertex list
00336             track->clearDecayVertices();
00337 
00338             // Add the remaining vertexes
00339             for (TrackingVertexRefVector::const_iterator idecay = decayVertices.begin(); idecay != decayVertices.end(); ++idecay)
00340                 if ( (*idecay).key() != index ) track->addDecayVertex(*idecay);
00341 
00342             // Redirect all the decay source vertexes to those in the electron daughter
00343             for (TrackingParticle::tv_iterator idecay = daughter->decayVertices_begin(); idecay != daughter->decayVertices_end(); ++idecay)
00344             {
00345                 // Add the vertexes to the decay list of the source particles
00346                 track->addDecayVertex(*idecay);
00347                 // Get a reference to decay vertex
00348                 TrackingVertex * vertex = &trackingVertexes_->at( idecay->key() );
00349                 // Copy all the source tracks from of the decay vertex
00350                 TrackingParticleRefVector sources( vertex->sourceTracks() );
00351                 // Clear the source track references
00352                 vertex->clearParentTracks();
00353                 // Add the new source tracks by excluding the one with the segment merged
00354                 for (TrackingVertex::tp_iterator isource = sources.begin(); isource != sources.end(); ++isource)
00355                     if (*isource != daughterRef)
00356                         vertex->addParentTrack(*isource);
00357                 // Add the track reference to the list of sources
00358                 vertex->addParentTrack(trackRef);
00359             }
00360 
00361             // Adding the vertex to the exlusion list
00362             excludedTV.insert(index);
00363 
00364             // Adding the electron segment tp into the exlusion list
00365             excludedTP.insert( daughterRef.key() );
00366         }
00367     }
00368 
00369     edm::LogInfo(MessageCategory_) << "Generating the merged collection." << std::endl;
00370 
00371     // Reserved the same amount of memory for the merged collections
00372     mergedTrackingParticles_->reserve(trackingParticles_->size());
00373     mergedTrackingVertexes_->reserve(trackingVertexes_->size());
00374 
00375     index = 0;
00376     std::map<unsigned int, unsigned int> vertexMap;
00377 
00378     // Copy non-excluded vertices discarding parent & child tracks
00379     for (TrackingVertexCollection::const_iterator iVC = trackingVertexes_->begin(); iVC != trackingVertexes_->end(); ++iVC, ++index)
00380     {
00381         if ( excludedTV.find(index) != excludedTV.end() ) continue;
00382         // Save the new location of the non excluded vertexes (take in consideration those were removed)
00383         vertexMap.insert( std::make_pair(index, mergedTrackingVertexes_->size()) );
00384         // Copy those vertexes are not excluded
00385         TrackingVertex newVertex = (*iVC);
00386         newVertex.clearDaughterTracks();
00387         newVertex.clearParentTracks();
00388         mergedTrackingVertexes_->push_back(newVertex);
00389     }
00390 
00391     index = 0;
00392 
00393     // Copy and cross reference the non-excluded tp to the merged collection
00394     for (TrackingParticleCollection::const_iterator iTP = trackingParticles_->begin(); iTP != trackingParticles_->end(); ++iTP, ++index)
00395     {
00396         if ( excludedTP.find(index) != excludedTP.end() ) continue;
00397 
00398         TrackingVertexRef       sourceV = iTP->parentVertex();
00399         TrackingVertexRefVector decayVs = iTP->decayVertices();
00400         TrackingParticle newTrack = *iTP;
00401 
00402         newTrack.clearParentVertex();
00403         newTrack.clearDecayVertices();
00404 
00405         // Set vertex indices for new vertex product and track references in those vertices
00406 
00407         // Index of parent vertex in vertex container
00408         unsigned int parentIndex = vertexMap[sourceV.key()];
00409         // Index of this track in track container
00410         unsigned int tIndex = mergedTrackingParticles_->size();
00411 
00412         // Add vertex to track
00413         newTrack.setParentVertex( TrackingVertexRef(refMergedTrackingVertexes_, parentIndex) );
00414         // Add track to vertex
00415         (mergedTrackingVertexes_->at(parentIndex)).addDaughterTrack(TrackingParticleRef(refMergedTrackingParticles_, tIndex));
00416 
00417         for (TrackingVertexRefVector::const_iterator iDecayV = decayVs.begin(); iDecayV != decayVs.end(); ++iDecayV)
00418         {
00419             // Index of decay vertex in vertex container
00420             unsigned int daughterIndex = vertexMap[iDecayV->key()];
00421             // Add vertex to track
00422             newTrack.addDecayVertex( TrackingVertexRef(refMergedTrackingVertexes_, daughterIndex) );
00423             // Add track to vertex
00424             (mergedTrackingVertexes_->at(daughterIndex)).addParentTrack( TrackingParticleRef(refMergedTrackingParticles_, tIndex) );
00425         }
00426 
00427         mergedTrackingParticles_->push_back(newTrack);
00428     }
00429 }
00430 
00431 
00432 bool TrackingTruthProducer::isBremsstrahlungVertex(
00433     TrackingVertex const & vertex,
00434     std::auto_ptr<TrackingParticleCollection> & tPC
00435 )
00436 {
00437     const TrackingParticleRefVector parents(vertex.sourceTracks());
00438 
00439     // Check for the basic parent conditions
00440     if ( parents.size() != 1)
00441         return false;
00442 
00443     // Check for the parent particle is a |electron| (electron or positron)
00444     if ( std::abs( tPC->at(parents.begin()->key()).pdgId() ) != 11 )
00445         return false;
00446 
00447     unsigned int nElectrons = 0;
00448     unsigned int nOthers = 0;
00449 
00450     // Loop over the daughter particles and counts the number of |electrons|, others (non photons)
00451     for ( TrackingVertex::tp_iterator it = vertex.daughterTracks_begin(); it != vertex.daughterTracks_end(); ++it )
00452     {
00453         // Stronger rejection for looping particles
00454         if ( parents[0] == *it )
00455             return false;
00456 
00457         if ( std::abs( tPC->at(it->key()).pdgId() ) == 11 )
00458             nElectrons++;
00459         else if ( tPC->at(it->key()).pdgId() != 22 )
00460             nOthers++;
00461     }
00462 
00463     // Condition to be a Bremsstrahlung Vertex
00464     if (nElectrons == 1 && nOthers == 0)
00465         return true;
00466 
00467     return false;
00468 }
00469 
00470 
00471 void TrackingTruthProducer::createTrackingTruth()
00472 {
00473     // Reset the event counter (use for define vertexId)
00474     eventIdCounter_.clear();
00475 
00476     // Define a container of vetoed traks
00477     std::map<int,std::size_t> vetoedTracks;
00478 
00479     // Define map between parent simtrack and tv indexes
00480     std::map<int,std::size_t> vetoedSimVertexes;
00481 
00482     for (int simTrackIndex = 0; simTrackIndex != simTracks_->size(); ++simTrackIndex)
00483     {
00484         // Check if the simTrack is excluded (includes non traceable and recovered by history)
00485         if ( vetoedTracks.find(simTrackIndex) != vetoedTracks.end() ) continue;
00486 
00487         SimTrack const & simTrack = simTracks_->getObject(simTrackIndex);
00488 
00489         TrackingParticle trackingParticle;
00490 
00491         // Set a bare tp (only with the psimhit) with a given simtrack
00492         // the function return true if it is tracable
00493         if ( setTrackingParticle(simTrack, trackingParticle) )
00494         {
00495             // Follows the path upward recovering the history of the particle
00496             SimTrack const * currentSimTrack = & simTrack;
00497 
00498             // Initial condition for the tp and tv indexes
00499             int trackingParticleIndex = -1;
00500             int trackingVertexIndex = -1;
00501 
00502             do
00503             {
00504                 // Set a new tracking particle for the current simtrack
00505                 // and add it to the list of parent tracks of previous vertex
00506                 if (trackingParticleIndex >= 0)
00507                 {
00508                     setTrackingParticle(*currentSimTrack, trackingParticle);
00509 
00510                     // Set the tp index to its new value
00511                     trackingParticleIndex = trackingParticles_->size();
00512                     // Push the tp in to the collection
00513                     trackingParticles_->push_back(trackingParticle);
00514 
00515                     // Add the previous track to the list of decay vertexes of the new tp
00516                     trackingParticles_->at(trackingParticleIndex).addDecayVertex(
00517                         TrackingVertexRef(refTrackingVertexes_, trackingVertexIndex)
00518                     );
00519 
00520                     // Add the new tp to the list of parent tracks of the previous tv
00521                     trackingVertexes_->at(trackingVertexIndex).addParentTrack(
00522                         TrackingParticleRef(refTrackingParticles_, trackingParticleIndex)
00523                     );
00524                 }
00525                 else
00526                 {
00527                     // Set the tp index to its new value
00528                     trackingParticleIndex = trackingParticles_->size();
00529                     // Push the tp in to the collection
00530                     trackingParticles_->push_back(trackingParticle);
00531                     // Vetoed the simTrack
00532                     vetoedTracks.insert( std::make_pair(simTrackIndex, trackingParticleIndex) );
00533                 }
00534 
00535                 // Verify if the parent simVertex has a simTrack or if the source is a vetoSimVertex
00536                 if (currentSimTrack->noVertex()) break;
00537 
00538                 // Get the simTrack parent index (it is implicit should be in the same event as current)
00539                 unsigned int parentSimVertexIndex = vertexIdToIndex_[
00540                                                         EncodedTruthId(
00541                                                             currentSimTrack->eventId(),
00542                                                             currentSimTrack->vertIndex()
00543                                                         )
00544                                                     ];
00545                 // Create a new tv
00546                 TrackingVertex trackingVertex;
00547                 // Get the parent simVertex associated to the current simTrack
00548                 SimVertex const * parentSimVertex = & simVertexes_->getObject(parentSimVertexIndex);
00549 
00550                 bool vetoSimVertex = vetoedSimVertexes.find(parentSimVertexIndex) != vetoedSimVertexes.end();
00551 
00552                 // Check for a already visited parent simTrack
00553                 if ( !vetoSimVertex )
00554                 {
00555                     // Set the tv by using simvertex
00556                     trackingVertexIndex = setTrackingVertex(*parentSimVertex, trackingVertex);
00557 
00558                     // Check if a new vertex needs to be created
00559                     if (trackingVertexIndex < 0)
00560                     {
00561                         // Set the tv index ot its new value
00562                         trackingVertexIndex = trackingVertexes_->size();
00563                         // Push the new tv in to the collection
00564                         trackingVertexes_->push_back(trackingVertex);
00565                     }
00566                     else
00567                     {
00568                         // Get the postion and time of the vertex
00569                         const LorentzVector & position = trackingVertexes_->at(trackingVertexIndex).position();
00570                         Vector xyz = Vector(position.x(), position.y(), position.z());
00571                         double t = position.t();
00572                         // Set the vertex postion of the tp to the closest vertex
00573                         trackingParticles_->at(trackingParticleIndex).setVertex(xyz, t);
00574                     }
00575 
00576                     vetoedSimVertexes.insert( std::make_pair(parentSimVertexIndex, trackingVertexIndex) );
00577                 }
00578                 else
00579                     trackingVertexIndex = vetoedSimVertexes[parentSimVertexIndex];
00580 
00581                 // Set the newly created tv as parent vertex
00582                 trackingParticles_->at(trackingParticleIndex).setParentVertex(
00583                     TrackingVertexRef(refTrackingVertexes_, trackingVertexIndex)
00584                 );
00585 
00586                 // Add the newly created tp to the tv daughter list
00587                 trackingVertexes_->at(trackingVertexIndex).addDaughterTrack(
00588                     TrackingParticleRef(refTrackingParticles_, trackingParticleIndex)
00589                 );
00590 
00591                 // Verify if the parent simVertex has a simTrack or if the source is a vetoSimVertex
00592                 if (parentSimVertex->noParent() || vetoSimVertex) break;
00593 
00594                 // Get the next simTrack index (it is implicit should be in the same event as current).
00595                 unsigned int nextSimTrackIndex = trackIdToIndex_[
00596                                                      EncodedTruthId(
00597                                                          currentSimTrack->eventId(),
00598                                                          parentSimVertex->parentIndex()
00599                                                      )
00600                                                  ];
00601 
00602                 // Check if the next track exist
00603                 if ( vetoedTracks.find(nextSimTrackIndex) != vetoedTracks.end() )
00604                 {
00605                     // Add to the newly created tv the existent next simtrack in to parent list.
00606                     trackingVertexes_->at(trackingVertexIndex).addParentTrack(
00607                         TrackingParticleRef(refTrackingParticles_, vetoedTracks[nextSimTrackIndex])
00608                     );
00609                     // Add the vertex to list of decay vertexes of the new tp
00610                     trackingParticles_->at(vetoedTracks[nextSimTrackIndex]).addDecayVertex(
00611                         TrackingVertexRef(refTrackingVertexes_, trackingVertexIndex)
00612                     );
00613                     break;
00614                 }
00615 
00616                 // Vetoed the next simTrack
00617                 vetoedTracks.insert( std::make_pair(nextSimTrackIndex, trackingParticleIndex) );
00618 
00619                 // Set the current simTrack as the next simTrack
00620                 currentSimTrack = & simTracks_->getObject(nextSimTrackIndex);
00621             }
00622             while (!currentSimTrack->noVertex());
00623         }
00624     }
00625 }
00626 
00627 
00628 bool TrackingTruthProducer::setTrackingParticle(
00629     SimTrack const & simTrack,
00630     TrackingParticle & trackingParticle
00631 )
00632 {
00633     // Get the eventid associated to the track
00634     EncodedEventId trackEventId = simTrack.eventId();
00635     // Get the simtrack id
00636     EncodedTruthId simTrackId = EncodedTruthId( trackEventId, simTrack.trackId() );
00637 
00638     // Location of the parent vertex
00639     LorentzVector position;
00640     // If not parent then location is (0,0,0,0)
00641     if (simTrack.noVertex())
00642         position = LorentzVector(0, 0, 0, 0);
00643     else
00644         position = simVertexes_->getObject(simTrack.vertIndex()). position();
00645 
00646     // Define the default status and pdgid
00647     int status = -99;
00648     int pdgId = simTrack.type();
00649 
00650     int genParticleIndex = simTrack.genpartIndex();
00651     bool signalEvent = (trackEventId.event() == 0 && trackEventId.bunchCrossing() == 0);
00652 
00653     // In the case of a existing generated particle and track
00654     // event is signal redefine status a pdgId
00655 
00656     edm::Handle<edm::HepMCProduct> hepmc;
00657 
00658     if (genParticleIndex >= 0 && (signalEvent || useMultipleHepMCLabels_) )
00659     {
00660         // Get the generated particle
00661         hepmc = (useMultipleHepMCLabels_) ? hepMCProducts_.at(trackEventId.rawId()) : hepmc = hepMCProducts_.at(0);
00662 
00663         const HepMC::GenParticle * genParticle = hepmc->GetEvent()->barcode_to_particle(genParticleIndex);
00664 
00665         if (genParticle)
00666         {
00667             status = genParticle->status();
00668             pdgId  = genParticle->pdg_id();
00669         }
00670     }
00671 
00672     // Create a tp from the simtrack
00673     trackingParticle = TrackingParticle(
00674                            (char) simTrack.charge(),
00675                            simTrack.momentum(),
00676                            Vector(position.x(), position.y(), position.z()),
00677                            position.t(),
00678                            pdgId,
00679                            status,
00680                            trackEventId
00681                        );
00682 
00683     bool init = true;
00684 
00685     int processType = 0;
00686     int particleType = 0;
00687 
00688     // Counting the TP hits using the layers (as in ORCA).
00689     // Does seem to find less hits. maybe b/c layer is a number now, not a pointer
00690     int totalSimHits = 0;
00691     int oldLayer = 0;
00692     int newLayer = 0;
00693     int oldDetector = 0;
00694     int newDetector = 0;
00695 
00696     // Loop over the associated hits per track
00697     for (
00698         EncodedTruthIdToIndexes::const_iterator iEntry = trackIdToHits_.lower_bound(simTrackId);
00699         iEntry != trackIdToHits_.upper_bound(simTrackId);
00700         ++iEntry
00701     )
00702     {
00703         // Get a constant reference to the simhit
00704         PSimHit const & pSimHit = pSimHits_.at(iEntry->second);
00705 
00706         // Initial condition for consistent simhit selection
00707         if (init)
00708         {
00709             processType = pSimHit.processType();
00710             particleType = pSimHit.particleType();
00711             init = false;
00712         }
00713 
00714         // Check for delta and interaction products discards
00715         if (processType == pSimHit.processType() && particleType == pSimHit.particleType() && pdgId == pSimHit.particleType() )
00716         {
00717             trackingParticle.addPSimHit(pSimHit);
00718 
00719             unsigned int detectorIdIndex = pSimHit.detUnitId();
00720             DetId detectorId = DetId(detectorIdIndex);
00721             oldLayer = newLayer;
00722             oldDetector = newDetector;
00723             newLayer = LayerFromDetid(detectorIdIndex);
00724             newDetector = detectorId.subdetId();
00725 
00726             // Count hits using layers for glued detectors
00727             // newlayer !=0 excludes Muon layers set to 0 by LayerFromDetid
00728             if ( ( oldLayer != newLayer || (oldLayer==newLayer && oldDetector!=newDetector ) ) && newLayer != 0 ) totalSimHits++;
00729         }
00730     }
00731 
00732     // Set the number of matched simhits
00733     trackingParticle.setMatchedHit(totalSimHits);
00734 
00735     // Add the simtrack associated to the tp
00736     trackingParticle.addG4Track(simTrack);
00737 
00738     // Add the generator information
00739     if ( genParticleIndex >= 0 && (signalEvent || useMultipleHepMCLabels_) )
00740         trackingParticle.addGenParticle( GenParticleRef(hepmc, genParticleIndex) );
00741 
00742     if (selectorFlag_) return selector_(trackingParticle);
00743 
00744     return true;
00745 }
00746 
00747 
00748 int TrackingTruthProducer::setTrackingVertex(
00749     SimVertex const & simVertex,
00750     TrackingVertex & trackingVertex
00751 )
00752 {
00753     LorentzVector const & position = simVertex.position();
00754 
00755     // Look for close by vertexes
00756     for (std::size_t trackingVertexIndex = 0; trackingVertexIndex < trackingVertexes_->size(); ++trackingVertexIndex)
00757     {
00758         // Calculate the distance
00759         double distance = (position - trackingVertexes_->at(trackingVertexIndex).position()).P();
00760         // If the distance is under a given cut return the trackingVertex index (vertex merging)
00761         if (distance <= distanceCut_)
00762         {
00763             // Add simvertex to the pre existent tv
00764             trackingVertexes_->at(trackingVertexIndex).addG4Vertex(simVertex);
00765             // return tv index
00766             return trackingVertexIndex;
00767         }
00768     }
00769 
00770     // Get the event if from the simvertex
00771     EncodedEventId simVertexEventId = simVertex.eventId();
00772 
00773     // Initialize the event counter
00774     if (eventIdCounter_.find(simVertexEventId) == eventIdCounter_.end())
00775         eventIdCounter_[simVertexEventId] = 0;
00776 
00777     // Get the simVertex id
00778     EncodedTruthId simVertexId = EncodedTruthId(simVertexEventId, eventIdCounter_[simVertexEventId]);
00779 
00780     // Calculate if the vertex is in the tracker volume (it needs to be review for other detectors)
00781     bool inVolume = (position.Pt() < volumeRadius_ && std::abs(position.z()) < volumeZ_); // In or out of Tracker
00782 
00783     // Initialize the new vertex
00784     trackingVertex = TrackingVertex(position, inVolume, simVertexId);
00785 
00786     // Find the the closest GenVertexes to the created tv
00787     addCloseGenVertexes(trackingVertex);
00788 
00789     // Add the g4 vertex to the tv
00790     trackingVertex.addG4Vertex(simVertex);
00791 
00792     // Initialize the event counter
00793     eventIdCounter_[simVertexEventId]++;
00794 
00795     return -1;
00796 }
00797 
00798 
00799 void TrackingTruthProducer::addCloseGenVertexes(TrackingVertex & trackingVertex)
00800 {
00801     // Get the generated particle
00802     edm::Handle<edm::HepMCProduct> hepmc = (useMultipleHepMCLabels_) ? hepMCProducts_.at(trackingVertex.eventId().rawId()) : hepMCProducts_.at(0);
00803     const HepMC::GenEvent * genEvent = hepmc->GetEvent();
00804 
00805     // Get the postion of the tv
00806     Vector tvPosition(trackingVertex.position().x(), trackingVertex.position().y(), trackingVertex.position().z());
00807 
00808     // Find HepMC vertices, put them in a close TrackingVertex (this could conceivably add the same GenVertex to multiple TrackingVertices)
00809     for (
00810         HepMC::GenEvent::vertex_const_iterator iGenVertex = genEvent->vertices_begin();
00811         iGenVertex != genEvent->vertices_end();
00812         ++iGenVertex
00813     )
00814     {
00815         // Get the position of the genVertex
00816         HepMC::ThreeVector rawPosition = (*iGenVertex)->position();
00817 
00818         // Convert to cm
00819         Vector genPosition(rawPosition.x()/10.0, rawPosition.y()/10.0, rawPosition.z()/10.0);
00820 
00821         // Calculate the dis
00822         double distance = sqrt( (tvPosition - genPosition).mag2() );
00823 
00824         if (distance <= distanceCut_)
00825             trackingVertex.addGenVertex( GenVertexRef(hepmc, (*iGenVertex)->barcode()) );
00826     }
00827 }
00828 
00829 
00830 int TrackingTruthProducer::LayerFromDetid(const unsigned int & detid)
00831 {
00832     DetId detId = DetId(detid);
00833 
00834     if ( detId.det() != DetId::Tracker ) return 0;
00835 
00836     int layerNumber=0;
00837     unsigned int subdetId = static_cast<unsigned int>(detId.subdetId());
00838 
00839     if ( subdetId == StripSubdetector::TIB)
00840     {
00841         TIBDetId tibid(detId.rawId());
00842         layerNumber = tibid.layer();
00843     }
00844     else if ( subdetId ==  StripSubdetector::TOB )
00845     {
00846         TOBDetId tobid(detId.rawId());
00847         layerNumber = tobid.layer();
00848     }
00849     else if ( subdetId ==  StripSubdetector::TID)
00850     {
00851         TIDDetId tidid(detId.rawId());
00852         layerNumber = tidid.wheel();
00853     }
00854     else if ( subdetId ==  StripSubdetector::TEC )
00855     {
00856         TECDetId tecid(detId.rawId());
00857         layerNumber = tecid.wheel();
00858     }
00859     else if ( subdetId ==  PixelSubdetector::PixelBarrel )
00860     {
00861         PXBDetId pxbid(detId.rawId());
00862         layerNumber = pxbid.layer();
00863     }
00864     else if ( subdetId ==  PixelSubdetector::PixelEndcap )
00865     {
00866         PXFDetId pxfid(detId.rawId());
00867         layerNumber = pxfid.disk();
00868     }
00869     else
00870         edm::LogVerbatim("TrackingTruthProducer") << "Unknown subdetid: " <<  subdetId;
00871 
00872     return layerNumber;
00873 }
00874 
00875 DEFINE_FWK_MODULE(TrackingTruthProducer);