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 <map>
31 #include <string>
32 
33 //----------------------------------------------------------------------------------------------------
34 
36 public:
38 
39 private:
40  void analyze(const edm::Event &, const edm::EventSetup &) override;
41  void endJob() override;
42 
43  void fillPlots(unsigned int meth_idx,
44  unsigned int idx,
45  const reco::ForwardProton &rec_pr,
46  const HepMC::FourVector &vtx,
47  const HepMC::FourVector &mom,
48  const LHCInfo &lhcInfo);
49 
52 
55 
57 
59 
60  struct PlotGroup {
61  std::unique_ptr<TH1D> h_de_xi;
62  std::unique_ptr<TProfile> p_de_xi_vs_xi_simu;
63  std::unique_ptr<TH2D> h_xi_reco_vs_xi_simu;
64 
65  std::unique_ptr<TH1D> h_de_th_x;
66  std::unique_ptr<TProfile> p_de_th_x_vs_xi_simu;
67 
68  std::unique_ptr<TH1D> h_de_th_y;
69  std::unique_ptr<TProfile> p_de_th_y_vs_xi_simu;
70 
71  std::unique_ptr<TH1D> h_de_vtx_y;
72  std::unique_ptr<TProfile> p_de_vtx_y_vs_xi_simu;
73 
74  std::unique_ptr<TH1D> h_de_t;
75  std::unique_ptr<TProfile> p_de_t_vs_xi_simu, p_de_t_vs_t_simu;
76 
77  PlotGroup()
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)),
80  h_xi_reco_vs_xi_simu(new TH2D("", ";#xi_{simu};#xi_{reco}", 100, 0., 0.30, 100, 0., 0.30)),
81 
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)),
84 
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)),
87 
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)),
90 
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.)) {}
94 
95  static TGraphErrors profileToRMSGraph(TProfile *p, const char *name = "") {
96  TGraphErrors gr_err;
97  gr_err.SetName(name);
98 
99  for (int bi = 1; bi <= p->GetNbinsX(); ++bi) {
100  double c = p->GetBinCenter(bi);
101  double w = p->GetBinWidth(bi);
102 
103  double N = p->GetBinEntries(bi);
104  double Sy = p->GetBinContent(bi) * N;
105  double Syy = p->GetSumw2()->At(bi);
106 
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; // Gaussian approximation
110  double si_unc = (si_unc_sq >= 0.) ? sqrt(si_unc_sq) : 0.;
111 
112  int idx = gr_err.GetN();
113  gr_err.SetPoint(idx, c, si);
114  gr_err.SetPointError(idx, w / 2., si_unc);
115  }
116 
117  return gr_err;
118  }
119 
120  void write() const {
121  h_xi_reco_vs_xi_simu->Write("h_xi_reco_vs_xi_simu");
122  h_de_xi->Write("h_de_xi");
123  p_de_xi_vs_xi_simu->Write("p_de_xi_vs_xi_simu");
124  profileToRMSGraph(p_de_xi_vs_xi_simu.get(), "g_rms_de_xi_vs_xi_simu").Write();
125 
126  h_de_th_x->Write("h_de_th_x");
127  p_de_th_x_vs_xi_simu->Write("p_de_th_x_vs_xi_simu");
128  profileToRMSGraph(p_de_th_x_vs_xi_simu.get(), "g_rms_de_th_x_vs_xi_simu").Write();
129 
130  h_de_th_y->Write("h_de_th_y");
131  p_de_th_y_vs_xi_simu->Write("p_de_th_y_vs_xi_simu");
132  profileToRMSGraph(p_de_th_y_vs_xi_simu.get(), "g_rms_de_th_y_vs_xi_simu").Write();
133 
134  h_de_vtx_y->Write("h_de_vtx_y");
135  p_de_vtx_y_vs_xi_simu->Write("p_de_vtx_y_vs_xi_simu");
136  profileToRMSGraph(p_de_vtx_y_vs_xi_simu.get(), "g_rms_de_vtx_y_vs_xi_simu").Write();
137 
138  h_de_t->Write("h_de_t");
139  p_de_t_vs_xi_simu->Write("p_de_t_vs_xi_simu");
140  profileToRMSGraph(p_de_t_vs_xi_simu.get(), "g_rms_de_t_vs_xi_simu").Write();
141  p_de_t_vs_t_simu->Write("p_de_t_vs_t_simu");
142  profileToRMSGraph(p_de_t_vs_t_simu.get(), "g_rms_de_t_vs_t_simu").Write();
143  }
144  };
145 
146  std::map<unsigned int, std::map<unsigned int, PlotGroup>> plots_;
147 };
148 
149 //----------------------------------------------------------------------------------------------------
150 
151 using namespace std;
152 using namespace edm;
153 using namespace HepMC;
154 
155 //----------------------------------------------------------------------------------------------------
156 
158  const edm::ParameterSet &iConfig)
159  : tokenHepMCBeforeSmearing_(
160  consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("tagHepMCBeforeSmearing"))),
161  tokenHepMCAfterSmearing_(
162  consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("tagHepMCAfterSmearing"))),
163  tokenRecoProtonsSingleRP_(
164  consumes<reco::ForwardProtonCollection>(iConfig.getParameter<InputTag>("tagRecoProtonsSingleRP"))),
165  tokenRecoProtonsMultiRP_(
166  consumes<reco::ForwardProtonCollection>(iConfig.getParameter<InputTag>("tagRecoProtonsMultiRP"))),
167  lhcInfoLabel_(iConfig.getParameter<std::string>("lhcInfoLabel")),
168  outputFile_(iConfig.getParameter<string>("outputFile")) {}
169 
170 //----------------------------------------------------------------------------------------------------
171 
173  // get conditions
174  edm::ESHandle<LHCInfo> hLHCInfo;
175  iSetup.get<LHCInfoRcd>().get(lhcInfoLabel_, hLHCInfo);
176 
177  // get input
178  edm::Handle<edm::HepMCProduct> hHepMCBeforeSmearing;
179  iEvent.getByToken(tokenHepMCBeforeSmearing_, hHepMCBeforeSmearing);
180  HepMC::GenEvent *hepMCEventBeforeSmearing = (HepMC::GenEvent *)hHepMCBeforeSmearing->GetEvent();
181 
182  edm::Handle<edm::HepMCProduct> hHepMCAfterSmearing;
183  iEvent.getByToken(tokenHepMCAfterSmearing_, hHepMCAfterSmearing);
184  HepMC::GenEvent *hepMCEventAfterSmearing = (HepMC::GenEvent *)hHepMCAfterSmearing->GetEvent();
185 
186  Handle<reco::ForwardProtonCollection> hRecoProtonsSingleRP;
187  iEvent.getByToken(tokenRecoProtonsSingleRP_, hRecoProtonsSingleRP);
188 
189  Handle<reco::ForwardProtonCollection> hRecoProtonsMultiRP;
190  iEvent.getByToken(tokenRecoProtonsMultiRP_, hRecoProtonsMultiRP);
191 
192  // extract vertex position
193  bool vertex_set = false;
194  FourVector vtx;
195  for (auto it = hepMCEventAfterSmearing->vertices_begin(); it != hepMCEventAfterSmearing->vertices_end(); ++it) {
196  if (vertex_set) {
197  LogError("CTPPSProtonReconstructionSimulationValidator") << "Multiple vertices found.";
198  return;
199  }
200 
201  vertex_set = true;
202  vtx = (*it)->position();
203  }
204 
205  // extract forward protons
206  bool proton_45_set = false;
207  bool proton_56_set = false;
208  FourVector mom_45, mom_56;
209 
210  for (auto it = hepMCEventBeforeSmearing->particles_begin(); it != hepMCEventBeforeSmearing->particles_end(); ++it) {
211  const auto &part = *it;
212 
213  // accept only stable non-beam protons
214  if (part->pdg_id() != 2212)
215  continue;
216 
217  if (part->status() != 1)
218  continue;
219 
220  if (part->is_beam())
221  continue;
222 
223  const auto &mom = part->momentum();
224 
225  if (mom.e() < 4500.)
226  continue;
227 
228  if (mom.z() > 0) {
229  // 45
230  if (proton_45_set) {
231  LogError("CTPPSProtonReconstructionSimulationValidator") << "Found multiple protons in sector 45.";
232  return;
233  }
234 
235  proton_45_set = true;
236  mom_45 = mom;
237  } else {
238  // 56
239  if (proton_56_set) {
240  LogError("CTPPSProtonReconstructionSimulationValidator") << "Found multiple protons in sector 56.";
241  return;
242  }
243 
244  proton_56_set = true;
245  mom_56 = mom;
246  }
247  }
248 
249  // do comparison
250  for (const auto &handle : {hRecoProtonsSingleRP, hRecoProtonsMultiRP}) {
251  for (const auto &rec_pr : *handle) {
252  if (!rec_pr.validFit())
253  continue;
254 
255  unsigned int idx;
256 
257  bool mom_set = false;
258  FourVector mom;
259 
261  idx = 0;
262  mom_set = proton_45_set;
263  mom = mom_45;
264  } else {
265  idx = 1;
266  mom_set = proton_56_set;
267  mom = mom_56;
268  }
269 
270  if (!mom_set)
271  continue;
272 
273  unsigned int meth_idx = 1234;
274 
276  meth_idx = 0;
277 
279  idx = 100 * rpId.arm() + 10 * rpId.station() + rpId.rp();
280  }
281 
283  meth_idx = 1;
284 
285  fillPlots(meth_idx, idx, rec_pr, vtx, mom, *hLHCInfo);
286  }
287  }
288 }
289 
290 //----------------------------------------------------------------------------------------------------
291 
293  unsigned int idx,
294  const reco::ForwardProton &rec_pr,
295  const HepMC::FourVector &vtx,
296  const HepMC::FourVector &mom,
297  const LHCInfo &lhcInfo) {
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);
304  const double t_simu = -reco::ForwardProton::calculateT(p_nom, mom.rho(), th_simu);
305 
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.; // conversion: cm --> mm
310  const double t_reco = -rec_pr.t();
311 
312  auto &plt = plots_[meth_idx][idx];
313 
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);
317 
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);
320 
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);
323 
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);
326 
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);
330 }
331 
332 //----------------------------------------------------------------------------------------------------
333 
335  auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
336 
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);
340 
341  for (const auto &eit : mit.second) {
342  gDirectory = d_method->mkdir(Form("%i", eit.first));
343  eit.second.write();
344  }
345  }
346 }
347 
348 //----------------------------------------------------------------------------------------------------
349 
LHCInfo::energy
const float energy() const
Definition: LHCInfo.cc:190
CTPPSProtonReconstructionSimulationValidator::endJob
void endJob() override
Definition: CTPPSProtonReconstructionSimulationValidator.cc:333
CTPPSProtonReconstructionSimulationValidator::lhcInfoLabel_
std::string lhcInfoLabel_
Definition: CTPPSProtonReconstructionSimulationValidator.cc:57
EDAnalyzer.h
CTPPSProtonReconstructionSimulationValidator::tokenHepMCAfterSmearing_
edm::EDGetTokenT< edm::HepMCProduct > tokenHepMCAfterSmearing_
Definition: CTPPSProtonReconstructionSimulationValidator.cc:52
reco::ForwardProton::ReconstructionMethod::singleRP
CTPPSProtonReconstructionSimulationValidator::PlotGroup::h_de_t
std::unique_ptr< TH1D > h_de_t
Definition: CTPPSProtonReconstructionSimulationValidator.cc:75
CTPPSProtonReconstructionSimulationValidator::PlotGroup
Definition: CTPPSProtonReconstructionSimulationValidator.cc:61
ESHandle.h
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:67
edm::EDGetTokenT< edm::HepMCProduct >
edm
HLT enums.
Definition: AlignableModifier.h:19
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
reco::ForwardProton::thetaY
float thetaY() const
horizontal scattering angle, in rad
Definition: ForwardProton.h:89
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:291
edm::RefVector::begin
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:223
LHCInfo.h
reco::ForwardProton
Definition: ForwardProton.h:22
LHCInfo
Definition: LHCInfo.h:12
CTPPSLocalTrackLite.h
reco::ForwardProton::vertex
const Point & vertex() const
fitted vertex position
Definition: ForwardProton.h:57
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 >
training_settings.idx
idx
Definition: training_settings.py:16
CTPPSProtonReconstructionSimulationValidator::plots_
std::map< unsigned int, std::map< unsigned int, PlotGroup > > plots_
Definition: CTPPSProtonReconstructionSimulationValidator.cc:147
CTPPSProtonReconstructionSimulationValidator::CTPPSProtonReconstructionSimulationValidator
CTPPSProtonReconstructionSimulationValidator(const edm::ParameterSet &)
Definition: CTPPSProtonReconstructionSimulationValidator.cc:156
CTPPSProtonReconstructionSimulationValidator::PlotGroup::write
void write() const
Definition: CTPPSProtonReconstructionSimulationValidator.cc:121
HepMC::GenEvent
Definition: hepmc_rootio.cc:9
CTPPSProtonReconstructionSimulationValidator::analyze
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: CTPPSProtonReconstructionSimulationValidator.cc:171
MakerMacros.h
part
part
Definition: HCALResponse.h:20
edm::EventSetup::get
T get() const
Definition: EventSetup.h:73
CTPPSProtonReconstructionSimulationValidator::tokenRecoProtonsMultiRP_
edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsMultiRP_
Definition: CTPPSProtonReconstructionSimulationValidator.cc:55
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
ForwardProtonFwd.h
LHCInfoRcd
Definition: LHCInfoRcd.h:24
CTPPSProtonReconstructionSimulationValidator::tokenHepMCBeforeSmearing_
edm::EDGetTokenT< edm::HepMCProduct > tokenHepMCBeforeSmearing_
Definition: CTPPSProtonReconstructionSimulationValidator.cc:51
reco::ForwardProton::contributingLocalTracks
const CTPPSLocalTrackLiteRefVector & contributingLocalTracks() const
list of RP tracks that contributed to this global track
Definition: ForwardProton.h:137
CTPPSProtonReconstructionSimulationValidator::PlotGroup::h_de_th_y
std::unique_ptr< TH1D > h_de_th_y
Definition: CTPPSProtonReconstructionSimulationValidator.cc:69
badGlobalMuonTaggersAOD_cff.vtx
vtx
Definition: badGlobalMuonTaggersAOD_cff.py:5
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
reco::ForwardProtonCollection
std::vector< ForwardProton > ForwardProtonCollection
Collection of ForwardProton objects.
Definition: ForwardProtonFwd.h:25
edm::ParameterSet
Definition: ParameterSet.h:36
edm::LogError
Definition: MessageLogger.h:183
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:70
ForwardProton.h
CTPPSDetId
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:31
reco::ForwardProton::xi
float xi() const
longitudinal fractional momentum loss
Definition: ForwardProton.h:85
iEvent
int iEvent
Definition: GenABIO.cc:224
CTPPSProtonReconstructionSimulationValidator
Definition: CTPPSProtonReconstructionSimulationValidator.cc:34
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
reco::ForwardProton::method
ReconstructionMethod method() const
reconstruction method for this track
Definition: ForwardProton.h:132
reco::ForwardProton::thetaX
float thetaX() const
vertical scattering angle, in rad
Definition: ForwardProton.h:87
CTPPSProtonReconstructionSimulationValidator::PlotGroup::p_de_t_vs_xi_simu
std::unique_ptr< TProfile > p_de_t_vs_xi_simu
Definition: CTPPSProtonReconstructionSimulationValidator.cc:76
get
#define get
CTPPSProtonReconstructionSimulationValidator::PlotGroup::profileToRMSGraph
static TGraphErrors profileToRMSGraph(TProfile *p, const char *name="")
Definition: CTPPSProtonReconstructionSimulationValidator.cc:96
CTPPSProtonReconstructionSimulationValidator::PlotGroup::h_de_xi
std::unique_ptr< TH1D > h_de_xi
Definition: CTPPSProtonReconstructionSimulationValidator.cc:62
CTPPSProtonReconstructionSimulationValidator::PlotGroup::h_de_th_x
std::unique_ptr< TH1D > h_de_th_x
Definition: CTPPSProtonReconstructionSimulationValidator.cc:66
reco::ForwardProton::validFit
bool validFit() const
flag for the fit validity
Definition: ForwardProton.h:127
CTPPSProtonReconstructionSimulationValidator::PlotGroup::h_de_vtx_y
std::unique_ptr< TH1D > h_de_vtx_y
Definition: CTPPSProtonReconstructionSimulationValidator.cc:72
std
Definition: JetResolutionObject.h:76
reco::ForwardProton::ReconstructionMethod::multiRP
CTPPSProtonReconstructionSimulationValidator::PlotGroup::p_de_t_vs_t_simu
std::unique_ptr< TProfile > p_de_t_vs_t_simu
Definition: CTPPSProtonReconstructionSimulationValidator.cc:76
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:64
CTPPSProtonReconstructionSimulationValidator::outputFile_
std::string outputFile_
Definition: CTPPSProtonReconstructionSimulationValidator.cc:59
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
EventSetup.h
CTPPSDetId.h
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:78
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:73
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:63
CTPPSProtonReconstructionSimulationValidator::tokenRecoProtonsSingleRP_
edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsSingleRP_
Definition: CTPPSProtonReconstructionSimulationValidator.cc:54
reco::ForwardProton::lhcSector
LHCSector lhcSector() const
Definition: ForwardProton.h:141