CMS 3D CMS Logo

ProtonProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: PhysicsTools/NanoAOD
4 // Class: ProtonProducer
5 //
10 //
11 // Original Author: Justin Williams
12 // Created: 04 Jul 2019 15:27:53 GMT
13 //
14 //
15 
16 // system include files
17 #include <memory>
18 #include <map>
19 #include <string>
20 #include <vector>
21 #include <iostream>
22 
23 // user include files
26 
29 
32 
35 
41 
43 public:
46  mayConsume<reco::ForwardProtonCollection>(ps.getParameter<edm::InputTag>("tagRecoProtonsMulti"))),
48  mayConsume<reco::ForwardProtonCollection>(ps.getParameter<edm::InputTag>("tagRecoProtonsSingle"))),
49  tokenTracksLite_(mayConsume<std::vector<CTPPSLocalTrackLite>>(ps.getParameter<edm::InputTag>("tagTrackLite"))),
50  storeSingleRPProtons_(ps.getParameter<bool>("storeSingleRPProtons")) {
51  produces<edm::ValueMap<int>>("arm");
52  produces<nanoaod::FlatTable>("ppsTrackTable");
54  produces<edm::ValueMap<int>>("protonRPId");
55  }
56  ~ProtonProducer() override {}
57 
58  // ------------ method called to produce the data ------------
59  void produce(edm::StreamID id, edm::Event &iEvent, const edm::EventSetup &iSetup) const override {
60  // Get Forward Proton handle
62  iEvent.getByToken(tokenRecoProtonsMultiRP_, hRecoProtonsMultiRP);
63 
64  edm::Handle<reco::ForwardProtonCollection> hRecoProtonsSingleRP;
65  iEvent.getByToken(tokenRecoProtonsSingleRP_, hRecoProtonsSingleRP);
66 
67  // Get PPS Local Track handle
69  iEvent.getByToken(tokenTracksLite_, ppsTracksLite);
70 
71  // book output variables for protons
72  std::vector<int> multiRP_arm;
73 
74  // book output variables for tracks
75  std::vector<float> trackX, trackY, trackTime, trackTimeUnc;
76  std::vector<int> multiRPProtonIdx, decRPId, rpType;
77  std::vector<int> singleRPProtonIdx, singleRP_protonRPId;
78 
80  // process single-RP protons
81  {
82  const auto &num_proton = hRecoProtonsSingleRP->size();
83  singleRP_protonRPId.reserve(num_proton);
84 
85  for (const auto &proton : *hRecoProtonsSingleRP) {
86  CTPPSDetId rpId((*proton.contributingLocalTracks().begin())->getRPId());
87  unsigned int rpDecId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
88  singleRP_protonRPId.push_back(rpDecId);
89  }
90  }
91 
92  // process multi-RP protons
93  {
94  const auto &num_proton = hRecoProtonsMultiRP->size();
95  multiRP_arm.reserve(num_proton);
96 
97  for (const auto &proton : *hRecoProtonsMultiRP) {
98  multiRP_arm.push_back((proton.pz() < 0.) ? 1 : 0);
99  }
100  }
101 
102  // process local tracks
103  for (unsigned int tr_idx = 0; tr_idx < ppsTracksLite->size(); ++tr_idx) {
104  const auto &tr = ppsTracksLite->at(tr_idx);
105 
106  bool found = false;
107 
108  CTPPSDetId rpId(tr.getRPId());
109  unsigned int rpDecId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
110 
111  signed int singleRP_idx = -1;
112  for (unsigned int p_idx = 0; p_idx < hRecoProtonsSingleRP->size(); ++p_idx) {
113  const auto &proton = hRecoProtonsSingleRP->at(p_idx);
114 
115  for (const auto &ref : proton.contributingLocalTracks()) {
116  if (ref.key() == tr_idx) {
117  singleRP_idx = p_idx;
118  found = true;
119  }
120  }
121  }
122 
123  signed int multiRP_idx = -1;
124  for (unsigned int p_idx = 0; p_idx < hRecoProtonsMultiRP->size(); ++p_idx) {
125  const auto &proton = hRecoProtonsMultiRP->at(p_idx);
126 
127  for (const auto &ref : proton.contributingLocalTracks()) {
128  if (ref.key() == tr_idx) {
129  multiRP_idx = p_idx;
130  found = true;
131  }
132  }
133  }
134 
135  if (found) {
136  singleRPProtonIdx.push_back(singleRP_idx);
137  multiRPProtonIdx.push_back(multiRP_idx);
138  decRPId.push_back(rpDecId);
139  rpType.push_back(rpId.subdetId());
140  trackX.push_back(tr.getX());
141  trackY.push_back(tr.getY());
142  trackTime.push_back(tr.getTime());
143  trackTimeUnc.push_back(tr.getTimeUnc());
144  }
145  }
146  }
147 
148  else {
149  for (unsigned int p_idx = 0; p_idx < hRecoProtonsMultiRP->size(); ++p_idx) {
150  const auto &proton = hRecoProtonsMultiRP->at(p_idx);
151  multiRP_arm.push_back((proton.pz() < 0.) ? 1 : 0);
152 
153  for (const auto &ref : proton.contributingLocalTracks()) {
154  multiRPProtonIdx.push_back(p_idx);
155  CTPPSDetId rpId(ref->getRPId());
156  unsigned int rpDecId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
157  decRPId.push_back(rpDecId);
158  rpType.push_back(rpId.subdetId());
159  trackX.push_back(ref->getX());
160  trackY.push_back(ref->getY());
161  trackTime.push_back(ref->getTime());
162  trackTimeUnc.push_back(ref->getTimeUnc());
163  }
164  }
165  }
166 
167  // update proton tables
168  std::unique_ptr<edm::ValueMap<int>> multiRP_armV(new edm::ValueMap<int>());
169  edm::ValueMap<int>::Filler fillermultiArm(*multiRP_armV);
170  fillermultiArm.insert(hRecoProtonsMultiRP, multiRP_arm.begin(), multiRP_arm.end());
171  fillermultiArm.fill();
172 
173  std::unique_ptr<edm::ValueMap<int>> protonRPIdV(new edm::ValueMap<int>());
174  edm::ValueMap<int>::Filler fillerID(*protonRPIdV);
175  if (storeSingleRPProtons_) {
176  fillerID.insert(hRecoProtonsSingleRP, singleRP_protonRPId.begin(), singleRP_protonRPId.end());
177  fillerID.fill();
178  }
179 
180  // build track table
181  auto ppsTab = std::make_unique<nanoaod::FlatTable>(trackX.size(), "PPSLocalTrack", false);
182  ppsTab->addColumn<int>("multiRPProtonIdx", multiRPProtonIdx, "local track - proton correspondence", nanoaod::FlatTable::IntColumn);
184  ppsTab->addColumn<int>("singleRPProtonIdx", singleRPProtonIdx, "local track - proton correspondence", nanoaod::FlatTable::IntColumn);
185  ppsTab->addColumn<float>("x", trackX, "local track x", nanoaod::FlatTable::FloatColumn, 16);
186  ppsTab->addColumn<float>("y", trackY, "local track y", nanoaod::FlatTable::FloatColumn, 13);
187  ppsTab->addColumn<float>("time", trackTime, "local track time", nanoaod::FlatTable::FloatColumn, 16);
188  ppsTab->addColumn<float>("timeUnc", trackTimeUnc, "local track time uncertainty", nanoaod::FlatTable::FloatColumn, 13);
189  ppsTab->addColumn<int>("decRPId", decRPId, "local track detector dec id", nanoaod::FlatTable::IntColumn);
190  ppsTab->addColumn<int>("rpType", rpType, "strip=3, pixel=4, diamond=5, timing=6", nanoaod::FlatTable::IntColumn);
191  ppsTab->setDoc("ppsLocalTrack variables");
192 
193  // save output
194  iEvent.put(std::move(multiRP_armV), "arm");
195  iEvent.put(std::move(ppsTab), "ppsTrackTable");
197  iEvent.put(std::move(protonRPIdV), "protonRPId");
198  }
199 
200  // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
203  desc.add<edm::InputTag>("tagRecoProtonsMulti")->setComment("multiRP proton collection");
204  desc.add<edm::InputTag>("tagRecoProtonsSingle")->setComment("singleRP proton collection");
205  desc.add<edm::InputTag>("tagTrackLite")->setComment("pps local tracks lite collection");
206  desc.add<bool>("storeSingleRPProtons")->setComment("flag to store singleRP protons and associated tracks");
207  descriptions.add("ProtonProducer", desc);
208  }
209 
210 protected:
215 };
216 
const edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsMultiRP_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
Local (=single RP) track with essential information only.
~ProtonProducer() override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:53
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsSingleRP_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void produce(edm::StreamID id, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
fixed size matrix
HLT enums.
ProtonProducer(edm::ParameterSet const &ps)
const edm::EDGetTokenT< std::vector< CTPPSLocalTrackLite > > tokenTracksLite_
std::vector< ForwardProton > ForwardProtonCollection
Collection of ForwardProton objects.
EDGetTokenT< ProductType > mayConsume(edm::InputTag const &tag)
const bool storeSingleRPProtons_
def move(src, dest)
Definition: eostools.py:511