27 #include "TGraphErrors.h" 29 #include "CLHEP/Units/GlobalPhysicalConstants.h" 49 const HepMC::FourVector &
vtx,
50 const HepMC::FourVector &mom,
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)),
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)),
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)),
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)),
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.)) {}
105 gr_err.SetName(
name);
107 for (
int bi = 1; bi <=
p->GetNbinsX(); ++bi) {
108 double c =
p->GetBinCenter(bi);
109 double w =
p->GetBinWidth(bi);
111 double N =
p->GetBinEntries(bi);
112 double Sy =
p->GetBinContent(bi) *
N;
113 double Syy =
p->GetSumw2()->At(bi);
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;
118 double si_unc = (si_unc_sq >= 0.) ?
sqrt(si_unc_sq) : 0.;
120 int idx = gr_err.GetN();
121 gr_err.SetPoint(
idx,
c, si);
122 gr_err.SetPointError(
idx,
w / 2., si_unc);
155 std::map<unsigned int, std::map<unsigned int, PlotGroup>>
plots_;
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.)),
165 h_t_dh_minus_vtx_z(new TH1D(
"",
";(t_56 - t_45)/2 - vtx_z (mm)", 100, -100., +100.)) {}
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;
194 using namespace HepMC;
200 : tokenHepMCBeforeSmearing_(
202 tokenHepMCAfterSmearing_(
204 tokenRecoProtonsSingleRP_(
206 tokenRecoProtonsMultiRP_(
211 useNewLHCInfo_(iConfig.getParameter<
bool>(
"useNewLHCInfo")),
212 outputFile_(iConfig.getParameter<
string>(
"outputFile")) {}
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");
224 desc.add<
std::string>(
"outputFile",
"output.root")->setComment(
"output file name");
226 desc.addUntracked<
unsigned int>(
"verbosity", 0)->setComment(
"verbosity level");
228 descriptions.
add(
"ctppsProtonReconstructionSimulationValidatorDefault",
desc);
254 bool vertex_set =
false;
256 for (
auto it = hepMCEventAfterSmearing->vertices_begin();
it != hepMCEventAfterSmearing->vertices_end(); ++
it) {
258 LogError(
"CTPPSProtonReconstructionSimulationValidator") <<
"Multiple vertices found.";
263 vtx = (*it)->position();
267 bool proton_45_set =
false;
268 bool proton_56_set =
false;
269 FourVector mom_45, mom_56;
271 for (
auto it = hepMCEventBeforeSmearing->particles_begin();
it != hepMCEventBeforeSmearing->particles_end(); ++
it) {
275 if (
part->pdg_id() != 2212)
278 if (
part->status() != 1)
284 const auto &mom =
part->momentum();
292 LogError(
"CTPPSProtonReconstructionSimulationValidator") <<
"Found multiple protons in sector 45.";
296 proton_45_set =
true;
301 LogError(
"CTPPSProtonReconstructionSimulationValidator") <<
"Found multiple protons in sector 56.";
305 proton_56_set =
true;
311 for (
const auto &
handle : {hRecoProtonsSingleRP, hRecoProtonsMultiRP}) {
312 for (
const auto &rec_pr : *
handle) {
313 if (!rec_pr.validFit())
316 unsigned int idx = 12345;
318 bool mom_set =
false;
323 mom_set = proton_45_set;
328 mom_set = proton_56_set;
335 unsigned int meth_idx = 1234;
340 CTPPSDetId rpId((*rec_pr.contributingLocalTracks().begin())->rpId());
341 idx = 100 * rpId.arm() + 10 * rpId.station() + rpId.rp();
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())
359 if (rec_pr.timeError() < 1E-3)
364 time_45 = rec_pr.time();
368 time_56 = rec_pr.time();
372 if (time_set_45 && time_set_56)
381 const HepMC::FourVector &
vtx,
382 const HepMC::FourVector &mom,
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);
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.;
396 const double t_reco = -rec_pr.
t();
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);
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);
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);
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);
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);
422 auto f_out = std::make_unique<TFile>(
outputFile_.c_str(),
"recreate");
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);
428 for (
const auto &eit : mit.second) {
429 gDirectory = d_method->mkdir(Form(
"%i", eit.first));
434 gDirectory = f_out->mkdir(
"double arm");
std::unique_ptr< TH1D > h_de_t
CTPPSProtonReconstructionSimulationValidator(const edm::ParameterSet &)
std::unique_ptr< TProfile > p_de_th_x_vs_xi_simu
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...
float thetaY() const
horizontal scattering angle, in rad
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)
std::unique_ptr< TH1D > h_t_dh_minus_vtx_z
edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsSingleRP_
const edm::ESGetToken< LHCInfoPerFill, LHCInfoPerFillRcd > lhcInfoPerFillToken_
float thetaX() const
vertical scattering angle, in rad
Log< level::Error, false > LogError
std::unique_ptr< TProfile > p_de_th_y_vs_xi_simu
std::unique_ptr< TH2D > h2_t_dh_vs_vtx_z
std::unique_ptr< TH1D > h_de_th_y
std::unique_ptr< TProfile > p_de_t_vs_xi_simu
edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsMultiRP_
const bool useNewLHCInfo_
#define DEFINE_FWK_MODULE(type)
float xi() const
longitudinal fractional momentum loss
std::unique_ptr< TH1D > h_de_vtx_y
const HepMC::GenEvent * GetEvent() const
std::unique_ptr< TH1D > h_de_th_x
std::unique_ptr< TH2D > h2_t_sh_vs_vtx_t
std::unique_ptr< TProfile > p_de_t_vs_t_simu
static TGraphErrors profileToRMSGraph(TProfile *p, const char *name="")
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const Point & vertex() const
fitted vertex position
std::unique_ptr< TH1D > h_t_sh_minus_vtx_t
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.
std::unique_ptr< TH2D > h_xi_reco_vs_xi_simu
edm::EDGetTokenT< edm::HepMCProduct > tokenHepMCBeforeSmearing_
DoubleArmPlotGroup double_arm_plots_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
const edm::ESGetToken< LHCInfo, LHCInfoRcd > lhcInfoToken_
std::unique_ptr< TH1D > h_de_xi
std::unique_ptr< TProfile > p_de_xi_vs_xi_simu
std::vector< ForwardProton > ForwardProtonCollection
Collection of ForwardProton objects.
edm::EDGetTokenT< edm::HepMCProduct > tokenHepMCAfterSmearing_
float t() const
four-momentum transfer squared, in GeV^2
std::unique_ptr< TH2D > h2_de_t_vs_t_simu
const edm::ESGetToken< LHCInfoPerLS, LHCInfoPerLSRcd > lhcInfoPerLSToken_
std::map< unsigned int, std::map< unsigned int, PlotGroup > > plots_
std::unique_ptr< TProfile > p_de_vtx_y_vs_xi_simu