CMS 3D CMS Logo

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