28 #include "TGraphErrors.h"
46 const HepMC::FourVector &
vtx,
47 const HepMC::FourVector &mom,
74 std::unique_ptr<TH1D>
h_de_t;
78 :
h_de_xi(new TH1D(
"",
";#xi_{reco} - #xi_{simu}", 100, 0., 0.)),
79 p_de_xi_vs_xi_simu(new TProfile(
"",
";#xi_{simu};#xi_{reco} - #xi_{simu}", 19, 0.015, 0.205)),
82 h_de_th_x(new TH1D(
"",
";#theta_{x,reco} - #theta_{x,simu}", 100, 0., 0.)),
83 p_de_th_x_vs_xi_simu(new TProfile(
"",
";#xi_{simu};#theta_{x,reco} - #theta_{x,simu}", 19, 0.015, 0.205)),
85 h_de_th_y(new TH1D(
"",
";#theta_{y,reco} - #theta_{y,simu}", 100, 0., 0.)),
86 p_de_th_y_vs_xi_simu(new TProfile(
"",
";#xi_{simu};#theta_{y,reco} - #theta_{y,simu}", 19, 0.015, 0.205)),
88 h_de_vtx_y(new TH1D(
"",
";vtx_{y,reco} - vtx_{y,simu} (mm)", 100, 0., 0.)),
89 p_de_vtx_y_vs_xi_simu(new TProfile(
"",
";#xi_{simu};vtx_{y,reco} - vtx_{y,simu} (mm)", 19, 0.015, 0.205)),
91 h_de_t(new TH1D(
"",
";t_{reco} - t_{simu}", 100, -1., +1.)),
92 p_de_t_vs_xi_simu(new TProfile(
"",
";xi_{simu};t_{reco} - t_{simu}", 19, 0.015, 0.205)),
93 p_de_t_vs_t_simu(new TProfile(
"",
";t_{simu};t_{reco} - t_{simu}", 20, 0., 5.)) {}
99 for (
int bi = 1; bi <=
p->GetNbinsX(); ++bi) {
100 double c =
p->GetBinCenter(bi);
101 double w =
p->GetBinWidth(bi);
103 double N =
p->GetBinEntries(bi);
104 double Sy =
p->GetBinContent(bi) *
N;
105 double Syy =
p->GetSumw2()->At(bi);
107 double si_sq = Syy /
N - Sy * Sy /
N /
N;
108 double si = (si_sq >= 0.) ?
sqrt(si_sq) : 0.;
109 double si_unc_sq = si_sq / 2. /
N;
110 double si_unc = (si_unc_sq >= 0.) ?
sqrt(si_unc_sq) : 0.;
112 int idx = gr_err.GetN();
113 gr_err.SetPoint(
idx,
c, si);
114 gr_err.SetPointError(
idx,
w / 2., si_unc);
146 std::map<unsigned int, std::map<unsigned int, PlotGroup>>
plots_;
153 using namespace HepMC;
159 : tokenHepMCBeforeSmearing_(
161 tokenHepMCAfterSmearing_(
163 tokenRecoProtonsSingleRP_(
165 tokenRecoProtonsMultiRP_(
167 lhcInfoLabel_(iConfig.getParameter<
std::
string>(
"lhcInfoLabel")),
168 outputFile_(iConfig.getParameter<
string>(
"outputFile")) {}
193 bool vertex_set =
false;
195 for (
auto it = hepMCEventAfterSmearing->vertices_begin(); it != hepMCEventAfterSmearing->vertices_end(); ++it) {
197 LogError(
"CTPPSProtonReconstructionSimulationValidator") <<
"Multiple vertices found.";
202 vtx = (*it)->position();
206 bool proton_45_set =
false;
207 bool proton_56_set =
false;
208 FourVector mom_45, mom_56;
210 for (
auto it = hepMCEventBeforeSmearing->particles_begin(); it != hepMCEventBeforeSmearing->particles_end(); ++it) {
211 const auto &
part = *it;
214 if (
part->pdg_id() != 2212)
217 if (
part->status() != 1)
223 const auto &mom =
part->momentum();
231 LogError(
"CTPPSProtonReconstructionSimulationValidator") <<
"Found multiple protons in sector 45.";
235 proton_45_set =
true;
240 LogError(
"CTPPSProtonReconstructionSimulationValidator") <<
"Found multiple protons in sector 56.";
244 proton_56_set =
true;
250 for (
const auto &
handle : {hRecoProtonsSingleRP, hRecoProtonsMultiRP}) {
251 for (
const auto &rec_pr : *
handle) {
257 bool mom_set =
false;
262 mom_set = proton_45_set;
266 mom_set = proton_56_set;
273 unsigned int meth_idx = 1234;
295 const HepMC::FourVector &
vtx,
296 const HepMC::FourVector &mom,
298 const double p_nom = lhcInfo.
energy();
299 const double xi_simu = (p_nom - mom.rho()) / p_nom;
300 const double th_x_simu = mom.x() / mom.rho();
301 const double th_y_simu = mom.y() / mom.rho();
302 const double vtx_y_simu =
vtx.y();
303 const double th_simu =
sqrt(th_x_simu * th_x_simu + th_y_simu * th_y_simu);
306 const double xi_reco = rec_pr.
xi();
307 const double th_x_reco = rec_pr.
thetaX();
308 const double th_y_reco = rec_pr.
thetaY();
309 const double vtx_y_reco = rec_pr.
vertex().y() * 10.;
310 const double t_reco = -rec_pr.
t();
314 plt.h_xi_reco_vs_xi_simu->Fill(xi_simu, xi_reco);
315 plt.h_de_xi->Fill(xi_reco - xi_simu);
316 plt.p_de_xi_vs_xi_simu->Fill(xi_simu, xi_reco - xi_simu);
318 plt.h_de_th_x->Fill(th_x_reco - th_x_simu);
319 plt.p_de_th_x_vs_xi_simu->Fill(xi_simu, th_x_reco - th_x_simu);
321 plt.h_de_th_y->Fill(th_y_reco - th_y_simu);
322 plt.p_de_th_y_vs_xi_simu->Fill(xi_simu, th_y_reco - th_y_simu);
324 plt.h_de_vtx_y->Fill(vtx_y_reco - vtx_y_simu);
325 plt.p_de_vtx_y_vs_xi_simu->Fill(xi_simu, vtx_y_reco - vtx_y_simu);
327 plt.h_de_t->Fill(t_reco - t_simu);
328 plt.p_de_t_vs_xi_simu->Fill(xi_simu, t_reco - t_simu);
329 plt.p_de_t_vs_t_simu->Fill(t_simu, t_reco - t_simu);
335 auto f_out = std::make_unique<TFile>(
outputFile_.c_str(),
"recreate");
337 for (
const auto &mit :
plots_) {
338 const char *
method = (mit.first == 0) ?
"single rp" :
"multi rp";
339 TDirectory *d_method = f_out->mkdir(
method);
341 for (
const auto &eit : mit.second) {
342 gDirectory = d_method->mkdir(Form(
"%i", eit.first));