CMS 3D CMS Logo

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, edm::EventSetup const& iSetup, std::vector<PileupSummaryInfo> const& ps, int bunchSpacing) override;
24 
25 private:
26  void add(const std::vector<TrackingParticle>& particles, const std::vector<TrackingVertex>& vertices);
27 
29  edm::EDGetTokenT<std::vector<TrackingVertex> > VtxSigToken_; // Token to retrieve information
30 
31  edm::InputTag TrackingParticlePileInputTag_ ; // InputTag for pileup tracks
32 
33  std::string TrackingParticleCollectionDM_; // secondary name to be given to new TrackingParticle
34 
35 
36  std::unique_ptr<std::vector<TrackingParticle>> NewTrackList_;
37  std::unique_ptr<std::vector<TrackingVertex>> NewVertexList_;
40 };
41 
43  TrackSigToken_(iC.consumes<std::vector<TrackingParticle> >(ps.getParameter<edm::InputTag>("labelSig"))),
44  VtxSigToken_ (iC.consumes<std::vector<TrackingVertex> >(ps.getParameter<edm::InputTag>("labelSig"))),
45  TrackingParticlePileInputTag_(ps.getParameter<edm::InputTag>("pileInputTag")),
46  TrackingParticleCollectionDM_(ps.getParameter<std::string>("collectionDM"))
47 {
48  producer.produces< std::vector<TrackingParticle> >(TrackingParticleCollectionDM_);
49  producer.produces< std::vector<TrackingVertex> >(TrackingParticleCollectionDM_);
50 }
51 
52 
54  NewTrackList_ = std::make_unique<std::vector<TrackingParticle>>();
55  NewVertexList_ = std::make_unique<std::vector<TrackingVertex>>();
56 
57  // need RefProds in order to re-key the particle<->vertex refs
58  // TODO: try to remove const_cast, requires making Event non-const in BMixingModule::initializeEvent
59  TrackListRef_ = const_cast<edm::Event&>(iEvent).getRefBeforePut< std::vector<TrackingParticle> >(TrackingParticleCollectionDM_);
60  VertexListRef_ = const_cast<edm::Event&>(iEvent).getRefBeforePut< std::vector<TrackingVertex> >(TrackingParticleCollectionDM_);
61 }
62 
65  iEvent.getByToken(TrackSigToken_, tracks);
66 
68  iEvent.getByToken(VtxSigToken_, vtxs);
69 
70  if(tracks.isValid() && vtxs.isValid()) {
71  add(*tracks, *vtxs);
72  }
73 }
74 
76  LogDebug("PreMixingTrackingParticleWorker") <<"\n===============> adding pileups from event "<<pep.principal().id()<<" for bunchcrossing "<<pep.bunchCrossing();
77 
79  pep.getByLabel(TrackingParticlePileInputTag_, inputHandle);
80 
82  pep.getByLabel(TrackingParticlePileInputTag_, inputVHandle);
83 
84  if(inputHandle.isValid() && inputVHandle.isValid()) {
85  add(*inputHandle, *inputVHandle);
86  }
87 }
88 
89 void PreMixingTrackingParticleWorker::add(const std::vector<TrackingParticle>& particles, const std::vector<TrackingVertex>& vertices) {
90  const size_t StartingIndexV = NewVertexList_->size();
91  const size_t StartingIndexT = NewTrackList_->size();
92 
93  // grab Vertices, store copy, preserving indices. Easier to loop over vertices first - fewer links
94  for(const auto& vtx: vertices) {
95  NewVertexList_->push_back(vtx);
96  }
97 
98  // grab tracks, store copy
99  for(const auto& track: particles) {
100  const auto& oldRef=track.parentVertex();
101  auto newRef=TrackingVertexRef( VertexListRef_, oldRef.index()+StartingIndexV );
102  NewTrackList_->push_back(track);
103 
104  auto & Ntrack = NewTrackList_->back(); //modify copy
105 
106  Ntrack.setParentVertex( newRef );
107  Ntrack.clearDecayVertices();
108 
109  // next, loop over daughter vertices, same strategy
110  for( auto const& vertexRef : track.decayVertices() ) {
111  auto newRef=TrackingVertexRef( VertexListRef_, vertexRef.index()+StartingIndexV );
112  Ntrack.addDecayVertex(newRef);
113  }
114  }
115 
116  // Now that tracks are handled, go back and put correct Refs in vertices
117  // Operate only on the added pileup vertices, and leave the already-existing vertices untouched
118  std::vector<decltype(TrackingParticleRef().index())> sourceTrackIndices;
119  std::vector<decltype(TrackingParticleRef().index())> daughterTrackIndices;
120  for(size_t iVertex = StartingIndexV; iVertex != NewVertexList_->size(); ++iVertex) {
121  auto& vertex = (*NewVertexList_)[iVertex];
122 
123  // Need to copy the indices before clearing the vectors
124  sourceTrackIndices.reserve(vertex.sourceTracks().size());
125  daughterTrackIndices.reserve(vertex.daughterTracks().size());
126  for(auto const& ref: vertex.sourceTracks()) sourceTrackIndices.push_back(ref.index());
127  for(auto const& ref: vertex.daughterTracks()) daughterTrackIndices.push_back(ref.index());
128 
129  vertex.clearParentTracks();
130  vertex.clearDaughterTracks();
131 
132  for( auto index : sourceTrackIndices ) {
133  auto newRef=TrackingParticleRef( TrackListRef_, index+StartingIndexT );
134  vertex.addParentTrack(newRef);
135  }
136 
137  // next, loop over daughter tracks, same strategy
138  for( auto index : daughterTrackIndices ) {
139  auto newRef=TrackingParticleRef( TrackListRef_, index+StartingIndexT );
140  vertex.addDaughterTrack(newRef);
141  }
142 
143  sourceTrackIndices.clear();
144  daughterTrackIndices.clear();
145  }
146 }
147 
148 void PreMixingTrackingParticleWorker::put(edm::Event& iEvent, edm::EventSetup const& iSetup, std::vector<PileupSummaryInfo> const& ps, int bunchSpacing) {
149  edm::LogInfo("PreMixingTrackingParticleWorker") << "total # Merged Tracks: " << NewTrackList_->size() ;
152 }
153 
#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:137
edm::EDGetTokenT< std::vector< TrackingParticle > > TrackSigToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
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:230
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:510