CMS 3D CMS Logo

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