CMS 3D CMS Logo

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