CMS 3D CMS Logo

CTPPSProtonReconstructionEfficiencyEstimatorMC.cc
Go to the documentation of this file.
1 /****************************************************************************
2  * Authors:
3  * Jan Kašpar
4  ****************************************************************************/
5 
12 
14 
17 
19 
27 
28 #include "TFile.h"
29 #include "TH1D.h"
30 #include "TH2D.h"
31 #include "TProfile.h"
32 #include "TGraphErrors.h"
33 
34 #include <map>
35 #include <string>
36 
37 //----------------------------------------------------------------------------------------------------
38 
40 public:
42 
43 private:
44  void analyze(const edm::Event &, const edm::EventSetup &) override;
45  void endJob() override;
46 
48 
51 
53 
55 
57 
58  unsigned int rpId_45_N_, rpId_45_F_;
59  unsigned int rpId_56_N_, rpId_56_F_;
60 
61  std::map<unsigned int, unsigned int> rpDecId_near_, rpDecId_far_;
62 
64 
65  unsigned int verbosity_;
66 
67  struct PlotGroup {
68  std::unique_ptr<TProfile> p_eff_vs_xi;
69 
70  std::unique_ptr<TH1D> h_n_part_acc_nr, h_n_part_acc_fr;
71 
72  PlotGroup()
73  : p_eff_vs_xi(new TProfile("", ";#xi_{simu};efficiency", 19, 0.015, 0.205)),
74  h_n_part_acc_nr(new TH1D("", ";n particles in acceptance", 6, -0.5, +5.5)),
75  h_n_part_acc_fr(new TH1D("", ";n particles in acceptance", 6, -0.5, +5.5)) {}
76 
77  void write() const {
78  p_eff_vs_xi->Write("p_eff_vs_xi");
79  h_n_part_acc_nr->Write("h_n_part_acc_nr");
80  h_n_part_acc_fr->Write("h_n_part_acc_fr");
81  }
82  };
83 
84  std::map<unsigned int, std::map<unsigned int, PlotGroup>> plots_; // map: arm, n(particles in acceptance) --> plots
85 };
86 
87 //----------------------------------------------------------------------------------------------------
88 
89 using namespace std;
90 using namespace edm;
91 using namespace HepMC;
92 
93 //----------------------------------------------------------------------------------------------------
94 
96  const edm::ParameterSet &iConfig)
97  : tokenHepMCAfterSmearing_(
98  consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("tagHepMCAfterSmearing"))),
99  tokenStripRecHitsPerParticle_(consumes<std::map<int, edm::DetSetVector<TotemRPRecHit>>>(
100  iConfig.getParameter<edm::InputTag>("tagStripRecHitsPerParticle"))),
101  tokenPixelRecHitsPerParticle_(consumes<std::map<int, edm::DetSetVector<CTPPSPixelRecHit>>>(
102  iConfig.getParameter<edm::InputTag>("tagPixelRecHitsPerParticle"))),
103  tracksToken_(consumes<CTPPSLocalTrackLiteCollection>(iConfig.getParameter<edm::InputTag>("tagTracks"))),
104  tokenRecoProtonsMultiRP_(
105  consumes<reco::ForwardProtonCollection>(iConfig.getParameter<InputTag>("tagRecoProtonsMultiRP"))),
106  lhcInfoLabel_(iConfig.getParameter<std::string>("lhcInfoLabel")),
107 
108  rpId_45_N_(iConfig.getParameter<unsigned int>("rpId_45_N")),
109  rpId_45_F_(iConfig.getParameter<unsigned int>("rpId_45_F")),
110  rpId_56_N_(iConfig.getParameter<unsigned int>("rpId_56_N")),
111  rpId_56_F_(iConfig.getParameter<unsigned int>("rpId_56_F")),
112 
113  outputFile_(iConfig.getParameter<string>("outputFile")),
114 
115  verbosity_(iConfig.getUntrackedParameter<unsigned int>("verbosity", 0)) {
118 
121 
122  // book plots
123  for (unsigned int arm = 0; arm < 2; ++arm) {
124  for (unsigned int np = 1; np <= 5; ++np)
125  plots_[arm][np] = PlotGroup();
126  }
127 }
128 
129 //----------------------------------------------------------------------------------------------------
130 
132  std::ostringstream os;
133 
134  // get conditions
135  edm::ESHandle<LHCInfo> hLHCInfo;
136  iSetup.get<LHCInfoRcd>().get(lhcInfoLabel_, hLHCInfo);
137 
138  // get input
139  edm::Handle<edm::HepMCProduct> hHepMCAfterSmearing;
140  iEvent.getByToken(tokenHepMCAfterSmearing_, hHepMCAfterSmearing);
141  HepMC::GenEvent *hepMCEventAfterSmearing = (HepMC::GenEvent *)hHepMCAfterSmearing->GetEvent();
142 
144  iEvent.getByToken(tokenStripRecHitsPerParticle_, hStripRecHitsPerParticle);
145 
147  iEvent.getByToken(tokenPixelRecHitsPerParticle_, hPixelRecHitsPerParticle);
148 
150  iEvent.getByToken(tracksToken_, tracks);
151 
152  Handle<reco::ForwardProtonCollection> hRecoProtonsMultiRP;
153  iEvent.getByToken(tokenRecoProtonsMultiRP_, hRecoProtonsMultiRP);
154 
155  // buffer for particle information
156  struct ParticleInfo {
157  unsigned int arm = 2;
158  double xi = 0.;
159  std::map<unsigned int, unsigned int> recHitsPerRP;
160  bool inAcceptanceNear = false, inAcceptanceFar = false, inAcceptance = false;
161  };
162 
163  std::map<int, ParticleInfo> particleInfo; // barcode --> info
164 
165  // process HepMC
166  for (auto it = hepMCEventAfterSmearing->particles_begin(); it != hepMCEventAfterSmearing->particles_end(); ++it) {
167  const auto &part = *it;
168 
169  // accept only stable non-beam protons
170  if (part->pdg_id() != 2212)
171  continue;
172 
173  if (part->status() != 1)
174  continue;
175 
176  if (part->is_beam())
177  continue;
178 
179  const auto &mom = part->momentum();
180 
181  if (mom.e() < 4500.)
182  continue;
183 
185 
186  info.arm = (mom.z() > 0.) ? 0 : 1;
187 
188  const double p_nom = hLHCInfo->energy();
189  info.xi = (p_nom - mom.rho()) / p_nom;
190 
191  particleInfo[part->barcode()] = std::move(info);
192  }
193 
194  // check acceptance
195  for (const auto &pp : *hStripRecHitsPerParticle) {
196  const auto barcode = pp.first;
197 
198  for (const auto &ds : pp.second) {
199  CTPPSDetId detId(ds.detId());
200  CTPPSDetId rpId = detId.rpId();
201  particleInfo[barcode].recHitsPerRP[rpId] += ds.size();
202  }
203  }
204 
205  for (const auto &pp : *hPixelRecHitsPerParticle) {
206  const auto barcode = pp.first;
207 
208  for (const auto &ds : pp.second) {
209  CTPPSDetId detId(ds.detId());
210  CTPPSDetId rpId = detId.rpId();
211  particleInfo[barcode].recHitsPerRP[rpId] += ds.size();
212  }
213  }
214 
215  std::map<unsigned int, bool> isStripRPNear, isStripRPFar;
216 
217  for (auto &pp : particleInfo) {
218  if (verbosity_)
219  os << "* barcode=" << pp.first << ", arm=" << pp.second.arm << ", xi=" << pp.second.xi << std::endl;
220 
221  for (const auto &rpp : pp.second.recHitsPerRP) {
222  CTPPSDetId rpId(rpp.first);
223  unsigned int needed_rec_hits = 1000;
224  if (rpId.subdetId() == CTPPSDetId::sdTrackingStrip)
225  needed_rec_hits = 6;
226  if (rpId.subdetId() == CTPPSDetId::sdTrackingPixel)
227  needed_rec_hits = 3;
228 
229  unsigned int rpDecId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
230 
231  if (rpId.subdetId() == CTPPSDetId::sdTrackingStrip) {
232  if (rpDecId == rpDecId_near_[rpId.arm()])
233  isStripRPNear[rpId.arm()] = true;
234  if (rpDecId == rpDecId_far_[rpId.arm()])
235  isStripRPFar[rpId.arm()] = true;
236  }
237 
238  if (rpp.second >= needed_rec_hits) {
239  if (rpDecId == rpDecId_near_[rpId.arm()])
240  pp.second.inAcceptanceNear = true;
241  if (rpDecId == rpDecId_far_[rpId.arm()])
242  pp.second.inAcceptanceFar = true;
243  }
244 
245  if (verbosity_)
246  os << " RP " << rpDecId << ": " << rpp.second << " hits" << std::endl;
247  }
248 
249  pp.second.inAcceptance = pp.second.inAcceptanceNear && pp.second.inAcceptanceFar;
250 
251  if (verbosity_)
252  os << " inAcceptance: near=" << pp.second.inAcceptanceNear << ", far=" << pp.second.inAcceptanceFar
253  << ", global=" << pp.second.inAcceptance << std::endl;
254  }
255 
256  // count particles in acceptance
257  struct ArmCounter {
258  unsigned int near = 0, far = 0, global = 0;
259  };
260  std::map<unsigned int, ArmCounter> nParticlesInAcceptance;
261  for (auto &pp : particleInfo) {
262  auto &counter = nParticlesInAcceptance[pp.second.arm];
263  if (pp.second.inAcceptanceNear)
264  counter.near++;
265  if (pp.second.inAcceptanceFar)
266  counter.far++;
267  if (pp.second.inAcceptance)
268  counter.global++;
269  }
270 
271  // count reconstructed tracks
272  std::map<unsigned int, ArmCounter> nReconstructedTracks;
273  for (const auto &tr : *tracks) {
274  CTPPSDetId rpId(tr.rpId());
275  unsigned int rpDecId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
276 
277  if (rpDecId == rpDecId_near_[rpId.arm()])
278  nReconstructedTracks[rpId.arm()].near++;
279  if (rpDecId == rpDecId_far_[rpId.arm()])
280  nReconstructedTracks[rpId.arm()].far++;
281  }
282 
283  // count reconstructed protons
284  std::map<unsigned int, unsigned int> nReconstructedProtons;
285  for (const auto &pr : *hRecoProtonsMultiRP) {
286  if (!pr.validFit())
287  continue;
288 
289  unsigned int arm = 2;
290  if (pr.lhcSector() == reco::ForwardProton::LHCSector::sector45)
291  arm = 0;
292  if (pr.lhcSector() == reco::ForwardProton::LHCSector::sector56)
293  arm = 1;
294 
295  nReconstructedProtons[arm]++;
296  }
297 
298  // fill plots
299  for (unsigned int arm = 0; arm < 2; arm++) {
300  const auto &npa = nParticlesInAcceptance[arm];
301  const auto &nrt = nReconstructedTracks[arm];
302 
303  if (verbosity_)
304  os << "* arm " << arm << ": nRecoProtons=" << nReconstructedProtons[arm]
305  << " (tracks near=" << nReconstructedTracks[arm].near << ", far=" << nReconstructedTracks[arm].far
306  << "), nAcc=" << npa.global << " (near=" << npa.near << ", far=" << npa.far << ")" << std::endl;
307 
308  // skip event if no track in global acceptance
309  if (npa.global < 1)
310  continue;
311 
312  const auto &p = plots_[arm][npa.global];
313 
314  p.h_n_part_acc_nr->Fill(npa.near);
315  p.h_n_part_acc_fr->Fill(npa.far);
316 
317  // skip events with some local reconstruction inefficiency
318  if (nrt.near != npa.near || nrt.far != npa.far)
319  continue;
320 
321  const double eff = double(nReconstructedProtons[arm]) / npa.global;
322 
323  if (verbosity_)
324  os << " eff=" << eff << std::endl;
325 
326  for (auto &pp : particleInfo) {
327  if (pp.second.arm != arm || !pp.second.inAcceptance)
328  continue;
329 
330  p.p_eff_vs_xi->Fill(pp.second.xi, eff);
331  }
332  }
333 
334  if (verbosity_)
335  edm::LogInfo("CTPPSProtonReconstructionEfficiencyEstimatorMC") << os.str();
336 }
337 
338 //----------------------------------------------------------------------------------------------------
339 
341  auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
342 
343  for (const auto &ait : plots_) {
344  char buf[100];
345  sprintf(buf, "arm%u", ait.first);
346  TDirectory *d_arm = f_out->mkdir(buf);
347 
348  for (const auto &npit : ait.second) {
349  sprintf(buf, "%u", npit.first);
350  TDirectory *d_np = d_arm->mkdir(buf);
351  gDirectory = d_np;
352 
353  npit.second.write();
354  }
355  }
356 }
357 
358 //----------------------------------------------------------------------------------------------------
359 
LHCInfo::energy
const float energy() const
Definition: LHCInfo.cc:190
TotemRPRecHit
Reconstructed hit in TOTEM RP.
Definition: TotemRPRecHit.h:17
edm::DetSetVector
Definition: DetSetVector.h:61
counter
Definition: counter.py:1
CTPPSLocalTrackLiteCollection
std::vector< CTPPSLocalTrackLite > CTPPSLocalTrackLiteCollection
Collection of CTPPSLocalTrackLite objects.
Definition: CTPPSLocalTrackLiteFwd.h:18
PDWG_EXOHSCP_cff.tracks
tracks
Definition: PDWG_EXOHSCP_cff.py:28
EDAnalyzer.h
CTPPSPixelRecHit
Definition: CTPPSPixelRecHit.h:17
CTPPSProtonReconstructionEfficiencyEstimatorMC::tokenHepMCAfterSmearing_
edm::EDGetTokenT< edm::HepMCProduct > tokenHepMCAfterSmearing_
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:48
CTPPSProtonReconstructionEfficiencyEstimatorMC::PlotGroup::h_n_part_acc_fr
std::unique_ptr< TH1D > h_n_part_acc_fr
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:71
ESHandle.h
edm::EDGetTokenT< edm::HepMCProduct >
edm
HLT enums.
Definition: AlignableModifier.h:19
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
np
int np
Definition: AMPTWrapper.h:43
hybridSuperClusters_cfi.xi
xi
Definition: hybridSuperClusters_cfi.py:10
CTPPSProtonReconstructionEfficiencyEstimatorMC::lhcInfoLabel_
std::string lhcInfoLabel_
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:57
LHCInfo.h
CTPPSProtonReconstructionEfficiencyEstimatorMC::PlotGroup
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:68
info
static const TGPicture * info(bool iBackgroundIsBlack)
Definition: FWCollectionSummaryWidget.cc:153
CTPPSLocalTrackLite.h
CTPPSProtonReconstructionEfficiencyEstimatorMC::endJob
void endJob() override
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:339
TotemRPRecHit.h
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
edm::LogInfo
Log< level::Info, false > LogInfo
Definition: MessageLogger.h:125
edm::one::EDAnalyzer
Definition: EDAnalyzer.h:30
CTPPSProtonReconstructionEfficiencyEstimatorMC::rpDecId_far_
std::map< unsigned int, unsigned int > rpDecId_far_
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:62
year_2016_postTS2_cff.rpId
rpId
Definition: year_2016_postTS2_cff.py:23
edm::Handle< edm::HepMCProduct >
CTPPSPixelRecHit.h
CTPPSProtonReconstructionEfficiencyEstimatorMC::rpId_56_N_
unsigned int rpId_56_N_
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:60
HepMC::GenEvent
Definition: hepmc_rootio.cc:9
CTPPSDetId::sdTrackingStrip
Definition: CTPPSDetId.h:44
MakerMacros.h
part
part
Definition: HCALResponse.h:20
edm::EventSetup::get
T get() const
Definition: EventSetup.h:80
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
CTPPSProtonReconstructionEfficiencyEstimatorMC::PlotGroup::p_eff_vs_xi
std::unique_ptr< TProfile > p_eff_vs_xi
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:69
edm::ESHandle< LHCInfo >
ForwardProtonFwd.h
LHCInfoRcd
Definition: LHCInfoRcd.h:24
CTPPSDetId::sdTrackingPixel
Definition: CTPPSDetId.h:44
CTPPSProtonReconstructionEfficiencyEstimatorMC::rpId_45_F_
unsigned int rpId_45_F_
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:59
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
CTPPSProtonReconstructionEfficiencyEstimatorMC::verbosity_
unsigned int verbosity_
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:66
CTPPSProtonReconstructionEfficiencyEstimatorMC::tokenStripRecHitsPerParticle_
edm::EDGetTokenT< std::map< int, edm::DetSetVector< TotemRPRecHit > > > tokenStripRecHitsPerParticle_
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:50
reco::ForwardProtonCollection
std::vector< ForwardProton > ForwardProtonCollection
Collection of ForwardProton objects.
Definition: ForwardProtonFwd.h:25
edm::ParameterSet
Definition: ParameterSet.h:47
CTPPSProtonReconstructionEfficiencyEstimatorMC::rpId_56_F_
unsigned int rpId_56_F_
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:60
sipixeldigitoraw
Definition: SiPixelDigiToRaw.cc:39
Event.h
CTPPSProtonReconstructionEfficiencyEstimatorMC::rpId_45_N_
unsigned int rpId_45_N_
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:59
CTPPSProtonReconstructionEfficiencyEstimatorMC::PlotGroup::h_n_part_acc_nr
std::unique_ptr< TH1D > h_n_part_acc_nr
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:71
ForwardProton.h
reco::ForwardProton::LHCSector::sector56
CTPPSDetId
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:31
createfilelist.int
int
Definition: createfilelist.py:10
iEvent
int iEvent
Definition: GenABIO.cc:224
CTPPSProtonReconstructionEfficiencyEstimatorMC
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:38
CTPPSProtonReconstructionEfficiencyEstimatorMC::outputFile_
std::string outputFile_
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:64
CTPPSProtonReconstructionEfficiencyEstimatorMC::CTPPSProtonReconstructionEfficiencyEstimatorMC
CTPPSProtonReconstructionEfficiencyEstimatorMC(const edm::ParameterSet &)
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:94
edm::EventSetup
Definition: EventSetup.h:57
CTPPSProtonReconstructionEfficiencyEstimatorMC::tracksToken_
edm::EDGetTokenT< CTPPSLocalTrackLiteCollection > tracksToken_
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:53
edm::HepMCProduct::GetEvent
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:34
DetSetVector.h
get
#define get
visDQMUpload.buf
buf
Definition: visDQMUpload.py:154
CTPPSProtonReconstructionEfficiencyEstimatorMC::rpDecId_near_
std::map< unsigned int, unsigned int > rpDecId_near_
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:62
CTPPSProtonReconstructionEfficiencyEstimatorMC::PlotGroup::PlotGroup
PlotGroup()
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:73
CTPPSProtonReconstructionEfficiencyEstimatorMC::tokenPixelRecHitsPerParticle_
edm::EDGetTokenT< std::map< int, edm::DetSetVector< CTPPSPixelRecHit > > > tokenPixelRecHitsPerParticle_
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:51
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
CTPPSProtonReconstructionEfficiencyEstimatorMC::PlotGroup::write
void write() const
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:78
Frameworkfwd.h
reco::ForwardProton::LHCSector::sector45
HepMC
Definition: GenParticle.h:15
CTPPSProtonReconstructionEfficiencyEstimatorMC::tokenRecoProtonsMultiRP_
edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsMultiRP_
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:55
CTPPSLocalTrackLiteFwd.h
EventSetup.h
CTPPSDetId.h
CTPPSProtonReconstructionEfficiencyEstimatorMC::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:130
CTPPSProtonReconstructionEfficiencyEstimatorMC::plots_
std::map< unsigned int, std::map< unsigned int, PlotGroup > > plots_
Definition: CTPPSProtonReconstructionEfficiencyEstimatorMC.cc:85
genParticles_cff.map
map
Definition: genParticles_cff.py:11
edm::HepMCProduct
Definition: HepMCProduct.h:18
LHCInfoRcd.h
createTree.pp
pp
Definition: createTree.py:17
HepMCProduct.h
edm::Event
Definition: Event.h:73
ParticleInfo
Definition: ParticleTable.h:15
edm::InputTag
Definition: InputTag.h:15