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
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
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
00090 hepMCProducts_.clear();
00091
00092
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
00119 edm::Handle<CrossingFrame<SimTrack> > cfSimTracks;
00120 event.getByLabel("mix", simHitLabel_, cfSimTracks);
00121
00122
00123 simTracks_ = std::auto_ptr<MixCollection<SimTrack> >( new MixCollection<SimTrack>(cfSimTracks.product()) );
00124
00125
00126 edm::Handle<CrossingFrame<SimVertex> > cfSimVertexes;
00127 event.getByLabel("mix", simHitLabel_, cfSimVertexes);
00128
00129
00130 simVertexes_ = std::auto_ptr<MixCollection<SimVertex> >( new MixCollection<SimVertex>(cfSimVertexes.product()) );
00131
00132
00133 trackingParticles_ = std::auto_ptr<TrackingParticleCollection>( new TrackingParticleCollection );
00134 trackingVertexes_ = std::auto_ptr<TrackingVertexCollection>( new TrackingVertexCollection );
00135
00136
00137 refTrackingParticles_ = event.getRefBeforePut<TrackingParticleCollection>();
00138 refTrackingVertexes_ = event.getRefBeforePut<TrackingVertexCollection>();
00139
00140
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
00155 associator(pSimHits_, trackIdToHits_);
00156
00157
00158 associator(simTracks_, trackIdToIndex_);
00159
00160
00161 associator(simVertexes_, vertexIdToIndex_);
00162
00163 createTrackingTruth();
00164
00165 if (mergedBremsstrahlung_)
00166 {
00167
00168 mergedTrackingParticles_ = std::auto_ptr<TrackingParticleCollection>( new TrackingParticleCollection );
00169 mergedTrackingVertexes_ = std::auto_ptr<TrackingVertexCollection>( new TrackingVertexCollection );
00170
00171
00172 refMergedTrackingParticles_ = event.getRefBeforePut<TrackingParticleCollection>("MergedTrackTruth");
00173 refMergedTrackingVertexes_ = event.getRefBeforePut<TrackingVertexCollection>("MergedTrackTruth");
00174
00175
00176 mergeBremsstrahlung();
00177
00178
00179 event.put(mergedTrackingParticles_, "MergedTrackTruth");
00180 event.put(mergedTrackingVertexes_, "MergedTrackTruth");
00181 event.put(trackingParticles_);
00182 event.put(trackingVertexes_);
00183 }
00184 else
00185 {
00186
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
00199 association.clear();
00200
00201
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
00217 association.clear();
00218
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
00235 bool useVertexId = true;
00236 EncodedEventIdToIndex vertexId;
00237 EncodedEventId oldEventId;
00238 unsigned int oldVertexId = 0;
00239
00240
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
00259 index = 0;
00260
00261
00262 association.clear();
00263
00264
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
00284 for (TrackingVertexCollection::iterator iVC = trackingVertexes_->begin(); iVC != trackingVertexes_->end(); ++iVC, ++index)
00285 {
00286
00287 if ( isBremsstrahlungVertex(*iVC, trackingParticles_) )
00288 {
00289
00290 TrackingParticle * track = &trackingParticles_->at(iVC->sourceTracks_begin()->key());
00291
00292 TrackingParticleRef trackRef = *iVC->sourceTracks_begin();
00293
00294 TrackingParticle * daughter = 0;
00295
00296 TrackingParticleRef daughterRef;
00297
00298
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
00305 daughter = pointer;
00306
00307 daughterRef = *idaughter;
00308 }
00309 else if ( pointer->pdgId() == 22 )
00310 {
00311
00312 pointer->clearParentVertex();
00313
00314 pointer->setParentVertex( track->parentVertex() );
00315
00316 TrackingVertex * vertex = &trackingVertexes_->at( track->parentVertex().key() );
00317
00318 vertex->addDaughterTrack( *idaughter );
00319 }
00320 }
00321
00322
00323
00324
00325 for (TrackingParticle::g4t_iterator isegment = daughter->g4Track_begin(); isegment != daughter->g4Track_end(); ++isegment)
00326 track->addG4Track(*isegment);
00327
00328
00329 for (std::vector<PSimHit>::const_iterator ihit = daughter->pSimHit_begin(); ihit != daughter->pSimHit_end(); ++ihit)
00330 track->addPSimHit(*ihit);
00331
00332
00333 TrackingVertexRefVector decayVertices( track->decayVertices() );
00334
00335
00336 track->clearDecayVertices();
00337
00338
00339 for (TrackingVertexRefVector::const_iterator idecay = decayVertices.begin(); idecay != decayVertices.end(); ++idecay)
00340 if ( (*idecay).key() != index ) track->addDecayVertex(*idecay);
00341
00342
00343 for (TrackingParticle::tv_iterator idecay = daughter->decayVertices_begin(); idecay != daughter->decayVertices_end(); ++idecay)
00344 {
00345
00346 track->addDecayVertex(*idecay);
00347
00348 TrackingVertex * vertex = &trackingVertexes_->at( idecay->key() );
00349
00350 TrackingParticleRefVector sources( vertex->sourceTracks() );
00351
00352 vertex->clearParentTracks();
00353
00354 for (TrackingVertex::tp_iterator isource = sources.begin(); isource != sources.end(); ++isource)
00355 if (*isource != daughterRef)
00356 vertex->addParentTrack(*isource);
00357
00358 vertex->addParentTrack(trackRef);
00359 }
00360
00361
00362 excludedTV.insert(index);
00363
00364
00365 excludedTP.insert( daughterRef.key() );
00366 }
00367 }
00368
00369 edm::LogInfo(MessageCategory_) << "Generating the merged collection." << std::endl;
00370
00371
00372 mergedTrackingParticles_->reserve(trackingParticles_->size());
00373 mergedTrackingVertexes_->reserve(trackingVertexes_->size());
00374
00375 index = 0;
00376 std::map<unsigned int, unsigned int> vertexMap;
00377
00378
00379 for (TrackingVertexCollection::const_iterator iVC = trackingVertexes_->begin(); iVC != trackingVertexes_->end(); ++iVC, ++index)
00380 {
00381 if ( excludedTV.find(index) != excludedTV.end() ) continue;
00382
00383 vertexMap.insert( std::make_pair(index, mergedTrackingVertexes_->size()) );
00384
00385 TrackingVertex newVertex = (*iVC);
00386 newVertex.clearDaughterTracks();
00387 newVertex.clearParentTracks();
00388 mergedTrackingVertexes_->push_back(newVertex);
00389 }
00390
00391 index = 0;
00392
00393
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
00406
00407
00408 unsigned int parentIndex = vertexMap[sourceV.key()];
00409
00410 unsigned int tIndex = mergedTrackingParticles_->size();
00411
00412
00413 newTrack.setParentVertex( TrackingVertexRef(refMergedTrackingVertexes_, parentIndex) );
00414
00415 (mergedTrackingVertexes_->at(parentIndex)).addDaughterTrack(TrackingParticleRef(refMergedTrackingParticles_, tIndex));
00416
00417 for (TrackingVertexRefVector::const_iterator iDecayV = decayVs.begin(); iDecayV != decayVs.end(); ++iDecayV)
00418 {
00419
00420 unsigned int daughterIndex = vertexMap[iDecayV->key()];
00421
00422 newTrack.addDecayVertex( TrackingVertexRef(refMergedTrackingVertexes_, daughterIndex) );
00423
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
00440 if ( parents.size() != 1)
00441 return false;
00442
00443
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
00451 for ( TrackingVertex::tp_iterator it = vertex.daughterTracks_begin(); it != vertex.daughterTracks_end(); ++it )
00452 {
00453
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
00464 if (nElectrons == 1 && nOthers == 0)
00465 return true;
00466
00467 return false;
00468 }
00469
00470
00471 void TrackingTruthProducer::createTrackingTruth()
00472 {
00473
00474 eventIdCounter_.clear();
00475
00476
00477 std::map<int,std::size_t> vetoedTracks;
00478
00479
00480 std::map<int,std::size_t> vetoedSimVertexes;
00481
00482 for (int simTrackIndex = 0; simTrackIndex != simTracks_->size(); ++simTrackIndex)
00483 {
00484
00485 if ( vetoedTracks.find(simTrackIndex) != vetoedTracks.end() ) continue;
00486
00487 SimTrack const & simTrack = simTracks_->getObject(simTrackIndex);
00488
00489 TrackingParticle trackingParticle;
00490
00491
00492
00493 if ( setTrackingParticle(simTrack, trackingParticle) )
00494 {
00495
00496 SimTrack const * currentSimTrack = & simTrack;
00497
00498
00499 int trackingParticleIndex = -1;
00500 int trackingVertexIndex = -1;
00501
00502 do
00503 {
00504
00505
00506 if (trackingParticleIndex >= 0)
00507 {
00508 setTrackingParticle(*currentSimTrack, trackingParticle);
00509
00510
00511 trackingParticleIndex = trackingParticles_->size();
00512
00513 trackingParticles_->push_back(trackingParticle);
00514
00515
00516 trackingParticles_->at(trackingParticleIndex).addDecayVertex(
00517 TrackingVertexRef(refTrackingVertexes_, trackingVertexIndex)
00518 );
00519
00520
00521 trackingVertexes_->at(trackingVertexIndex).addParentTrack(
00522 TrackingParticleRef(refTrackingParticles_, trackingParticleIndex)
00523 );
00524 }
00525 else
00526 {
00527
00528 trackingParticleIndex = trackingParticles_->size();
00529
00530 trackingParticles_->push_back(trackingParticle);
00531
00532 vetoedTracks.insert( std::make_pair(simTrackIndex, trackingParticleIndex) );
00533 }
00534
00535
00536 if (currentSimTrack->noVertex()) break;
00537
00538
00539 unsigned int parentSimVertexIndex = vertexIdToIndex_[
00540 EncodedTruthId(
00541 currentSimTrack->eventId(),
00542 currentSimTrack->vertIndex()
00543 )
00544 ];
00545
00546 TrackingVertex trackingVertex;
00547
00548 SimVertex const * parentSimVertex = & simVertexes_->getObject(parentSimVertexIndex);
00549
00550 bool vetoSimVertex = vetoedSimVertexes.find(parentSimVertexIndex) != vetoedSimVertexes.end();
00551
00552
00553 if ( !vetoSimVertex )
00554 {
00555
00556 trackingVertexIndex = setTrackingVertex(*parentSimVertex, trackingVertex);
00557
00558
00559 if (trackingVertexIndex < 0)
00560 {
00561
00562 trackingVertexIndex = trackingVertexes_->size();
00563
00564 trackingVertexes_->push_back(trackingVertex);
00565 }
00566 else
00567 {
00568
00569
00570
00571
00572
00573
00574
00575 }
00576
00577 vetoedSimVertexes.insert( std::make_pair(parentSimVertexIndex, trackingVertexIndex) );
00578 }
00579 else
00580 trackingVertexIndex = vetoedSimVertexes[parentSimVertexIndex];
00581
00582
00583 trackingParticles_->at(trackingParticleIndex).setParentVertex(
00584 TrackingVertexRef(refTrackingVertexes_, trackingVertexIndex)
00585 );
00586
00587
00588 trackingVertexes_->at(trackingVertexIndex).addDaughterTrack(
00589 TrackingParticleRef(refTrackingParticles_, trackingParticleIndex)
00590 );
00591
00592
00593 if (parentSimVertex->noParent() || vetoSimVertex) break;
00594
00595
00596 unsigned int nextSimTrackIndex = trackIdToIndex_[
00597 EncodedTruthId(
00598 currentSimTrack->eventId(),
00599 parentSimVertex->parentIndex()
00600 )
00601 ];
00602
00603
00604 if ( vetoedTracks.find(nextSimTrackIndex) != vetoedTracks.end() )
00605 {
00606
00607 trackingVertexes_->at(trackingVertexIndex).addParentTrack(
00608 TrackingParticleRef(refTrackingParticles_, vetoedTracks[nextSimTrackIndex])
00609 );
00610
00611 trackingParticles_->at(vetoedTracks[nextSimTrackIndex]).addDecayVertex(
00612 TrackingVertexRef(refTrackingVertexes_, trackingVertexIndex)
00613 );
00614 break;
00615 }
00616
00617
00618 vetoedTracks.insert( std::make_pair(nextSimTrackIndex, trackingParticleIndex) );
00619
00620
00621 currentSimTrack = & simTracks_->getObject(nextSimTrackIndex);
00622 }
00623 while (!currentSimTrack->noVertex());
00624 }
00625 }
00626 }
00627
00628
00629 bool TrackingTruthProducer::setTrackingParticle(
00630 SimTrack const & simTrack,
00631 TrackingParticle & trackingParticle
00632 )
00633 {
00634
00635 EncodedEventId trackEventId = simTrack.eventId();
00636
00637 EncodedTruthId simTrackId = EncodedTruthId( trackEventId, simTrack.trackId() );
00638
00639
00640 LorentzVector position;
00641
00642 if (simTrack.noVertex())
00643 position = LorentzVector(0, 0, 0, 0);
00644 else
00645 {
00646 unsigned int parentSimVertexIndex = vertexIdToIndex_[ EncodedTruthId( simTrack.eventId(), simTrack.vertIndex() ) ];
00647 position = simVertexes_->getObject(parentSimVertexIndex).position();
00648 }
00649
00650
00651 int status = -99;
00652 int pdgId = simTrack.type();
00653
00654 int genParticleIndex = simTrack.genpartIndex();
00655 bool signalEvent = (trackEventId.event() == 0 && trackEventId.bunchCrossing() == 0);
00656
00657
00658
00659
00660 edm::Handle<edm::HepMCProduct> hepmc;
00661
00662 if (genParticleIndex >= 0 && (signalEvent || useMultipleHepMCLabels_) )
00663 {
00664
00665 hepmc = (useMultipleHepMCLabels_) ? hepMCProducts_.at(trackEventId.rawId()) : hepmc = hepMCProducts_.at(0);
00666
00667 const HepMC::GenParticle * genParticle = hepmc->GetEvent()->barcode_to_particle(genParticleIndex);
00668
00669 if (genParticle)
00670 {
00671 status = genParticle->status();
00672 pdgId = genParticle->pdg_id();
00673 }
00674 }
00675
00676
00677 trackingParticle = TrackingParticle(
00678 (char) simTrack.charge(),
00679 simTrack.momentum(),
00680 Vector(position.x(), position.y(), position.z()),
00681 position.t(),
00682 pdgId,
00683 status,
00684 trackEventId
00685 );
00686
00687 bool init = true;
00688
00689 int processType = 0;
00690 int particleType = 0;
00691
00692
00693
00694 int totalSimHits = 0;
00695 int oldLayer = 0;
00696 int newLayer = 0;
00697 int oldDetector = 0;
00698 int newDetector = 0;
00699
00700
00701 for (
00702 EncodedTruthIdToIndexes::const_iterator iEntry = trackIdToHits_.lower_bound(simTrackId);
00703 iEntry != trackIdToHits_.upper_bound(simTrackId);
00704 ++iEntry
00705 )
00706 {
00707
00708 PSimHit const & pSimHit = pSimHits_.at(iEntry->second);
00709
00710
00711 if (init)
00712 {
00713 processType = pSimHit.processType();
00714 particleType = pSimHit.particleType();
00715 init = false;
00716 }
00717
00718
00719 if (processType == pSimHit.processType() && particleType == pSimHit.particleType() && pdgId == pSimHit.particleType() )
00720 {
00721 trackingParticle.addPSimHit(pSimHit);
00722
00723 unsigned int detectorIdIndex = pSimHit.detUnitId();
00724 DetId detectorId = DetId(detectorIdIndex);
00725 oldLayer = newLayer;
00726 oldDetector = newDetector;
00727 newLayer = LayerFromDetid(detectorIdIndex);
00728 newDetector = detectorId.subdetId();
00729
00730
00731
00732 if ( ( oldLayer != newLayer || (oldLayer==newLayer && oldDetector!=newDetector ) ) && newLayer != 0 ) totalSimHits++;
00733 }
00734 }
00735
00736
00737 trackingParticle.setMatchedHit(totalSimHits);
00738
00739
00740 trackingParticle.addG4Track(simTrack);
00741
00742
00743 if ( genParticleIndex >= 0 && (signalEvent || useMultipleHepMCLabels_) )
00744 trackingParticle.addGenParticle( GenParticleRef(hepmc, genParticleIndex) );
00745
00746 if (selectorFlag_) return selector_(trackingParticle);
00747
00748 return true;
00749 }
00750
00751
00752 int TrackingTruthProducer::setTrackingVertex(
00753 SimVertex const & simVertex,
00754 TrackingVertex & trackingVertex
00755 )
00756 {
00757 LorentzVector const & position = simVertex.position();
00758
00759
00760 for (std::size_t trackingVertexIndex = 0; trackingVertexIndex < trackingVertexes_->size(); ++trackingVertexIndex)
00761 {
00762
00763 double distance = (position - trackingVertexes_->at(trackingVertexIndex).position()).P();
00764
00765 if (distance <= distanceCut_)
00766 {
00767
00768 trackingVertexes_->at(trackingVertexIndex).addG4Vertex(simVertex);
00769
00770 return trackingVertexIndex;
00771 }
00772 }
00773
00774
00775 EncodedEventId simVertexEventId = simVertex.eventId();
00776
00777
00778 if (eventIdCounter_.find(simVertexEventId) == eventIdCounter_.end())
00779 eventIdCounter_[simVertexEventId] = 0;
00780
00781
00782 EncodedTruthId simVertexId = EncodedTruthId(simVertexEventId, eventIdCounter_[simVertexEventId]);
00783
00784
00785 bool inVolume = (position.Pt() < volumeRadius_ && std::abs(position.z()) < volumeZ_);
00786
00787
00788 trackingVertex = TrackingVertex(position, inVolume, simVertexId);
00789
00790
00791 addCloseGenVertexes(trackingVertex);
00792
00793
00794 trackingVertex.addG4Vertex(simVertex);
00795
00796
00797 eventIdCounter_[simVertexEventId]++;
00798
00799 return -1;
00800 }
00801
00802
00803 void TrackingTruthProducer::addCloseGenVertexes(TrackingVertex & trackingVertex)
00804 {
00805
00806 edm::Handle<edm::HepMCProduct> hepmc = (useMultipleHepMCLabels_) ? hepMCProducts_.at(trackingVertex.eventId().rawId()) : hepMCProducts_.at(0);
00807 const HepMC::GenEvent * genEvent = hepmc->GetEvent();
00808
00809
00810 Vector tvPosition(trackingVertex.position().x(), trackingVertex.position().y(), trackingVertex.position().z());
00811
00812
00813 for (
00814 HepMC::GenEvent::vertex_const_iterator iGenVertex = genEvent->vertices_begin();
00815 iGenVertex != genEvent->vertices_end();
00816 ++iGenVertex
00817 )
00818 {
00819
00820 HepMC::ThreeVector rawPosition = (*iGenVertex)->position();
00821
00822
00823 Vector genPosition(rawPosition.x()/10.0, rawPosition.y()/10.0, rawPosition.z()/10.0);
00824
00825
00826 double distance = sqrt( (tvPosition - genPosition).mag2() );
00827
00828 if (distance <= distanceCut_)
00829 trackingVertex.addGenVertex( GenVertexRef(hepmc, (*iGenVertex)->barcode()) );
00830 }
00831 }
00832
00833
00834 int TrackingTruthProducer::LayerFromDetid(const unsigned int & detid)
00835 {
00836 DetId detId = DetId(detid);
00837
00838 if ( detId.det() != DetId::Tracker ) return 0;
00839
00840 int layerNumber=0;
00841 unsigned int subdetId = static_cast<unsigned int>(detId.subdetId());
00842
00843 if ( subdetId == StripSubdetector::TIB)
00844 {
00845 TIBDetId tibid(detId.rawId());
00846 layerNumber = tibid.layer();
00847 }
00848 else if ( subdetId == StripSubdetector::TOB )
00849 {
00850 TOBDetId tobid(detId.rawId());
00851 layerNumber = tobid.layer();
00852 }
00853 else if ( subdetId == StripSubdetector::TID)
00854 {
00855 TIDDetId tidid(detId.rawId());
00856 layerNumber = tidid.wheel();
00857 }
00858 else if ( subdetId == StripSubdetector::TEC )
00859 {
00860 TECDetId tecid(detId.rawId());
00861 layerNumber = tecid.wheel();
00862 }
00863 else if ( subdetId == PixelSubdetector::PixelBarrel )
00864 {
00865 PXBDetId pxbid(detId.rawId());
00866 layerNumber = pxbid.layer();
00867 }
00868 else if ( subdetId == PixelSubdetector::PixelEndcap )
00869 {
00870 PXFDetId pxfid(detId.rawId());
00871 layerNumber = pxfid.disk();
00872 }
00873 else
00874 edm::LogVerbatim("TrackingTruthProducer") << "Unknown subdetid: " << subdetId;
00875
00876 return layerNumber;
00877 }
00878
00879 DEFINE_FWK_MODULE(TrackingTruthProducer);