CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackingTruthAccumulator.cc
Go to the documentation of this file.
1 
29 
44 
51 
52 
53 
54 
55 //---------------------------------------------------------------------------------
56 //---------------------------------------------------------------------------------
57 //---------------------------------------------------------------------------------
58 //------ ------
59 //------ Declarations for utility classes and functions in the unnamed ------
60 //------ namespace. The definitions are right at the bottom of the file. ------
61 //------ ------
62 //---------------------------------------------------------------------------------
63 //---------------------------------------------------------------------------------
64 //---------------------------------------------------------------------------------
65 
66 namespace
67 {
72  struct DecayChainTrack
73  {
74  int simTrackIndex;
75  struct DecayChainVertex* pParentVertex;
76  // Some things have multiple decay vertices. Not sure how that works but it seems to be mostly electrons and some photons.
77  std::vector<struct DecayChainVertex*> daughterVertices;
78  DecayChainTrack* pMergedBremSource;
79  DecayChainTrack() : simTrackIndex(-1), pParentVertex(NULL), pMergedBremSource(NULL) {}
80  DecayChainTrack( int newSimTrackIndex ) : simTrackIndex(newSimTrackIndex), pParentVertex(NULL), pMergedBremSource() {}
81  };
82 
87  struct DecayChainVertex
88  {
89  int simVertexIndex;
90  DecayChainTrack* pParentTrack;
91  std::vector<DecayChainTrack*> daughterTracks;
92  DecayChainVertex* pMergedBremSource;
93  DecayChainVertex() : simVertexIndex(-1), pParentTrack(NULL), pMergedBremSource(NULL) {}
94  DecayChainVertex( int newIndex ) : simVertexIndex(newIndex), pParentTrack(NULL), pMergedBremSource(NULL) {}
95  };
96 
107  struct DecayChain
108  {
109  public:
110  DecayChain( const std::vector<SimTrack>& trackCollection, const std::vector<SimVertex>& vertexCollection );
111  const size_t decayTracksSize;
112  const size_t decayVerticesSize;
113 
118  void integrityCheck();
119  const SimTrack& getSimTrack( const ::DecayChainTrack* pDecayTrack ) const { return simTrackCollection_.at( pDecayTrack->simTrackIndex ); }
120  const SimVertex& getSimVertex( const ::DecayChainVertex* pDecayVertex ) const { return simVertexCollection_.at( pDecayVertex->simVertexIndex ); }
121  private:
122  void findBrem( const std::vector<SimTrack>& trackCollection, const std::vector<SimVertex>& vertexCollection );
123  std::unique_ptr< ::DecayChainTrack[]> decayTracks_;
124  std::unique_ptr< ::DecayChainVertex[]> decayVertices_;
125 
127  std::vector< ::DecayChainVertex*> rootVertices_;
128 
129  // Keep a record of the constructor parameters
130  const std::vector<SimTrack>& simTrackCollection_;
131  const std::vector<SimVertex>& simVertexCollection_;
132  public:
133  const std::unique_ptr< ::DecayChainTrack[]>& decayTracks;
134  const std::unique_ptr< ::DecayChainVertex[]>& decayVertices;
135  const std::vector< ::DecayChainVertex*>& rootVertices;
136  };
137 
142  class TrackingParticleFactory
143  {
144  public:
145  TrackingParticleFactory( const ::DecayChain& decayChain, const edm::Handle< std::vector<reco::GenParticle> >& hGenParticles,
146  const edm::Handle< edm::HepMCProduct >& hepMCproduct, const edm::Handle< std::vector<int> >& hHepMCGenParticleIndices,
147  const std::vector<const PSimHit*>& simHits, double volumeRadius, double volumeZ, double vertexDistanceCut, bool allowDifferentProcessTypes );
148  TrackingParticle createTrackingParticle( const DecayChainTrack* pTrack, const TrackerTopology *tTopo ) const;
149  TrackingVertex createTrackingVertex( const DecayChainVertex* pVertex) const;
150  bool vectorIsInsideVolume( const math::XYZTLorentzVectorD& vector ) const;
151  private:
152  const ::DecayChain& decayChain_;
153  const edm::Handle< std::vector<reco::GenParticle> >& hGenParticles_;
154  const edm::Handle< edm::HepMCProduct >& hepMCproduct_;
155  std::vector<int> genParticleIndices_;
156  const std::vector<const PSimHit*>& simHits_;
157  const double volumeRadius_;
158  const double volumeZ_;
159  const double vertexDistanceCut2_; // distance based on which HepMC::GenVertexs are added to SimVertexs
160  std::multimap<unsigned int, size_t> trackIdToHitIndex_;
161  bool allowDifferentProcessTypeForDifferentDetectors_;
162  };
163 
168  class OutputCollectionWrapper
169  {
170  public:
171  OutputCollectionWrapper( const DecayChain& decayChain, TrackingTruthAccumulator::OutputCollections& outputCollections );
172  TrackingParticle* addTrackingParticle( const ::DecayChainTrack* pDecayTrack, const TrackingParticle& trackingParticle );
173  TrackingVertex* addTrackingVertex( const ::DecayChainVertex* pDecayVertex, const TrackingVertex& trackingVertex );
174  TrackingParticle* getTrackingParticle( const ::DecayChainTrack* pDecayTrack );
176  void setProxy( const ::DecayChainTrack* pOriginalTrack, const ::DecayChainTrack* pProxyTrack );
177  void setProxy( const ::DecayChainVertex* pOriginalVertex, const ::DecayChainVertex* pProxyVertex );
178  TrackingVertex* getTrackingVertex( const ::DecayChainVertex* pDecayVertex );
179  TrackingParticleRef getRef( const ::DecayChainTrack* pDecayTrack );
180  TrackingVertexRef getRef( const ::DecayChainVertex* pDecayVertex );
181  private:
182  void associateToExistingObjects( const ::DecayChainVertex* pChainVertex );
183  void associateToExistingObjects( const ::DecayChainTrack* pChainTrack );
185  std::vector<int> trackingParticleIndices_;
186  std::vector<int> trackingVertexIndices_;
187  };
188 
193  TrackingParticle* addTrackAndParentVertex( ::DecayChainTrack* pDecayTrack, const TrackingParticle& trackingParticle, ::OutputCollectionWrapper* pOutput );
194 
199  void addTrack( ::DecayChainTrack* pDecayChainTrack, const TrackingParticleSelector* pSelector, ::OutputCollectionWrapper* pUnmergedOutput, ::OutputCollectionWrapper* pMergedOutput, const ::TrackingParticleFactory& objectFactory, bool addAncestors, const TrackerTopology *tTopo);
200 
201 } // end of the unnamed namespace
202 
203 
204 
205 
206 //---------------------------------------------------------------------------------
207 //---------------------------------------------------------------------------------
208 //---------------------------------------------------------------------------------
209 //------ ------
210 //------ Definitions for the TrackingTruthAccumulator methods ------
211 //------ are below. ------
212 //------ ------
213 //---------------------------------------------------------------------------------
214 //---------------------------------------------------------------------------------
215 //---------------------------------------------------------------------------------
216 
218  messageCategory_("TrackingTruthAccumulator"),
219  volumeRadius_( config.getParameter<double>("volumeRadius") ),
220  volumeZ_( config.getParameter<double>("volumeZ") ),
221  vertexDistanceCut_( config.getParameter<double>("vertexDistanceCut") ),
222  ignoreTracksOutsideVolume_( config.getParameter<bool>("ignoreTracksOutsideVolume") ),
223  maximumPreviousBunchCrossing_( config.getParameter<unsigned int>("maximumPreviousBunchCrossing") ),
224  maximumSubsequentBunchCrossing_( config.getParameter<unsigned int>("maximumSubsequentBunchCrossing") ),
225  createUnmergedCollection_( config.getParameter<bool>("createUnmergedCollection") ),
226  createMergedCollection_(config.getParameter<bool>("createMergedBremsstrahlung") ),
227  addAncestors_( config.getParameter<bool>("alwaysAddAncestors") ),
228  removeDeadModules_( config.getParameter<bool>("removeDeadModules") ),
229  simTrackLabel_( config.getParameter<edm::InputTag>("simTrackCollection") ),
230  simVertexLabel_( config.getParameter<edm::InputTag>("simVertexCollection") ),
231  collectionTags_( ),
232  genParticleLabel_( config.getParameter<edm::InputTag>("genParticleCollection") ),
233  hepMCproductLabel_( config.getParameter<edm::InputTag>("HepMCProductLabel") ),
234  allowDifferentProcessTypeForDifferentDetectors_( config.getParameter<bool>("allowDifferentSimHitProcesses") )
235 {
236  //
237  // Make sure at least one of the merged and unmerged collections have been set
238  // to be created.
239  //
241  edm::LogError(messageCategory_) << "Both \"createUnmergedCollection\" and \"createMergedBremsstrahlung\" have been"
242  << "set to false, which means no collections will be created";
243 
244  // Initialize selection for building TrackingParticles
245  //
246  if( config.exists( "select" ) )
247  {
248  edm::ParameterSet param=config.getParameter<edm::ParameterSet>("select");
249  selector_=TrackingParticleSelector( param.getParameter<double>( "ptMinTP" ),
250  param.getParameter<double>( "minRapidityTP" ),
251  param.getParameter<double>( "maxRapidityTP" ),
252  param.getParameter<double>( "tipTP" ),
253  param.getParameter<double>( "lipTP" ),
254  param.getParameter<int>( "minHitTP" ),
255  param.getParameter<bool>( "signalOnlyTP" ),
256  param.getParameter<bool>( "chargedOnlyTP" ),
257  param.getParameter<bool>( "stableOnlyTP" ),
258  param.getParameter<std::vector<int> >("pdgIdTP") );
259  selectorFlag_=true;
260 
261  // Also set these two variables, which are used to drop out early if the SimTrack doesn't conform.
262  // The selector_ requires a full TrackingParticle object, but these two variables can veto things early.
263  chargedOnly_=param.getParameter<bool>( "chargedOnlyTP" );
264  signalOnly_=param.getParameter<bool>( "signalOnlyTP" );
265  }
266  else
267  {
268  selectorFlag_=false;
269  chargedOnly_=false;
270  signalOnly_=false;
271  }
272 
273  //
274  // Need to state what collections are going to be added to the event. This
275  // depends on which of the merged and unmerged collections have been configured
276  // to be created.
277  //
279  {
280  mixMod.produces<TrackingVertexCollection>();
281  mixMod.produces<TrackingParticleCollection>();
282  }
283 
285  {
286  mixMod.produces<TrackingParticleCollection>("MergedTrackTruth");
287  mixMod.produces<TrackingVertexCollection>("MergedTrackTruth");
288  }
289 
290  iC.consumes<std::vector<SimTrack> >(simTrackLabel_);
291  iC.consumes<std::vector<SimVertex> >(simVertexLabel_);
292  iC.consumes<std::vector<reco::GenParticle> >(genParticleLabel_);
293  iC.consumes<std::vector<int> >(genParticleLabel_);
294  iC.consumes<std::vector<int> >(hepMCproductLabel_);
295 
296  // Fill the collection tags
297  const edm::ParameterSet& simHitCollectionConfig=config.getParameterSet("simHitCollections");
298  std::vector<std::string> parameterNames=simHitCollectionConfig.getParameterNames();
299 
300  for( const auto& parameterName : parameterNames )
301  {
302  std::vector<edm::InputTag> tags=simHitCollectionConfig.getParameter<std::vector<edm::InputTag> >(parameterName);
303  collectionTags_.insert(collectionTags_.end(), tags.begin(), tags.end());
304  }
305 
306  for( const auto& collectionTag : collectionTags_ ) {
307  iC.consumes<std::vector<PSimHit> >(collectionTag);
308  }
309 
310 }
311 
313 {
315  {
318  unmergedOutput_.refTrackingParticles=const_cast<edm::Event&>( event ).getRefBeforePut<TrackingParticleCollection>();
319  unmergedOutput_.refTrackingVertexes=const_cast<edm::Event&>( event ).getRefBeforePut<TrackingVertexCollection>();
320  }
321 
323  {
326  mergedOutput_.refTrackingParticles=const_cast<edm::Event&>( event ).getRefBeforePut<TrackingParticleCollection>("MergedTrackTruth");
327  mergedOutput_.refTrackingVertexes=const_cast<edm::Event&>( event ).getRefBeforePut<TrackingVertexCollection>("MergedTrackTruth");
328  }
329 }
330 
335 
337 {
338  // Call the templated version that does the same for both signal and pileup events
339 
341  event.getByLabel(hepMCproductLabel_, hepmc);
342 
343  accumulateEvent( event, setup, hepmc );
344 }
345 
347 {
348  // If this bunch crossing is outside the user configured limit, don't do anything.
349  if( event.bunchCrossing()>=-static_cast<int>(maximumPreviousBunchCrossing_) && event.bunchCrossing()<=static_cast<int>(maximumSubsequentBunchCrossing_) )
350  {
351  //edm::LogInfo(messageCategory_) << "Analysing pileup event for bunch crossing " << event.bunchCrossing();
352 
353  //simply create empty handle as we do not have a HepMCProduct in PU anyway
355  accumulateEvent( event, setup, hepmc );
356  }
357  else edm::LogInfo(messageCategory_) << "Skipping pileup event for bunch crossing " << event.bunchCrossing();
358 }
359 
361 {
362 
364  {
365  edm::LogInfo("TrackingTruthAccumulator") << "Adding " << unmergedOutput_.pTrackingParticles->size() << " TrackingParticles and " << unmergedOutput_.pTrackingVertices->size()
366  << " TrackingVertexs to the event.";
367 
368  event.put( unmergedOutput_.pTrackingParticles );
369  event.put( unmergedOutput_.pTrackingVertices );
370  }
371 
373  {
374  edm::LogInfo("TrackingTruthAccumulator") << "Adding " << mergedOutput_.pTrackingParticles->size() << " merged TrackingParticles and " << mergedOutput_.pTrackingVertices->size()
375  << " merged TrackingVertexs to the event.";
376 
377  event.put( mergedOutput_.pTrackingParticles, "MergedTrackTruth" );
378  event.put( mergedOutput_.pTrackingVertices, "MergedTrackTruth" );
379  }
380 
381 }
382 
384 {
385  //
386  // Get the collections
387  //
389  edm::Handle<std::vector<SimVertex> > hSimVertices;
391  edm::Handle< std::vector<int> > hGenParticleIndices;
392 
393  event.getByLabel( simTrackLabel_, hSimTracks );
394  event.getByLabel( simVertexLabel_, hSimVertices );
395 
396  try
397  {
398  event.getByLabel( genParticleLabel_, hGenParticles );
399  event.getByLabel( genParticleLabel_, hGenParticleIndices );
400  }
401  catch( cms::Exception& exception )
402  {
403  //
404  // The Monte Carlo is not always available, e.g. for pileup events. The information
405  // is only used if it's available, but for some reason the PileUpEventPrincipal
406  // wrapper throws an exception here rather than waiting to see if the handle is
407  // used (as is the case for edm::Event). So I just want to catch this exception
408  // and use the normal handle checking later on.
409  //
410  }
411 
412  //Retrieve tracker topology from geometry
413  edm::ESHandle<TrackerTopology> tTopoHandle;
414  setup.get<IdealGeometryRecord>().get(tTopoHandle);
415  const TrackerTopology* const tTopo = tTopoHandle.product();
416 
417 
418  // Run through the collections and work out the decay chain of each track/vertex. The
419  // information in SimTrack and SimVertex only allows traversing upwards, but this will
420  // allow traversal in both directions. This is required for things like grouping electrons
421  // that bremsstrahlung as one TrackingParticle if "mergedBremsstrahlung" is set in the
422  // config file.
423  DecayChain decayChain( *hSimTracks, *hSimVertices );
424 
425  // I only want to create these collections if they're actually required
426  std::auto_ptr< ::OutputCollectionWrapper> pUnmergedCollectionWrapper;
427  std::auto_ptr< ::OutputCollectionWrapper> pMergedCollectionWrapper;
428  if( createUnmergedCollection_ ) pUnmergedCollectionWrapper.reset( new ::OutputCollectionWrapper( decayChain, unmergedOutput_ ) );
429  if( createMergedCollection_ ) pMergedCollectionWrapper.reset( new ::OutputCollectionWrapper( decayChain, mergedOutput_ ) );
430 
431  std::vector<const PSimHit*> simHitPointers;
432  fillSimHits( simHitPointers, event, setup );
433  TrackingParticleFactory objectFactory( decayChain, hGenParticles, hepMCproduct, hGenParticleIndices, simHitPointers, volumeRadius_, volumeZ_, vertexDistanceCut_, allowDifferentProcessTypeForDifferentDetectors_ );
434 
435  // While I'm testing, perform some checks.
436  // TODO - drop this call once I'm happy it works in all situations.
437  //decayChain.integrityCheck();
438 
439  TrackingParticleSelector* pSelector=NULL;
440  if( selectorFlag_ ) pSelector=&selector_;
441 
442  // Run over all of the SimTracks, but because I'm interested in the decay hierarchy
443  // do it through the DecayChainTrack objects. These are looped over in sequence here
444  // but they have the hierarchy information for the functions called to traverse the
445  // decay chain.
446 
447  for( size_t index=0; index<decayChain.decayTracksSize; ++index )
448  {
449  ::DecayChainTrack* pDecayTrack=&decayChain.decayTracks[index];
450  const SimTrack& simTrack=hSimTracks->at(pDecayTrack->simTrackIndex);
451 
452 
453  // Perform some quick checks to see if we can drop out early. Note that these are
454  // a subset of the cuts in the selector_ so the created TrackingParticle could still
455  // fail. The selector_ requires the full TrackingParticle to be made however, which
456  // can be computationally expensive.
457  if( chargedOnly_ && simTrack.charge()==0 ) continue;
458  if( signalOnly_ && (simTrack.eventId().bunchCrossing()!=0 || simTrack.eventId().event()!=0) ) continue;
459 
460  // Also perform a check to see if the production vertex is inside the tracker volume (if required).
462  {
463  const SimVertex& simVertex=hSimVertices->at( pDecayTrack->pParentVertex->simVertexIndex );
464  if( !objectFactory.vectorIsInsideVolume( simVertex.position() ) ) continue;
465  }
466 
467 
468  // This function creates the TrackinParticle and adds it to the collection if it
469  // passes the selection criteria specified in the configuration. If the config
470  // specifies adding ancestors, the function is called recursively to do that.
471  ::addTrack( pDecayTrack, pSelector, pUnmergedCollectionWrapper.get(), pMergedCollectionWrapper.get(), objectFactory, addAncestors_, tTopo );
472  }
473 }
474 
475 template<class T> void TrackingTruthAccumulator::fillSimHits( std::vector<const PSimHit*>& returnValue, const T& event, const edm::EventSetup& setup )
476 {
477  // loop over the collections
478  for( const auto& collectionTag : collectionTags_ )
479  {
481  event.getByLabel( collectionTag, hSimHits );
482 
483  // TODO - implement removing the dead modules
484  for( const auto& simHit : *hSimHits )
485  {
486  returnValue.push_back( &simHit );
487  }
488 
489  } // end of loop over InputTags
490 }
491 
492 
493 //---------------------------------------------------------------------------------
494 //---------------------------------------------------------------------------------
495 //---------------------------------------------------------------------------------
496 //------ ------
497 //------ End of the definitions for the TrackingTruthAccumulator methods. ------
498 //------ Definitions for the functions and classes declared in the unnamed ------
499 //------ namespace are below. Documentation for those is above, where ------
500 //------ they're declared. ------
501 //------ ------
502 //---------------------------------------------------------------------------------
503 //---------------------------------------------------------------------------------
504 //---------------------------------------------------------------------------------
505 
506 
507 namespace // Unnamed namespace for things only used in this file
508 {
509  //---------------------------------------------------------------------------------
510  //---------------------------------------------------------------------------------
511  //---- TrackingParticleFactory methods ----------------------------------------
512  //---------------------------------------------------------------------------------
513  //---------------------------------------------------------------------------------
514 
515  ::TrackingParticleFactory::TrackingParticleFactory( const ::DecayChain& decayChain, const edm::Handle< std::vector<reco::GenParticle> >& hGenParticles,
516  const edm::Handle< edm::HepMCProduct >& hepMCproduct, const edm::Handle< std::vector<int> >& hHepMCGenParticleIndices, const std::vector<const PSimHit*>& simHits,
517  double volumeRadius, double volumeZ, double vertexDistanceCut, bool allowDifferentProcessTypes )
518  : decayChain_(decayChain), hGenParticles_(hGenParticles), hepMCproduct_(hepMCproduct), simHits_(simHits), volumeRadius_(volumeRadius),
519  volumeZ_(volumeZ), vertexDistanceCut2_(vertexDistanceCut*vertexDistanceCut), allowDifferentProcessTypeForDifferentDetectors_(allowDifferentProcessTypes)
520  {
521  // Need to create a multimap to get from a SimTrackId to all of the hits in it. The SimTrackId
522  // is an unsigned int.
523  for( size_t index=0; index<simHits_.size(); ++index )
524  {
525  trackIdToHitIndex_.insert( std::make_pair( simHits_[index]->trackId(), index ) );
526  }
527 
528  if( hHepMCGenParticleIndices.isValid() ) // Monte Carlo might not be available for the pileup events
529  {
530  genParticleIndices_.resize( hHepMCGenParticleIndices->size()+1 );
531 
532  // What I need is the reverse mapping of this vector. The sizes are already equivalent because I set
533  // the size in the initialiser list.
534  for( size_t recoGenParticleIndex=0; recoGenParticleIndex<hHepMCGenParticleIndices->size(); ++recoGenParticleIndex )
535  {
536  size_t hepMCGenParticleIndex=(*hHepMCGenParticleIndices)[recoGenParticleIndex];
537 
538  // They should be the same size, give or take a fencepost error, so this should never happen - but just in case
539  if( genParticleIndices_.size()<=hepMCGenParticleIndex ) genParticleIndices_.resize(hepMCGenParticleIndex+1);
540 
541  genParticleIndices_[ hepMCGenParticleIndex ]=recoGenParticleIndex;
542  }
543  }
544  }
545 
546  TrackingParticle TrackingParticleFactory::createTrackingParticle( const ::DecayChainTrack* pChainTrack, const TrackerTopology *tTopo ) const
547  {
549  typedef math::XYZPoint Vector;
550 
551  const SimTrack& simTrack=decayChain_.getSimTrack( pChainTrack );
552  const SimVertex& parentSimVertex=decayChain_.getSimVertex( pChainTrack->pParentVertex );
553 
554  LorentzVector position( 0, 0, 0, 0 );
555  if( !simTrack.noVertex() ) position=parentSimVertex.position();
556 
557  int pdgId=simTrack.type();
558 
559  TrackingParticle returnValue;
560  // N.B. The sim track is added a few lines below, the parent and decay vertices are added in
561  // another function later on.
562 
563  //
564  // If there is some valid Monte Carlo for this track, take some information from that.
565  // Only do so if it is from the signal event however. Not sure why but that's what the
566  // old code did.
567  //
568  if( simTrack.eventId().event()==0 && simTrack.eventId().bunchCrossing()==0 ) // if this is a track in the signal event
569  {
570  int hepMCGenParticleIndex=simTrack.genpartIndex();
571  if( hepMCGenParticleIndex>=0 && hepMCGenParticleIndex<static_cast<int>(genParticleIndices_.size()) && hGenParticles_.isValid() )
572  {
573  int recoGenParticleIndex=genParticleIndices_[hepMCGenParticleIndex];
574  reco::GenParticleRef generatorParticleRef( hGenParticles_, recoGenParticleIndex );
575  pdgId=generatorParticleRef->pdgId();
576  returnValue.addGenParticle( generatorParticleRef );
577  }
578  }
579 
580  returnValue.addG4Track( simTrack );
581 
582  // I need to run over the sim hits and count up the different types of hits. To be honest
583  // I don't fully understand this next section of code, I copied it from the old TrackingTruthProducer
584  // to provide the same results. The different types of hits that I need to count are:
585  size_t matchedHits=0; // As far as I can tell, this is the number of tracker layers with hits on them, i.e. taking account of overlaps.
586  size_t numberOfHits=0; // The total number of hits not taking account of overlaps.
587  size_t numberOfTrackerHits=0; // The number of tracker hits not taking account of overlaps.
588 
589  bool init=true;
590  int processType=0;
591  int particleType=0;
592  int oldLayer = 0;
593  int newLayer = 0;
594  DetId oldDetector;
595  DetId newDetector;
596 
597  for( std::multimap<unsigned int,size_t>::const_iterator iHitIndex=trackIdToHitIndex_.lower_bound( simTrack.trackId() );
598  iHitIndex!=trackIdToHitIndex_.upper_bound( simTrack.trackId() );
599  ++iHitIndex )
600  {
601  const auto& pSimHit=simHits_[ iHitIndex->second ];
602 
603  // Initial condition for consistent simhit selection
604  if( init )
605  {
606  processType=pSimHit->processType();
607  particleType=pSimHit->particleType();
608  newDetector=DetId( pSimHit->detUnitId() );
609  init=false;
610  }
611 
612  oldDetector=newDetector;
613  newDetector=DetId( pSimHit->detUnitId() );
614 
615  // Fast sim seems to have different processType between tracker and muons, so if this flag
616  // has been set allow the processType to change if the detector changes.
617  if( allowDifferentProcessTypeForDifferentDetectors_ && newDetector.det()!=oldDetector.det() ) processType=pSimHit->processType();
618 
619  // Check for delta and interaction products discards
620  if( processType==pSimHit->processType() && particleType==pSimHit->particleType() && pdgId==pSimHit->particleType() )
621  {
622  ++numberOfHits;
623  oldLayer=newLayer;
624  newLayer=0;
625  if( newDetector.det() == DetId::Tracker ) {
626  ++numberOfTrackerHits;
627 
628  newLayer=tTopo->layer( newDetector );
629 
630  // Count hits using layers for glued detectors
631  if( (oldLayer!=newLayer || (oldLayer==newLayer && oldDetector.subdetId()!=newDetector.subdetId())) ) ++matchedHits;
632  }
633  }
634  } // end of loop over the sim hits for this sim track
635 
636  returnValue.setNumberOfTrackerLayers( matchedHits );
637  returnValue.setNumberOfHits( numberOfHits );
638  returnValue.setNumberOfTrackerHits( numberOfTrackerHits );
639 
640  return returnValue;
641  }
642 
643  TrackingVertex TrackingParticleFactory::createTrackingVertex( const ::DecayChainVertex* pChainVertex) const
644  {
645 
647  typedef math::XYZPoint Vector;
649 
650  const SimVertex& simVertex=decayChain_.getSimVertex( pChainVertex );
651 
652  bool isInVolume=this->vectorIsInsideVolume( simVertex.position() );
653 
654  TrackingVertex returnValue( simVertex.position(), isInVolume, simVertex.eventId() );
655 
656  // add the SimVertex to the TrackingVertex
657  returnValue.addG4Vertex(simVertex);
658 
659  // also add refs to nearby HepMC::GenVertexs; the way this is done (i.e. based on the position) is transcribed over from the old TrackingTruthProducer
660  if( simVertex.eventId().event()==0 && simVertex.eventId().bunchCrossing()==0 && hepMCproduct_.isValid()) // if this is a track in the signal event
661  {
662  const HepMC::GenEvent* genEvent = hepMCproduct_->GetEvent();
663 
664  if (genEvent != NULL)
665  {
666  Vector tvPosition(returnValue.position().x(), returnValue.position().y(), returnValue.position().z());
667 
668  for (HepMC::GenEvent::vertex_const_iterator iGenVertex = genEvent->vertices_begin(); iGenVertex != genEvent->vertices_end(); ++iGenVertex)
669  {
670  HepMC::ThreeVector rawPosition = (*iGenVertex)->position();
671 
672  Vector genPosition(rawPosition.x()*0.1, rawPosition.y()*0.1, rawPosition.z()*0.1);
673 
674  auto distance2 = (tvPosition - genPosition).mag2();
675 
676  if (distance2 < vertexDistanceCut2_)
677  returnValue.addGenVertex( GenVertexRef(hepMCproduct_, (*iGenVertex)->barcode()) );
678  }
679  }
680  }
681 
682  return returnValue;
683  }
684 
685  bool ::TrackingParticleFactory::vectorIsInsideVolume( const math::XYZTLorentzVectorD& vector ) const
686  {
687  return ( vector.Pt()<volumeRadius_ && vector.z()<volumeZ_ );
688  }
689 
690  //---------------------------------------------------------------------------------
691  //---------------------------------------------------------------------------------
692  //---- DecayChain methods -----------------------------------------------------
693  //---------------------------------------------------------------------------------
694  //---------------------------------------------------------------------------------
695 
696  ::DecayChain::DecayChain( const std::vector<SimTrack>& trackCollection, const std::vector<SimVertex>& vertexCollection )
697  : decayTracksSize( trackCollection.size() ),
698  decayVerticesSize( vertexCollection.size() ),
699  decayTracks_( new DecayChainTrack[decayTracksSize] ),
700  decayVertices_( new DecayChainVertex[decayVerticesSize] ),
701  simTrackCollection_( trackCollection ),
702  simVertexCollection_( vertexCollection ),
703  decayTracks( decayTracks_ ), // These const references map onto the actual members for public const access
704  decayVertices( decayVertices_ ),
705  rootVertices( rootVertices_ )
706  {
707  // I need some maps to be able to get object pointers from the track/vertex ID
708  std::map<int,::DecayChainTrack*> trackIdToDecayTrack;
709  std::map<int,::DecayChainVertex*> vertexIdToDecayVertex;
710 
711  // First create a DecayChainTrack for every SimTrack and make a note of the
712  // trackIds in the map. Also add a pointer to the daughter list of the parent
713  // DecayChainVertex, which might include creating the vertex object if it
714  // doesn't already exist.
715  size_t decayVertexIndex=0; // The index of the next free slot in the DecayChainVertex array.
716  for( size_t index=0; index<trackCollection.size(); ++index )
717  {
718  ::DecayChainTrack* pDecayTrack=&decayTracks_[index]; // Get a pointer for ease of use
719 
720  // This is the array index of the SimTrack that this DecayChainTrack corresponds to. At
721  // the moment this is a 1 to 1 relationship with the DecayChainTrack index, but they won't
722  // necessarily be accessed through the array later so it's still required to store it.
723  pDecayTrack->simTrackIndex=index;
724 
725  trackIdToDecayTrack[ trackCollection[index].trackId() ]=pDecayTrack;
726 
727  int parentVertexIndex=trackCollection[index].vertIndex();
728  if( parentVertexIndex>=0 )
729  {
730  // Get the DecayChainVertex corresponding to this SimVertex, or initialise it if it hasn't been done already.
731  ::DecayChainVertex*& pParentVertex=vertexIdToDecayVertex[parentVertexIndex];
732  if( pParentVertex==NULL )
733  {
734  // Note that I'm using a reference, so changing pParentVertex will change the entry in the map too.
735  pParentVertex=&decayVertices_[decayVertexIndex];
736  ++decayVertexIndex;
737  pParentVertex->simVertexIndex=parentVertexIndex;
738  }
739  pParentVertex->daughterTracks.push_back(pDecayTrack);
740  pDecayTrack->pParentVertex=pParentVertex;
741  }
742  else throw std::runtime_error( "TrackingTruthAccumulator: Found a track with an invalid parent vertex index." );
743  }
744 
745  // This assert was originally in to check the internal consistency of the decay chain. Fast sim
746  // pileup seems to have a load of vertices with no tracks pointing to them though, so fast sim fails
747  // this assert if pileup is added. I don't think the problem is vital however, so if the assert is
748  // taken out these extra vertices are ignored.
749  //assert( vertexIdToDecayVertex.size()==vertexCollection.size() && vertexCollection.size()==decayVertexIndex );
750 
751  // I still need to set DecayChainTrack::daughterVertices and DecayChainVertex::pParentTrack.
752  // The information to do this comes from SimVertex::parentIndex. I couldn't do this before
753  // because I need all of the DecayChainTracks initialised.
754  for( auto& decayVertexMapPair : vertexIdToDecayVertex )
755  {
756  ::DecayChainVertex* pDecayVertex=decayVertexMapPair.second;
757  int parentTrackIndex=vertexCollection[pDecayVertex->simVertexIndex].parentIndex();
758  if( parentTrackIndex!=-1 )
759  {
760  std::map<int,::DecayChainTrack*>::iterator iParentTrackMapPair=trackIdToDecayTrack.find(parentTrackIndex);
761  if( iParentTrackMapPair==trackIdToDecayTrack.end() )
762  {
763  std::stringstream errorStream;
764  errorStream << "TrackingTruthAccumulator: Something has gone wrong with the indexing. Parent track index is " << parentTrackIndex << ".";
765  throw std::runtime_error( errorStream.str() );
766  }
767 
768  ::DecayChainTrack* pParentTrackHierarchy=iParentTrackMapPair->second;
769 
770  pParentTrackHierarchy->daughterVertices.push_back( pDecayVertex );
771  pDecayVertex->pParentTrack=pParentTrackHierarchy;
772  }
773  else rootVertices_.push_back(pDecayVertex); // Has no parent so is at the top of the decay chain.
774  } // end of loop over the vertexIdToDecayVertex map
775 
776  findBrem( trackCollection, vertexCollection );
777 
778  } // end of ::DecayChain constructor
779 
780  // Function documentation is with the declaration above. This function is only used for testing.
781  void ::DecayChain::integrityCheck()
782  {
783  //
784  // Start off checking the tracks
785  //
786  for( size_t index=0; index<decayTracksSize; ++index )
787  {
788  const auto& decayTrack=decayTracks[index];
789  //
790  // Tracks should always have a production vertex
791  //
792  if( decayTrack.pParentVertex==NULL ) throw std::runtime_error( "TrackingTruthAccumulator.cc integrityCheck(): Found DecayChainTrack with no parent vertex." );
793 
794  //
795  // Make sure the parent has this as a child. Also make sure it's only listed once.
796  //
797  size_t numberOfTimesListed=0;
798  for( const auto pSiblingTrack : decayTrack.pParentVertex->daughterTracks )
799  {
800  if( pSiblingTrack==&decayTrack ) ++numberOfTimesListed;
801  }
802  if( numberOfTimesListed!=1 ) throw std::runtime_error( "TrackingTruthAccumulator.cc integrityCheck(): Found DecayChainTrack whose parent does not have it listed once and only once as a child." );
803 
804  //
805  // Make sure all of the children have this listed as a parent.
806  //
807  for( const auto pDaughterVertex : decayTrack.daughterVertices )
808  {
809  if( pDaughterVertex->pParentTrack!=&decayTrack ) throw std::runtime_error( "TrackingTruthAccumulator.cc integrityCheck(): Found DecayChainTrack whose child does not have it listed as the parent." );
810  }
811 
812  //
813  // Follow the chain up to the root vertex, and make sure that is listed in rootVertices
814  //
815  const DecayChainVertex* pAncestorVertex=decayTrack.pParentVertex;
816  while( pAncestorVertex->pParentTrack!=NULL )
817  {
818  if( pAncestorVertex->pParentTrack->pParentVertex==NULL ) throw std::runtime_error( "TrackingTruthAccumulator.cc integrityCheck(): Found DecayChainTrack with no parent vertex higher in the decay chain." );
819  pAncestorVertex=pAncestorVertex->pParentTrack->pParentVertex;
820  }
821  if( std::find( rootVertices.begin(), rootVertices.end(), pAncestorVertex )==rootVertices.end() )
822  {
823  throw std::runtime_error( "TrackingTruthAccumulator.cc integrityCheck(): Found DecayChainTrack whose root vertex is not recorded anywhere." );
824  }
825  } // end of loop over decayTracks
826 
827  //
828  // Now check the vertices
829  //
830  for( size_t index=0; index<decayVerticesSize; ++index )
831  {
832  const auto& decayVertex=decayVertices[index];
833 
834  //
835  // Make sure this, or a vertex higher in the chain, is in the list of root vertices.
836  //
837  const DecayChainVertex* pAncestorVertex=&decayVertex;
838  while( pAncestorVertex->pParentTrack!=NULL )
839  {
840  if( pAncestorVertex->pParentTrack->pParentVertex==NULL ) throw std::runtime_error( "TrackingTruthAccumulator.cc integrityCheck(): Found DecayChainTrack with no parent vertex higher in the vertex decay chain." );
841  pAncestorVertex=pAncestorVertex->pParentTrack->pParentVertex;
842  }
843  if( std::find( rootVertices.begin(), rootVertices.end(), pAncestorVertex )==rootVertices.end() )
844  {
845  throw std::runtime_error( "TrackingTruthAccumulator.cc integrityCheck(): Found DecayChainTrack whose root vertex is not recorded anywhere." );
846  }
847 
848  //
849  // Make sure the parent has this listed in its daughters once and only once.
850  //
851  if( decayVertex.pParentTrack!=NULL )
852  {
853  size_t numberOfTimesListed=0;
854  for( const auto pSibling : decayVertex.pParentTrack->daughterVertices )
855  {
856  if( pSibling==&decayVertex ) ++numberOfTimesListed;
857  }
858  if( numberOfTimesListed!=1 ) throw std::runtime_error( "TrackingTruthAccumulator.cc integrityCheck(): Found DecayChainVertex whose parent does not have it listed once and only once as a child." );
859  }
860 
861  //
862  // Make sure all of the children have this listed as the parent
863  //
864  for( const auto pDaughter : decayVertex.daughterTracks )
865  {
866  if( pDaughter->pParentVertex!=&decayVertex ) throw std::runtime_error( "TrackingTruthAccumulator.cc integrityCheck(): Found DecayChainVertex whose child does not have it listed as the parent." );
867  }
868  } // end of loop over decay vertices
869 
870  std::cout << "TrackingTruthAccumulator.cc integrityCheck() completed successfully" << std::endl;
871  } // end of ::DecayChain::integrityCheck()
872 
873  void ::DecayChain::findBrem( const std::vector<SimTrack>& trackCollection, const std::vector<SimVertex>& vertexCollection )
874  {
875  for( size_t index=0; index<decayVerticesSize; ++index )
876  {
877  auto& vertex=decayVertices_[index];
878 
879  // Make sure parent is an electron
880  if( vertex.pParentTrack==NULL ) continue;
881  int parentTrackPDG=trackCollection[vertex.pParentTrack->simTrackIndex].type();
882  if( std::abs( parentTrackPDG )!=11 ) continue;
883 
884  size_t numberOfElectrons=0;
885  size_t numberOfNonElectronsOrPhotons=0;
886  for( auto& pDaughterTrack : vertex.daughterTracks )
887  {
888  const auto& simTrack=trackCollection[pDaughterTrack->simTrackIndex];
889  if( simTrack.type()==11 || simTrack.type()==-11 ) ++numberOfElectrons;
890  else if( simTrack.type()!=22 ) ++numberOfNonElectronsOrPhotons;
891  }
892  if( numberOfElectrons==1 && numberOfNonElectronsOrPhotons==0 )
893  {
894  // This is a valid brem, run through and tell all daughters to use the parent
895  // as the brem
896  for( auto& pDaughterTrack : vertex.daughterTracks ) pDaughterTrack->pMergedBremSource=vertex.pParentTrack;
897  vertex.pMergedBremSource=vertex.pParentTrack->pParentVertex;
898  }
899  }
900  } // end of ::DecayChain::findBrem()
901 
902 
903  //---------------------------------------------------------------------------------
904  //---------------------------------------------------------------------------------
905  //---- OutputCollectionWrapper methods ----------------------------------------
906  //---------------------------------------------------------------------------------
907  //---------------------------------------------------------------------------------
908 
909  ::OutputCollectionWrapper::OutputCollectionWrapper( const ::DecayChain& decayChain, TrackingTruthAccumulator::OutputCollections& outputCollections )
910  : output_(outputCollections),
911  trackingParticleIndices_(decayChain.decayTracksSize,-1),
912  trackingVertexIndices_(decayChain.decayVerticesSize,-1)
913  {
914  // No operation besides the initialiser list
915  }
916 
917  TrackingParticle* ::OutputCollectionWrapper::addTrackingParticle( const ::DecayChainTrack* pDecayTrack, const TrackingParticle& trackingParticle )
918  {
919  if( trackingParticleIndices_[pDecayTrack->simTrackIndex]!=-1 ) throw std::runtime_error( "OutputCollectionWrapper::addTrackingParticle - trying to add a particle twice" );
920 
921  trackingParticleIndices_[pDecayTrack->simTrackIndex]=output_.pTrackingParticles->size();
922  output_.pTrackingParticles->push_back( trackingParticle );
923 
924  // Clear any associations in case there were any beforehand
925  output_.pTrackingParticles->back().clearDecayVertices();
926  output_.pTrackingParticles->back().clearParentVertex();
927 
928  // Associate to the objects that are already in the output collections
929  associateToExistingObjects( pDecayTrack );
930 
931  return &output_.pTrackingParticles->back();
932  }
933 
934  TrackingVertex* ::OutputCollectionWrapper::addTrackingVertex( const ::DecayChainVertex* pDecayVertex, const TrackingVertex& trackingVertex )
935  {
936  if( trackingVertexIndices_[pDecayVertex->simVertexIndex]!=-1 ) throw std::runtime_error( "OutputCollectionWrapper::addTrackingVertex - trying to add a vertex twice" );
937 
938  trackingVertexIndices_[pDecayVertex->simVertexIndex]=output_.pTrackingVertices->size();
939  output_.pTrackingVertices->push_back( trackingVertex );
940 
941  // Associate the new vertex to any parents or children that are already in the output collections
942  associateToExistingObjects( pDecayVertex );
943 
944  return &output_.pTrackingVertices->back();
945  }
946 
947  TrackingParticle* ::OutputCollectionWrapper::getTrackingParticle( const ::DecayChainTrack* pDecayTrack )
948  {
949  const int index=trackingParticleIndices_[pDecayTrack->simTrackIndex];
950  if( index==-1 ) return NULL;
951  else return &(*output_.pTrackingParticles)[index];
952  }
953 
954  TrackingVertex* ::OutputCollectionWrapper::getTrackingVertex( const ::DecayChainVertex* pDecayVertex )
955  {
956  const int index=trackingVertexIndices_[pDecayVertex->simVertexIndex];
957  if( index==-1 ) return NULL;
958  else return &(*output_.pTrackingVertices)[index];
959  }
960 
961  TrackingParticleRef OutputCollectionWrapper::getRef( const ::DecayChainTrack* pDecayTrack )
962  {
963  const int index=trackingParticleIndices_[pDecayTrack->simTrackIndex];
964  if( index==-1 ) throw std::runtime_error( "OutputCollectionWrapper::getRefTrackingParticle - ref requested for a non existent TrackingParticle" );
965  else return TrackingParticleRef( output_.refTrackingParticles, index );
966  }
967 
968  TrackingVertexRef OutputCollectionWrapper::getRef( const ::DecayChainVertex* pDecayVertex )
969  {
970  const int index=trackingVertexIndices_[pDecayVertex->simVertexIndex];
971  if( index==-1 ) throw std::runtime_error( "OutputCollectionWrapper::getRefTrackingParticle - ref requested for a non existent TrackingVertex" );
972  else return TrackingVertexRef( output_.refTrackingVertexes, index );
973  }
974 
975  void ::OutputCollectionWrapper::setProxy( const ::DecayChainTrack* pOriginalTrack, const ::DecayChainTrack* pProxyTrack )
976  {
977  int& index=trackingParticleIndices_[pOriginalTrack->simTrackIndex];
978  if( index!=-1 ) throw std::runtime_error( "OutputCollectionWrapper::setProxy() was called for a TrackingParticle that has already been created" );
979  // Note that index is a reference so this call changes the value in the vector too
980  index=trackingParticleIndices_[pProxyTrack->simTrackIndex];
981  }
982 
983  void ::OutputCollectionWrapper::setProxy( const ::DecayChainVertex* pOriginalVertex, const ::DecayChainVertex* pProxyVertex )
984  {
985  int& index=trackingVertexIndices_[pOriginalVertex->simVertexIndex];
986  const int newIndex=trackingVertexIndices_[pProxyVertex->simVertexIndex];
987 
988  if( index!=-1 && index!=newIndex ) throw std::runtime_error( "OutputCollectionWrapper::setProxy() was called for a TrackingVertex that has already been created" );
989 
990  // Note that index is a reference so this call changes the value in the vector too
991  index=newIndex;
992  }
993 
994  void ::OutputCollectionWrapper::associateToExistingObjects( const ::DecayChainVertex* pChainVertex )
995  {
996  // First make sure the DecayChainVertex supplied has been turned into a TrackingVertex
997  TrackingVertex* pTrackingVertex=getTrackingVertex( pChainVertex );
998  if( pTrackingVertex==NULL ) throw std::runtime_error( "associateToExistingObjects was passed a non existent TrackingVertex" );
999 
1000  //
1001  // Associate to the parent track (if there is one)
1002  //
1003  ::DecayChainTrack* pParentChainTrack=pChainVertex->pParentTrack;
1004  if( pParentChainTrack!=NULL ) // Make sure there is a parent track first
1005  {
1006  // There is a parent track, but it might not have been turned into a TrackingParticle yet
1007  TrackingParticle* pParentTrackingParticle=getTrackingParticle(pParentChainTrack);
1008  if( pParentTrackingParticle!=NULL )
1009  {
1010  pParentTrackingParticle->addDecayVertex( getRef(pChainVertex) );
1011  pTrackingVertex->addParentTrack( getRef(pParentChainTrack) );
1012  }
1013  // Don't worry about the "else" case - the parent track might not necessarily get added
1014  // to the output collection at all. A check is done on daughter vertices when tracks
1015  // are added, so the association will be picked up then.
1016  }
1017 
1018  //
1019  // A parent TrackingVertex is always associated to a daughter TrackingParticle when
1020  // the TrackingParticle is created. Hence there is no need to check the list of
1021  // daughter tracks.
1022  //
1023  }
1024 
1025  void ::OutputCollectionWrapper::associateToExistingObjects( const ::DecayChainTrack* pChainTrack )
1026  {
1027  //
1028  // First make sure this DecayChainTrack has been turned into a TrackingParticle
1029  //
1030  TrackingParticle* pTrackingParticle=getTrackingParticle( pChainTrack );
1031  if( pTrackingParticle==NULL ) throw std::runtime_error( "associateToExistingObjects was passed a non existent TrackingParticle" );
1032 
1033  // Get the parent vertex. This should always already have been turned into a TrackingVertex, and
1034  // there should always be a parent DecayChainVertex.
1035  ::DecayChainVertex* pParentChainVertex=pChainTrack->pParentVertex;
1036  TrackingVertex* pParentTrackingVertex=getTrackingVertex( pParentChainVertex );
1037 
1038  //
1039  // First associate to the parent vertex.
1040  //
1041  pTrackingParticle->setParentVertex( getRef(pParentChainVertex) );
1042  pParentTrackingVertex->addDaughterTrack( getRef(pChainTrack) );
1043 
1044  // If there are any daughter vertices that have already been made into a TrackingVertex
1045  // make sure they know about each other. If the daughter vertices haven't been made into
1046  // TrackingVertexs yet, I won't do anything and create the association when the vertex is
1047  // made.
1048  for( auto pDaughterChainVertex : pChainTrack->daughterVertices )
1049  {
1050  TrackingVertex* pDaughterTrackingVertex=getTrackingVertex( pDaughterChainVertex );
1051  if( pDaughterTrackingVertex!=NULL )
1052  {
1053  pTrackingParticle->addDecayVertex( getRef(pDaughterChainVertex) );
1054  pDaughterTrackingVertex->addParentTrack( getRef(pChainTrack) );
1055  }
1056  }
1057  }
1058 
1059 
1060  TrackingParticle* addTrackAndParentVertex( ::DecayChainTrack* pDecayTrack, const TrackingParticle& trackingParticle, ::OutputCollectionWrapper* pOutput )
1061  {
1062  // See if this TrackingParticle has already been created (could be if the DecayChainTracks are
1063  // looped over in a funny order). If it has then there's no need to do anything.
1064  TrackingParticle* pTrackingParticle=pOutput->getTrackingParticle( pDecayTrack );
1065  if( pTrackingParticle==NULL )
1066  {
1067  // Need to make sure the production vertex has been created first
1068  TrackingVertex* pProductionVertex=pOutput->getTrackingVertex( pDecayTrack->pParentVertex );
1069  if( pProductionVertex==NULL )
1070  {
1071  // TrackingVertex doesn't exist in the output collection yet. However, it's already been
1072  // created in the addTrack() function and a temporary reference to it set in the TrackingParticle.
1073  // I'll use that reference to create a copy in the output collection. When the TrackingParticle
1074  // is added to the output collection a few lines below the temporary reference to the parent
1075  // vertex will be cleared, and the correct one referring to the output collection will be set.
1076  pProductionVertex=pOutput->addTrackingVertex( pDecayTrack->pParentVertex, *trackingParticle.parentVertex() );
1077  }
1078 
1079 
1080  pTrackingParticle=pOutput->addTrackingParticle( pDecayTrack, trackingParticle );
1081  }
1082 
1083  return pTrackingParticle;
1084  }
1085 
1086  void addTrack( ::DecayChainTrack* pDecayChainTrack, const TrackingParticleSelector* pSelector, ::OutputCollectionWrapper* pUnmergedOutput, ::OutputCollectionWrapper* pMergedOutput, const ::TrackingParticleFactory& objectFactory, bool addAncestors, const TrackerTopology *tTopo )
1087  {
1088  if( pDecayChainTrack==NULL ) return; // This is required for when the addAncestors_ recursive call reaches the top of the chain
1089 
1090  // Check to see if this particle has already been processed while traversing up the parents
1091  // of another split in the decay chain. The check in the line above only stops when the top
1092  // of the chain is reached, whereas this will stop when a previously traversed split is reached.
1093  { // block to limit the scope of temporary variables
1094  bool alreadyProcessed=true;
1095  if( pUnmergedOutput!=NULL )
1096  {
1097  if( pUnmergedOutput->getTrackingParticle( pDecayChainTrack )==NULL ) alreadyProcessed=false;
1098  }
1099  if( pMergedOutput!=NULL )
1100  {
1101  if( pMergedOutput->getTrackingParticle( pDecayChainTrack )==NULL ) alreadyProcessed=false;
1102  }
1103  if( alreadyProcessed ) return;
1104  }
1105 
1106  // Create a TrackingParticle.
1107  TrackingParticle newTrackingParticle=objectFactory.createTrackingParticle( pDecayChainTrack, tTopo );
1108 
1109  // The selector checks the impact parameters from the vertex, so I need to have a valid reference
1110  // to the parent vertex in the TrackingParticle before that can be called. TrackingParticle needs
1111  // an edm::Ref for the parent TrackingVertex though. I still don't know if this is going to be
1112  // added to the collection so I can't take it from there, so I need to create a temporary one.
1113  // When the addTrackAndParentVertex() is called (assuming it passes selection) it will use the
1114  // temporary reference to create a copy of the parent vertex, put that in the output collection,
1115  // and then set the reference in the TrackingParticle properly.
1116  TrackingVertexCollection dummyCollection; // Only needed to create an edm::Ref
1117  dummyCollection.push_back( objectFactory.createTrackingVertex( pDecayChainTrack->pParentVertex) );
1118  TrackingVertexRef temporaryRef( &dummyCollection, 0 );
1119  newTrackingParticle.setParentVertex( temporaryRef );
1120 
1121  // If a selector has been supplied apply it on the new TrackingParticle and return if it fails.
1122  if( pSelector )
1123  {
1124  if( !(*pSelector)( newTrackingParticle ) ) return; // Return if the TrackingParticle fails selection
1125  }
1126 
1127  // Add the ancestors first (if required) so that the collection has some kind of chronological
1128  // order. I don't know how important that is but other code might assume chronological order.
1129  // If adding ancestors, no selection is applied. Note that I've already checked that all
1130  // DecayChainTracks have a pParentVertex.
1131  if( addAncestors ) addTrack( pDecayChainTrack->pParentVertex->pParentTrack, NULL, pUnmergedOutput, pMergedOutput, objectFactory, addAncestors, tTopo );
1132 
1133  // If creation of the unmerged collection has been turned off in the config this pointer
1134  // will be null.
1135  if( pUnmergedOutput!=NULL ) addTrackAndParentVertex( pDecayChainTrack, newTrackingParticle, pUnmergedOutput );
1136 
1137  // If creation of the merged collection has been turned off in the config this pointer
1138  // will be null.
1139  if( pMergedOutput!=NULL )
1140  {
1141  ::DecayChainTrack* pBremParentChainTrack=pDecayChainTrack;
1142  while( pBremParentChainTrack->pMergedBremSource!=NULL ) pBremParentChainTrack=pBremParentChainTrack->pMergedBremSource;
1143 
1144  if( pBremParentChainTrack!=pDecayChainTrack )
1145  {
1146  TrackingParticle* pBremParentTrackingParticle=addTrackAndParentVertex( pBremParentChainTrack, newTrackingParticle, pMergedOutput );
1147  // The original particle in the bremsstrahlung decay chain has been
1148  // created (or retrieved if it already existed), now copy in the
1149  // extra information.
1150  // TODO - copy extra information.
1151 
1152  if( std::abs( newTrackingParticle.pdgId() ) == 22 )
1153  {
1154  // Photons are added as separate TrackingParticles, but with the production
1155  // vertex changed to be the production vertex of the original electron.
1156 
1157  // Set up a proxy, so that all requests for the parent TrackingVertex get
1158  // redirected to the brem parent's TrackingVertex.
1159  pMergedOutput->setProxy( pDecayChainTrack->pParentVertex, pBremParentChainTrack->pParentVertex );
1160 
1161  // Now that pMergedOutput thinks the parent vertex is the brem parent's vertex I can
1162  // call this and it will set the TrackingParticle parent vertex correctly to the brem
1163  // parent vertex.
1164  addTrackAndParentVertex( pDecayChainTrack, newTrackingParticle, pMergedOutput );
1165  }
1166  else if( std::abs( newTrackingParticle.pdgId() ) == 11 )
1167  {
1168  // Electrons have their G4 tracks and SimHits copied to the parent TrackingParticle.
1169  for( const auto& trackSegment : newTrackingParticle.g4Tracks() )
1170  {
1171  pBremParentTrackingParticle->addG4Track( trackSegment );
1172  }
1173 
1174  // Also copy the generator particle references
1175  for( const auto& genParticleRef : newTrackingParticle.genParticles() )
1176  {
1177  pBremParentTrackingParticle->addGenParticle( genParticleRef );
1178  }
1179 
1180  pBremParentTrackingParticle->setNumberOfHits(pBremParentTrackingParticle->numberOfHits()+newTrackingParticle.numberOfHits());
1181  pBremParentTrackingParticle->setNumberOfTrackerHits(pBremParentTrackingParticle->numberOfTrackerHits()+newTrackingParticle.numberOfTrackerHits());
1182  pBremParentTrackingParticle->setNumberOfTrackerLayers(pBremParentTrackingParticle->numberOfTrackerLayers()+newTrackingParticle.numberOfTrackerLayers());
1183 
1184  // Set a proxy in the output collection wrapper so that any attempt to get objects for
1185  // this DecayChainTrack again get redirected to the brem parent.
1186  pMergedOutput->setProxy( pDecayChainTrack, pBremParentChainTrack );
1187  }
1188  }
1189  else
1190  {
1191  // This is not the result of bremsstrahlung so add to the collection as normal.
1192  addTrackAndParentVertex( pDecayChainTrack, newTrackingParticle, pMergedOutput );
1193  }
1194  } // end of "if( pMergedOutput!=NULL )", i.e. end of "if bremsstrahlung merging is turned on"
1195 
1196  } // end of addTrack function
1197 
1198 } // end of the unnamed namespace
1199 
1200 // Register with the framework
helper::MatcherGetRef< C >::ref_type getRef(const Handle< C > &c, size_t k)
Definition: getRef.h:28
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
EncodedEventId eventId() const
Definition: CoreSimTrack.h:31
int event() const
get the contents of the subdetector field (should be protected?)
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
virtual void initializeEvent(const edm::Event &event, const edm::EventSetup &setup)
std::vector< TrackingParticle > TrackingParticleCollection
ROOT::Math::Plane3D::Vector Vector
Definition: EcalHitMaker.cc:29
int pdgId() const
PDG ID.
const bool addAncestors_
Whether or not to add the full parentage of any TrackingParticle that is inserted in the collection...
int init
Definition: HydjetWrapper.h:62
bool exists(std::string const &parameterName) const
checks if a parameter exists
void addDaughterTrack(const TrackingParticleRef &)
const double vertexDistanceCut_
maximum distance for HepMC::GenVertex to be added to SimVertex
#define NULL
Definition: scimark2.h:8
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
void addParentTrack(const TrackingParticleRef &)
float charge() const
charge
Definition: CoreSimTrack.cc:18
tuple vertexCollection
bool signalOnly_
Uses the same config as selector_, but can be used to drop out early since selector_ requires the Tra...
void fillSimHits(std::vector< const PSimHit * > &returnValue, const T &event, const edm::EventSetup &setup)
Fills the supplied vector with pointers to the SimHits, checking for bad modules if required...
const unsigned int maximumSubsequentBunchCrossing_
void setParentVertex(const TrackingVertexRef &ref)
void addGenParticle(const reco::GenParticleRef &ref)
SingleObjectSelector< TrackingParticleCollection,::TrackingParticleSelector > TrackingParticleSelector
const reco::GenParticleRefVector & genParticles() const
void accumulateEvent(const T &event, const edm::EventSetup &setup, const edm::Handle< edm::HepMCProduct > &hepMCproduct)
Both forms of accumulate() delegate to this templated method.
double y() const
Same as rapidity().
int bunchCrossing() const
get the detector field from this detid
void addG4Track(const SimTrack &t)
int genpartIndex() const
index of the corresponding Generator particle in the Event container (-1 if no Genpart) ...
Definition: SimTrack.h:33
int numberOfTrackerLayers() const
The number of tracker layers with a hit.
virtual void finalizeEvent(edm::Event &event, const edm::EventSetup &setup)
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
const std::string messageCategory_
The message category used to send messages to MessageLogger.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool noVertex() const
Definition: SimTrack.h:30
tuple genEvent
Definition: MCTruth.py:33
const math::XYZTLorentzVectorD & position() const
Definition: CoreSimVertex.h:26
std::vector< std::string > getParameterNames() const
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
bool allowDifferentProcessTypeForDifferentDetectors_
When counting hits, allows hits in different detectors to have a different process type...
std::auto_ptr< TrackingParticleCollection > pTrackingParticles
virtual void accumulate(const edm::Event &event, const edm::EventSetup &setup)
tuple tags
Definition: o2o.py:248
edm::Ref< TrackingVertexCollection > TrackingVertexRef
Definition: DetId.h:18
edm::Ref< edm::HepMCProduct, HepMC::GenVertex > GenVertexRef
unsigned int trackId() const
Definition: CoreSimTrack.h:34
ParameterSet const & getParameterSet(std::string const &) const
std::auto_ptr< TrackingVertexCollection > pTrackingVertices
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
void addDecayVertex(const TrackingVertexRef &ref)
std::vector< TrackingVertex > TrackingVertexCollection
const T & get() const
Definition: EventSetup.h:55
const std::vector< SimTrack > & g4Tracks() const
T const * product() const
Definition: ESHandle.h:86
tuple simHits
Definition: trackerHits.py:16
int numberOfHits() const
Gives the total number of hits, including muon hits. Hits on overlaps in the same layer count separat...
const unsigned int maximumPreviousBunchCrossing_
unsigned int layer(const DetId &id) const
void setNumberOfTrackerHits(int numberOfTrackerHits)
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:25
const TrackingVertexRef & parentVertex() const
EncodedEventId eventId() const
Definition: CoreSimVertex.h:30
#define DEFINE_DIGI_ACCUMULATOR(type)
static int position[264][3]
Definition: ReadPGInfo.cc:509
bool chargedOnly_
Uses the same config as selector_, but can be used to drop out early since selector_ requires the Tra...
const bool createUnmergedCollection_
If bremsstrahlung merging, whether to also add the unmerged collection to the event or not...
TrackingParticleSelector selector_
Monte Carlo truth information used for tracking validation.
tuple cout
Definition: gather_cfg.py:121
int numberOfTrackerHits() const
The number of hits in the tracker. Hits on overlaps in the same layer count separately.
void setNumberOfTrackerLayers(const int numberOfTrackerLayers)
TrackingTruthAccumulator(const edm::ParameterSet &config, edm::one::EDProducerBase &mixMod, edm::ConsumesCollector &iC)
Replacement for TrackingTruthProducer in the new pileup mixing setup.
Detector det() const
get the detector field from this detid
Definition: DetId.h:35
long double T
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
tuple size
Write out results.
edm::Ref< TrackingParticleCollection > TrackingParticleRef
edm::InputTag hepMCproductLabel_
Needed to add HepMC::GenVertex to SimVertex.
math::PtEtaPhiELorentzVectorF LorentzVector
void setNumberOfHits(int numberOfHits)
std::vector< edm::InputTag > collectionTags_