CMS 3D CMS Logo

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