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