CMS 3D CMS Logo

CTPPSProtonReconstructionPlotter.cc
Go to the documentation of this file.
1 /****************************************************************************
2 * Authors:
3 * Jan Kašpar (jan.kaspar@gmail.com)
4 ****************************************************************************/
5 
14 
16 
20 
23 
24 #include "TFile.h"
25 #include "TGraphErrors.h"
26 #include "TH1D.h"
27 #include "TH2D.h"
28 #include "TProfile.h"
29 
30 //----------------------------------------------------------------------------------------------------
31 
33 public:
36 
37 private:
38  void analyze(const edm::Event &, const edm::EventSetup &) override;
39 
40  void endJob() override;
41 
45 
46  unsigned int rpId_45_N_, rpId_45_F_;
47  unsigned int rpId_56_N_, rpId_56_F_;
48 
49  struct AssociationCuts {
50  double ti_tr_min;
51  double ti_tr_max;
52 
53  void load(const edm::ParameterSet &ps) {
54  ti_tr_min = ps.getParameter<double>("ti_tr_min");
55  ti_tr_max = ps.getParameter<double>("ti_tr_max");
56  }
57  };
58 
59  std::map<unsigned int, AssociationCuts> association_cuts_;
60 
62 
63  signed int maxNonEmptyEvents_;
64 
65  static void profileToRMSGraph(TProfile *p, TGraphErrors *g) {
66  for (int bi = 1; bi <= p->GetNbinsX(); ++bi) {
67  double c = p->GetBinCenter(bi);
68 
69  double N = p->GetBinEntries(bi);
70  double Sy = p->GetBinContent(bi) * N;
71  double Syy = p->GetSumw2()->At(bi);
72 
73  double si_sq = Syy / N - Sy * Sy / N / N;
74  double si = (si_sq >= 0.) ? sqrt(si_sq) : 0.;
75  double si_unc_sq = si_sq / 2. / N; // Gaussian approximation
76  double si_unc = (si_unc_sq >= 0.) ? sqrt(si_unc_sq) : 0.;
77 
78  int idx = g->GetN();
79  g->SetPoint(idx, c, si);
80  g->SetPointError(idx, 0., si_unc);
81  }
82  }
83 
85  const CTPPSLocalTrackLite &tr,
86  const CTPPSGeometry &geometry,
87  double &de_x,
88  double &de_x_unc);
89 
90  struct SingleRPPlots {
91  std::unique_ptr<TH1D> h_multiplicity;
92  std::unique_ptr<TH1D> h_xi;
93  std::unique_ptr<TH2D> h2_th_y_vs_xi;
94  std::unique_ptr<TProfile> p_th_y_vs_xi;
95 
96  std::map<unsigned int, TH1D *> m_h_xi_nTracks;
97  std::unique_ptr<TH1D> h_xi_n1f1;
98 
100  : h_multiplicity(new TH1D("", ";reconstructed protons", 11, -0.5, 10.5)),
101  h_xi(new TH1D("", ";#xi", 100, 0., 0.3)),
102  h2_th_y_vs_xi(new TH2D("", ";#xi;#theta_{y} (rad)", 100, 0., 0.3, 100, -500E-6, +500E-6)),
103  p_th_y_vs_xi(new TProfile("", ";#xi;#theta_{y} (rad)", 100, 0., 0.3)),
104  h_xi_n1f1(new TH1D("", ";#xi", 100, 0., 0.3)) {
105  for (unsigned int n = 2; n <= 10; ++n)
106  m_h_xi_nTracks[n] = new TH1D(*h_xi);
107  }
108 
109  void fill(const reco::ForwardProton &p, unsigned int nTracks, bool n1f1) {
110  if (p.validFit()) {
111  h_xi->Fill(p.xi());
112 
113  const double th_y = p.thetaY();
114  h2_th_y_vs_xi->Fill(p.xi(), th_y);
115  p_th_y_vs_xi->Fill(p.xi(), th_y);
116 
117  auto it = m_h_xi_nTracks.find(nTracks);
118  if (it != m_h_xi_nTracks.end())
119  it->second->Fill(p.xi());
120 
121  if (n1f1)
122  h_xi_n1f1->Fill(p.xi());
123  }
124  }
125 
126  void write() const {
127  h_multiplicity->Write("h_multiplicity");
128  h_xi->Write("h_xi");
129 
130  h2_th_y_vs_xi->Write("h2_th_y_vs_xi");
131  p_th_y_vs_xi->Write("p_th_y_vs_xi");
132 
133  TDirectory *d_top = gDirectory;
134 
135  gDirectory = d_top->mkdir("h_xi_nTracks");
136  for (const auto p : m_h_xi_nTracks) {
137  char buf[100];
138  sprintf(buf, "h_xi_nTracks_%u", p.first);
139  p.second->Write(buf);
140  }
141 
142  gDirectory = d_top;
143 
144  h_xi_n1f1->Write("h_xi_n1f1");
145  }
146  };
147 
148  std::map<unsigned int, SingleRPPlots> singleRPPlots_;
149 
150  struct MultiRPPlots {
151  std::unique_ptr<TH1D> h_multiplicity;
152  std::unique_ptr<TH1D> h_xi, h_th_x, h_th_y, h_vtx_y, h_t_unif, h_t, h_chi_sq, h_log_chi_sq, h_chi_sq_norm;
153  std::unique_ptr<TH1D> h_t_xi_range1, h_t_xi_range2, h_t_xi_range3;
154  std::unique_ptr<TH1D> h_time;
155  std::unique_ptr<TH1D> h_n_contrib_tracking_tracks, h_n_contrib_timing_tracks;
156  std::unique_ptr<TH2D> h2_th_x_vs_xi, h2_th_y_vs_xi, h2_vtx_y_vs_xi, h2_t_vs_xi;
157  std::unique_ptr<TProfile> p_th_x_vs_xi, p_th_y_vs_xi, p_vtx_y_vs_xi;
158 
159  std::unique_ptr<TH2D> h2_timing_tracks_vs_prot_mult;
160 
161  std::map<unsigned int, TH1D *> m_h_xi_nTracks;
162  std::unique_ptr<TH1D> h_xi_n1f1;
163 
164  std::unique_ptr<TH1D> h_de_x_timing_vs_tracking, h_de_x_rel_timing_vs_tracking, h_de_x_match_timing_vs_tracking;
165  std::unique_ptr<TH1D> h_de_x_timing_vs_tracking_ClCo, h_de_x_rel_timing_vs_tracking_ClCo,
166  h_de_x_match_timing_vs_tracking_ClCo;
167 
168  std::unique_ptr<TH2D> h2_y_vs_x_tt0_ClCo, h2_y_vs_x_tt1_ClCo, h2_y_vs_x_ttm_ClCo;
169 
171  : h_multiplicity(new TH1D("", ";reconstructed protons per event", 11, -0.5, 10.5)),
172  h_xi(new TH1D("", ";#xi", 100, 0., 0.3)),
173  h_th_x(new TH1D("", ";#theta_{x} (rad)", 250, -500E-6, +500E-6)),
174  h_th_y(new TH1D("", ";#theta_{y} (rad)", 500, -1000E-6, +1000E-6)),
175  h_vtx_y(new TH1D("", ";vtx_{y} (cm)", 100, -100E-3, +100E-3)),
176  h_chi_sq(new TH1D("", ";#chi^{2}", 100, 0., 10.)),
177  h_log_chi_sq(new TH1D("", ";log_{10} #chi^{2}", 100, -20., 5.)),
178  h_chi_sq_norm(new TH1D("", ";#chi^{2}/ndf", 100, 0., 5.)),
179  h_time(new TH1D("", ";time", 100, -2., +2.)),
180  h_n_contrib_tracking_tracks(new TH1D("", ";n of contrib. tracking tracks per reco proton", 4, -0.5, +3.5)),
181  h_n_contrib_timing_tracks(new TH1D("", ";n of contrib. timing tracks per reco proton", 4, -0.5, +3.5)),
182  h2_th_x_vs_xi(new TH2D("", ";#xi;#theta_{x} (rad)", 100, 0., 0.3, 100, -500E-6, +500E-6)),
183  h2_th_y_vs_xi(new TH2D("", ";#xi;#theta_{y} (rad)", 100, 0., 0.3, 200, -1000E-6, +1000E-6)),
184  h2_vtx_y_vs_xi(new TH2D("", ";#xi;vtx_{y} (cm)", 100, 0., 0.3, 100, -100E-3, +100E-3)),
185  p_th_x_vs_xi(new TProfile("", ";#xi;#theta_{x} (rad)", 100, 0., 0.3)),
186  p_th_y_vs_xi(new TProfile("", ";#xi;#theta_{y} (rad)", 100, 0., 0.3)),
187  p_vtx_y_vs_xi(new TProfile("", ";#xi;vtx_{y} (cm)", 100, 0., 0.3)),
188  h2_timing_tracks_vs_prot_mult(
189  new TH2D("", ";reco protons per event;timing tracks per event", 11, -0.5, 10.5, 11, -0.5, 10.5)),
190  h_xi_n1f1(new TH1D("", ";#xi", 100, 0., 0.3)),
191 
192  h_de_x_timing_vs_tracking(new TH1D("", ";#Delta x (mm)", 200, -1., +1.)),
193  h_de_x_rel_timing_vs_tracking(new TH1D("", ";#Delta x / #sigma(x)", 200, -20., +20.)),
194  h_de_x_match_timing_vs_tracking(new TH1D("", ";match between tracking and timing tracks", 2, -0.5, +1.5)),
195 
196  h_de_x_timing_vs_tracking_ClCo(new TH1D("", ";#Delta x (mm)", 200, -1., +1.)),
197  h_de_x_rel_timing_vs_tracking_ClCo(new TH1D("", ";#Delta x / #sigma(x)", 200, -20., +20.)),
198  h_de_x_match_timing_vs_tracking_ClCo(
199  new TH1D("", ";match between tracking and timing tracks", 2, -0.5, +1.5)),
200 
201  h2_y_vs_x_tt0_ClCo(new TH2D("", ";x (mm);y (mm)", 100, -5., 25., 100, -15., +15.)),
202  h2_y_vs_x_tt1_ClCo(new TH2D("", ";x (mm);y (mm)", 100, -5., 25., 100, -15., +15.)),
203  h2_y_vs_x_ttm_ClCo(new TH2D("", ";x (mm);y (mm)", 100, -5., 25., 100, -15., +15.)) {
204  std::vector<double> v_t_bin_edges;
205  for (double t = 0; t <= 5.;) {
206  v_t_bin_edges.push_back(t);
207  const double de_t = 0.05 + 0.09 * t + 0.02 * t * t;
208  t += de_t;
209  }
210  h_t_unif.reset(new TH1D("", ";|t| (GeV^2)", 100, 0., 5.));
211  h_t.reset(new TH1D("", ";|t| (GeV^2)", v_t_bin_edges.size() - 1, v_t_bin_edges.data()));
212  h_t_xi_range1.reset(new TH1D("", ";|t| (GeV^2)", v_t_bin_edges.size() - 1, v_t_bin_edges.data()));
213  h_t_xi_range2.reset(new TH1D("", ";|t| (GeV^2)", v_t_bin_edges.size() - 1, v_t_bin_edges.data()));
214  h_t_xi_range3.reset(new TH1D("", ";|t| (GeV^2)", v_t_bin_edges.size() - 1, v_t_bin_edges.data()));
215  h2_t_vs_xi.reset(
216  new TH2D("", ";#xi;|t| (GeV^2)", 100, 0., 0.3, v_t_bin_edges.size() - 1, v_t_bin_edges.data()));
217 
218  for (unsigned int n = 2; n <= 10; ++n)
219  m_h_xi_nTracks[n] = new TH1D(*h_xi);
220  }
221 
222  void fill(const reco::ForwardProton &p, unsigned int nTracks, bool n1f1) {
223  if (!p.validFit())
224  return;
225 
226  unsigned int n_contrib_tracking_tracks = 0, n_contrib_timing_tracks = 0;
227  for (const auto &tr : p.contributingLocalTracks()) {
228  CTPPSDetId detId(tr->getRPId());
229  if (detId.subdetId() == CTPPSDetId::sdTrackingStrip || detId.subdetId() == CTPPSDetId::sdTrackingPixel)
230  n_contrib_tracking_tracks++;
231  if (detId.subdetId() == CTPPSDetId::sdTimingDiamond || detId.subdetId() == CTPPSDetId::sdTimingFastSilicon)
232  n_contrib_timing_tracks++;
233  }
234 
235  const double th_x = p.thetaX();
236  const double th_y = p.thetaY();
237  const double mt = -p.t();
238 
239  h_chi_sq->Fill(p.chi2());
240  h_log_chi_sq->Fill(log10(p.chi2()));
241  if (p.ndof() > 0)
242  h_chi_sq_norm->Fill(p.normalizedChi2());
243 
244  h_n_contrib_tracking_tracks->Fill(n_contrib_tracking_tracks);
245  h_n_contrib_timing_tracks->Fill(n_contrib_timing_tracks);
246 
247  h_xi->Fill(p.xi());
248 
249  h_th_x->Fill(th_x);
250  h_th_y->Fill(th_y);
251 
252  h_vtx_y->Fill(p.vertex().y());
253 
254  h_t_unif->Fill(mt);
255  h_t->Fill(mt);
256  if (p.xi() > 0.04 && p.xi() < 0.07)
257  h_t_xi_range1->Fill(mt);
258  if (p.xi() > 0.07 && p.xi() < 0.10)
259  h_t_xi_range2->Fill(mt);
260  if (p.xi() > 0.10 && p.xi() < 0.13)
261  h_t_xi_range3->Fill(mt);
262 
263  h_time->Fill(p.time());
264 
265  h2_th_x_vs_xi->Fill(p.xi(), th_x);
266  h2_th_y_vs_xi->Fill(p.xi(), th_y);
267  h2_vtx_y_vs_xi->Fill(p.xi(), p.vertex().y());
268  h2_t_vs_xi->Fill(p.xi(), mt);
269 
270  p_th_x_vs_xi->Fill(p.xi(), th_x);
271  p_th_y_vs_xi->Fill(p.xi(), th_y);
272  p_vtx_y_vs_xi->Fill(p.xi(), p.vertex().y());
273 
274  auto it = m_h_xi_nTracks.find(nTracks);
275  if (it != m_h_xi_nTracks.end())
276  it->second->Fill(p.xi());
277 
278  if (n1f1)
279  h_xi_n1f1->Fill(p.xi());
280  }
281 
282  void write() const {
283  h_multiplicity->Write("h_multiplicity");
284 
285  h_chi_sq->Write("h_chi_sq");
286  h_log_chi_sq->Write("h_log_chi_sq");
287  h_chi_sq_norm->Write("h_chi_sq_norm");
288 
289  h_n_contrib_tracking_tracks->Write("h_n_contrib_tracking_tracks");
290  h_n_contrib_timing_tracks->Write("h_n_contrib_timing_tracks");
291 
292  h2_timing_tracks_vs_prot_mult->Write("h2_timing_tracks_vs_prot_mult");
293 
294  h_xi->Write("h_xi");
295 
296  h_th_x->Write("h_th_x");
297  h2_th_x_vs_xi->Write("h2_th_x_vs_xi");
298  p_th_x_vs_xi->Write("p_th_x_vs_xi");
299  auto g_th_x_RMS_vs_xi = std::make_unique<TGraphErrors>();
300  profileToRMSGraph(p_th_x_vs_xi.get(), g_th_x_RMS_vs_xi.get());
301  g_th_x_RMS_vs_xi->Write("g_th_x_RMS_vs_xi");
302 
303  h_th_y->Write("h_th_y");
304  h2_th_y_vs_xi->Write("h2_th_y_vs_xi");
305  p_th_y_vs_xi->Write("p_th_y_vs_xi");
306  auto g_th_y_RMS_vs_xi = std::make_unique<TGraphErrors>();
307  profileToRMSGraph(p_th_y_vs_xi.get(), g_th_y_RMS_vs_xi.get());
308  g_th_y_RMS_vs_xi->Write("g_th_y_RMS_vs_xi");
309 
310  h_vtx_y->Write("h_vtx_y");
311  h2_vtx_y_vs_xi->Write("h2_vtx_y_vs_xi");
312  p_vtx_y_vs_xi->Write("p_vtx_y_vs_xi");
313  auto g_vtx_y_RMS_vs_xi = std::make_unique<TGraphErrors>();
314  profileToRMSGraph(p_vtx_y_vs_xi.get(), g_vtx_y_RMS_vs_xi.get());
315  g_vtx_y_RMS_vs_xi->Write("g_vtx_y_RMS_vs_xi");
316 
317  h_t->Scale(1., "width");
318 
319  h_t_unif->Write("h_t_unif");
320  h_t->Write("h_t");
321  h_t_xi_range1->Write("h_t_xi_range1");
322  h_t_xi_range2->Write("h_t_xi_range2");
323  h_t_xi_range3->Write("h_t_xi_range3");
324 
325  h2_t_vs_xi->Write("h2_t_vs_xi");
326 
327  h_time->Write("h_time");
328 
329  TDirectory *d_top = gDirectory;
330 
331  gDirectory = d_top->mkdir("h_xi_nTracks");
332  for (const auto p : m_h_xi_nTracks) {
333  char buf[100];
334  sprintf(buf, "h_xi_nTracks_%u", p.first);
335  p.second->Write(buf);
336  }
337 
338  gDirectory = d_top;
339 
340  h_xi_n1f1->Write("h_xi_n1f1");
341 
342  h_de_x_timing_vs_tracking->Write("h_de_x_timing_vs_tracking");
343  h_de_x_rel_timing_vs_tracking->Write("h_de_x_rel_timing_vs_tracking");
344  h_de_x_match_timing_vs_tracking->Write("h_de_x_match_timing_vs_tracking");
345 
346  h_de_x_timing_vs_tracking_ClCo->Write("h_de_x_timing_vs_tracking_ClCo");
347  h_de_x_rel_timing_vs_tracking_ClCo->Write("h_de_x_rel_timing_vs_tracking_ClCo");
348  h_de_x_match_timing_vs_tracking_ClCo->Write("h_de_x_match_timing_vs_tracking_ClCo");
349 
350  h2_y_vs_x_tt0_ClCo->Write("h2_y_vs_x_tt0_ClCo");
351  h2_y_vs_x_tt1_ClCo->Write("h2_y_vs_x_tt1_ClCo");
352  h2_y_vs_x_ttm_ClCo->Write("h2_y_vs_x_ttm_ClCo");
353  }
354  };
355 
356  std::map<unsigned int, MultiRPPlots> multiRPPlots_;
357 
359  std::unique_ptr<TH2D> h2_xi_mu_vs_xi_si;
360  std::unique_ptr<TH1D> h_xi_diff_mu_si, h_xi_diff_si_mu;
361 
362  std::unique_ptr<TH2D> h2_xi_diff_si_mu_vs_xi_mu;
363  std::unique_ptr<TProfile> p_xi_diff_si_mu_vs_xi_mu;
364 
365  std::unique_ptr<TH2D> h2_th_y_mu_vs_th_y_si;
366 
368  : h2_xi_mu_vs_xi_si(new TH2D("", ";#xi_{single};#xi_{multi}", 100, 0., 0.3, 100, 0., 0.3)),
369  h_xi_diff_mu_si(new TH1D("", ";#xi_{multi} - #xi_{single}", 100, -0.1, +0.1)),
370  h_xi_diff_si_mu(new TH1D("", ";#xi_{single} - #xi_{multi}", 100, -0.1, +0.1)),
371  h2_xi_diff_si_mu_vs_xi_mu(
372  new TH2D("", ";#xi_{multi};#xi_{single} - #xi_{multi}", 100, 0., 0.3, 100, -0.10, 0.10)),
373  p_xi_diff_si_mu_vs_xi_mu(new TProfile("", ";#xi_{multi};#xi_{single} - #xi_{multi}", 100, 0., 0.3)),
374  h2_th_y_mu_vs_th_y_si(
375  new TH2D("", ";#theta^{*}_{y,si};#theta^{*}_{y,mu}", 100, -500E-6, +500E-6, 100, -500E-6, +500E-6)) {}
376 
377  void fill(const reco::ForwardProton &p_single, const reco::ForwardProton &p_multi) {
378  if (p_single.validFit() && p_multi.validFit()) {
379  h2_xi_mu_vs_xi_si->Fill(p_single.xi(), p_multi.xi());
380  h_xi_diff_mu_si->Fill(p_multi.xi() - p_single.xi());
381  h_xi_diff_si_mu->Fill(p_single.xi() - p_multi.xi());
382 
383  h2_xi_diff_si_mu_vs_xi_mu->Fill(p_multi.xi(), p_single.xi() - p_multi.xi());
384  p_xi_diff_si_mu_vs_xi_mu->Fill(p_multi.xi(), p_single.xi() - p_multi.xi());
385 
386  const double th_y_si = p_single.thetaY();
387  const double th_y_mu = p_multi.thetaY();
388 
389  h2_th_y_mu_vs_th_y_si->Fill(th_y_si, th_y_mu);
390  }
391  }
392 
393  void write() const {
394  h2_xi_mu_vs_xi_si->Write("h2_xi_mu_vs_xi_si");
395  h_xi_diff_mu_si->Write("h_xi_diff_mu_si");
396  h_xi_diff_si_mu->Write("h_xi_diff_si_mu");
397 
398  h2_xi_diff_si_mu_vs_xi_mu->Write("h2_xi_diff_si_mu_vs_xi_mu");
399  p_xi_diff_si_mu_vs_xi_mu->Write("p_xi_diff_si_mu_vs_xi_mu");
400 
401  h2_th_y_mu_vs_th_y_si->Write("h2_th_y_mu_vs_th_y_si");
402  }
403  };
404 
405  std::map<unsigned int, SingleMultiCorrelationPlots> singleMultiCorrelationPlots_;
406 
408  std::unique_ptr<TH1D> h_xi_si_diffNF;
409  std::unique_ptr<TH2D> h2_xi_si_diffNF_vs_xi_mu;
410  std::unique_ptr<TProfile> p_xi_si_diffNF_vs_xi_mu;
411 
412  std::unique_ptr<TH1D> h_th_y_si_diffNF;
413  std::unique_ptr<TProfile> p_th_y_si_diffNF_vs_xi_mu;
414 
416  : h_xi_si_diffNF(new TH1D("", ";#xi_{sF} - #xi_{sN}", 100, -0.02, +0.02)),
417  h2_xi_si_diffNF_vs_xi_mu(new TH2D("", ";#xi_{m};#xi_{sF} - #xi_{sN}", 100, 0., 0.3, 100, -0.02, +0.02)),
418  p_xi_si_diffNF_vs_xi_mu(new TProfile("", ";#xi_{m};#xi_{sF} - #xi_{sN}", 100, 0., 0.3)),
419  h_th_y_si_diffNF(new TH1D("", ";#theta_{y,sF} - #theta_{y,sN}", 100, -100E-6, +100E-6)),
420  p_th_y_si_diffNF_vs_xi_mu(new TProfile("", ";#xi_{m};#theta_{y,sF} - #theta_{y,sN}", 100, 0., 0.3)) {}
421 
422  void fill(const reco::ForwardProton &p_s_N, const reco::ForwardProton &p_s_F, const reco::ForwardProton &p_m) {
423  if (!p_s_N.validFit() || !p_s_F.validFit() || !p_m.validFit())
424  return;
425 
426  const double th_y_s_N = p_s_N.thetaY();
427  const double th_y_s_F = p_s_F.thetaY();
428 
429  h_xi_si_diffNF->Fill(p_s_F.xi() - p_s_N.xi());
430  h2_xi_si_diffNF_vs_xi_mu->Fill(p_m.xi(), p_s_F.xi() - p_s_N.xi());
431  p_xi_si_diffNF_vs_xi_mu->Fill(p_m.xi(), p_s_F.xi() - p_s_N.xi());
432 
433  h_th_y_si_diffNF->Fill(th_y_s_F - th_y_s_N);
434  p_th_y_si_diffNF_vs_xi_mu->Fill(p_m.xi(), th_y_s_F - th_y_s_N);
435  }
436 
437  void write() const {
438  h_xi_si_diffNF->Write("h_xi_si_diffNF");
439  h2_xi_si_diffNF_vs_xi_mu->Write("h2_xi_si_diffNF_vs_xi_mu");
440  p_xi_si_diffNF_vs_xi_mu->Write("p_xi_si_diffNF_vs_xi_mu");
441 
442  h_th_y_si_diffNF->Write("h_th_y_si_diffNF");
443  p_th_y_si_diffNF_vs_xi_mu->Write("p_th_y_si_diffNF_vs_xi_mu");
444  }
445  };
446 
447  std::map<unsigned int, ArmCorrelationPlots> armCorrelationPlots_;
448 
449  std::unique_ptr<TProfile> p_x_L_diffNF_vs_x_L_N_, p_x_R_diffNF_vs_x_R_N_;
450  std::unique_ptr<TProfile> p_y_L_diffNF_vs_y_L_N_, p_y_R_diffNF_vs_y_R_N_;
451 
453 };
454 
455 //----------------------------------------------------------------------------------------------------
456 
457 using namespace std;
458 using namespace edm;
459 
460 //----------------------------------------------------------------------------------------------------
461 
463  : tokenTracks_(consumes<CTPPSLocalTrackLiteCollection>(ps.getParameter<edm::InputTag>("tagTracks"))),
465  consumes<reco::ForwardProtonCollection>(ps.getParameter<InputTag>("tagRecoProtonsSingleRP"))),
467  consumes<reco::ForwardProtonCollection>(ps.getParameter<InputTag>("tagRecoProtonsMultiRP"))),
468 
469  rpId_45_N_(ps.getParameter<unsigned int>("rpId_45_N")),
470  rpId_45_F_(ps.getParameter<unsigned int>("rpId_45_F")),
471  rpId_56_N_(ps.getParameter<unsigned int>("rpId_56_N")),
472  rpId_56_F_(ps.getParameter<unsigned int>("rpId_56_F")),
473 
474  outputFile_(ps.getParameter<string>("outputFile")),
475  maxNonEmptyEvents_(ps.getUntrackedParameter<signed int>("maxNonEmptyEvents", -1)),
476 
477  p_x_L_diffNF_vs_x_L_N_(new TProfile("p_x_L_diffNF_vs_x_L_N", ";x_{LN};x_{LF} - x_{LN}", 100, 0., +20.)),
478  p_x_R_diffNF_vs_x_R_N_(new TProfile("p_x_R_diffNF_vs_x_R_N", ";x_{RN};x_{RF} - x_{RN}", 100, 0., +20.)),
479 
480  p_y_L_diffNF_vs_y_L_N_(new TProfile("p_y_L_diffNF_vs_y_L_N", ";y_{LN};y_{LF} - y_{LN}", 100, -20., +20.)),
481  p_y_R_diffNF_vs_y_R_N_(new TProfile("p_y_R_diffNF_vs_y_R_N", ";y_{RN};y_{RF} - y_{RN}", 100, -20., +20.)),
482 
484  for (const std::string &sector : {"45", "56"}) {
485  const unsigned int arm = (sector == "45") ? 0 : 1;
486  association_cuts_[arm].load(ps.getParameterSet("association_cuts_" + sector));
487  }
488 }
489 
490 //----------------------------------------------------------------------------------------------------
491 
493  const CTPPSLocalTrackLite &tr_ti,
494  const CTPPSGeometry &geometry,
495  double &de_x,
496  double &de_x_unc) {
497  // identify tracking-RP tracks
498  const auto &tr_i = *proton.contributingLocalTracks().at(0);
499  const auto &tr_j = *proton.contributingLocalTracks().at(1);
500 
501  const double z_i = geometry.getRPTranslation(tr_i.getRPId()).z();
502  const double z_j = geometry.getRPTranslation(tr_j.getRPId()).z();
503 
504  // interpolation from tracking RPs
505  const double z_ti =
506  -geometry.getRPTranslation(tr_ti.getRPId()).z(); // the minus sign fixes a bug in the diamond geometry
507  const double f_i = (z_ti - z_j) / (z_i - z_j), f_j = (z_i - z_ti) / (z_i - z_j);
508  const double x_inter = f_i * tr_i.getX() + f_j * tr_j.getX();
509  const double x_inter_unc_sq =
510  f_i * f_i * tr_i.getXUnc() * tr_i.getXUnc() + f_j * f_j * tr_j.getXUnc() * tr_j.getXUnc();
511 
512  // save distance
513  de_x = tr_ti.getX() - x_inter;
514  de_x_unc = sqrt(tr_ti.getXUnc() * tr_ti.getXUnc() + x_inter_unc_sq);
515 }
516 
517 //----------------------------------------------------------------------------------------------------
518 
520  // get input
522  event.getByToken(tokenTracks_, hTracks);
523 
524  Handle<reco::ForwardProtonCollection> hRecoProtonsSingleRP;
525  event.getByToken(tokenRecoProtonsSingleRP_, hRecoProtonsSingleRP);
526 
527  Handle<reco::ForwardProtonCollection> hRecoProtonsMultiRP;
528  event.getByToken(tokenRecoProtonsMultiRP_, hRecoProtonsMultiRP);
529 
530  if (!hRecoProtonsSingleRP->empty())
532 
534  throw cms::Exception("CTPPSProtonReconstructionPlotter") << "Number of non empty events reached maximum.";
535 
536  // get conditions
538  iSetup.get<VeryForwardRealGeometryRecord>().get(hGeometry);
539 
540  // track plots
541  const CTPPSLocalTrackLite *tr_L_N = nullptr;
542  const CTPPSLocalTrackLite *tr_L_F = nullptr;
543  const CTPPSLocalTrackLite *tr_R_N = nullptr;
544  const CTPPSLocalTrackLite *tr_R_F = nullptr;
545  std::map<unsigned int, unsigned int> armTrackCounter, armTimingTrackCounter;
546  std::map<unsigned int, unsigned int> armTrackCounter_N, armTrackCounter_F;
547 
548  for (const auto &tr : *hTracks) {
549  CTPPSDetId rpId(tr.getRPId());
550  unsigned int decRPId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
551 
552  if (decRPId == rpId_45_N_) {
553  tr_L_N = &tr;
554  armTrackCounter_N[0]++;
555  }
556  if (decRPId == rpId_45_F_) {
557  tr_L_F = &tr;
558  armTrackCounter_F[0]++;
559  }
560  if (decRPId == rpId_56_N_) {
561  tr_R_N = &tr;
562  armTrackCounter_N[1]++;
563  }
564  if (decRPId == rpId_56_F_) {
565  tr_R_F = &tr;
566  armTrackCounter_F[1]++;
567  }
568 
569  armTrackCounter[rpId.arm()]++;
570 
571  const bool trackerRP =
572  (rpId.subdetId() == CTPPSDetId::sdTrackingStrip || rpId.subdetId() == CTPPSDetId::sdTrackingPixel);
573  if (!trackerRP)
574  armTimingTrackCounter[rpId.arm()]++;
575  }
576 
577  if (tr_L_N && tr_L_F) {
578  p_x_L_diffNF_vs_x_L_N_->Fill(tr_L_N->getX(), tr_L_F->getX() - tr_L_N->getX());
579  p_y_L_diffNF_vs_y_L_N_->Fill(tr_L_N->getY(), tr_L_F->getY() - tr_L_N->getY());
580  }
581 
582  if (tr_R_N && tr_R_F) {
583  p_x_R_diffNF_vs_x_R_N_->Fill(tr_R_N->getX(), tr_R_F->getX() - tr_R_N->getX());
584  p_y_R_diffNF_vs_y_R_N_->Fill(tr_R_N->getY(), tr_R_F->getY() - tr_R_N->getY());
585  }
586 
587  // initialise multiplicity counters
588  std::map<unsigned int, unsigned int> singleRPMultiplicity, multiRPMultiplicity;
589  singleRPMultiplicity[rpId_45_N_] = singleRPMultiplicity[rpId_45_F_] = singleRPMultiplicity[rpId_56_N_] =
590  singleRPMultiplicity[rpId_56_F_] = 0;
591  multiRPMultiplicity[0] = multiRPMultiplicity[1] = 0;
592 
593  // make single-RP-reco plots
594  for (const auto &proton : *hRecoProtonsSingleRP) {
595  CTPPSDetId rpId((*proton.contributingLocalTracks().begin())->getRPId());
596  unsigned int decRPId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
597 
598  const bool n1f1 = (armTrackCounter_N[rpId.arm()] == 1 && armTrackCounter_F[rpId.arm()] == 1);
599 
600  singleRPPlots_[decRPId].fill(proton, armTrackCounter[rpId.arm()], n1f1);
601 
602  if (proton.validFit())
603  singleRPMultiplicity[decRPId]++;
604  }
605 
606  for (const auto it : singleRPMultiplicity)
607  singleRPPlots_[it.first].h_multiplicity->Fill(it.second);
608 
609  // make multi-RP-reco plots
610  for (const auto &proton : *hRecoProtonsMultiRP) {
611  CTPPSDetId rpId((*proton.contributingLocalTracks().begin())->getRPId());
612  unsigned int armId = rpId.arm();
613 
614  const bool n1f1 = (armTrackCounter_N[armId] == 1 && armTrackCounter_F[armId] == 1);
615 
616  multiRPPlots_[armId].fill(proton, armTrackCounter[armId], n1f1);
617 
618  if (proton.validFit())
619  multiRPMultiplicity[armId]++;
620  }
621 
622  for (const auto it : multiRPMultiplicity) {
623  const auto &pl = multiRPPlots_[it.first];
624  pl.h_multiplicity->Fill(it.second);
625  pl.h2_timing_tracks_vs_prot_mult->Fill(it.second, armTimingTrackCounter[it.first]);
626  }
627 
628  // define "clean condition" for each arm
629  map<unsigned int, bool> clCo;
630  clCo[0] = (singleRPMultiplicity[rpId_45_N_] && singleRPMultiplicity[rpId_45_F_] && multiRPMultiplicity[0] == 1);
631  clCo[1] = (singleRPMultiplicity[rpId_56_N_] && singleRPMultiplicity[rpId_56_F_] && multiRPMultiplicity[1] == 1);
632 
633  // plot distances between multi-RP protons and timing tracks in the same arm
634  for (const auto &proton : *hRecoProtonsMultiRP) {
635  if (!proton.validFit())
636  continue;
637 
638  CTPPSDetId rpId_proton((*proton.contributingLocalTracks().begin())->getRPId());
639  unsigned int armId = rpId_proton.arm();
640  const auto &pl = multiRPPlots_[armId];
641 
642  for (const auto &tr : *hTracks) {
643  CTPPSDetId rpId_tr(tr.getRPId());
644  if (rpId_tr.arm() != armId)
645  continue;
646 
647  const bool trackerRP =
648  (rpId_tr.subdetId() == CTPPSDetId::sdTrackingStrip || rpId_tr.subdetId() == CTPPSDetId::sdTrackingPixel);
649  if (trackerRP)
650  continue;
651 
652  double de_x = 0., de_x_unc = 0.;
653  CalculateTimingTrackingDistance(proton, tr, *hGeometry, de_x, de_x_unc);
654 
655  const double rd = (de_x_unc > 0.) ? de_x / de_x_unc : -1E10;
656  const auto &ac = association_cuts_[armId];
657  const bool match = (ac.ti_tr_min <= fabs(rd) && fabs(rd) <= ac.ti_tr_max);
658 
659  pl.h_de_x_timing_vs_tracking->Fill(de_x);
660  pl.h_de_x_rel_timing_vs_tracking->Fill(rd);
661  pl.h_de_x_match_timing_vs_tracking->Fill(match ? 1. : 0.);
662 
663  if (clCo[armId] && armTimingTrackCounter[armId] == 1) {
664  pl.h_de_x_timing_vs_tracking_ClCo->Fill(de_x);
665  pl.h_de_x_rel_timing_vs_tracking_ClCo->Fill(rd);
666  pl.h_de_x_match_timing_vs_tracking_ClCo->Fill(match ? 1. : 0.);
667  }
668  }
669  }
670 
671  // plot xy maps
672  for (const auto &proton : *hRecoProtonsMultiRP) {
673  if (!proton.validFit())
674  continue;
675 
676  CTPPSDetId rpId((*proton.contributingLocalTracks().begin())->getRPId());
677  unsigned int armId = rpId.arm();
678  const auto &pl = multiRPPlots_[armId];
679  const auto &nTimingTracks = armTimingTrackCounter[armId];
680 
681  if (!clCo[armId])
682  continue;
683 
684  double x_ref = 0., y_ref = 0.;
685  if (armId == 0) {
686  x_ref = tr_L_N->getX();
687  y_ref = tr_L_N->getY();
688  }
689  if (armId == 1) {
690  x_ref = tr_R_N->getX();
691  y_ref = tr_R_N->getY();
692  }
693 
694  if (nTimingTracks == 0)
695  pl.h2_y_vs_x_tt0_ClCo->Fill(x_ref, y_ref);
696  if (nTimingTracks == 1)
697  pl.h2_y_vs_x_tt1_ClCo->Fill(x_ref, y_ref);
698  if (nTimingTracks > 1)
699  pl.h2_y_vs_x_ttm_ClCo->Fill(x_ref, y_ref);
700  }
701 
702  // make correlation plots
703  for (const auto &proton_m : *hRecoProtonsMultiRP) {
704  CTPPSDetId rpId_m((*proton_m.contributingLocalTracks().begin())->getRPId());
705  unsigned int arm = rpId_m.arm();
706 
707  const reco::ForwardProton *p_s_N = nullptr, *p_s_F = nullptr;
708 
709  for (const auto &proton_s : *hRecoProtonsSingleRP) {
710  // skip if source of single-RP reco not included in multi-RP reco
711  const auto key_s = proton_s.contributingLocalTracks()[0].key();
712  bool compatible = false;
713  for (const auto &tr_m : proton_m.contributingLocalTracks()) {
714  if (tr_m.key() == key_s) {
715  compatible = true;
716  break;
717  }
718  }
719 
720  if (!compatible)
721  continue;
722 
723  // fill single-multi plots
724  CTPPSDetId rpId_s((*proton_s.contributingLocalTracks().begin())->getRPId());
725  const unsigned int idx = rpId_s.arm() * 1000 + rpId_s.station() * 100 + rpId_s.rp() * 10 + rpId_s.arm();
726  singleMultiCorrelationPlots_[idx].fill(proton_s, proton_m);
727 
728  // select protons for arm-correlation plots
729  const unsigned int rpDecId_s = rpId_s.arm() * 100 + rpId_s.station() * 10 + rpId_s.rp();
730  if (rpDecId_s == rpId_45_N_ || rpDecId_s == rpId_56_N_)
731  p_s_N = &proton_s;
732  if (rpDecId_s == rpId_45_F_ || rpDecId_s == rpId_56_F_)
733  p_s_F = &proton_s;
734  }
735 
736  // fill arm-correlation plots
737  if (p_s_N && p_s_F)
738  armCorrelationPlots_[arm].fill(*p_s_N, *p_s_F, proton_m);
739  }
740 }
741 
742 //----------------------------------------------------------------------------------------------------
743 
745  LogInfo("CTPPSProtonReconstructionPlotter") << "n_non_empty_events = " << n_non_empty_events_;
746 
747  auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
748 
749  p_x_L_diffNF_vs_x_L_N_->Write();
750  p_x_R_diffNF_vs_x_R_N_->Write();
751 
752  p_y_L_diffNF_vs_y_L_N_->Write();
753  p_y_R_diffNF_vs_y_R_N_->Write();
754 
755  TDirectory *d_singleRPPlots = f_out->mkdir("singleRPPlots");
756  for (const auto &it : singleRPPlots_) {
757  gDirectory = d_singleRPPlots->mkdir(Form("rp%u", it.first));
758  it.second.write();
759  }
760 
761  TDirectory *d_multiRPPlots = f_out->mkdir("multiRPPlots");
762  for (const auto &it : multiRPPlots_) {
763  gDirectory = d_multiRPPlots->mkdir(Form("arm%u", it.first));
764  it.second.write();
765  }
766 
767  TDirectory *d_singleMultiCorrelationPlots = f_out->mkdir("singleMultiCorrelationPlots");
768  for (const auto &it : singleMultiCorrelationPlots_) {
769  unsigned int si_rp = it.first / 10;
770  unsigned int mu_arm = it.first % 10;
771 
772  gDirectory = d_singleMultiCorrelationPlots->mkdir(Form("si_rp%u_mu_arm%u", si_rp, mu_arm));
773  it.second.write();
774  }
775 
776  TDirectory *d_armCorrelationPlots = f_out->mkdir("armCorrelationPlots");
777  for (const auto &it : armCorrelationPlots_) {
778  gDirectory = d_armCorrelationPlots->mkdir(Form("arm%u", it.first));
779  it.second.write();
780  }
781 }
782 
783 //----------------------------------------------------------------------------------------------------
784 
const Point & vertex() const
fitted vertex position
Definition: ForwardProton.h:46
T getParameter(std::string const &) const
static void profileToRMSGraph(TProfile *p, TGraphErrors *g)
std::map< unsigned int, MultiRPPlots > multiRPPlots_
const unsigned int nTracks(const reco::Vertex &sv)
edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsSingleRP_
float getXUnc() const
returns the horizontal track position uncertainty
std::map< unsigned int, SingleRPPlots > singleRPPlots_
Local (=single RP) track with essential information only.
void analyze(const edm::Event &, const edm::EventSetup &) override
float t() const
four-momentum transfer squared, in GeV^2
edm::EDGetTokenT< CTPPSLocalTrackLiteCollection > tokenTracks_
void fill(const reco::ForwardProton &p_single, const reco::ForwardProton &p_multi)
float getX() const
returns the horizontal track position
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
Event setup record containing the real (actual) geometry information.
float chi2() const
chi-squared of the fit
Definition: ForwardProton.h:67
float thetaX() const
vertical scattering angle, in rad
Definition: ForwardProton.h:78
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
uint32_t getRPId() const
returns the RP id
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
T sqrt(T t)
Definition: SSEVec.h:18
float normalizedChi2() const
chi-squared divided by ndof (or chi-squared * 1e6 if ndof is zero)
Definition: ForwardProton.h:71
float xi() const
longitudinal fractional momentum loss
Definition: ForwardProton.h:76
float getY() const
returns the vertical track position
uint32_t arm() const
Definition: CTPPSDetId.h:51
const CTPPSLocalTrackLiteRefVector & contributingLocalTracks() const
list of RP tracks that contributed to this global track
std::map< unsigned int, ArmCorrelationPlots > armCorrelationPlots_
The manager class for TOTEM RP geometry.
Definition: CTPPSGeometry.h:33
void CalculateTimingTrackingDistance(const reco::ForwardProton &proton, const CTPPSLocalTrackLite &tr, const CTPPSGeometry &geometry, double &de_x, double &de_x_unc)
CLHEP::Hep3Vector getRPTranslation(unsigned int id) const
bool validFit() const
flag for the fit validity
#define N
Definition: blowfish.cc:9
std::vector< CTPPSLocalTrackLite > CTPPSLocalTrackLiteCollection
Collection of CTPPSLocalTrackLite objects.
ParameterSet const & getParameterSet(std::string const &) const
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
fixed size matrix
value_type const at(size_type idx) const
Retrieve an element of the RefVector.
Definition: RefVector.h:88
HLT enums.
CTPPSProtonReconstructionPlotter(const edm::ParameterSet &)
std::map< unsigned int, SingleMultiCorrelationPlots > singleMultiCorrelationPlots_
T get() const
Definition: EventSetup.h:71
float thetaY() const
horizontal scattering angle, in rad
Definition: ForwardProton.h:80
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.
void fill(const reco::ForwardProton &p_s_N, const reco::ForwardProton &p_s_F, const reco::ForwardProton &p_m)
float time() const
time of proton arrival at forward stations
edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsMultiRP_
std::map< unsigned int, AssociationCuts > association_cuts_
Definition: event.py:1
void fill(const reco::ForwardProton &p, unsigned int nTracks, bool n1f1)
void fill(const reco::ForwardProton &p, unsigned int nTracks, bool n1f1)
unsigned int ndof() const
number of degrees of freedom for the track fit
Definition: ForwardProton.h:69