CMS 3D CMS Logo

PreMixingTrackingParticleWorker.cc
Go to the documentation of this file.
7 
11 
14 
16 public:
20  ~PreMixingTrackingParticleWorker() override = default;
21 
22  void initializeEvent(edm::Event const &iEvent, edm::EventSetup const &iSetup) override;
23  void addSignals(edm::Event const &iEvent, edm::EventSetup const &iSetup) override;
24  void addPileups(PileUpEventPrincipal const &pep, edm::EventSetup const &iSetup) override;
25  void put(edm::Event &iEvent,
26  edm::EventSetup const &iSetup,
27  std::vector<PileupSummaryInfo> const &ps,
28  int bunchSpacing) override;
29 
30 private:
31  void add(const std::vector<TrackingParticle> &particles, const std::vector<TrackingVertex> &vertices);
32 
35 
36  edm::InputTag TrackingParticlePileInputTag_; // InputTag for pileup tracks
37 
38  std::string TrackingParticleCollectionDM_; // secondary name to be given to
39  // new TrackingParticle
40 
41  std::unique_ptr<std::vector<TrackingParticle>> NewTrackList_;
42  std::unique_ptr<std::vector<TrackingVertex>> NewVertexList_;
45 };
46 
50  : TrackSigToken_(iC.consumes<std::vector<TrackingParticle>>(ps.getParameter<edm::InputTag>("labelSig"))),
51  VtxSigToken_(iC.consumes<std::vector<TrackingVertex>>(ps.getParameter<edm::InputTag>("labelSig"))),
52  TrackingParticlePileInputTag_(ps.getParameter<edm::InputTag>("pileInputTag")),
53  TrackingParticleCollectionDM_(ps.getParameter<std::string>("collectionDM")) {
54  producer.produces<std::vector<TrackingParticle>>(TrackingParticleCollectionDM_);
55  producer.produces<std::vector<TrackingVertex>>(TrackingParticleCollectionDM_);
56 }
57 
59  NewTrackList_ = std::make_unique<std::vector<TrackingParticle>>();
60  NewVertexList_ = std::make_unique<std::vector<TrackingVertex>>();
61 
62  // need RefProds in order to re-key the particle<->vertex refs
63  // TODO: try to remove const_cast, requires making Event non-const in
64  // BMixingModule::initializeEvent
66  const_cast<edm::Event &>(iEvent).getRefBeforePut<std::vector<TrackingParticle>>(TrackingParticleCollectionDM_);
68  const_cast<edm::Event &>(iEvent).getRefBeforePut<std::vector<TrackingVertex>>(TrackingParticleCollectionDM_);
69 }
70 
73  iEvent.getByToken(TrackSigToken_, tracks);
74 
76  iEvent.getByToken(VtxSigToken_, vtxs);
77 
78  if (tracks.isValid() && vtxs.isValid()) {
79  add(*tracks, *vtxs);
80  }
81 }
82 
84  LogDebug("PreMixingTrackingParticleWorker") << "\n===============> adding pileups from event "
85  << pep.principal().id() << " for bunchcrossing " << pep.bunchCrossing();
86 
88  pep.getByLabel(TrackingParticlePileInputTag_, inputHandle);
89 
91  pep.getByLabel(TrackingParticlePileInputTag_, inputVHandle);
92 
93  if (inputHandle.isValid() && inputVHandle.isValid()) {
94  add(*inputHandle, *inputVHandle);
95  }
96 }
97 
98 void PreMixingTrackingParticleWorker::add(const std::vector<TrackingParticle> &particles,
99  const std::vector<TrackingVertex> &vertices) {
100  const size_t StartingIndexV = NewVertexList_->size();
101  const size_t StartingIndexT = NewTrackList_->size();
102 
103  // grab Vertices, store copy, preserving indices. Easier to loop over
104  // vertices first - fewer links
105  for (const auto &vtx : vertices) {
106  NewVertexList_->push_back(vtx);
107  }
108 
109  // grab tracks, store copy
110  for (const auto &track : particles) {
111  const auto &oldRef = track.parentVertex();
112  auto newRef = TrackingVertexRef(VertexListRef_, oldRef.index() + StartingIndexV);
113  NewTrackList_->push_back(track);
114 
115  auto &Ntrack = NewTrackList_->back(); // modify copy
116 
117  Ntrack.setParentVertex(newRef);
118  Ntrack.clearDecayVertices();
119 
120  // next, loop over daughter vertices, same strategy
121  for (auto const &vertexRef : track.decayVertices()) {
122  auto newRef = TrackingVertexRef(VertexListRef_, vertexRef.index() + StartingIndexV);
123  Ntrack.addDecayVertex(newRef);
124  }
125  }
126 
127  // Now that tracks are handled, go back and put correct Refs in vertices
128  // Operate only on the added pileup vertices, and leave the already-existing
129  // vertices untouched
130  std::vector<decltype(TrackingParticleRef().index())> sourceTrackIndices;
131  std::vector<decltype(TrackingParticleRef().index())> daughterTrackIndices;
132  for (size_t iVertex = StartingIndexV; iVertex != NewVertexList_->size(); ++iVertex) {
133  auto &vertex = (*NewVertexList_)[iVertex];
134 
135  // Need to copy the indices before clearing the vectors
136  sourceTrackIndices.reserve(vertex.sourceTracks().size());
137  daughterTrackIndices.reserve(vertex.daughterTracks().size());
138  for (auto const &ref : vertex.sourceTracks())
139  sourceTrackIndices.push_back(ref.index());
140  for (auto const &ref : vertex.daughterTracks())
141  daughterTrackIndices.push_back(ref.index());
142 
143  vertex.clearParentTracks();
144  vertex.clearDaughterTracks();
145 
146  for (auto index : sourceTrackIndices) {
147  auto newRef = TrackingParticleRef(TrackListRef_, index + StartingIndexT);
148  vertex.addParentTrack(newRef);
149  }
150 
151  // next, loop over daughter tracks, same strategy
152  for (auto index : daughterTrackIndices) {
153  auto newRef = TrackingParticleRef(TrackListRef_, index + StartingIndexT);
154  vertex.addDaughterTrack(newRef);
155  }
156 
157  sourceTrackIndices.clear();
158  daughterTrackIndices.clear();
159  }
160 }
161 
163  edm::EventSetup const &iSetup,
164  std::vector<PileupSummaryInfo> const &ps,
165  int bunchSpacing) {
166  edm::LogInfo("PreMixingTrackingParticleWorker") << "total # Merged Tracks: " << NewTrackList_->size();
169 }
170 
#define LogDebug(id)
BranchAliasSetterT< ProductType > produces()
declare what type of product will make and with which optional label
~PreMixingTrackingParticleWorker() override=default
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
edm::EDGetTokenT< std::vector< TrackingParticle > > TrackSigToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
void add(const std::vector< TrackingParticle > &particles, const std::vector< TrackingVertex > &vertices)
EventID const & id() const
std::unique_ptr< std::vector< TrackingParticle > > NewTrackList_
std::unique_ptr< std::vector< TrackingVertex > > NewVertexList_
edm::EventPrincipal const & principal()
int iEvent
Definition: GenABIO.cc:224
PreMixingTrackingParticleWorker(const edm::ParameterSet &ps, edm::ProducerBase &producer, edm::ConsumesCollector &&iC)
bool isValid() const
Definition: HandleBase.h:74
edm::Ref< TrackingVertexCollection > TrackingVertexRef
void addPileups(PileUpEventPrincipal const &pep, edm::EventSetup const &iSetup) override
HLT enums.
edm::EDGetTokenT< std::vector< TrackingVertex > > VtxSigToken_
bool getByLabel(edm::InputTag const &tag, edm::Handle< T > &result) const
void put(edm::Event &iEvent, edm::EventSetup const &iSetup, std::vector< PileupSummaryInfo > const &ps, int bunchSpacing) override
Monte Carlo truth information used for tracking validation.
void initializeEvent(edm::Event const &iEvent, edm::EventSetup const &iSetup) override
void addSignals(edm::Event const &iEvent, edm::EventSetup const &iSetup) override
#define DEFINE_PREMIXING_WORKER(TYPE)
edm::Ref< TrackingParticleCollection > TrackingParticleRef
def move(src, dest)
Definition: eostools.py:511