CMS 3D CMS Logo

CTPPSProtonReconstructionEfficiencyEstimatorData.cc
Go to the documentation of this file.
1 /****************************************************************************
2  * Authors:
3  * Jan Kašpar
4  * Mauricio Thiel
5  ****************************************************************************/
6 
14 
17 
19 
24 
25 #include "TFile.h"
26 #include "TH1D.h"
27 #include "TH2D.h"
28 #include "TProfile.h"
29 #include "TF1.h"
30 #include "TGraph.h"
31 
32 #include <map>
33 #include <string>
34 #include <sstream>
35 
36 //----------------------------------------------------------------------------------------------------
37 
39 public:
41 
42  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
43 
44 private:
45  void analyze(const edm::Event &, const edm::EventSetup &) override;
46  void endJob() override;
47 
49 
51 
53 
55 
57 
58  unsigned int n_prep_events_;
59 
60  unsigned int n_exp_prot_max_;
61 
62  std::vector<double> n_sigmas_;
63 
65 
66  unsigned int verbosity_;
67 
69 
70  struct ArmData {
71  unsigned int rpId_N, rpId_F;
72 
73  std::shared_ptr<const TSpline3> s_x_to_xi_N;
74 
75  unsigned int n_events_with_tracks;
76 
77  std::unique_ptr<TH1D> h_de_x, h_de_y;
78  std::unique_ptr<TH2D> h2_de_y_vs_de_x;
79 
82 
83  std::vector<double> n_sigmas;
84 
85  unsigned int n_exp_prot_max;
86 
87  struct EffPlots {
88  std::unique_ptr<TProfile> p_eff1_vs_x_N;
89  std::unique_ptr<TProfile> p_eff1_vs_xi_N;
90 
91  std::unique_ptr<TProfile> p_eff2_vs_x_N;
92  std::unique_ptr<TProfile> p_eff2_vs_xi_N;
93 
95  : p_eff1_vs_x_N(new TProfile("", ";x_{N} (mm);efficiency", 50, 0., 25.)),
96  p_eff1_vs_xi_N(new TProfile("", ";#xi_{si,N};efficiency", 50, 0., 0.25)),
97 
98  p_eff2_vs_x_N(new TProfile("", ";x_{N} (mm);efficiency", 50, 0., 25.)),
99  p_eff2_vs_xi_N(new TProfile("", ";#xi_{si,N};efficiency", 50, 0., 0.25)) {}
100 
101  void write() const {
102  p_eff1_vs_x_N->Write("p_eff1_vs_x_N");
103  p_eff1_vs_xi_N->Write("p_eff1_vs_xi_N");
104 
105  p_eff2_vs_x_N->Write("p_eff2_vs_x_N");
106  p_eff2_vs_xi_N->Write("p_eff2_vs_xi_N");
107  }
108  };
109 
110  // (n exp protons, index in n_sigmas) --> plots
111  // n exp protons = 0 --> no condition (summary case)
112  std::map<unsigned int, std::map<unsigned int, EffPlots>> effPlots;
113 
115  : n_events_with_tracks(0),
116 
117  h_de_x(new TH1D("", ";x_{F} - x_{N} (mm)", 200, -1., +1.)),
118  h_de_y(new TH1D("", ";y_{F} - y_{N} (mm)", 200, -1., +1.)),
119  h2_de_y_vs_de_x(new TH2D("", ";x_{F} - x_{N} (mm);y_{F} - y_{N} (mm)", 100, -1., +1., 100, -1., +1.)),
120 
121  de_x_mean(0.),
122  de_x_sigma(0.),
123  de_y_mean(0.),
124  de_y_sigma(0.) {
125  for (unsigned int nep = 0; nep <= n_exp_prot_max; ++nep) {
126  for (unsigned int nsi = 0; nsi < n_sigmas.size(); ++nsi) {
127  effPlots[nep][nsi] = EffPlots();
128  }
129  }
130  }
131 
132  void write() const {
133  h_de_x->Write("h_de_x");
134  h_de_y->Write("h_de_y");
135  h2_de_y_vs_de_x->Write("h2_de_y_vs_de_x");
136 
137  char buf[100];
138 
139  for (const auto &n_si : n_sigmas) {
140  sprintf(buf, "g_de_y_vs_de_x_n_si=%.1f", n_si);
141  TGraph *g = new TGraph();
142  g->SetName(buf);
143 
144  g->SetPoint(0, de_x_mean - n_si * de_x_sigma, de_y_mean - n_si * de_y_sigma);
145  g->SetPoint(1, de_x_mean + n_si * de_x_sigma, de_y_mean - n_si * de_y_sigma);
146  g->SetPoint(2, de_x_mean + n_si * de_x_sigma, de_y_mean + n_si * de_y_sigma);
147  g->SetPoint(3, de_x_mean - n_si * de_x_sigma, de_y_mean + n_si * de_y_sigma);
148  g->SetPoint(4, de_x_mean - n_si * de_x_sigma, de_y_mean - n_si * de_y_sigma);
149 
150  g->Write();
151  }
152 
153  TDirectory *d_top = gDirectory;
154 
155  for (const auto &nep_p : effPlots) {
156  if (nep_p.first == 0)
157  sprintf(buf, "exp prot all");
158  else
159  sprintf(buf, "exp prot %u", nep_p.first);
160 
161  TDirectory *d_nep = d_top->mkdir(buf);
162 
163  for (const auto &nsi_p : nep_p.second) {
164  sprintf(buf, "nsi = %.1f", n_sigmas[nsi_p.first]);
165 
166  TDirectory *d_nsi = d_nep->mkdir(buf);
167 
168  gDirectory = d_nsi;
169 
170  nsi_p.second.write();
171  }
172  }
173 
174  gDirectory = d_top;
175  }
176 
179 
180  for (const auto &p : ofc) {
181  CTPPSDetId rpId(p.first);
182  unsigned int decRPId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
183 
184  if (decRPId == rpId_N) {
185  of = &p.second;
186  break;
187  }
188  }
189 
190  if (!of) {
191  edm::LogError("ArmData::UpdateOptics") << "Cannot find optical functions for RP " << rpId_N;
192  return;
193  }
194 
195  std::vector<double> xiValues =
196  of->getXiValues(); // local copy made since the TSpline constructor needs non-const parameters
197  std::vector<double> xDValues = of->getFcnValues()[LHCOpticalFunctionsSet::exd];
198  s_x_to_xi_N = std::make_shared<TSpline3>("", xDValues.data(), xiValues.data(), xiValues.size());
199  }
200  };
201 
202  std::map<unsigned int, ArmData> data_;
203 
204  std::unique_ptr<TF1> ff_;
205 };
206 
207 //----------------------------------------------------------------------------------------------------
208 
209 using namespace std;
210 using namespace edm;
211 
212 //----------------------------------------------------------------------------------------------------
213 
215  const edm::ParameterSet &iConfig)
216  : tokenTracks_(consumes<CTPPSLocalTrackLiteCollection>(iConfig.getParameter<edm::InputTag>("tagTracks"))),
217 
219  consumes<reco::ForwardProtonCollection>(iConfig.getParameter<InputTag>("tagRecoProtonsMultiRP"))),
220 
221  pixelDiscardBXShiftedTracks_(iConfig.getParameter<bool>("pixelDiscardBXShiftedTracks")),
222 
223  localAngleXMin_(iConfig.getParameter<double>("localAngleXMin")),
224  localAngleXMax_(iConfig.getParameter<double>("localAngleXMax")),
225  localAngleYMin_(iConfig.getParameter<double>("localAngleYMin")),
226  localAngleYMax_(iConfig.getParameter<double>("localAngleYMax")),
227 
228  opticsLabel_(iConfig.getParameter<std::string>("opticsLabel")),
229  n_prep_events_(iConfig.getParameter<unsigned int>("n_prep_events")),
230  n_exp_prot_max_(iConfig.getParameter<unsigned int>("n_exp_prot_max")),
231  n_sigmas_(iConfig.getParameter<std::vector<double>>("n_sigmas")),
232 
233  outputFile_(iConfig.getParameter<string>("outputFile")),
234 
235  verbosity_(iConfig.getUntrackedParameter<unsigned int>("verbosity")),
236 
237  ff_(new TF1("ff", "[0] * exp(- (x-[1])*(x-[1]) / 2 / [2]/[2]) + [4]")) {
238  data_[0].n_exp_prot_max = n_exp_prot_max_;
239  data_[0].n_sigmas = n_sigmas_;
240  data_[0].rpId_N = iConfig.getParameter<unsigned int>("rpId_45_N");
241  data_[0].rpId_F = iConfig.getParameter<unsigned int>("rpId_45_F");
242 
243  data_[1].n_exp_prot_max = n_exp_prot_max_;
244  data_[1].n_sigmas = n_sigmas_;
245  data_[1].rpId_N = iConfig.getParameter<unsigned int>("rpId_56_N");
246  data_[1].rpId_F = iConfig.getParameter<unsigned int>("rpId_56_F");
247 }
248 
249 //----------------------------------------------------------------------------------------------------
250 
253 
254  desc.add<edm::InputTag>("tagTracks", edm::InputTag())->setComment("input tag for local lite tracks");
255  desc.add<edm::InputTag>("tagRecoProtonsMultiRP", edm::InputTag())->setComment("input tag for multi-RP reco protons");
256 
257  desc.add<bool>("pixelDiscardBXShiftedTracks", false)
258  ->setComment("whether to discard pixel tracks built from BX-shifted planes");
259 
260  desc.add<double>("localAngleXMin", -0.03)->setComment("minimal accepted value of local horizontal angle (rad)");
261  desc.add<double>("localAngleXMax", +0.03)->setComment("maximal accepted value of local horizontal angle (rad)");
262  desc.add<double>("localAngleYMin", -0.04)->setComment("minimal accepted value of local vertical angle (rad)");
263  desc.add<double>("localAngleYMax", +0.04)->setComment("maximal accepted value of local vertical angle (rad)");
264 
265  desc.add<std::string>("opticsLabel", "")->setComment("label of the optics records");
266 
267  desc.add<unsigned int>("n_prep_events", 1000)
268  ->setComment("number of preparatory events (to determine de x and de y window)");
269 
270  desc.add<unsigned int>("n_exp_prot_max", 5)->setComment("maximum number of expected protons per event and per arm");
271 
272  desc.add<std::vector<double>>("n_sigmas", {3., 5., 7.})->setComment("list of n_sigma values");
273 
274  desc.add<unsigned int>("rpId_45_N", 0)->setComment("decimal RP id for 45 near");
275  desc.add<unsigned int>("rpId_45_F", 0)->setComment("decimal RP id for 45 far");
276  desc.add<unsigned int>("rpId_56_N", 0)->setComment("decimal RP id for 56 near");
277  desc.add<unsigned int>("rpId_56_F", 0)->setComment("decimal RP id for 56 far");
278 
279  desc.add<std::string>("outputFile", "output.root")->setComment("output file name");
280 
281  desc.addUntracked<unsigned int>("verbosity", 0)->setComment("verbosity level");
282 
283  descriptions.add("ctppsProtonReconstructionEfficiencyEstimatorData", desc);
284 }
285 
286 //----------------------------------------------------------------------------------------------------
287 
289  const edm::EventSetup &iSetup) {
290  std::ostringstream os;
291 
292  // get conditions
294  iSetup.get<CTPPSInterpolatedOpticsRcd>().get(opticsLabel_, hOpticalFunctions);
295 
296  // check optics change
297  if (opticsWatcher_.check(iSetup)) {
298  data_[0].UpdateOptics(*hOpticalFunctions);
299  data_[1].UpdateOptics(*hOpticalFunctions);
300  }
301 
302  // get input
304  iEvent.getByToken(tokenTracks_, hTracks);
305 
306  Handle<reco::ForwardProtonCollection> hRecoProtonsMultiRP;
307  iEvent.getByToken(tokenRecoProtonsMultiRP_, hRecoProtonsMultiRP);
308 
309  // pre-selection of tracks
310  std::vector<unsigned int> sel_track_idc;
311  for (unsigned int idx = 0; idx < hTracks->size(); ++idx) {
312  const auto &tr = hTracks->at(idx);
313 
314  if (tr.getTx() < localAngleXMin_ || tr.getTx() > localAngleXMax_ || tr.getTy() < localAngleYMin_ ||
315  tr.getTy() > localAngleYMax_)
316  continue;
317 
319  if (tr.getPixelTrackRecoInfo() == CTPPSpixelLocalTrackReconstructionInfo::allShiftedPlanes ||
320  tr.getPixelTrackRecoInfo() == CTPPSpixelLocalTrackReconstructionInfo::mixedPlanes)
321  continue;
322  }
323 
324  sel_track_idc.push_back(idx);
325  }
326 
327  // debug print
328  if (verbosity_ > 1) {
329  os << "* tracks (pre-selected):" << std::endl;
330 
331  for (const auto idx : sel_track_idc) {
332  const auto &tr = hTracks->at(idx);
333  CTPPSDetId rpId(tr.getRPId());
334  unsigned int decRPId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
335 
336  os << " [" << idx << "] RP=" << decRPId << ", x=" << tr.getX() << ", y=" << tr.getY() << std::endl;
337  }
338 
339  os << "* protons:" << std::endl;
340  for (unsigned int i = 0; i < hRecoProtonsMultiRP->size(); ++i) {
341  const auto &pr = (*hRecoProtonsMultiRP)[i];
342  os << " [" << i << "] track indices: ";
343  for (const auto &trr : pr.contributingLocalTracks())
344  os << trr.key() << ", ";
345  os << std::endl;
346  }
347  }
348 
349  // make de_x and de_y plots
350  map<unsigned int, bool> hasTracksInArm;
351 
352  for (const auto idx_i : sel_track_idc) {
353  const auto &tr_i = hTracks->at(idx_i);
354  CTPPSDetId rpId_i(tr_i.getRPId());
355  unsigned int decRPId_i = rpId_i.arm() * 100 + rpId_i.station() * 10 + rpId_i.rp();
356 
357  for (const auto idx_j : sel_track_idc) {
358  const auto &tr_j = hTracks->at(idx_j);
359  CTPPSDetId rpId_j(tr_j.getRPId());
360  unsigned int decRPId_j = rpId_j.arm() * 100 + rpId_j.station() * 10 + rpId_j.rp();
361 
362  // check whether desired RP combination
363  unsigned int arm = 123;
364  for (const auto &ap : data_) {
365  if (ap.second.rpId_N == decRPId_i && ap.second.rpId_F == decRPId_j)
366  arm = ap.first;
367  }
368 
369  if (arm > 1)
370  continue;
371 
372  // update flag
373  hasTracksInArm[arm] = true;
374 
375  // update plots
376  auto &ad = data_[arm];
377  const double de_x = tr_j.getX() - tr_i.getX();
378  const double de_y = tr_j.getY() - tr_i.getY();
379 
380  if (ad.n_events_with_tracks < n_prep_events_) {
381  ad.h_de_x->Fill(de_x);
382  ad.h_de_y->Fill(de_y);
383  }
384 
385  ad.h2_de_y_vs_de_x->Fill(de_x, de_y);
386  }
387  }
388 
389  // update event counter
390  for (auto &ap : data_) {
391  if (hasTracksInArm[ap.first])
392  ap.second.n_events_with_tracks++;
393  }
394 
395  // if threshold reached do fits
396  for (auto &ap : data_) {
397  auto &ad = ap.second;
398 
399  if (ad.n_events_with_tracks == n_prep_events_) {
400  if (ad.de_x_sigma > 0. && ad.de_y_sigma > 0.)
401  continue;
402 
403  std::vector<std::pair<unsigned int, TH1D *>> m;
404  m.emplace_back(0, ad.h_de_x.get());
405  m.emplace_back(1, ad.h_de_y.get());
406 
407  for (const auto &p : m) {
408  double max_pos = -1E100, max_val = -1.;
409  for (int bi = 1; bi < p.second->GetNbinsX(); ++bi) {
410  const double pos = p.second->GetBinCenter(bi);
411  const double val = p.second->GetBinContent(bi);
412 
413  if (val > max_val) {
414  max_val = val;
415  max_pos = pos;
416  }
417  }
418 
419  const double sig = 0.2;
420 
421  ff_->SetParameters(max_val, max_pos, sig, 0.);
422  p.second->Fit(ff_.get(), "Q", "", max_pos - 3. * sig, max_pos + 3. * sig);
423  p.second->Fit(ff_.get(), "Q", "", max_pos - 3. * sig, max_pos + 3. * sig);
424 
425  if (p.first == 0) {
426  ad.de_x_mean = ff_->GetParameter(1);
427  ad.de_x_sigma = fabs(ff_->GetParameter(2));
428  }
429  if (p.first == 1) {
430  ad.de_y_mean = ff_->GetParameter(1);
431  ad.de_y_sigma = fabs(ff_->GetParameter(2));
432  }
433  }
434 
435  if (verbosity_) {
436  os << "* fitting arm " << ap.first << std::endl
437  << " de_x: mean = " << ad.de_x_mean << ", sigma = " << ad.de_x_sigma << std::endl
438  << " de_y: mean = " << ad.de_y_mean << ", sigma = " << ad.de_y_sigma;
439  }
440  }
441  }
442 
443  // data structures for efficiency analysis
444  struct ArmEventData {
445  std::map<unsigned int, std::set<unsigned int>> matched_track_idc;
446 
447  std::set<unsigned int> reco_proton_idc;
448 
449  std::map<unsigned int, std::set<unsigned int>> matched_track_with_prot_idc, matched_track_without_prot_idc;
450  };
451 
452  std::map<unsigned int, ArmEventData> armEventData;
453 
454  // determine the number of expected protons
455  for (const auto idx_i : sel_track_idc) {
456  const auto &tr_i = hTracks->at(idx_i);
457  CTPPSDetId rpId_i(tr_i.getRPId());
458  unsigned int decRPId_i = rpId_i.arm() * 100 + rpId_i.station() * 10 + rpId_i.rp();
459 
460  for (const auto idx_j : sel_track_idc) {
461  const auto &tr_j = hTracks->at(idx_j);
462  CTPPSDetId rpId_j(tr_j.getRPId());
463  unsigned int decRPId_j = rpId_j.arm() * 100 + rpId_j.station() * 10 + rpId_j.rp();
464 
465  // check whether desired RP combination
466  unsigned int arm = 123;
467  for (const auto &ap : data_) {
468  if (ap.second.rpId_N == decRPId_i && ap.second.rpId_F == decRPId_j)
469  arm = ap.first;
470  }
471 
472  if (arm > 1)
473  continue;
474 
475  // check near-far matching
476  auto &ad = data_[arm];
477  const double de_x = tr_j.getX() - tr_i.getX();
478  const double de_y = tr_j.getY() - tr_i.getY();
479 
480  for (unsigned int nsi = 0; nsi < ad.n_sigmas.size(); ++nsi) {
481  const double n_si = ad.n_sigmas[nsi];
482  const bool match_x = fabs(de_x - ad.de_x_mean) < n_si * ad.de_x_sigma;
483  const bool match_y = fabs(de_y - ad.de_y_mean) < n_si * ad.de_y_sigma;
484  if (match_x && match_y)
485  armEventData[arm].matched_track_idc[nsi].insert(idx_i);
486  }
487  }
488  }
489 
490  // determine the number of reconstructed protons
491  for (unsigned int i = 0; i < hRecoProtonsMultiRP->size(); ++i) {
492  const auto &proton = (*hRecoProtonsMultiRP)[i];
493 
494  CTPPSDetId rpId((*proton.contributingLocalTracks().begin())->getRPId());
495  unsigned int arm = rpId.arm();
496 
497  if (proton.validFit())
498  armEventData[arm].reco_proton_idc.insert(i);
499  }
500 
501  // compare matched tracks with reco protons
502  if (verbosity_ > 1)
503  os << "* cmp matched tracks vs. reco protons" << std::endl;
504 
505  for (auto &ap : armEventData) {
506  auto &ad = data_[ap.first];
507 
508  if (verbosity_ > 1)
509  os << " arm " << ap.first << std::endl;
510 
511  for (unsigned int nsi = 0; nsi < ad.n_sigmas.size(); ++nsi) {
512  if (verbosity_ > 1)
513  os << " nsi = " << nsi << std::endl;
514 
515  for (const auto &tri : ap.second.matched_track_idc[nsi]) {
516  const auto &track = hTracks->at(tri);
517 
518  bool some_proton_matching = false;
519 
520  if (verbosity_ > 1)
521  os << " tri = " << tri << std::endl;
522 
523  for (const auto &pri : ap.second.reco_proton_idc) {
524  const auto &proton = (*hRecoProtonsMultiRP)[pri];
525 
526  bool proton_matching = false;
527 
528  for (const auto &pr_tr : proton.contributingLocalTracks()) {
529  CTPPSDetId rpId(pr_tr->getRPId());
530  unsigned int decRPId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
531 
532  if (decRPId != ad.rpId_N)
533  continue;
534 
535  const double x = pr_tr->getX();
536  const double y = pr_tr->getY();
537  const double th = 1E-3; // 1 um
538 
539  const bool match = (fabs(x - track.getX()) < th) && (fabs(y - track.getY()) < th);
540 
541  if (verbosity_ > 1)
542  os << " pri = " << pri << ": x_tr = " << track.getX() << ", x_pr = " << x
543  << ", match = " << match << std::endl;
544 
545  if (match) {
546  proton_matching = true;
547  break;
548  }
549  }
550 
551  if (proton_matching) {
552  some_proton_matching = true;
553  break;
554  }
555  }
556 
557  if (verbosity_ > 1)
558  os << " --> some_proton_matching " << some_proton_matching << std::endl;
559 
560  if (some_proton_matching)
561  ap.second.matched_track_with_prot_idc[nsi].insert(tri);
562  else
563  ap.second.matched_track_without_prot_idc[nsi].insert(tri);
564  }
565  }
566  }
567 
568  // debug print
569  if (verbosity_ > 1) {
570  for (auto &ap : armEventData) {
571  auto &ad = data_[ap.first];
572 
573  if (ad.de_x_sigma <= 0. && ad.de_y_sigma <= 0.)
574  continue;
575 
576  os << "* results for arm " << ap.first << std::endl;
577 
578  os << " reco_proton_idc: ";
579  for (const auto &idx : ap.second.reco_proton_idc)
580  os << idx << ", ";
581  os << std::endl;
582 
583  for (unsigned int nsi = 0; nsi < ad.n_sigmas.size(); ++nsi) {
584  os << " n_si = " << ad.n_sigmas[nsi] << std::endl;
585 
586  os << " matched_track_idc: ";
587  for (const auto &idx : ap.second.matched_track_idc[nsi])
588  os << idx << ", ";
589  os << std::endl;
590 
591  os << " matched_track_with_prot_idc: ";
592  for (const auto &idx : ap.second.matched_track_with_prot_idc[nsi])
593  os << idx << ", ";
594  os << std::endl;
595 
596  os << " matched_track_without_prot_idc: ";
597  for (const auto &idx : ap.second.matched_track_without_prot_idc[nsi])
598  os << idx << ", ";
599  os << std::endl;
600  }
601  }
602  }
603 
604  // update efficiency plots
605  for (auto &ap : armEventData) {
606  auto &ad = data_[ap.first];
607 
608  // stop if sigmas not yet determined
609  if (ad.de_x_sigma <= 0. && ad.de_y_sigma <= 0.)
610  continue;
611 
612  // loop over n_sigma choices
613  for (unsigned int nsi = 0; nsi < ad.n_sigmas.size(); ++nsi) {
614  const unsigned int n_exp_prot = ap.second.matched_track_idc[nsi].size();
615  const unsigned int n_rec_prot = ap.second.reco_proton_idc.size();
616 
617  // stop if N(expected protons) out of range
618  if (n_exp_prot < 1 || n_exp_prot > ad.n_exp_prot_max)
619  continue;
620 
621  // update method 1 plots
622  const double eff = double(n_rec_prot) / n_exp_prot;
623 
624  for (unsigned int tri : ap.second.matched_track_idc[nsi]) {
625  const double x_N = hTracks->at(tri).getX();
626  const double xi_N = ad.s_x_to_xi_N->Eval(x_N * 1E-1); // conversion mm to cm
627 
628  ad.effPlots[0][nsi].p_eff1_vs_x_N->Fill(x_N, eff);
629  ad.effPlots[0][nsi].p_eff1_vs_xi_N->Fill(xi_N, eff);
630 
631  ad.effPlots[n_exp_prot][nsi].p_eff1_vs_x_N->Fill(x_N, eff);
632  ad.effPlots[n_exp_prot][nsi].p_eff1_vs_xi_N->Fill(xi_N, eff);
633  }
634 
635  // update method 2 plots
636  for (const auto &tri : ap.second.matched_track_with_prot_idc[nsi]) {
637  const double x_N = hTracks->at(tri).getX();
638  const double xi_N = ad.s_x_to_xi_N->Eval(x_N * 1E-1); // conversion mm to cm
639 
640  ad.effPlots[0][nsi].p_eff2_vs_x_N->Fill(x_N, 1.);
641  ad.effPlots[0][nsi].p_eff2_vs_xi_N->Fill(xi_N, 1.);
642 
643  ad.effPlots[n_exp_prot][nsi].p_eff2_vs_x_N->Fill(x_N, 1.);
644  ad.effPlots[n_exp_prot][nsi].p_eff2_vs_xi_N->Fill(xi_N, 1.);
645  }
646 
647  for (const auto &tri : ap.second.matched_track_without_prot_idc[nsi]) {
648  const double x_N = hTracks->at(tri).getX();
649  const double xi_N = ad.s_x_to_xi_N->Eval(x_N * 1E-1); // conversion mm to cm
650 
651  ad.effPlots[0][nsi].p_eff2_vs_x_N->Fill(x_N, 0.);
652  ad.effPlots[0][nsi].p_eff2_vs_xi_N->Fill(xi_N, 0.);
653 
654  ad.effPlots[n_exp_prot][nsi].p_eff2_vs_x_N->Fill(x_N, 0.);
655  ad.effPlots[n_exp_prot][nsi].p_eff2_vs_xi_N->Fill(xi_N, 0.);
656  }
657  }
658  }
659 
660  if (verbosity_)
661  edm::LogInfo("CTPPSProtonReconstructionEfficiencyEstimatorData") << os.str();
662 }
663 
664 //----------------------------------------------------------------------------------------------------
665 
667  auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
668 
669  for (const auto &ait : data_) {
670  char buf[100];
671  sprintf(buf, "arm %u", ait.first);
672  TDirectory *d_arm = f_out->mkdir(buf);
673  gDirectory = d_arm;
674 
675  ait.second.write();
676  }
677 }
678 
679 //----------------------------------------------------------------------------------------------------
680 
uint32_t station() const
Definition: CTPPSDetId.h:58
T getParameter(std::string const &) const
void setComment(std::string const &value)
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
#define NULL
Definition: scimark2.h:8
const std::vector< double > & getXiValues() const
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
void UpdateOptics(const LHCInterpolatedOpticalFunctionsSetCollection &ofc)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsMultiRP_
uint32_t arm() const
Definition: CTPPSDetId.h:51
ParameterDescriptionBase * add(U const &iLabel, T const &value)
std::vector< CTPPSLocalTrackLite > CTPPSLocalTrackLiteCollection
Collection of CTPPSLocalTrackLite objects.
Set of optical functions corresponding to one scoring plane along LHC, including splines for interpol...
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:52
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
fixed size matrix
HLT enums.
T get() const
Definition: EventSetup.h:71
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
std::vector< ForwardProton > ForwardProtonCollection
Collection of ForwardProton objects.
const std::vector< std::vector< double > > & getFcnValues() const
uint32_t rp() const
Definition: CTPPSDetId.h:65
void analyze(const edm::Event &, const edm::EventSetup &) override