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 
34 
37 
43 
45 public:
48  mayConsume<reco::ForwardProtonCollection>(ps.getParameter<edm::InputTag>("tagRecoProtonsMulti"))),
50  mayConsume<reco::ForwardProtonCollection>(ps.getParameter<edm::InputTag>("tagRecoProtonsSingle"))),
51  tokenTracksLite_(mayConsume<std::vector<CTPPSLocalTrackLite>>(ps.getParameter<edm::InputTag>("tagTrackLite"))),
52  storeSingleRPProtons_(ps.getParameter<bool>("storeSingleRPProtons")) {
53  produces<edm::ValueMap<int>>("arm");
54  produces<nanoaod::FlatTable>("ppsTrackTable");
56  produces<edm::ValueMap<int>>("protonRPId");
57  }
58  ~ProtonProducer() override {}
59 
60  // ------------ method called to produce the data ------------
61  void produce(edm::StreamID id, edm::Event &iEvent, const edm::EventSetup &iSetup) const override {
62  // Get Forward Proton handle
63  auto hRecoProtonsMultiRP = iEvent.getHandle(tokenRecoProtonsMultiRP_);
64  auto hRecoProtonsSingleRP = iEvent.getHandle(tokenRecoProtonsSingleRP_);
65 
66  // Get PPS Local Track handle
67  auto ppsTracksLite = iEvent.getHandle(tokenTracksLite_);
68 
69  // book output variables for protons
70  std::vector<int> multiRP_arm;
71 
72  // book output variables for tracks
73  std::vector<float> trackX, trackY, trackTime, trackTimeUnc;
74  std::vector<int> multiRPProtonIdx, decRPId, rpType;
75  std::vector<int> singleRPProtonIdx, singleRP_protonRPId;
76 
78  // process single-RP protons
79  {
80  const auto &num_proton = hRecoProtonsSingleRP->size();
81  singleRP_protonRPId.reserve(num_proton);
82 
83  for (const auto &proton : *hRecoProtonsSingleRP) {
84  CTPPSDetId rpId((*proton.contributingLocalTracks().begin())->rpId());
85  unsigned int rpDecId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
86  singleRP_protonRPId.push_back(rpDecId);
87  }
88  }
89 
90  // process multi-RP protons
91  {
92  const auto &num_proton = hRecoProtonsMultiRP->size();
93  multiRP_arm.reserve(num_proton);
94 
95  for (const auto &proton : *hRecoProtonsMultiRP) {
96  multiRP_arm.push_back((proton.pz() < 0.) ? 1 : 0);
97  }
98  }
99 
100  // process local tracks
101  for (unsigned int tr_idx = 0; tr_idx < ppsTracksLite->size(); ++tr_idx) {
102  const auto &tr = ppsTracksLite->at(tr_idx);
103 
104  bool found = false;
105 
106  CTPPSDetId rpId(tr.rpId());
107  unsigned int rpDecId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
108 
109  signed int singleRP_idx = -1;
110  for (unsigned int p_idx = 0; p_idx < hRecoProtonsSingleRP->size(); ++p_idx) {
111  const auto &proton = hRecoProtonsSingleRP->at(p_idx);
112 
113  for (const auto &ref : proton.contributingLocalTracks()) {
114  if (ref.key() == tr_idx) {
115  singleRP_idx = p_idx;
116  found = true;
117  }
118  }
119  }
120 
121  signed int multiRP_idx = -1;
122  for (unsigned int p_idx = 0; p_idx < hRecoProtonsMultiRP->size(); ++p_idx) {
123  const auto &proton = hRecoProtonsMultiRP->at(p_idx);
124 
125  for (const auto &ref : proton.contributingLocalTracks()) {
126  if (ref.key() == tr_idx) {
127  multiRP_idx = p_idx;
128  found = true;
129  }
130  }
131  }
132 
133  if (found) {
134  singleRPProtonIdx.push_back(singleRP_idx);
135  multiRPProtonIdx.push_back(multiRP_idx);
136  decRPId.push_back(rpDecId);
137  rpType.push_back(rpId.subdetId());
138  trackX.push_back(tr.x());
139  trackY.push_back(tr.y());
140  trackTime.push_back(tr.time());
141  trackTimeUnc.push_back(tr.timeUnc());
142  }
143  }
144  }
145 
146  else {
147  for (unsigned int p_idx = 0; p_idx < hRecoProtonsMultiRP->size(); ++p_idx) {
148  const auto &proton = hRecoProtonsMultiRP->at(p_idx);
149  multiRP_arm.push_back((proton.pz() < 0.) ? 1 : 0);
150 
151  for (const auto &ref : proton.contributingLocalTracks()) {
152  multiRPProtonIdx.push_back(p_idx);
153  CTPPSDetId rpId(ref->rpId());
154  unsigned int rpDecId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
155  decRPId.push_back(rpDecId);
156  rpType.push_back(rpId.subdetId());
157  trackX.push_back(ref->x());
158  trackY.push_back(ref->y());
159  trackTime.push_back(ref->time());
160  trackTimeUnc.push_back(ref->timeUnc());
161  }
162  }
163  }
164 
165  // update proton tables
166  std::unique_ptr<edm::ValueMap<int>> multiRP_armV(new edm::ValueMap<int>());
167  edm::ValueMap<int>::Filler fillermultiArm(*multiRP_armV);
168  fillermultiArm.insert(hRecoProtonsMultiRP, multiRP_arm.begin(), multiRP_arm.end());
169  fillermultiArm.fill();
170 
171  std::unique_ptr<edm::ValueMap<int>> protonRPIdV(new edm::ValueMap<int>());
172  edm::ValueMap<int>::Filler fillerID(*protonRPIdV);
173  if (storeSingleRPProtons_) {
174  fillerID.insert(hRecoProtonsSingleRP, singleRP_protonRPId.begin(), singleRP_protonRPId.end());
175  fillerID.fill();
176  }
177 
178  // build track table
179  auto ppsTab = std::make_unique<nanoaod::FlatTable>(trackX.size(), "PPSLocalTrack", false);
180  ppsTab->addColumn<int>("multiRPProtonIdx", multiRPProtonIdx, "local track - proton correspondence");
182  ppsTab->addColumn<int>("singleRPProtonIdx", singleRPProtonIdx, "local track - proton correspondence");
183  ppsTab->addColumn<float>("x", trackX, "local track x", 16);
184  ppsTab->addColumn<float>("y", trackY, "local track y", 13);
185  ppsTab->addColumn<float>("time", trackTime, "local track time", 16);
186  ppsTab->addColumn<float>("timeUnc", trackTimeUnc, "local track time uncertainty", 13);
187  ppsTab->addColumn<int>("decRPId", decRPId, "local track detector dec id");
188  ppsTab->addColumn<int>("rpType", rpType, "strip=3, pixel=4, diamond=5, timing=6");
189  ppsTab->setDoc("ppsLocalTrack variables");
190 
191  // save output
192  iEvent.put(std::move(multiRP_armV), "arm");
193  iEvent.put(std::move(ppsTab), "ppsTrackTable");
195  iEvent.put(std::move(protonRPIdV), "protonRPId");
196  }
197 
198  // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
201  desc.add<edm::InputTag>("tagRecoProtonsMulti")->setComment("multiRP proton collection");
202  desc.add<edm::InputTag>("tagRecoProtonsSingle")->setComment("singleRP proton collection");
203  desc.add<edm::InputTag>("tagTrackLite")->setComment("pps local tracks lite collection");
204  desc.add<bool>("storeSingleRPProtons")->setComment("flag to store singleRP protons and associated tracks");
205  descriptions.add("ProtonProducer", desc);
206  }
207 
208 protected:
213 };
214 
const edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsMultiRP_
Local (=single RP) track with essential information only.
~ProtonProducer() override
uint32_t arm() const
Definition: CTPPSDetId.h:57
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
const edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsSingleRP_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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
void produce(edm::StreamID id, edm::Event &iEvent, const edm::EventSetup &iSetup) const override