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