CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
CTPPSProtonReconstructionEfficiencyEstimatorData.cc
Go to the documentation of this file.
1 /****************************************************************************
2  * Authors:
3  * Jan Kašpar
4  * Mauricio Thiel
5  ****************************************************************************/
6 
14 
21 
23 
28 
29 #include "TFile.h"
30 #include "TH1D.h"
31 #include "TH2D.h"
32 #include "TProfile.h"
33 #include "TF1.h"
34 #include "TGraph.h"
35 
36 #include <map>
37 #include <string>
38 #include <sstream>
39 #include <utility>
40 
41 //----------------------------------------------------------------------------------------------------
42 
44 public:
46 
47  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
48 
49 private:
50  void analyze(const edm::Event &, const edm::EventSetup &) override;
51  void endJob() override;
52 
55 
59 
61 
63 
64  unsigned int n_prep_events_;
65 
66  unsigned int n_exp_prot_max_;
67 
68  std::vector<double> n_sigmas_;
69 
71 
72  unsigned int verbosity_;
73 
75 
76  struct ArmData {
77  unsigned int rpId_N, rpId_F;
78 
79  std::shared_ptr<const TSpline3> s_x_to_xi_N;
80 
81  unsigned int n_events_with_tracks;
82 
83  std::unique_ptr<TH1D> h_de_x, h_de_y;
84  std::unique_ptr<TH2D> h2_de_y_vs_de_x;
85 
88 
89  std::vector<double> n_sigmas;
90 
91  unsigned int n_exp_prot_max;
92 
93  struct EffPlots {
94  std::unique_ptr<TProfile> p_eff1_vs_x_N;
95  std::unique_ptr<TProfile> p_eff1_vs_xi_N;
96 
97  std::unique_ptr<TProfile> p_eff2_vs_x_N;
98  std::unique_ptr<TProfile> p_eff2_vs_xi_N;
99 
100  std::unique_ptr<TProfile> p_fr_match_unique, p_fr_match_non_unique;
102 
104  : p_eff1_vs_x_N(new TProfile("", ";x_{N} (mm);efficiency", 50, 0., 25.)),
105  p_eff1_vs_xi_N(new TProfile("", ";#xi_{si,N};efficiency", 50, 0., 0.25)),
106 
107  p_eff2_vs_x_N(new TProfile("", ";x_{N} (mm);efficiency", 50, 0., 25.)),
108  p_eff2_vs_xi_N(new TProfile("", ";#xi_{si,N};efficiency", 50, 0., 0.25)),
109 
110  p_fr_match_unique(new TProfile("", ";x_{N} (mm);fraction", 50, 0., 25.)),
111  p_fr_match_non_unique(new TProfile("", ";x_{N} (mm);fraction", 50, 0., 25.)),
112 
113  p_fr_signature_0(new TProfile("", ";x_{N} (mm);fraction", 50, 0., 25.)),
114  p_fr_signature_1(new TProfile("", ";x_{N} (mm);fraction", 50, 0., 25.)),
115  p_fr_signature_2(new TProfile("", ";x_{N} (mm);fraction", 50, 0., 25.)),
116  p_fr_signature_12(new TProfile("", ";x_{N} (mm);fraction", 50, 0., 25.)) {}
117 
118  void write() const {
119  p_eff1_vs_x_N->Write("p_eff1_vs_x_N");
120  p_eff1_vs_xi_N->Write("p_eff1_vs_xi_N");
121 
122  p_eff2_vs_x_N->Write("p_eff2_vs_x_N");
123  p_eff2_vs_xi_N->Write("p_eff2_vs_xi_N");
124 
125  p_fr_match_unique->Write("p_fr_match_unique");
126  p_fr_match_non_unique->Write("p_fr_match_non_unique");
127 
128  p_fr_signature_0->Write("p_fr_signature_0");
129  p_fr_signature_1->Write("p_fr_signature_1");
130  p_fr_signature_2->Write("p_fr_signature_2");
131  p_fr_signature_12->Write("p_fr_signature_12");
132  }
133  };
134 
135  // (n exp protons, index in n_sigmas) --> plots
136  // n exp protons = 0 --> no condition (summary case)
137  std::map<unsigned int, std::map<unsigned int, EffPlots>> effPlots;
138 
141 
142  h_de_x(new TH1D("", ";x_{F} - x_{N} (mm)", 200, -1., +1.)),
143  h_de_y(new TH1D("", ";y_{F} - y_{N} (mm)", 200, -1., +1.)),
144  h2_de_y_vs_de_x(new TH2D("", ";x_{F} - x_{N} (mm);y_{F} - y_{N} (mm)", 100, -1., +1., 100, -1., +1.)),
145 
146  de_x_mean(0.),
147  de_x_sigma(0.),
148  de_y_mean(0.),
149  de_y_sigma(0.) {
150  for (unsigned int nep = 0; nep <= n_exp_prot_max; ++nep) {
151  for (unsigned int nsi = 0; nsi < n_sigmas.size(); ++nsi) {
152  effPlots[nep][nsi] = EffPlots();
153  }
154  }
155  }
156 
157  void write() const {
158  h_de_x->Write("h_de_x");
159  h_de_y->Write("h_de_y");
160  h2_de_y_vs_de_x->Write("h2_de_y_vs_de_x");
161 
162  char buf[100];
163 
164  for (const auto &n_si : n_sigmas) {
165  sprintf(buf, "g_de_y_vs_de_x_n_si=%.1f", n_si);
166  TGraph *g = new TGraph();
167  g->SetName(buf);
168 
169  g->SetPoint(0, de_x_mean - n_si * de_x_sigma, de_y_mean - n_si * de_y_sigma);
170  g->SetPoint(1, de_x_mean + n_si * de_x_sigma, de_y_mean - n_si * de_y_sigma);
171  g->SetPoint(2, de_x_mean + n_si * de_x_sigma, de_y_mean + n_si * de_y_sigma);
172  g->SetPoint(3, de_x_mean - n_si * de_x_sigma, de_y_mean + n_si * de_y_sigma);
173  g->SetPoint(4, de_x_mean - n_si * de_x_sigma, de_y_mean - n_si * de_y_sigma);
174 
175  g->Write();
176  }
177 
178  TDirectory *d_top = gDirectory;
179 
180  for (const auto &nep_p : effPlots) {
181  if (nep_p.first == 0)
182  sprintf(buf, "exp prot all");
183  else
184  sprintf(buf, "exp prot %u", nep_p.first);
185 
186  TDirectory *d_nep = d_top->mkdir(buf);
187 
188  for (const auto &nsi_p : nep_p.second) {
189  sprintf(buf, "nsi = %.1f", n_sigmas[nsi_p.first]);
190 
191  TDirectory *d_nsi = d_nep->mkdir(buf);
192 
193  gDirectory = d_nsi;
194 
195  nsi_p.second.write();
196  }
197  }
198 
199  gDirectory = d_top;
200  }
201 
203  const LHCInterpolatedOpticalFunctionsSet *of = nullptr;
204 
205  for (const auto &p : ofc) {
206  CTPPSDetId rpId(p.first);
207  unsigned int decRPId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
208 
209  if (decRPId == rpId_N) {
210  of = &p.second;
211  break;
212  }
213  }
214 
215  if (!of) {
216  edm::LogError("ArmData::UpdateOptics") << "Cannot find optical functions for RP " << rpId_N;
217  return;
218  }
219 
220  std::vector<double> xiValues =
221  of->getXiValues(); // local copy made since the TSpline constructor needs non-const parameters
222  std::vector<double> xDValues = of->getFcnValues()[LHCOpticalFunctionsSet::exd];
223  s_x_to_xi_N = std::make_shared<TSpline3>("", xDValues.data(), xiValues.data(), xiValues.size());
224  }
225  };
226 
227  std::map<unsigned int, ArmData> data_;
228 
229  std::unique_ptr<TF1> ff_;
230 };
231 
232 //----------------------------------------------------------------------------------------------------
233 
234 using namespace std;
235 using namespace edm;
236 
237 //----------------------------------------------------------------------------------------------------
238 
240  const edm::ParameterSet &iConfig)
241  : tokenTracks_(consumes<CTPPSLocalTrackLiteCollection>(iConfig.getParameter<edm::InputTag>("tagTracks"))),
242 
243  tokenRecoProtonsMultiRP_(
244  consumes<reco::ForwardProtonCollection>(iConfig.getParameter<InputTag>("tagRecoProtonsMultiRP"))),
245 
246  lhcInfoESToken_(esConsumes(ESInputTag("", iConfig.getParameter<std::string>("lhcInfoLabel")))),
247  opticsESToken_(esConsumes(ESInputTag("", iConfig.getParameter<std::string>("opticsLabel")))),
248  ppsAssociationCutsToken_(
249  esConsumes(ESInputTag("", iConfig.getParameter<std::string>("ppsAssociationCutsLabel")))),
250 
251  pixelDiscardBXShiftedTracks_(iConfig.getParameter<bool>("pixelDiscardBXShiftedTracks")),
252 
253  localAngleXMin_(iConfig.getParameter<double>("localAngleXMin")),
254  localAngleXMax_(iConfig.getParameter<double>("localAngleXMax")),
255  localAngleYMin_(iConfig.getParameter<double>("localAngleYMin")),
256  localAngleYMax_(iConfig.getParameter<double>("localAngleYMax")),
257 
258  n_prep_events_(iConfig.getParameter<unsigned int>("n_prep_events")),
259  n_exp_prot_max_(iConfig.getParameter<unsigned int>("n_exp_prot_max")),
260  n_sigmas_(iConfig.getParameter<std::vector<double>>("n_sigmas")),
261 
262  outputFile_(iConfig.getParameter<string>("outputFile")),
263 
264  verbosity_(iConfig.getUntrackedParameter<unsigned int>("verbosity")),
265 
266  ff_(new TF1("ff", "[0] * exp(- (x-[1])*(x-[1]) / 2 / [2]/[2]) + [4]")) {
267  data_[0].n_exp_prot_max = n_exp_prot_max_;
268  data_[0].n_sigmas = n_sigmas_;
269  data_[0].rpId_N = iConfig.getParameter<unsigned int>("rpId_45_N");
270  data_[0].rpId_F = iConfig.getParameter<unsigned int>("rpId_45_F");
271 
272  data_[1].n_exp_prot_max = n_exp_prot_max_;
273  data_[1].n_sigmas = n_sigmas_;
274  data_[1].rpId_N = iConfig.getParameter<unsigned int>("rpId_56_N");
275  data_[1].rpId_F = iConfig.getParameter<unsigned int>("rpId_56_F");
276 }
277 
278 //----------------------------------------------------------------------------------------------------
279 
282 
283  desc.add<edm::InputTag>("tagTracks", edm::InputTag())->setComment("input tag for local lite tracks");
284  desc.add<edm::InputTag>("tagRecoProtonsMultiRP", edm::InputTag())->setComment("input tag for multi-RP reco protons");
285 
286  desc.add<string>("lhcInfoLabel", "")->setComment("label of LHCInfo data");
287  desc.add<string>("opticsLabel", "")->setComment("label of optics data");
288  desc.add<string>("ppsAssociationCutsLabel", "")->setComment("label of PPSAssociationCuts data");
289 
290  desc.add<bool>("pixelDiscardBXShiftedTracks", false)
291  ->setComment("whether to discard pixel tracks built from BX-shifted planes");
292 
293  desc.add<double>("localAngleXMin", -0.03)->setComment("minimal accepted value of local horizontal angle (rad)");
294  desc.add<double>("localAngleXMax", +0.03)->setComment("maximal accepted value of local horizontal angle (rad)");
295  desc.add<double>("localAngleYMin", -0.04)->setComment("minimal accepted value of local vertical angle (rad)");
296  desc.add<double>("localAngleYMax", +0.04)->setComment("maximal accepted value of local vertical angle (rad)");
297 
298  desc.add<unsigned int>("n_prep_events", 1000)
299  ->setComment("number of preparatory events (to determine de x and de y window)");
300 
301  desc.add<unsigned int>("n_exp_prot_max", 5)->setComment("maximum number of expected protons per event and per arm");
302 
303  desc.add<std::vector<double>>("n_sigmas", {3., 5., 7.})->setComment("list of n_sigma values");
304 
305  desc.add<unsigned int>("rpId_45_N", 0)->setComment("decimal RP id for 45 near");
306  desc.add<unsigned int>("rpId_45_F", 0)->setComment("decimal RP id for 45 far");
307  desc.add<unsigned int>("rpId_56_N", 0)->setComment("decimal RP id for 56 near");
308  desc.add<unsigned int>("rpId_56_F", 0)->setComment("decimal RP id for 56 far");
309 
310  desc.add<std::string>("outputFile", "output.root")->setComment("output file name");
311 
312  desc.addUntracked<unsigned int>("verbosity", 0)->setComment("verbosity level");
313 
314  descriptions.add("ctppsProtonReconstructionEfficiencyEstimatorData", desc);
315 }
316 
317 //----------------------------------------------------------------------------------------------------
318 
320  const edm::EventSetup &iSetup) {
321  std::ostringstream os;
322 
323  // get conditions
324  const auto &lhcInfo = iSetup.getData(lhcInfoESToken_);
325  const auto &opticalFunctions = iSetup.getData(opticsESToken_);
326  const auto &ppsAssociationCuts = iSetup.getData(ppsAssociationCutsToken_);
327 
328  // check optics change
329  if (opticsWatcher_.check(iSetup)) {
330  data_[0].UpdateOptics(opticalFunctions);
331  data_[1].UpdateOptics(opticalFunctions);
332  }
333 
334  // get input
336  iEvent.getByToken(tokenTracks_, hTracks);
337 
338  Handle<reco::ForwardProtonCollection> hRecoProtonsMultiRP;
339  iEvent.getByToken(tokenRecoProtonsMultiRP_, hRecoProtonsMultiRP);
340 
341  // pre-selection of tracks
342  std::vector<unsigned int> sel_track_idc;
343 
344  std::map<unsigned int, std::vector<unsigned int>> trackingSelection;
345 
346  for (unsigned int idx = 0; idx < hTracks->size(); ++idx) {
347  const auto &tr = hTracks->at(idx);
348 
349  if (tr.tx() < localAngleXMin_ || tr.tx() > localAngleXMax_ || tr.ty() < localAngleYMin_ ||
350  tr.ty() > localAngleYMax_)
351  continue;
352 
354  if (tr.pixelTrackRecoInfo() == CTPPSpixelLocalTrackReconstructionInfo::allShiftedPlanes ||
355  tr.pixelTrackRecoInfo() == CTPPSpixelLocalTrackReconstructionInfo::mixedPlanes)
356  continue;
357  }
358 
359  sel_track_idc.push_back(idx);
360 
361  const CTPPSDetId rpId(tr.rpId());
362 
363  const bool trackerRP =
364  (rpId.subdetId() == CTPPSDetId::sdTrackingStrip || rpId.subdetId() == CTPPSDetId::sdTrackingPixel);
365 
366  if (trackerRP)
367  trackingSelection[rpId.arm()].push_back(idx);
368  }
369 
370  // debug print
371  if (verbosity_ > 1) {
372  os << "* tracks (pre-selected):" << std::endl;
373 
374  for (const auto idx : sel_track_idc) {
375  const auto &tr = hTracks->at(idx);
376  CTPPSDetId rpId(tr.rpId());
377  unsigned int decRPId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
378 
379  os << " [" << idx << "] RP=" << decRPId << ", x=" << tr.x() << ", y=" << tr.y() << std::endl;
380  }
381 
382  os << "* protons:" << std::endl;
383  for (unsigned int i = 0; i < hRecoProtonsMultiRP->size(); ++i) {
384  const auto &pr = (*hRecoProtonsMultiRP)[i];
385  os << " [" << i << "] track indices: ";
386  for (const auto &trr : pr.contributingLocalTracks())
387  os << trr.key() << ", ";
388  os << std::endl;
389  }
390  }
391 
392  // make de_x and de_y plots
393  map<unsigned int, bool> hasTracksInArm;
394 
395  for (const auto idx_i : sel_track_idc) {
396  const auto &tr_i = hTracks->at(idx_i);
397  CTPPSDetId rpId_i(tr_i.rpId());
398  unsigned int decRPId_i = rpId_i.arm() * 100 + rpId_i.station() * 10 + rpId_i.rp();
399 
400  for (const auto idx_j : sel_track_idc) {
401  const auto &tr_j = hTracks->at(idx_j);
402  CTPPSDetId rpId_j(tr_j.rpId());
403  unsigned int decRPId_j = rpId_j.arm() * 100 + rpId_j.station() * 10 + rpId_j.rp();
404 
405  // check whether desired RP combination
406  unsigned int arm = 123;
407  for (const auto &ap : data_) {
408  if (ap.second.rpId_N == decRPId_i && ap.second.rpId_F == decRPId_j)
409  arm = ap.first;
410  }
411 
412  if (arm > 1)
413  continue;
414 
415  // update flag
416  hasTracksInArm[arm] = true;
417 
418  // update plots
419  auto &ad = data_[arm];
420  const double de_x = tr_j.x() - tr_i.x();
421  const double de_y = tr_j.y() - tr_i.y();
422 
423  if (ad.n_events_with_tracks < n_prep_events_) {
424  ad.h_de_x->Fill(de_x);
425  ad.h_de_y->Fill(de_y);
426  }
427 
428  ad.h2_de_y_vs_de_x->Fill(de_x, de_y);
429  }
430  }
431 
432  // update event counter
433  for (auto &ap : data_) {
434  if (hasTracksInArm[ap.first])
435  ap.second.n_events_with_tracks++;
436  }
437 
438  // if threshold reached do fits
439  for (auto &ap : data_) {
440  auto &ad = ap.second;
441 
442  if (ad.n_events_with_tracks == n_prep_events_) {
443  if (ad.de_x_sigma > 0. && ad.de_y_sigma > 0.)
444  continue;
445 
446  std::vector<std::pair<unsigned int, TH1D *>> m;
447  m.emplace_back(0, ad.h_de_x.get());
448  m.emplace_back(1, ad.h_de_y.get());
449 
450  for (const auto &p : m) {
451  double max_pos = -1E100, max_val = -1.;
452  for (int bi = 1; bi < p.second->GetNbinsX(); ++bi) {
453  const double pos = p.second->GetBinCenter(bi);
454  const double val = p.second->GetBinContent(bi);
455 
456  if (val > max_val) {
457  max_val = val;
458  max_pos = pos;
459  }
460  }
461 
462  const double sig = 0.2;
463 
464  ff_->SetParameters(max_val, max_pos, sig, 0.);
465  p.second->Fit(ff_.get(), "Q", "", max_pos - 3. * sig, max_pos + 3. * sig);
466  p.second->Fit(ff_.get(), "Q", "", max_pos - 3. * sig, max_pos + 3. * sig);
467 
468  if (p.first == 0) {
469  ad.de_x_mean = ff_->GetParameter(1);
470  ad.de_x_sigma = fabs(ff_->GetParameter(2));
471  }
472  if (p.first == 1) {
473  ad.de_y_mean = ff_->GetParameter(1);
474  ad.de_y_sigma = fabs(ff_->GetParameter(2));
475  }
476  }
477 
478  if (verbosity_) {
479  os << "* fitting arm " << ap.first << std::endl
480  << " de_x: mean = " << ad.de_x_mean << ", sigma = " << ad.de_x_sigma << std::endl
481  << " de_y: mean = " << ad.de_y_mean << ", sigma = " << ad.de_y_sigma;
482  }
483  }
484  }
485 
486  // logic below:
487  // 1) evaluate "proton candiates" or "expected protons" = local tracks in the near RP which have a compatible
488  // partner in the far RP; the compatibility is evaluated at multiple n_sigma choices
489  // 2) check how often each "proton candidate" can be matched to a reconstructed proton
490  // 3) re-evaluate association cuts --> check details for each "proton candidate"
491 
492  // data structures for efficiency analysis
493 
494  struct ArmEventData {
495  // n_sigma --> list of "proton candidates" (local track indices in the near RP)
496  std::map<unsigned int, std::set<unsigned int>> matched_track_idc;
497 
498  // list of valid reco-proton indices (in the chosen arm)
499  std::set<unsigned int> reco_proton_idc;
500 
501  // n_sigma --> list of "proton candicates" which have/don't have matching reco proton
502  std::map<unsigned int, std::set<unsigned int>> matched_track_with_prot_idc, matched_track_without_prot_idc;
503 
504  // for re-evaluation of association cuts
505  struct EvaluatedPair {
506  unsigned idx_N, idx_F;
507  bool x_cut, y_cut;
508  bool match = false;
509  bool unique = false;
510  };
511 
512  vector<EvaluatedPair> evaluatedPairs;
513  };
514 
515  std::map<unsigned int, ArmEventData> armEventData;
516 
517  // determine the number of proton candiates or expected protons
518  for (const auto idx_i : sel_track_idc) {
519  const auto &tr_i = hTracks->at(idx_i);
520  CTPPSDetId rpId_i(tr_i.rpId());
521  unsigned int decRPId_i = rpId_i.arm() * 100 + rpId_i.station() * 10 + rpId_i.rp();
522 
523  for (const auto idx_j : sel_track_idc) {
524  const auto &tr_j = hTracks->at(idx_j);
525  CTPPSDetId rpId_j(tr_j.rpId());
526  unsigned int decRPId_j = rpId_j.arm() * 100 + rpId_j.station() * 10 + rpId_j.rp();
527 
528  // check whether desired RP combination
529  unsigned int arm = 123;
530  for (const auto &ap : data_) {
531  if (ap.second.rpId_N == decRPId_i && ap.second.rpId_F == decRPId_j)
532  arm = ap.first;
533  }
534 
535  if (arm > 1)
536  continue;
537 
538  // check near-far matching
539  auto &ad = data_[arm];
540  const double de_x = tr_j.x() - tr_i.x();
541  const double de_y = tr_j.y() - tr_i.y();
542 
543  for (unsigned int nsi = 0; nsi < ad.n_sigmas.size(); ++nsi) {
544  const double n_si = ad.n_sigmas[nsi];
545  const bool match_x = fabs(de_x - ad.de_x_mean) < n_si * ad.de_x_sigma;
546  const bool match_y = fabs(de_y - ad.de_y_mean) < n_si * ad.de_y_sigma;
547  if (match_x && match_y)
548  armEventData[arm].matched_track_idc[nsi].insert(idx_i);
549  }
550  }
551  }
552 
553  // determine the number of reconstructed protons
554  for (unsigned int i = 0; i < hRecoProtonsMultiRP->size(); ++i) {
555  const auto &proton = (*hRecoProtonsMultiRP)[i];
556 
557  CTPPSDetId rpId((*proton.contributingLocalTracks().begin())->rpId());
558  unsigned int arm = rpId.arm();
559 
560  if (proton.validFit())
561  armEventData[arm].reco_proton_idc.insert(i);
562  }
563 
564  // compare matched tracks with reco protons
565  if (verbosity_ > 1)
566  os << "* cmp matched tracks vs. reco protons" << std::endl;
567 
568  for (auto &ap : armEventData) {
569  auto &ad = data_[ap.first];
570 
571  if (verbosity_ > 1)
572  os << " arm " << ap.first << std::endl;
573 
574  for (unsigned int nsi = 0; nsi < ad.n_sigmas.size(); ++nsi) {
575  if (verbosity_ > 1)
576  os << " nsi = " << nsi << std::endl;
577 
578  for (const auto &tri : ap.second.matched_track_idc[nsi]) {
579  const auto &track = hTracks->at(tri);
580 
581  bool some_proton_matching = false;
582 
583  if (verbosity_ > 1)
584  os << " tri = " << tri << std::endl;
585 
586  for (const auto &pri : ap.second.reco_proton_idc) {
587  const auto &proton = (*hRecoProtonsMultiRP)[pri];
588 
589  bool proton_matching = false;
590 
591  for (const auto &pr_tr : proton.contributingLocalTracks()) {
592  CTPPSDetId rpId(pr_tr->rpId());
593  unsigned int decRPId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
594 
595  if (decRPId != ad.rpId_N)
596  continue;
597 
598  const double x = pr_tr->x();
599  const double y = pr_tr->y();
600  const double th = 1E-3; // 1 um
601 
602  const bool match = (fabs(x - track.x()) < th) && (fabs(y - track.y()) < th);
603 
604  if (verbosity_ > 1)
605  os << " pri = " << pri << ": x_tr = " << track.x() << ", x_pr = " << x
606  << ", match = " << match << std::endl;
607 
608  if (match) {
609  proton_matching = true;
610  break;
611  }
612  }
613 
614  if (proton_matching) {
615  some_proton_matching = true;
616  break;
617  }
618  }
619 
620  if (verbosity_ > 1)
621  os << " --> some_proton_matching " << some_proton_matching << std::endl;
622 
623  if (some_proton_matching)
624  ap.second.matched_track_with_prot_idc[nsi].insert(tri);
625  else
626  ap.second.matched_track_without_prot_idc[nsi].insert(tri);
627  }
628  }
629  }
630 
631  // redo association cuts
632  for (auto &ap : armEventData) {
633  const auto &arm = ap.first;
634 
635  auto &ad = data_[arm];
636 
637  auto &aed = ap.second;
638 
639  auto &ass_cut = ppsAssociationCuts.getAssociationCuts(arm);
640 
641  const auto &indices = trackingSelection[arm];
642 
643  map<unsigned int, unsigned> matching_multiplicity;
644 
645  for (const auto &i : indices) {
646  for (const auto &j : indices) {
647  const auto &tr_i = hTracks->at(i);
648  const auto &tr_j = hTracks->at(j);
649 
650  const CTPPSDetId id_i(tr_i.rpId());
651  const CTPPSDetId id_j(tr_j.rpId());
652 
653  const unsigned rpDecId_i = id_i.arm() * 100 + id_i.station() * 10 + id_i.rp();
654  const unsigned rpDecId_j = id_j.arm() * 100 + id_j.station() * 10 + id_j.rp();
655 
656  if (rpDecId_i != ad.rpId_N || rpDecId_j != ad.rpId_F)
657  continue;
658 
659  ArmEventData::EvaluatedPair evp;
660  evp.idx_N = i;
661  evp.idx_F = j;
662 
663  evp.x_cut = ass_cut.isSatisfied(ass_cut.qX, tr_i.x(), tr_i.y(), lhcInfo.crossingAngle(), tr_i.x() - tr_j.x());
664  evp.y_cut = ass_cut.isSatisfied(ass_cut.qY, tr_i.x(), tr_i.y(), lhcInfo.crossingAngle(), tr_i.y() - tr_j.y());
665 
666  evp.match = evp.x_cut && evp.y_cut;
667 
668  aed.evaluatedPairs.push_back(evp);
669 
670  if (evp.match) {
671  matching_multiplicity[i]++;
672  matching_multiplicity[j]++;
673  }
674  }
675  }
676 
677  for (auto &evp : aed.evaluatedPairs)
678  evp.unique = (matching_multiplicity[evp.idx_N] == 1) && (matching_multiplicity[evp.idx_F] == 1);
679  }
680 
681  // debug print
682  if (verbosity_ > 1) {
683  for (auto &ap : armEventData) {
684  auto &ad = data_[ap.first];
685 
686  if (ad.de_x_sigma <= 0. && ad.de_y_sigma <= 0.)
687  continue;
688 
689  os << "* results for arm " << ap.first << std::endl;
690 
691  os << " reco_proton_idc: ";
692  for (const auto &idx : ap.second.reco_proton_idc)
693  os << idx << ", ";
694  os << std::endl;
695 
696  for (unsigned int nsi = 0; nsi < ad.n_sigmas.size(); ++nsi) {
697  os << " n_si = " << ad.n_sigmas[nsi] << std::endl;
698 
699  os << " matched_track_idc: ";
700  for (const auto &idx : ap.second.matched_track_idc[nsi])
701  os << idx << ", ";
702  os << std::endl;
703 
704  os << " matched_track_with_prot_idc: ";
705  for (const auto &idx : ap.second.matched_track_with_prot_idc[nsi])
706  os << idx << ", ";
707  os << std::endl;
708 
709  os << " matched_track_without_prot_idc: ";
710  for (const auto &idx : ap.second.matched_track_without_prot_idc[nsi])
711  os << idx << ", ";
712  os << std::endl;
713  }
714 
715  os << " evaluated pairs: ";
716  for (const auto &p : ap.second.evaluatedPairs)
717  os << "(" << p.idx_N << "-" << p.idx_F << ",M" << p.match << ",U" << p.unique << ")"
718  << ", ";
719  os << std::endl;
720  }
721  }
722 
723  // update efficiency plots
724  for (auto &ap : armEventData) {
725  const auto &arm = ap.first;
726  auto &ad = data_[arm];
727 
728  // stop if sigmas not yet determined
729  if (ad.de_x_sigma <= 0. && ad.de_y_sigma <= 0.)
730  continue;
731 
732  // loop over n_sigma choices
733  for (unsigned int nsi = 0; nsi < ad.n_sigmas.size(); ++nsi) {
734  const unsigned int n_exp_prot = ap.second.matched_track_idc[nsi].size();
735  const unsigned int n_rec_prot = ap.second.reco_proton_idc.size();
736 
737  // stop if N(expected protons) out of range
738  if (n_exp_prot < 1 || n_exp_prot > ad.n_exp_prot_max)
739  continue;
740 
741  // update method 1 plots
742  const double eff = double(n_rec_prot) / n_exp_prot;
743 
744  for (unsigned int tri : ap.second.matched_track_idc[nsi]) {
745  const double x_N = hTracks->at(tri).x();
746  const double xi_N = ad.s_x_to_xi_N->Eval(x_N * 1E-1); // conversion mm to cm
747 
748  ad.effPlots[0][nsi].p_eff1_vs_x_N->Fill(x_N, eff);
749  ad.effPlots[0][nsi].p_eff1_vs_xi_N->Fill(xi_N, eff);
750 
751  ad.effPlots[n_exp_prot][nsi].p_eff1_vs_x_N->Fill(x_N, eff);
752  ad.effPlots[n_exp_prot][nsi].p_eff1_vs_xi_N->Fill(xi_N, eff);
753  }
754 
755  // update method 2 plots
756  for (const auto &tri : ap.second.matched_track_with_prot_idc[nsi]) {
757  const double x_N = hTracks->at(tri).x();
758  const double xi_N = ad.s_x_to_xi_N->Eval(x_N * 1E-1); // conversion mm to cm
759 
760  ad.effPlots[0][nsi].p_eff2_vs_x_N->Fill(x_N, 1.);
761  ad.effPlots[0][nsi].p_eff2_vs_xi_N->Fill(xi_N, 1.);
762 
763  ad.effPlots[n_exp_prot][nsi].p_eff2_vs_x_N->Fill(x_N, 1.);
764  ad.effPlots[n_exp_prot][nsi].p_eff2_vs_xi_N->Fill(xi_N, 1.);
765  }
766 
767  for (const auto &tri : ap.second.matched_track_without_prot_idc[nsi]) {
768  const double x_N = hTracks->at(tri).x();
769  const double xi_N = ad.s_x_to_xi_N->Eval(x_N * 1E-1); // conversion mm to cm
770 
771  ad.effPlots[0][nsi].p_eff2_vs_x_N->Fill(x_N, 0.);
772  ad.effPlots[0][nsi].p_eff2_vs_xi_N->Fill(xi_N, 0.);
773 
774  ad.effPlots[n_exp_prot][nsi].p_eff2_vs_x_N->Fill(x_N, 0.);
775  ad.effPlots[n_exp_prot][nsi].p_eff2_vs_xi_N->Fill(xi_N, 0.);
776  }
777 
778  // update association cut plots
779  for (const auto &cand : ap.second.matched_track_idc[nsi]) // loop over candidates
780  {
781  const double x_N = hTracks->at(cand).x();
782 
783  auto &plots = ad.effPlots[n_exp_prot][nsi];
784 
785  bool hasMatchingAndUniquePair = false;
786  bool hasMatchingAndNonUniquePair = false;
787  for (const auto &pair : ap.second.evaluatedPairs) // loop over pairs
788  {
789  // stop if cand not compatible with pair
790  if (cand != pair.idx_N)
791  continue;
792 
793  if (pair.match && pair.unique)
794  hasMatchingAndUniquePair = true;
795 
796  if (pair.match && !pair.unique)
797  hasMatchingAndNonUniquePair = true;
798 
799  unsigned int signature = 999999;
800  if (!pair.x_cut && !pair.y_cut)
801  signature = 0;
802  if (pair.x_cut && !pair.y_cut)
803  signature = 1;
804  if (!pair.x_cut && pair.y_cut)
805  signature = 2;
806  if (pair.x_cut && pair.y_cut)
807  signature = 12;
808 
809  plots.p_fr_signature_0->Fill(x_N, (signature == 0) ? 1 : 0);
810  plots.p_fr_signature_1->Fill(x_N, (signature == 1) ? 1 : 0);
811  plots.p_fr_signature_2->Fill(x_N, (signature == 2) ? 1 : 0);
812  plots.p_fr_signature_12->Fill(x_N, (signature == 12) ? 1 : 0);
813  }
814 
815  plots.p_fr_match_unique->Fill(x_N, (hasMatchingAndUniquePair) ? 1 : 0);
816  plots.p_fr_match_non_unique->Fill(x_N, (hasMatchingAndNonUniquePair) ? 1 : 0);
817  }
818  }
819  }
820 
821  if (verbosity_)
822  edm::LogInfo("CTPPSProtonReconstructionEfficiencyEstimatorData") << os.str();
823 }
824 
825 //----------------------------------------------------------------------------------------------------
826 
828  auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
829 
830  for (const auto &ait : data_) {
831  char buf[100];
832  sprintf(buf, "arm %u", ait.first);
833  TDirectory *d_arm = f_out->mkdir(buf);
834  gDirectory = d_arm;
835 
836  ait.second.write();
837  }
838 }
839 
840 //----------------------------------------------------------------------------------------------------
841 
uint32_t station() const
Definition: CTPPSDetId.h:58
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:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
def unique
Definition: tier0.py:24
const std::vector< double > & getXiValues() const
Log< level::Error, false > LogError
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
void UpdateOptics(const LHCInterpolatedOpticalFunctionsSetCollection &ofc)
bool getData(T &iHolder) const
Definition: EventSetup.h:122
int iEvent
Definition: GenABIO.cc:224
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
tuple of
Definition: haddnano.py:40
edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsMultiRP_
uint32_t arm() const
Definition: CTPPSDetId.h:51
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Log< level::Info, false > LogInfo
std::vector< CTPPSLocalTrackLite > CTPPSLocalTrackLiteCollection
Collection of CTPPSLocalTrackLite objects.
edm::ESGetToken< LHCInterpolatedOpticalFunctionsSetCollection, CTPPSInterpolatedOpticsRcd > opticsESToken_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
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:57
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
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.
list indices
Definition: dqmdumpme.py:50
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 signature
Definition: Activities.doc:4
edm::ESGetToken< PPSAssociationCuts, PPSAssociationCutsRcd > ppsAssociationCutsToken_
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
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