CMS 3D CMS Logo

DeDxEstimatorRekeyer.cc
Go to the documentation of this file.
8 
9 //
10 // class declaration
11 //
12 
14 public:
15  explicit DeDxEstimatorRekeyer(const edm::ParameterSet&);
16  ~DeDxEstimatorRekeyer() override{};
17  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
18 
19 private:
20  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
21 
22  template <typename T>
23  std::map<std::string, edm::EDGetTokenT<T>> getTokens(const std::vector<edm::InputTag>& tags) {
24  std::map<std::string, edm::EDGetTokenT<T>> tokens;
25  for (const auto& tag : tags)
26  tokens.emplace(tag.label(), consumes<T>(tag));
27  return tokens;
28  };
29 
30  // ----------member data ---------------------------
34  const std::map<std::string, edm::EDGetTokenT<edm::ValueMap<reco::DeDxData>>> dedxEstimatorsTokens_;
35  const std::map<std::string, edm::EDGetTokenT<pat::PackedCandidateCollection>> packedCandidatesTokens_;
36  const std::map<std::string, edm::EDGetTokenT<edm::Association<pat::PackedCandidateCollection>>> trk2pcTokens_;
37 };
38 
41  desc.add<edm::InputTag>("tracks", {"generalTracks"});
42  desc.add<edm::InputTag>("dedxHits", {"dedxHitInfo"});
43  desc.add<edm::InputTag>("dedxMomentum", {"dedxHitInfo:momentumAtHit"});
44  desc.add<std::vector<edm::InputTag>>(
45  "packedCandidates",
46  {edm::InputTag("packedPFCandidates"), edm::InputTag("lostTracks"), edm::InputTag("lostTracks:eleTracks")});
47  desc.add<std::vector<edm::InputTag>>("dedxEstimators",
48  {edm::InputTag("dedxHarmonic2"), edm::InputTag("dedxPixelHarmonic2")});
49  descriptions.addWithDefaultLabel(desc);
50 }
51 
53  : tracksToken_(consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("tracks"))),
54  dedxHitAssToken_(consumes<reco::DeDxHitInfoAss>(iConfig.getParameter<edm::InputTag>("dedxHits"))),
55  dedxHitMomToken_(
56  consumes<edm::ValueMap<std::vector<float>>>(iConfig.getParameter<edm::InputTag>("dedxMomentum"))),
57  dedxEstimatorsTokens_(
58  getTokens<edm::ValueMap<reco::DeDxData>>(iConfig.getParameter<std::vector<edm::InputTag>>("dedxEstimators"))),
59  packedCandidatesTokens_(getTokens<pat::PackedCandidateCollection>(
60  iConfig.getParameter<std::vector<edm::InputTag>>("packedCandidates"))),
61  trk2pcTokens_(getTokens<edm::Association<pat::PackedCandidateCollection>>(
62  iConfig.getParameter<std::vector<edm::InputTag>>("packedCandidates"))) {
63  for (const auto& d : dedxEstimatorsTokens_)
64  produces<edm::ValueMap<reco::DeDxData>>(d.first);
65  produces<reco::DeDxHitInfoCollection>();
66  produces<reco::DeDxHitInfoAss>();
67  produces<edm::ValueMap<std::vector<float>>>("momentumAtHit");
68 }
69 
71  // Get input collections
72  const auto& tracks = iEvent.getHandle(tracksToken_);
73  std::vector<reco::TrackRef> trackRefs;
74  trackRefs.reserve(tracks->size());
75  for (auto track = tracks->begin(); track != tracks->end(); track++)
76  trackRefs.emplace_back(tracks, track - tracks->begin());
77 
78  typedef std::map<pat::PackedCandidateRef, reco::TrackRef> PCTrkMap;
79  std::vector<std::pair<edm::Handle<pat::PackedCandidateCollection>, PCTrkMap>> pcTrkMap;
80  pcTrkMap.reserve(packedCandidatesTokens_.size());
81  for (const auto& p : packedCandidatesTokens_) {
82  PCTrkMap map;
83  const auto& trk2pc = iEvent.get(trk2pcTokens_.at(p.first));
84  for (const auto& track : trackRefs) {
85  const auto& pc = trk2pc[track];
86  if (pc.isNonnull())
87  map.emplace(pc, track);
88  }
89  pcTrkMap.emplace_back(iEvent.getHandle(p.second), map);
90  }
91 
92  // Rekey dEdx estimators
93  for (const auto& d : dedxEstimatorsTokens_) {
94  const auto& dedxEstimators = iEvent.get(d.second);
95  auto trackDeDxValueMap = std::make_unique<edm::ValueMap<reco::DeDxData>>();
97  // Loop over packed candidates
98  for (const auto& h : pcTrkMap) {
99  std::vector<reco::DeDxData> dedxEstimate(h.first->size());
100  for (const auto& p : h.second)
101  dedxEstimate[p.first.key()] = dedxEstimators[p.second];
102  filler.insert(h.first, dedxEstimate.begin(), dedxEstimate.end());
103  }
104  // Fill the value map and put it into the event
105  filler.fill();
106  iEvent.put(std::move(trackDeDxValueMap), d.first);
107  }
108 
109  // Rekey dEdx hit info
110  const auto& dedxHitMom = iEvent.get(dedxHitMomToken_);
111  const auto& dedxHitAss = iEvent.get(dedxHitAssToken_);
112  const auto& dedxHitInfoHandle = iEvent.getRefBeforePut<reco::DeDxHitInfoCollection>();
113  auto dedxHitInfoAssociation = std::make_unique<reco::DeDxHitInfoAss>(dedxHitInfoHandle);
114  reco::DeDxHitInfoAss::Filler filler(*dedxHitInfoAssociation);
115  auto resultdedxHitColl = std::make_unique<reco::DeDxHitInfoCollection>();
116  resultdedxHitColl->reserve(pcTrkMap.size() > 0 ? pcTrkMap.size() * pcTrkMap[0].second.size() : 0);
117  std::vector<std::vector<float>> momenta;
118  momenta.reserve(resultdedxHitColl->capacity());
119  // Loop over packed candidates
120  for (const auto& h : pcTrkMap) {
121  std::vector<int> indices(h.first->size(), -1);
122  for (const auto& p : h.second) {
123  const auto& dedxHit = dedxHitAss[p.second];
124  if (dedxHit.isNull())
125  continue;
126  indices[p.first.key()] = resultdedxHitColl->size();
127  resultdedxHitColl->emplace_back(*dedxHit);
128  momenta.emplace_back(dedxHitMom[dedxHit]);
129  }
130  filler.insert(h.first, indices.begin(), indices.end());
131  }
132  const auto& dedxHitCollHandle = iEvent.put(std::move(resultdedxHitColl));
133  // Fill the association map and put it into the event
134  filler.fill();
135  iEvent.put(std::move(dedxHitInfoAssociation));
136  // Fill the value map and put it into the event
137  auto dedxMomenta = std::make_unique<edm::ValueMap<std::vector<float>>>();
138  edm::ValueMap<std::vector<float>>::Filler mfiller(*dedxMomenta);
139  mfiller.insert(dedxHitCollHandle, momenta.begin(), momenta.end());
140  mfiller.fill();
141  iEvent.put(std::move(dedxMomenta), "momentumAtHit");
142 }
143 
144 //define this as a plug-in
const edm::EDGetTokenT< edm::ValueMap< std::vector< float > > > dedxHitMomToken_
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
const std::map< std::string, edm::EDGetTokenT< edm::Association< pat::PackedCandidateCollection > > > trk2pcTokens_
std::vector< pat::PackedCandidate > PackedCandidateCollection
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
DeDxEstimatorRekeyer(const edm::ParameterSet &)
Definition: HeavyIon.h:7
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const edm::EDGetTokenT< reco::DeDxHitInfoAss > dedxHitAssToken_
int iEvent
Definition: GenABIO.cc:224
edm::Association< DeDxHitInfoCollection > DeDxHitInfoAss
Definition: DeDxHitInfo.h:123
const std::map< std::string, edm::EDGetTokenT< pat::PackedCandidateCollection > > packedCandidatesTokens_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const edm::EDGetTokenT< reco::TrackCollection > tracksToken_
std::map< std::string, edm::EDGetTokenT< T > > getTokens(const std::vector< edm::InputTag > &tags)
d
Definition: ztail.py:151
std::vector< DeDxHitInfo > DeDxHitInfoCollection
Definition: DeDxHitInfo.h:119
const std::map< std::string, edm::EDGetTokenT< edm::ValueMap< reco::DeDxData > > > dedxEstimatorsTokens_
fixed size matrix
HLT enums.
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
def move(src, dest)
Definition: eostools.py:511