CMS 3D CMS Logo

CTPPSProtonReconstructionSimulationValidator.cc
Go to the documentation of this file.
1 /****************************************************************************
2  * Authors:
3  * Jan Kašpar
4  ****************************************************************************/
5 
12 
14 
17 
19 
23 
24 #include "TFile.h"
25 #include "TH1D.h"
26 #include "TH2D.h"
27 #include "TProfile.h"
28 #include "TGraphErrors.h"
29 
30 #include "CLHEP/Units/GlobalPhysicalConstants.h"
31 
32 #include <map>
33 #include <string>
34 
35 //----------------------------------------------------------------------------------------------------
36 
38 public:
40 
41 private:
42  void analyze(const edm::Event &, const edm::EventSetup &) override;
43  void endJob() override;
44 
45  void fillPlots(unsigned int meth_idx,
46  unsigned int idx,
47  const reco::ForwardProton &rec_pr,
48  const HepMC::FourVector &vtx,
49  const HepMC::FourVector &mom,
50  const LHCInfo &lhcInfo);
51 
54 
57 
59 
61 
62  struct PlotGroup {
63  std::unique_ptr<TH1D> h_de_xi;
64  std::unique_ptr<TProfile> p_de_xi_vs_xi_simu;
65  std::unique_ptr<TH2D> h_xi_reco_vs_xi_simu;
66 
67  std::unique_ptr<TH1D> h_de_th_x;
68  std::unique_ptr<TProfile> p_de_th_x_vs_xi_simu;
69 
70  std::unique_ptr<TH1D> h_de_th_y;
71  std::unique_ptr<TProfile> p_de_th_y_vs_xi_simu;
72 
73  std::unique_ptr<TH1D> h_de_vtx_y;
74  std::unique_ptr<TProfile> p_de_vtx_y_vs_xi_simu;
75 
76  std::unique_ptr<TH1D> h_de_t;
77  std::unique_ptr<TProfile> p_de_t_vs_xi_simu, p_de_t_vs_t_simu;
78 
79  PlotGroup()
80  : h_de_xi(new TH1D("", ";#xi_{reco} - #xi_{simu}", 100, 0., 0.)),
81  p_de_xi_vs_xi_simu(new TProfile("", ";#xi_{simu};#xi_{reco} - #xi_{simu}", 50, 0., 0.25)),
82  h_xi_reco_vs_xi_simu(new TH2D("", ";#xi_{simu};#xi_{reco}", 100, 0., 0.30, 100, 0., 0.30)),
83 
84  h_de_th_x(new TH1D("", ";#theta_{x,reco} - #theta_{x,simu}", 100, 0., 0.)),
85  p_de_th_x_vs_xi_simu(new TProfile("", ";#xi_{simu};#theta_{x,reco} - #theta_{x,simu}", 50, 0., 0.25)),
86 
87  h_de_th_y(new TH1D("", ";#theta_{y,reco} - #theta_{y,simu}", 100, 0., 0.)),
88  p_de_th_y_vs_xi_simu(new TProfile("", ";#xi_{simu};#theta_{y,reco} - #theta_{y,simu}", 50, 0., 0.25)),
89 
90  h_de_vtx_y(new TH1D("", ";vtx_{y,reco} - vtx_{y,simu} (mm)", 100, 0., 0.)),
91  p_de_vtx_y_vs_xi_simu(new TProfile("", ";#xi_{simu};vtx_{y,reco} - vtx_{y,simu} (mm)", 50, 0., 0.25)),
92 
93  h_de_t(new TH1D("", ";t_{reco} - t_{simu}", 100, -1., +1.)),
94  p_de_t_vs_xi_simu(new TProfile("", ";xi_{simu};t_{reco} - t_{simu}", 50, 0., 0.25)),
95  p_de_t_vs_t_simu(new TProfile("", ";t_{simu};t_{reco} - t_{simu}", 20, 0., 5.)) {}
96 
97  static TGraphErrors profileToRMSGraph(TProfile *p, const char *name = "") {
98  TGraphErrors gr_err;
99  gr_err.SetName(name);
100 
101  for (int bi = 1; bi <= p->GetNbinsX(); ++bi) {
102  double c = p->GetBinCenter(bi);
103  double w = p->GetBinWidth(bi);
104 
105  double N = p->GetBinEntries(bi);
106  double Sy = p->GetBinContent(bi) * N;
107  double Syy = p->GetSumw2()->At(bi);
108 
109  double si_sq = Syy / N - Sy * Sy / N / N;
110  double si = (si_sq >= 0.) ? sqrt(si_sq) : 0.;
111  double si_unc_sq = si_sq / 2. / N; // Gaussian approximation
112  double si_unc = (si_unc_sq >= 0.) ? sqrt(si_unc_sq) : 0.;
113 
114  int idx = gr_err.GetN();
115  gr_err.SetPoint(idx, c, si);
116  gr_err.SetPointError(idx, w / 2., si_unc);
117  }
118 
119  return gr_err;
120  }
121 
122  void write() const {
123  h_xi_reco_vs_xi_simu->Write("h_xi_reco_vs_xi_simu");
124  h_de_xi->Write("h_de_xi");
125  p_de_xi_vs_xi_simu->Write("p_de_xi_vs_xi_simu");
126  profileToRMSGraph(p_de_xi_vs_xi_simu.get(), "g_rms_de_xi_vs_xi_simu").Write();
127 
128  h_de_th_x->Write("h_de_th_x");
129  p_de_th_x_vs_xi_simu->Write("p_de_th_x_vs_xi_simu");
130  profileToRMSGraph(p_de_th_x_vs_xi_simu.get(), "g_rms_de_th_x_vs_xi_simu").Write();
131 
132  h_de_th_y->Write("h_de_th_y");
133  p_de_th_y_vs_xi_simu->Write("p_de_th_y_vs_xi_simu");
134  profileToRMSGraph(p_de_th_y_vs_xi_simu.get(), "g_rms_de_th_y_vs_xi_simu").Write();
135 
136  h_de_vtx_y->Write("h_de_vtx_y");
137  p_de_vtx_y_vs_xi_simu->Write("p_de_vtx_y_vs_xi_simu");
138  profileToRMSGraph(p_de_vtx_y_vs_xi_simu.get(), "g_rms_de_vtx_y_vs_xi_simu").Write();
139 
140  h_de_t->Write("h_de_t");
141  p_de_t_vs_xi_simu->Write("p_de_t_vs_xi_simu");
142  profileToRMSGraph(p_de_t_vs_xi_simu.get(), "g_rms_de_t_vs_xi_simu").Write();
143  p_de_t_vs_t_simu->Write("p_de_t_vs_t_simu");
144  profileToRMSGraph(p_de_t_vs_t_simu.get(), "g_rms_de_t_vs_t_simu").Write();
145  }
146  };
147 
148  std::map<unsigned int, std::map<unsigned int, PlotGroup>> plots_;
149 
150  struct DoubleArmPlotGroup {
151  std::unique_ptr<TH2D> h2_t_sh_vs_vtx_t, h2_t_dh_vs_vtx_z;
152  std::unique_ptr<TH1D> h_t_sh_minus_vtx_t, h_t_dh_minus_vtx_z;
153 
155  : h2_t_sh_vs_vtx_t(new TH2D("", ";vtx_t (mm);(t_56 + t_45)/2 (mm)", 100, -250., -250., 100, +250., +250.)),
156  h2_t_dh_vs_vtx_z(new TH2D("", ";vtx_z (mm);(t_56 - t_45)/2 (mm)", 100, -250., -250., 100, +250., +250.)),
157  h_t_sh_minus_vtx_t(new TH1D("", ";(t_56 + t_45)/2 - vtx_t (mm)", 100, -100., +100.)),
158  h_t_dh_minus_vtx_z(new TH1D("", ";(t_56 - t_45)/2 - vtx_z (mm)", 100, -100., +100.)) {}
159 
160  void fill(double time_45, double time_56, double vtx_z, double vtx_t) {
161  const double t_sum_half = (time_56 + time_45) / 2. * CLHEP::c_light;
162  const double t_dif_half = (time_56 - time_45) / 2. * CLHEP::c_light;
163 
164  h2_t_sh_vs_vtx_t->Fill(t_sum_half, vtx_t);
165  h_t_sh_minus_vtx_t->Fill(t_sum_half - vtx_t);
166 
167  h2_t_dh_vs_vtx_z->Fill(t_dif_half, vtx_z);
168  h_t_dh_minus_vtx_z->Fill(t_dif_half - vtx_z);
169  }
170 
171  void write() const {
172  h2_t_sh_vs_vtx_t->Write("h2_t_sh_vs_vtx_t");
173  h_t_sh_minus_vtx_t->Write("h_t_sh_minus_vtx_t");
174 
175  h2_t_dh_vs_vtx_z->Write("h2_t_dh_vs_vtx_z");
176  h_t_dh_minus_vtx_z->Write("h_t_dh_minus_vtx_z");
177  }
178  };
179 
180  DoubleArmPlotGroup double_arm_plots_;
181 };
182 
183 //----------------------------------------------------------------------------------------------------
184 
185 using namespace std;
186 using namespace edm;
187 using namespace HepMC;
188 
189 //----------------------------------------------------------------------------------------------------
190 
192  const edm::ParameterSet &iConfig)
193  : tokenHepMCBeforeSmearing_(
194  consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("tagHepMCBeforeSmearing"))),
195  tokenHepMCAfterSmearing_(
196  consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("tagHepMCAfterSmearing"))),
197  tokenRecoProtonsSingleRP_(
198  consumes<reco::ForwardProtonCollection>(iConfig.getParameter<InputTag>("tagRecoProtonsSingleRP"))),
199  tokenRecoProtonsMultiRP_(
200  consumes<reco::ForwardProtonCollection>(iConfig.getParameter<InputTag>("tagRecoProtonsMultiRP"))),
201  lhcInfoLabel_(iConfig.getParameter<std::string>("lhcInfoLabel")),
202  outputFile_(iConfig.getParameter<string>("outputFile")) {}
203 
204 //----------------------------------------------------------------------------------------------------
205 
207  // get conditions
208  edm::ESHandle<LHCInfo> hLHCInfo;
209  iSetup.get<LHCInfoRcd>().get(lhcInfoLabel_, hLHCInfo);
210 
211  // get input
212  edm::Handle<edm::HepMCProduct> hHepMCBeforeSmearing;
213  iEvent.getByToken(tokenHepMCBeforeSmearing_, hHepMCBeforeSmearing);
214  HepMC::GenEvent *hepMCEventBeforeSmearing = (HepMC::GenEvent *)hHepMCBeforeSmearing->GetEvent();
215 
216  edm::Handle<edm::HepMCProduct> hHepMCAfterSmearing;
217  iEvent.getByToken(tokenHepMCAfterSmearing_, hHepMCAfterSmearing);
218  HepMC::GenEvent *hepMCEventAfterSmearing = (HepMC::GenEvent *)hHepMCAfterSmearing->GetEvent();
219 
220  Handle<reco::ForwardProtonCollection> hRecoProtonsSingleRP;
221  iEvent.getByToken(tokenRecoProtonsSingleRP_, hRecoProtonsSingleRP);
222 
223  Handle<reco::ForwardProtonCollection> hRecoProtonsMultiRP;
224  iEvent.getByToken(tokenRecoProtonsMultiRP_, hRecoProtonsMultiRP);
225 
226  // extract vertex position
227  bool vertex_set = false;
228  FourVector vtx;
229  for (auto it = hepMCEventAfterSmearing->vertices_begin(); it != hepMCEventAfterSmearing->vertices_end(); ++it) {
230  if (vertex_set) {
231  LogError("CTPPSProtonReconstructionSimulationValidator") << "Multiple vertices found.";
232  return;
233  }
234 
235  vertex_set = true;
236  vtx = (*it)->position();
237  }
238 
239  // extract forward protons
240  bool proton_45_set = false;
241  bool proton_56_set = false;
242  FourVector mom_45, mom_56;
243 
244  for (auto it = hepMCEventBeforeSmearing->particles_begin(); it != hepMCEventBeforeSmearing->particles_end(); ++it) {
245  const auto &part = *it;
246 
247  // accept only stable non-beam protons
248  if (part->pdg_id() != 2212)
249  continue;
250 
251  if (part->status() != 1)
252  continue;
253 
254  if (part->is_beam())
255  continue;
256 
257  const auto &mom = part->momentum();
258 
259  if (mom.e() < 4500.)
260  continue;
261 
262  if (mom.z() > 0) {
263  // 45
264  if (proton_45_set) {
265  LogError("CTPPSProtonReconstructionSimulationValidator") << "Found multiple protons in sector 45.";
266  return;
267  }
268 
269  proton_45_set = true;
270  mom_45 = mom;
271  } else {
272  // 56
273  if (proton_56_set) {
274  LogError("CTPPSProtonReconstructionSimulationValidator") << "Found multiple protons in sector 56.";
275  return;
276  }
277 
278  proton_56_set = true;
279  mom_56 = mom;
280  }
281  }
282 
283  // single-arm comparison
284  for (const auto &handle : {hRecoProtonsSingleRP, hRecoProtonsMultiRP}) {
285  for (const auto &rec_pr : *handle) {
286  if (!rec_pr.validFit())
287  continue;
288 
289  unsigned int idx = 12345;
290 
291  bool mom_set = false;
292  FourVector mom;
293 
295  idx = 0;
296  mom_set = proton_45_set;
297  mom = mom_45;
298  }
300  idx = 1;
301  mom_set = proton_56_set;
302  mom = mom_56;
303  }
304 
305  if (!mom_set)
306  continue;
307 
308  unsigned int meth_idx = 1234;
309 
311  meth_idx = 0;
312 
314  idx = 100 * rpId.arm() + 10 * rpId.station() + rpId.rp();
315  }
316 
318  meth_idx = 1;
319 
320  fillPlots(meth_idx, idx, rec_pr, vtx, mom, *hLHCInfo);
321  }
322  }
323 
324  // double-arm comparison
325  bool time_set_45 = false, time_set_56 = false;
326  double time_45 = 0., time_56 = 0.;
327  for (const auto &rec_pr : *hRecoProtonsMultiRP) {
328  if (!rec_pr.validFit())
329  continue;
330 
331  // skip protons with no timing information
332  if (rec_pr.timeError() < 1E-3)
333  continue;
334 
336  time_set_45 = true;
337  time_45 = rec_pr.time();
338  }
340  time_set_56 = true;
341  time_56 = rec_pr.time();
342  }
343  }
344 
345  if (time_set_45 && time_set_56)
346  double_arm_plots_.fill(time_45, time_56, vtx.z(), vtx.t());
347 }
348 
349 //----------------------------------------------------------------------------------------------------
350 
352  unsigned int idx,
353  const reco::ForwardProton &rec_pr,
354  const HepMC::FourVector &vtx,
355  const HepMC::FourVector &mom,
356  const LHCInfo &lhcInfo) {
357  const double p_nom = lhcInfo.energy();
358  const double xi_simu = (p_nom - mom.rho()) / p_nom;
359  const double th_x_simu = mom.x() / mom.rho();
360  const double th_y_simu = mom.y() / mom.rho();
361  const double vtx_y_simu = vtx.y();
362  const double th_simu = sqrt(th_x_simu * th_x_simu + th_y_simu * th_y_simu);
363  const double t_simu = -reco::ForwardProton::calculateT(p_nom, mom.rho(), th_simu);
364 
365  const double xi_reco = rec_pr.xi();
366  const double th_x_reco = rec_pr.thetaX();
367  const double th_y_reco = rec_pr.thetaY();
368  const double vtx_y_reco = rec_pr.vertex().y() * 10.; // conversion: cm --> mm
369  const double t_reco = -rec_pr.t();
370 
371  auto &plt = plots_[meth_idx][idx];
372 
373  plt.h_xi_reco_vs_xi_simu->Fill(xi_simu, xi_reco);
374  plt.h_de_xi->Fill(xi_reco - xi_simu);
375  plt.p_de_xi_vs_xi_simu->Fill(xi_simu, xi_reco - xi_simu);
376 
377  plt.h_de_th_x->Fill(th_x_reco - th_x_simu);
378  plt.p_de_th_x_vs_xi_simu->Fill(xi_simu, th_x_reco - th_x_simu);
379 
380  plt.h_de_th_y->Fill(th_y_reco - th_y_simu);
381  plt.p_de_th_y_vs_xi_simu->Fill(xi_simu, th_y_reco - th_y_simu);
382 
383  plt.h_de_vtx_y->Fill(vtx_y_reco - vtx_y_simu);
384  plt.p_de_vtx_y_vs_xi_simu->Fill(xi_simu, vtx_y_reco - vtx_y_simu);
385 
386  plt.h_de_t->Fill(t_reco - t_simu);
387  plt.p_de_t_vs_xi_simu->Fill(xi_simu, t_reco - t_simu);
388  plt.p_de_t_vs_t_simu->Fill(t_simu, t_reco - t_simu);
389 }
390 
391 //----------------------------------------------------------------------------------------------------
392 
394  auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
395 
396  for (const auto &mit : plots_) {
397  const char *method = (mit.first == 0) ? "single rp" : "multi rp";
398  TDirectory *d_method = f_out->mkdir(method);
399 
400  for (const auto &eit : mit.second) {
401  gDirectory = d_method->mkdir(Form("%i", eit.first));
402  eit.second.write();
403  }
404  }
405 
406  gDirectory = f_out->mkdir("double arm");
408 }
409 
410 //----------------------------------------------------------------------------------------------------
411 
LHCInfo::energy
const float energy() const
Definition: LHCInfo.cc:190
CTPPSProtonReconstructionSimulationValidator::endJob
void endJob() override
Definition: CTPPSProtonReconstructionSimulationValidator.cc:392
CTPPSProtonReconstructionSimulationValidator::lhcInfoLabel_
std::string lhcInfoLabel_
Definition: CTPPSProtonReconstructionSimulationValidator.cc:59
EDAnalyzer.h
CTPPSProtonReconstructionSimulationValidator::tokenHepMCAfterSmearing_
edm::EDGetTokenT< edm::HepMCProduct > tokenHepMCAfterSmearing_
Definition: CTPPSProtonReconstructionSimulationValidator.cc:54
reco::ForwardProton::ReconstructionMethod::singleRP
CTPPSProtonReconstructionSimulationValidator::PlotGroup::h_de_t
std::unique_ptr< TH1D > h_de_t
Definition: CTPPSProtonReconstructionSimulationValidator.cc:77
CTPPSProtonReconstructionSimulationValidator::PlotGroup
Definition: CTPPSProtonReconstructionSimulationValidator.cc:63
ESHandle.h
CTPPSProtonReconstructionSimulationValidator::DoubleArmPlotGroup::write
void write() const
Definition: CTPPSProtonReconstructionSimulationValidator.cc:172
patZpeak.handle
handle
Definition: patZpeak.py:23
CTPPSProtonReconstructionSimulationValidator::PlotGroup::p_de_th_x_vs_xi_simu
std::unique_ptr< TProfile > p_de_th_x_vs_xi_simu
Definition: CTPPSProtonReconstructionSimulationValidator.cc:69
edm::EDGetTokenT< edm::HepMCProduct >
edm
HLT enums.
Definition: AlignableModifier.h:19
CTPPSProtonReconstructionSimulationValidator::DoubleArmPlotGroup::DoubleArmPlotGroup
DoubleArmPlotGroup()
Definition: CTPPSProtonReconstructionSimulationValidator.cc:155
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
reco::ForwardProton::thetaY
float thetaY() const
horizontal scattering angle, in rad
Definition: ForwardProton.h:88
AlcaSiPixelAliHarvester0T_cff.method
method
Definition: AlcaSiPixelAliHarvester0T_cff.py:41
CTPPSProtonReconstructionSimulationValidator::fillPlots
void fillPlots(unsigned int meth_idx, unsigned int idx, const reco::ForwardProton &rec_pr, const HepMC::FourVector &vtx, const HepMC::FourVector &mom, const LHCInfo &lhcInfo)
Definition: CTPPSProtonReconstructionSimulationValidator.cc:350
edm::RefVector::begin
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:223
LHCInfo.h
reco::ForwardProton
Definition: ForwardProton.h:21
LHCInfo
Definition: LHCInfo.h:12
CTPPSLocalTrackLite.h
reco::ForwardProton::vertex
const Point & vertex() const
fitted vertex position
Definition: ForwardProton.h:56
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
edm::one::EDAnalyzer
Definition: EDAnalyzer.h:30
year_2016_postTS2_cff.rpId
rpId
Definition: year_2016_postTS2_cff.py:23
edm::Handle< edm::HepMCProduct >
CTPPSProtonReconstructionSimulationValidator::plots_
std::map< unsigned int, std::map< unsigned int, PlotGroup > > plots_
Definition: CTPPSProtonReconstructionSimulationValidator.cc:149
CTPPSProtonReconstructionSimulationValidator::CTPPSProtonReconstructionSimulationValidator
CTPPSProtonReconstructionSimulationValidator(const edm::ParameterSet &)
Definition: CTPPSProtonReconstructionSimulationValidator.cc:190
CTPPSProtonReconstructionSimulationValidator::PlotGroup::write
void write() const
Definition: CTPPSProtonReconstructionSimulationValidator.cc:123
HepMC::GenEvent
Definition: hepmc_rootio.cc:9
CTPPSProtonReconstructionSimulationValidator::DoubleArmPlotGroup
Definition: CTPPSProtonReconstructionSimulationValidator.cc:151
CTPPSProtonReconstructionSimulationValidator::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: CTPPSProtonReconstructionSimulationValidator.cc:205
heavyIonCSV_trainingSettings.idx
idx
Definition: heavyIonCSV_trainingSettings.py:5
MakerMacros.h
part
part
Definition: HCALResponse.h:20
edm::EventSetup::get
T get() const
Definition: EventSetup.h:80
CTPPSProtonReconstructionSimulationValidator::tokenRecoProtonsMultiRP_
edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsMultiRP_
Definition: CTPPSProtonReconstructionSimulationValidator.cc:57
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
w
const double w
Definition: UKUtility.cc:23
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
edm::ESHandle< LHCInfo >
N
#define N
Definition: blowfish.cc:9
reco::ForwardProton::time
float time() const
time of proton arrival at forward stations
Definition: ForwardProton.h:116
CTPPSProtonReconstructionSimulationValidator::DoubleArmPlotGroup::h_t_dh_minus_vtx_z
std::unique_ptr< TH1D > h_t_dh_minus_vtx_z
Definition: CTPPSProtonReconstructionSimulationValidator.cc:153
ForwardProtonFwd.h
LHCInfoRcd
Definition: LHCInfoRcd.h:24
CTPPSProtonReconstructionSimulationValidator::tokenHepMCBeforeSmearing_
edm::EDGetTokenT< edm::HepMCProduct > tokenHepMCBeforeSmearing_
Definition: CTPPSProtonReconstructionSimulationValidator.cc:53
reco::ForwardProton::contributingLocalTracks
const CTPPSLocalTrackLiteRefVector & contributingLocalTracks() const
list of RP tracks that contributed to this global track
Definition: ForwardProton.h:136
CTPPSProtonReconstructionSimulationValidator::PlotGroup::h_de_th_y
std::unique_ptr< TH1D > h_de_th_y
Definition: CTPPSProtonReconstructionSimulationValidator.cc:71
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
vtx_t
ap_int< 5 > vtx_t
Definition: data.h:8
reco::ForwardProtonCollection
std::vector< ForwardProton > ForwardProtonCollection
Collection of ForwardProton objects.
Definition: ForwardProtonFwd.h:25
edm::ParameterSet
Definition: ParameterSet.h:47
Event.h
CTPPSProtonReconstructionSimulationValidator::PlotGroup::p_de_th_y_vs_xi_simu
std::unique_ptr< TProfile > p_de_th_y_vs_xi_simu
Definition: CTPPSProtonReconstructionSimulationValidator.cc:72
ForwardProton.h
reco::ForwardProton::LHCSector::sector56
CTPPSDetId
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:31
reco::ForwardProton::xi
float xi() const
longitudinal fractional momentum loss
Definition: ForwardProton.h:84
iEvent
int iEvent
Definition: GenABIO.cc:224
CTPPSProtonReconstructionSimulationValidator
Definition: CTPPSProtonReconstructionSimulationValidator.cc:36
reco::ForwardProton::calculateT
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...
Definition: ForwardProton.cc:44
edm::EventSetup
Definition: EventSetup.h:57
edm::HepMCProduct::GetEvent
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:34
HltBtagPostValidation_cff.c
c
Definition: HltBtagPostValidation_cff.py:31
edm::LogError
Log< level::Error, false > LogError
Definition: MessageLogger.h:123
reco::ForwardProton::method
ReconstructionMethod method() const
reconstruction method for this track
Definition: ForwardProton.h:131
reco::ForwardProton::thetaX
float thetaX() const
vertical scattering angle, in rad
Definition: ForwardProton.h:86
CTPPSProtonReconstructionSimulationValidator::PlotGroup::p_de_t_vs_xi_simu
std::unique_ptr< TProfile > p_de_t_vs_xi_simu
Definition: CTPPSProtonReconstructionSimulationValidator.cc:78
get
#define get
CTPPSProtonReconstructionSimulationValidator::PlotGroup::profileToRMSGraph
static TGraphErrors profileToRMSGraph(TProfile *p, const char *name="")
Definition: CTPPSProtonReconstructionSimulationValidator.cc:98
CTPPSProtonReconstructionSimulationValidator::PlotGroup::h_de_xi
std::unique_ptr< TH1D > h_de_xi
Definition: CTPPSProtonReconstructionSimulationValidator.cc:64
CTPPSProtonReconstructionSimulationValidator::PlotGroup::h_de_th_x
std::unique_ptr< TH1D > h_de_th_x
Definition: CTPPSProtonReconstructionSimulationValidator.cc:68
reco::ForwardProton::validFit
bool validFit() const
flag for the fit validity
Definition: ForwardProton.h:126
CTPPSProtonReconstructionSimulationValidator::PlotGroup::h_de_vtx_y
std::unique_ptr< TH1D > h_de_vtx_y
Definition: CTPPSProtonReconstructionSimulationValidator.cc:74
std
Definition: JetResolutionObject.h:76
CTPPSProtonReconstructionSimulationValidator::DoubleArmPlotGroup::h2_t_dh_vs_vtx_z
std::unique_ptr< TH2D > h2_t_dh_vs_vtx_z
Definition: CTPPSProtonReconstructionSimulationValidator.cc:152
reco::ForwardProton::ReconstructionMethod::multiRP
CTPPSProtonReconstructionSimulationValidator::DoubleArmPlotGroup::h_t_sh_minus_vtx_t
std::unique_ptr< TH1D > h_t_sh_minus_vtx_t
Definition: CTPPSProtonReconstructionSimulationValidator.cc:153
extraflags_cff.vtx
vtx
Definition: extraflags_cff.py:18
CTPPSProtonReconstructionSimulationValidator::PlotGroup::p_de_t_vs_t_simu
std::unique_ptr< TProfile > p_de_t_vs_t_simu
Definition: CTPPSProtonReconstructionSimulationValidator.cc:78
Frameworkfwd.h
reco::ForwardProton::LHCSector::sector45
HepMC
Definition: GenParticle.h:15
CTPPSProtonReconstructionSimulationValidator::PlotGroup::h_xi_reco_vs_xi_simu
std::unique_ptr< TH2D > h_xi_reco_vs_xi_simu
Definition: CTPPSProtonReconstructionSimulationValidator.cc:66
CTPPSProtonReconstructionSimulationValidator::outputFile_
std::string outputFile_
Definition: CTPPSProtonReconstructionSimulationValidator.cc:61
CTPPSProtonReconstructionSimulationValidator::DoubleArmPlotGroup::fill
void fill(double time_45, double time_56, double vtx_z, double vtx_t)
Definition: CTPPSProtonReconstructionSimulationValidator.cc:161
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
EventSetup.h
CTPPSDetId.h
reco::ForwardProton::timeError
float timeError() const
uncertainty on time of proton arrival at forward stations
Definition: ForwardProton.h:120
edm::HepMCProduct
Definition: HepMCProduct.h:18
LHCInfoRcd.h
reco::ForwardProton::t
float t() const
four-momentum transfer squared, in GeV^2
Definition: ForwardProton.cc:51
HepMCProduct.h
CTPPSProtonReconstructionSimulationValidator::PlotGroup::PlotGroup
PlotGroup()
Definition: CTPPSProtonReconstructionSimulationValidator.cc:80
edm::Event
Definition: Event.h:73
CTPPSProtonReconstructionSimulationValidator::PlotGroup::p_de_vtx_y_vs_xi_simu
std::unique_ptr< TProfile > p_de_vtx_y_vs_xi_simu
Definition: CTPPSProtonReconstructionSimulationValidator.cc:75
CTPPSProtonReconstructionSimulationValidator::DoubleArmPlotGroup::h2_t_sh_vs_vtx_t
std::unique_ptr< TH2D > h2_t_sh_vs_vtx_t
Definition: CTPPSProtonReconstructionSimulationValidator.cc:152
CTPPSProtonReconstructionSimulationValidator::double_arm_plots_
DoubleArmPlotGroup double_arm_plots_
Definition: CTPPSProtonReconstructionSimulationValidator.cc:181
edm::InputTag
Definition: InputTag.h:15
CTPPSProtonReconstructionSimulationValidator::PlotGroup::p_de_xi_vs_xi_simu
std::unique_ptr< TProfile > p_de_xi_vs_xi_simu
Definition: CTPPSProtonReconstructionSimulationValidator.cc:65
CTPPSProtonReconstructionSimulationValidator::tokenRecoProtonsSingleRP_
edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsSingleRP_
Definition: CTPPSProtonReconstructionSimulationValidator.cc:56
reco::ForwardProton::lhcSector
LHCSector lhcSector() const
Definition: ForwardProton.h:140