CMS 3D CMS Logo

CTPPSProtonReconstructionSimulationValidator.cc
Go to the documentation of this file.
1 /****************************************************************************
2  * Authors:
3  * Jan Kašpar
4  ****************************************************************************/
5 
12 
14 
16 
18 
22 
23 #include "TFile.h"
24 #include "TH1D.h"
25 #include "TH2D.h"
26 #include "TProfile.h"
27 #include "TGraphErrors.h"
28 
29 #include "CLHEP/Units/GlobalPhysicalConstants.h"
30 
31 #include <map>
32 #include <string>
33 
34 //----------------------------------------------------------------------------------------------------
35 
37 public:
39 
40  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
41 
42 private:
43  void analyze(const edm::Event &, const edm::EventSetup &) override;
44  void endJob() override;
45 
46  void fillPlots(unsigned int meth_idx,
47  unsigned int idx,
48  const reco::ForwardProton &rec_pr,
49  const HepMC::FourVector &vtx,
50  const HepMC::FourVector &mom,
51  const double energy);
52 
55 
58 
62  const bool useNewLHCInfo_;
63 
65 
66  struct PlotGroup {
67  std::unique_ptr<TH1D> h_de_xi;
68  std::unique_ptr<TProfile> p_de_xi_vs_xi_simu;
69  std::unique_ptr<TH2D> h_xi_reco_vs_xi_simu;
70 
71  std::unique_ptr<TH1D> h_de_th_x;
72  std::unique_ptr<TProfile> p_de_th_x_vs_xi_simu;
73 
74  std::unique_ptr<TH1D> h_de_th_y;
75  std::unique_ptr<TProfile> p_de_th_y_vs_xi_simu;
76 
77  std::unique_ptr<TH1D> h_de_vtx_y;
78  std::unique_ptr<TProfile> p_de_vtx_y_vs_xi_simu;
79 
80  std::unique_ptr<TH1D> h_de_t;
81  std::unique_ptr<TProfile> p_de_t_vs_xi_simu, p_de_t_vs_t_simu;
82  std::unique_ptr<TH2D> h2_de_t_vs_t_simu;
83 
85  : h_de_xi(new TH1D("", ";#xi_{reco} - #xi_{simu}", 100, 0., 0.)),
86  p_de_xi_vs_xi_simu(new TProfile("", ";#xi_{simu};#xi_{reco} - #xi_{simu}", 50, 0., 0.25)),
87  h_xi_reco_vs_xi_simu(new TH2D("", ";#xi_{simu};#xi_{reco}", 100, 0., 0.30, 100, 0., 0.30)),
88 
89  h_de_th_x(new TH1D("", ";#theta_{x,reco} - #theta_{x,simu}", 100, 0., 0.)),
90  p_de_th_x_vs_xi_simu(new TProfile("", ";#xi_{simu};#theta_{x,reco} - #theta_{x,simu}", 50, 0., 0.25)),
91 
92  h_de_th_y(new TH1D("", ";#theta_{y,reco} - #theta_{y,simu}", 100, 0., 0.)),
93  p_de_th_y_vs_xi_simu(new TProfile("", ";#xi_{simu};#theta_{y,reco} - #theta_{y,simu}", 50, 0., 0.25)),
94 
95  h_de_vtx_y(new TH1D("", ";vtx_{y,reco} - vtx_{y,simu} (mm)", 100, 0., 0.)),
96  p_de_vtx_y_vs_xi_simu(new TProfile("", ";#xi_{simu};vtx_{y,reco} - vtx_{y,simu} (mm)", 50, 0., 0.25)),
97 
98  h_de_t(new TH1D("", ";t_{reco} - t_{simu}", 100, -1., +1.)),
99  p_de_t_vs_xi_simu(new TProfile("", ";xi_{simu};t_{reco} - t_{simu}", 50, 0., 0.25)),
100  p_de_t_vs_t_simu(new TProfile("", ";t_{simu};t_{reco} - t_{simu}", 20, 0., 5.)),
101  h2_de_t_vs_t_simu(new TH2D("", ";t_{simu};t_{reco} - t_{simu}", 150, 0., 5., 300, -2., +2.)) {}
102 
103  static TGraphErrors profileToRMSGraph(TProfile *p, const char *name = "") {
104  TGraphErrors gr_err;
105  gr_err.SetName(name);
106 
107  for (int bi = 1; bi <= p->GetNbinsX(); ++bi) {
108  double c = p->GetBinCenter(bi);
109  double w = p->GetBinWidth(bi);
110 
111  double N = p->GetBinEntries(bi);
112  double Sy = p->GetBinContent(bi) * N;
113  double Syy = p->GetSumw2()->At(bi);
114 
115  double si_sq = Syy / N - Sy * Sy / N / N;
116  double si = (si_sq >= 0.) ? sqrt(si_sq) : 0.;
117  double si_unc_sq = si_sq / 2. / N; // Gaussian approximation
118  double si_unc = (si_unc_sq >= 0.) ? sqrt(si_unc_sq) : 0.;
119 
120  int idx = gr_err.GetN();
121  gr_err.SetPoint(idx, c, si);
122  gr_err.SetPointError(idx, w / 2., si_unc);
123  }
124 
125  return gr_err;
126  }
127 
128  void write() const {
129  h_xi_reco_vs_xi_simu->Write("h_xi_reco_vs_xi_simu");
130  h_de_xi->Write("h_de_xi");
131  p_de_xi_vs_xi_simu->Write("p_de_xi_vs_xi_simu");
132  profileToRMSGraph(p_de_xi_vs_xi_simu.get(), "g_rms_de_xi_vs_xi_simu").Write();
133 
134  h_de_th_x->Write("h_de_th_x");
135  p_de_th_x_vs_xi_simu->Write("p_de_th_x_vs_xi_simu");
136  profileToRMSGraph(p_de_th_x_vs_xi_simu.get(), "g_rms_de_th_x_vs_xi_simu").Write();
137 
138  h_de_th_y->Write("h_de_th_y");
139  p_de_th_y_vs_xi_simu->Write("p_de_th_y_vs_xi_simu");
140  profileToRMSGraph(p_de_th_y_vs_xi_simu.get(), "g_rms_de_th_y_vs_xi_simu").Write();
141 
142  h_de_vtx_y->Write("h_de_vtx_y");
143  p_de_vtx_y_vs_xi_simu->Write("p_de_vtx_y_vs_xi_simu");
144  profileToRMSGraph(p_de_vtx_y_vs_xi_simu.get(), "g_rms_de_vtx_y_vs_xi_simu").Write();
145 
146  h_de_t->Write("h_de_t");
147  p_de_t_vs_xi_simu->Write("p_de_t_vs_xi_simu");
148  profileToRMSGraph(p_de_t_vs_xi_simu.get(), "g_rms_de_t_vs_xi_simu").Write();
149  p_de_t_vs_t_simu->Write("p_de_t_vs_t_simu");
150  profileToRMSGraph(p_de_t_vs_t_simu.get(), "g_rms_de_t_vs_t_simu").Write();
151  h2_de_t_vs_t_simu->Write("h2_de_t_vs_t_simu");
152  }
153  };
154 
155  std::map<unsigned int, std::map<unsigned int, PlotGroup>> plots_;
156 
158  std::unique_ptr<TH2D> h2_t_sh_vs_vtx_t, h2_t_dh_vs_vtx_z;
159  std::unique_ptr<TH1D> h_t_sh_minus_vtx_t, h_t_dh_minus_vtx_z;
160 
162  : h2_t_sh_vs_vtx_t(new TH2D("", ";vtx_t (mm);(t_56 + t_45)/2 (mm)", 100, -250., +250., 100, -250., +250.)),
163  h2_t_dh_vs_vtx_z(new TH2D("", ";vtx_z (mm);(t_56 - t_45)/2 (mm)", 100, -250., +250., 100, -250., +250.)),
164  h_t_sh_minus_vtx_t(new TH1D("", ";(t_56 + t_45)/2 - vtx_t (mm)", 100, -100., +100.)),
165  h_t_dh_minus_vtx_z(new TH1D("", ";(t_56 - t_45)/2 - vtx_z (mm)", 100, -100., +100.)) {}
166 
167  void fill(double time_45, double time_56, double vtx_z, double vtx_t) {
168  const double t_sum_half = (time_56 + time_45) / 2. * CLHEP::c_light;
169  const double t_dif_half = (time_56 - time_45) / 2. * CLHEP::c_light;
170 
171  h2_t_sh_vs_vtx_t->Fill(vtx_t, t_sum_half);
172  h_t_sh_minus_vtx_t->Fill(t_sum_half - vtx_t);
173 
174  h2_t_dh_vs_vtx_z->Fill(vtx_z, t_dif_half);
175  h_t_dh_minus_vtx_z->Fill(t_dif_half - vtx_z);
176  }
177 
178  void write() const {
179  h2_t_sh_vs_vtx_t->Write("h2_t_sh_vs_vtx_t");
180  h_t_sh_minus_vtx_t->Write("h_t_sh_minus_vtx_t");
181 
182  h2_t_dh_vs_vtx_z->Write("h2_t_dh_vs_vtx_z");
183  h_t_dh_minus_vtx_z->Write("h_t_dh_minus_vtx_z");
184  }
185  };
186 
188 };
189 
190 //----------------------------------------------------------------------------------------------------
191 
192 using namespace std;
193 using namespace edm;
194 using namespace HepMC;
195 
196 //----------------------------------------------------------------------------------------------------
197 
199  const edm::ParameterSet &iConfig)
200  : tokenHepMCBeforeSmearing_(
201  consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("tagHepMCBeforeSmearing"))),
202  tokenHepMCAfterSmearing_(
203  consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("tagHepMCAfterSmearing"))),
204  tokenRecoProtonsSingleRP_(
205  consumes<reco::ForwardProtonCollection>(iConfig.getParameter<InputTag>("tagRecoProtonsSingleRP"))),
206  tokenRecoProtonsMultiRP_(
207  consumes<reco::ForwardProtonCollection>(iConfig.getParameter<InputTag>("tagRecoProtonsMultiRP"))),
208  lhcInfoToken_(esConsumes(ESInputTag("", iConfig.getParameter<std::string>("lhcInfoLabel")))),
209  lhcInfoPerLSToken_(esConsumes(ESInputTag("", iConfig.getParameter<std::string>("lhcInfoPerLSLabel")))),
210  lhcInfoPerFillToken_(esConsumes(ESInputTag("", iConfig.getParameter<std::string>("lhcInfoPerFillLabel")))),
211  useNewLHCInfo_(iConfig.getParameter<bool>("useNewLHCInfo")),
212  outputFile_(iConfig.getParameter<string>("outputFile")) {}
213 
214 //----------------------------------------------------------------------------------------------------
215 
218 
219  desc.add<std::string>("lhcInfoLabel", "")->setComment("label of the LHCInfo record");
220  desc.add<std::string>("lhcInfoPerLSLabel", "")->setComment("label of the LHCInfoPerLS record");
221  desc.add<std::string>("lhcInfoPerFillLabel", "")->setComment("label of the LHCInfoPerFill record");
222  desc.add<bool>("useNewLHCInfo", false)->setComment("flag whether to use new LHCInfoPer* records or old LHCInfo");
223 
224  desc.add<std::string>("outputFile", "output.root")->setComment("output file name");
225 
226  desc.addUntracked<unsigned int>("verbosity", 0)->setComment("verbosity level");
227 
228  descriptions.add("ctppsProtonReconstructionSimulationValidatorDefault", desc);
229 }
230 
231 //----------------------------------------------------------------------------------------------------
232 
234  // get conditions
235  const LHCInfoCombined lhcInfoCombined(
237 
238  // get input
239  edm::Handle<edm::HepMCProduct> hHepMCBeforeSmearing;
240  iEvent.getByToken(tokenHepMCBeforeSmearing_, hHepMCBeforeSmearing);
241  HepMC::GenEvent *hepMCEventBeforeSmearing = (HepMC::GenEvent *)hHepMCBeforeSmearing->GetEvent();
242 
243  edm::Handle<edm::HepMCProduct> hHepMCAfterSmearing;
244  iEvent.getByToken(tokenHepMCAfterSmearing_, hHepMCAfterSmearing);
245  HepMC::GenEvent *hepMCEventAfterSmearing = (HepMC::GenEvent *)hHepMCAfterSmearing->GetEvent();
246 
247  Handle<reco::ForwardProtonCollection> hRecoProtonsSingleRP;
248  iEvent.getByToken(tokenRecoProtonsSingleRP_, hRecoProtonsSingleRP);
249 
250  Handle<reco::ForwardProtonCollection> hRecoProtonsMultiRP;
251  iEvent.getByToken(tokenRecoProtonsMultiRP_, hRecoProtonsMultiRP);
252 
253  // extract vertex position
254  bool vertex_set = false;
255  FourVector vtx;
256  for (auto it = hepMCEventAfterSmearing->vertices_begin(); it != hepMCEventAfterSmearing->vertices_end(); ++it) {
257  if (vertex_set) {
258  LogError("CTPPSProtonReconstructionSimulationValidator") << "Multiple vertices found.";
259  return;
260  }
261 
262  vertex_set = true;
263  vtx = (*it)->position();
264  }
265 
266  // extract forward protons
267  bool proton_45_set = false;
268  bool proton_56_set = false;
269  FourVector mom_45, mom_56;
270 
271  for (auto it = hepMCEventBeforeSmearing->particles_begin(); it != hepMCEventBeforeSmearing->particles_end(); ++it) {
272  const auto &part = *it;
273 
274  // accept only stable non-beam protons
275  if (part->pdg_id() != 2212)
276  continue;
277 
278  if (part->status() != 1)
279  continue;
280 
281  if (part->is_beam())
282  continue;
283 
284  const auto &mom = part->momentum();
285 
286  if (mom.e() < 4500.)
287  continue;
288 
289  if (mom.z() > 0) {
290  // 45
291  if (proton_45_set) {
292  LogError("CTPPSProtonReconstructionSimulationValidator") << "Found multiple protons in sector 45.";
293  return;
294  }
295 
296  proton_45_set = true;
297  mom_45 = mom;
298  } else {
299  // 56
300  if (proton_56_set) {
301  LogError("CTPPSProtonReconstructionSimulationValidator") << "Found multiple protons in sector 56.";
302  return;
303  }
304 
305  proton_56_set = true;
306  mom_56 = mom;
307  }
308  }
309 
310  // single-arm comparison
311  for (const auto &handle : {hRecoProtonsSingleRP, hRecoProtonsMultiRP}) {
312  for (const auto &rec_pr : *handle) {
313  if (!rec_pr.validFit())
314  continue;
315 
316  unsigned int idx = 12345;
317 
318  bool mom_set = false;
319  FourVector mom;
320 
321  if (rec_pr.lhcSector() == reco::ForwardProton::LHCSector::sector45) {
322  idx = 0;
323  mom_set = proton_45_set;
324  mom = mom_45;
325  }
326  if (rec_pr.lhcSector() == reco::ForwardProton::LHCSector::sector56) {
327  idx = 1;
328  mom_set = proton_56_set;
329  mom = mom_56;
330  }
331 
332  if (!mom_set)
333  continue;
334 
335  unsigned int meth_idx = 1234;
336 
338  meth_idx = 0;
339 
340  CTPPSDetId rpId((*rec_pr.contributingLocalTracks().begin())->rpId());
341  idx = 100 * rpId.arm() + 10 * rpId.station() + rpId.rp();
342  }
343 
345  meth_idx = 1;
346 
347  fillPlots(meth_idx, idx, rec_pr, vtx, mom, lhcInfoCombined.energy);
348  }
349  }
350 
351  // double-arm comparison
352  bool time_set_45 = false, time_set_56 = false;
353  double time_45 = 0., time_56 = 0.;
354  for (const auto &rec_pr : *hRecoProtonsMultiRP) {
355  if (!rec_pr.validFit())
356  continue;
357 
358  // skip protons with no timing information
359  if (rec_pr.timeError() < 1E-3)
360  continue;
361 
362  if (rec_pr.lhcSector() == reco::ForwardProton::LHCSector::sector45) {
363  time_set_45 = true;
364  time_45 = rec_pr.time();
365  }
366  if (rec_pr.lhcSector() == reco::ForwardProton::LHCSector::sector56) {
367  time_set_56 = true;
368  time_56 = rec_pr.time();
369  }
370  }
371 
372  if (time_set_45 && time_set_56)
373  double_arm_plots_.fill(time_45, time_56, vtx.z(), vtx.t());
374 }
375 
376 //----------------------------------------------------------------------------------------------------
377 
379  unsigned int idx,
380  const reco::ForwardProton &rec_pr,
381  const HepMC::FourVector &vtx,
382  const HepMC::FourVector &mom,
383  const double energy) {
384  const double p_nom = energy;
385  const double xi_simu = (p_nom - mom.rho()) / p_nom;
386  const double th_x_simu = mom.x() / mom.rho();
387  const double th_y_simu = mom.y() / mom.rho();
388  const double vtx_y_simu = vtx.y();
389  const double th_simu = sqrt(th_x_simu * th_x_simu + th_y_simu * th_y_simu);
390  const double t_simu = -reco::ForwardProton::calculateT(p_nom, mom.rho(), th_simu);
391 
392  const double xi_reco = rec_pr.xi();
393  const double th_x_reco = rec_pr.thetaX();
394  const double th_y_reco = rec_pr.thetaY();
395  const double vtx_y_reco = rec_pr.vertex().y() * 10.; // conversion: cm --> mm
396  const double t_reco = -rec_pr.t();
397 
398  auto &plt = plots_[meth_idx][idx];
399 
400  plt.h_xi_reco_vs_xi_simu->Fill(xi_simu, xi_reco);
401  plt.h_de_xi->Fill(xi_reco - xi_simu);
402  plt.p_de_xi_vs_xi_simu->Fill(xi_simu, xi_reco - xi_simu);
403 
404  plt.h_de_th_x->Fill(th_x_reco - th_x_simu);
405  plt.p_de_th_x_vs_xi_simu->Fill(xi_simu, th_x_reco - th_x_simu);
406 
407  plt.h_de_th_y->Fill(th_y_reco - th_y_simu);
408  plt.p_de_th_y_vs_xi_simu->Fill(xi_simu, th_y_reco - th_y_simu);
409 
410  plt.h_de_vtx_y->Fill(vtx_y_reco - vtx_y_simu);
411  plt.p_de_vtx_y_vs_xi_simu->Fill(xi_simu, vtx_y_reco - vtx_y_simu);
412 
413  plt.h_de_t->Fill(t_reco - t_simu);
414  plt.p_de_t_vs_xi_simu->Fill(xi_simu, t_reco - t_simu);
415  plt.p_de_t_vs_t_simu->Fill(t_simu, t_reco - t_simu);
416  plt.h2_de_t_vs_t_simu->Fill(t_simu, t_reco - t_simu);
417 }
418 
419 //----------------------------------------------------------------------------------------------------
420 
422  auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
423 
424  for (const auto &mit : plots_) {
425  const char *method = (mit.first == 0) ? "single rp" : "multi rp";
426  TDirectory *d_method = f_out->mkdir(method);
427 
428  for (const auto &eit : mit.second) {
429  gDirectory = d_method->mkdir(Form("%i", eit.first));
430  eit.second.write();
431  }
432  }
433 
434  gDirectory = f_out->mkdir("double arm");
436 }
437 
438 //----------------------------------------------------------------------------------------------------
439 
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
static float calculateT(double beam_mom, double proton_mom, double theta)
compute the squared four-momentum transfer from incident and scattered momenta, and angular informati...
T w() const
float thetaY() const
horizontal scattering angle, in rad
Definition: ForwardProton.h:83
void fillPlots(unsigned int meth_idx, unsigned int idx, const reco::ForwardProton &rec_pr, const HepMC::FourVector &vtx, const HepMC::FourVector &mom, const double energy)
edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsSingleRP_
const edm::ESGetToken< LHCInfoPerFill, LHCInfoPerFillRcd > lhcInfoPerFillToken_
float thetaX() const
vertical scattering angle, in rad
Definition: ForwardProton.h:81
Log< level::Error, false > LogError
int iEvent
Definition: GenABIO.cc:224
T sqrt(T t)
Definition: SSEVec.h:19
edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsMultiRP_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
float xi() const
longitudinal fractional momentum loss
Definition: ForwardProton.h:79
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:37
#define N
Definition: blowfish.cc:9
part
Definition: HCALResponse.h:20
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const Point & vertex() const
fitted vertex position
Definition: ForwardProton.h:51
void fill(double time_45, double time_56, double vtx_z, double vtx_t)
void analyze(const edm::Event &, const edm::EventSetup &) override
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
fixed size matrix
HLT enums.
ap_int< 5 > vtx_t
Definition: datatypes.h:20
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
std::vector< ForwardProton > ForwardProtonCollection
Collection of ForwardProton objects.
float t() const
four-momentum transfer squared, in GeV^2
const edm::ESGetToken< LHCInfoPerLS, LHCInfoPerLSRcd > lhcInfoPerLSToken_
std::map< unsigned int, std::map< unsigned int, PlotGroup > > plots_