CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PreMixingTrackingParticleWorker.cc
Go to the documentation of this file.
7 
11 
14 
16 public:
18  ~PreMixingTrackingParticleWorker() override = default;
19 
20  void initializeEvent(edm::Event const &iEvent, edm::EventSetup const &iSetup) override;
21  void addSignals(edm::Event const &iEvent, edm::EventSetup const &iSetup) override;
22  void addPileups(PileUpEventPrincipal const &pep, edm::EventSetup const &iSetup) override;
23  void put(edm::Event &iEvent,
24  edm::EventSetup const &iSetup,
25  std::vector<PileupSummaryInfo> const &ps,
26  int bunchSpacing) override;
27 
28 private:
29  void add(const std::vector<TrackingParticle> &particles, const std::vector<TrackingVertex> &vertices);
30 
33 
34  edm::InputTag TrackingParticlePileInputTag_; // InputTag for pileup tracks
35 
36  std::string TrackingParticleCollectionDM_; // secondary name to be given to
37  // new TrackingParticle
38 
39  std::unique_ptr<std::vector<TrackingParticle>> NewTrackList_;
40  std::unique_ptr<std::vector<TrackingVertex>> NewVertexList_;
43 };
44 
46  edm::ProducesCollector producesCollector,
48  : TrackSigToken_(iC.consumes<std::vector<TrackingParticle>>(ps.getParameter<edm::InputTag>("labelSig"))),
49  VtxSigToken_(iC.consumes<std::vector<TrackingVertex>>(ps.getParameter<edm::InputTag>("labelSig"))),
50  TrackingParticlePileInputTag_(ps.getParameter<edm::InputTag>("pileInputTag")),
51  TrackingParticleCollectionDM_(ps.getParameter<std::string>("collectionDM")) {
52  producesCollector.produces<std::vector<TrackingParticle>>(TrackingParticleCollectionDM_);
53  producesCollector.produces<std::vector<TrackingVertex>>(TrackingParticleCollectionDM_);
54 }
55 
57  NewTrackList_ = std::make_unique<std::vector<TrackingParticle>>();
58  NewVertexList_ = std::make_unique<std::vector<TrackingVertex>>();
59 
60  // need RefProds in order to re-key the particle<->vertex refs
61  // TODO: try to remove const_cast, requires making Event non-const in
62  // BMixingModule::initializeEvent
64  const_cast<edm::Event &>(iEvent).getRefBeforePut<std::vector<TrackingParticle>>(TrackingParticleCollectionDM_);
66  const_cast<edm::Event &>(iEvent).getRefBeforePut<std::vector<TrackingVertex>>(TrackingParticleCollectionDM_);
67 }
68 
71  iEvent.getByToken(TrackSigToken_, tracks);
72 
74  iEvent.getByToken(VtxSigToken_, vtxs);
75 
76  if (tracks.isValid() && vtxs.isValid()) {
77  add(*tracks, *vtxs);
78  }
79 }
80 
82  LogDebug("PreMixingTrackingParticleWorker") << "\n===============> adding pileups from event "
83  << pep.principal().id() << " for bunchcrossing " << pep.bunchCrossing();
84 
86  pep.getByLabel(TrackingParticlePileInputTag_, inputHandle);
87 
89  pep.getByLabel(TrackingParticlePileInputTag_, inputVHandle);
90 
91  if (inputHandle.isValid() && inputVHandle.isValid()) {
92  add(*inputHandle, *inputVHandle);
93  }
94 }
95 
96 void PreMixingTrackingParticleWorker::add(const std::vector<TrackingParticle> &particles,
97  const std::vector<TrackingVertex> &vertices) {
98  const size_t StartingIndexV = NewVertexList_->size();
99  const size_t StartingIndexT = NewTrackList_->size();
100 
101  // grab Vertices, store copy, preserving indices. Easier to loop over
102  // vertices first - fewer links
103  for (const auto &vtx : vertices) {
104  NewVertexList_->push_back(vtx);
105  }
106 
107  // grab tracks, store copy
108  for (const auto &track : particles) {
109  const auto &oldRef = track.parentVertex();
110  auto newRef = TrackingVertexRef(VertexListRef_, oldRef.index() + StartingIndexV);
111  NewTrackList_->push_back(track);
112 
113  auto &Ntrack = NewTrackList_->back(); // modify copy
114 
115  Ntrack.setParentVertex(newRef);
116  Ntrack.clearDecayVertices();
117 
118  // next, loop over daughter vertices, same strategy
119  for (auto const &vertexRef : track.decayVertices()) {
120  auto newRef = TrackingVertexRef(VertexListRef_, vertexRef.index() + StartingIndexV);
121  Ntrack.addDecayVertex(newRef);
122  }
123  }
124 
125  // Now that tracks are handled, go back and put correct Refs in vertices
126  // Operate only on the added pileup vertices, and leave the already-existing
127  // vertices untouched
128  std::vector<decltype(TrackingParticleRef().index())> sourceTrackIndices;
129  std::vector<decltype(TrackingParticleRef().index())> daughterTrackIndices;
130  for (size_t iVertex = StartingIndexV; iVertex != NewVertexList_->size(); ++iVertex) {
131  auto &vertex = (*NewVertexList_)[iVertex];
132 
133  // Need to copy the indices before clearing the vectors
134  sourceTrackIndices.reserve(vertex.sourceTracks().size());
135  daughterTrackIndices.reserve(vertex.daughterTracks().size());
136  for (auto const &ref : vertex.sourceTracks())
137  sourceTrackIndices.push_back(ref.index());
138  for (auto const &ref : vertex.daughterTracks())
139  daughterTrackIndices.push_back(ref.index());
140 
141  vertex.clearParentTracks();
142  vertex.clearDaughterTracks();
143 
144  for (auto index : sourceTrackIndices) {
145  auto newRef = TrackingParticleRef(TrackListRef_, index + StartingIndexT);
146  vertex.addParentTrack(newRef);
147  }
148 
149  // next, loop over daughter tracks, same strategy
150  for (auto index : daughterTrackIndices) {
151  auto newRef = TrackingParticleRef(TrackListRef_, index + StartingIndexT);
152  vertex.addDaughterTrack(newRef);
153  }
154 
155  sourceTrackIndices.clear();
156  daughterTrackIndices.clear();
157  }
158 }
159 
161  edm::EventSetup const &iSetup,
162  std::vector<PileupSummaryInfo> const &ps,
163  int bunchSpacing) {
164  edm::LogInfo("PreMixingTrackingParticleWorker") << "total # Merged Tracks: " << NewTrackList_->size();
167 }
168 
~PreMixingTrackingParticleWorker() override=default
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
edm::EDGetTokenT< std::vector< TrackingParticle > > TrackSigToken_
ProductRegistryHelper::BranchAliasSetterT< ProductType > produces()
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
void add(const std::vector< TrackingParticle > &particles, const std::vector< TrackingVertex > &vertices)
EventID const & id() const
auto const & tracks
cannot be loose
std::unique_ptr< std::vector< TrackingParticle > > NewTrackList_
std::unique_ptr< std::vector< TrackingVertex > > NewVertexList_
edm::EventPrincipal const & principal()
int iEvent
Definition: GenABIO.cc:224
def move
Definition: eostools.py:511
bool isValid() const
Definition: HandleBase.h:70
edm::Ref< TrackingVertexCollection > TrackingVertexRef
Log< level::Info, false > LogInfo
PreMixingTrackingParticleWorker(const edm::ParameterSet &ps, edm::ProducesCollector, edm::ConsumesCollector &&iC)
void addPileups(PileUpEventPrincipal const &pep, edm::EventSetup const &iSetup) override
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
#define LogDebug(id)