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