CMS 3D CMS Logo

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