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 {
37  public:
39 
40  private:
41  void analyze(const edm::Event&, const edm::EventSetup&) override;
42  void endJob() override;
43 
44  void fillPlots(unsigned int meth_idx, unsigned int idx, const reco::ForwardProton &rec_pr,
45  const HepMC::FourVector &vtx, const HepMC::FourVector &mom, const LHCInfo &lhcInfo);
46 
49 
52 
54 
56 
57  struct PlotGroup
58  {
59  std::unique_ptr<TH1D> h_de_xi;
60  std::unique_ptr<TProfile> p_de_xi_vs_xi_simu;
61  std::unique_ptr<TH2D> h_xi_reco_vs_xi_simu;
62 
63  std::unique_ptr<TH1D> h_de_th_x;
64  std::unique_ptr<TProfile> p_de_th_x_vs_xi_simu;
65 
66  std::unique_ptr<TH1D> h_de_th_y;
67  std::unique_ptr<TProfile> p_de_th_y_vs_xi_simu;
68 
69  std::unique_ptr<TH1D> h_de_vtx_y;
70  std::unique_ptr<TProfile> p_de_vtx_y_vs_xi_simu;
71 
72  std::unique_ptr<TH1D> h_de_t;
73  std::unique_ptr<TProfile> p_de_t_vs_xi_simu, p_de_t_vs_t_simu;
74 
76  h_de_xi(new TH1D("", ";#xi_{reco} - #xi_{simu}", 100, 0., 0.)),
77  p_de_xi_vs_xi_simu(new TProfile("", ";#xi_{simu};#xi_{reco} - #xi_{simu}", 19, 0.015, 0.205)),
78  h_xi_reco_vs_xi_simu(new TH2D("", ";#xi_{simu};#xi_{reco}", 100, 0., 0.30, 100, 0., 0.30)),
79 
80  h_de_th_x(new TH1D("", ";#theta_{x,reco} - #theta_{x,simu}", 100, 0., 0.)),
81  p_de_th_x_vs_xi_simu(new TProfile("", ";#xi_{simu};#theta_{x,reco} - #theta_{x,simu}", 19, 0.015, 0.205)),
82 
83  h_de_th_y(new TH1D("", ";#theta_{y,reco} - #theta_{y,simu}", 100, 0., 0.)),
84  p_de_th_y_vs_xi_simu(new TProfile("", ";#xi_{simu};#theta_{y,reco} - #theta_{y,simu}", 19, 0.015, 0.205)),
85 
86  h_de_vtx_y(new TH1D("", ";vtx_{y,reco} - vtx_{y,simu} (mm)", 100, 0., 0.)),
87  p_de_vtx_y_vs_xi_simu(new TProfile("", ";#xi_{simu};vtx_{y,reco} - vtx_{y,simu} (mm)", 19, 0.015, 0.205)),
88 
89  h_de_t(new TH1D("", ";t_{reco} - t_{simu}", 100, -1., +1.)),
90  p_de_t_vs_xi_simu(new TProfile("", ";xi_{simu};t_{reco} - t_{simu}", 19, 0.015, 0.205)),
91  p_de_t_vs_t_simu(new TProfile("", ";t_{simu};t_{reco} - t_{simu}", 20, 0., 5.))
92  {}
93 
94  static TGraphErrors profileToRMSGraph(TProfile *p, const char* name = "")
95  {
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  {
122  h_xi_reco_vs_xi_simu->Write("h_xi_reco_vs_xi_simu");
123  h_de_xi->Write("h_de_xi");
124  p_de_xi_vs_xi_simu->Write("p_de_xi_vs_xi_simu");
125  profileToRMSGraph(p_de_xi_vs_xi_simu.get(), "g_rms_de_xi_vs_xi_simu").Write();
126 
127  h_de_th_x->Write("h_de_th_x");
128  p_de_th_x_vs_xi_simu->Write("p_de_th_x_vs_xi_simu");
129  profileToRMSGraph(p_de_th_x_vs_xi_simu.get(), "g_rms_de_th_x_vs_xi_simu").Write();
130 
131  h_de_th_y->Write("h_de_th_y");
132  p_de_th_y_vs_xi_simu->Write("p_de_th_y_vs_xi_simu");
133  profileToRMSGraph(p_de_th_y_vs_xi_simu.get(), "g_rms_de_th_y_vs_xi_simu").Write();
134 
135  h_de_vtx_y->Write("h_de_vtx_y");
136  p_de_vtx_y_vs_xi_simu->Write("p_de_vtx_y_vs_xi_simu");
137  profileToRMSGraph(p_de_vtx_y_vs_xi_simu.get(), "g_rms_de_vtx_y_vs_xi_simu").Write();
138 
139  h_de_t->Write("h_de_t");
140  p_de_t_vs_xi_simu->Write("p_de_t_vs_xi_simu");
141  profileToRMSGraph(p_de_t_vs_xi_simu.get(), "g_rms_de_t_vs_xi_simu").Write();
142  p_de_t_vs_t_simu->Write("p_de_t_vs_t_simu");
143  profileToRMSGraph(p_de_t_vs_t_simu.get(), "g_rms_de_t_vs_t_simu").Write();
144  }
145  };
146 
147  std::map<unsigned int, std::map<unsigned int, PlotGroup>> plots_;
148 };
149 
150 //----------------------------------------------------------------------------------------------------
151 
152 using namespace std;
153 using namespace edm;
154 using namespace HepMC;
155 
156 //----------------------------------------------------------------------------------------------------
157 
159  tokenHepMCBeforeSmearing_(consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("tagHepMCBeforeSmearing"))),
160  tokenHepMCAfterSmearing_ (consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("tagHepMCAfterSmearing"))),
161  tokenRecoProtonsSingleRP_(consumes<reco::ForwardProtonCollection>(iConfig.getParameter<InputTag>("tagRecoProtonsSingleRP"))),
162  tokenRecoProtonsMultiRP_ (consumes<reco::ForwardProtonCollection>(iConfig.getParameter<InputTag>("tagRecoProtonsMultiRP"))),
163  lhcInfoLabel_(iConfig.getParameter<std::string>("lhcInfoLabel")),
164  outputFile_(iConfig.getParameter<string>("outputFile"))
165 {}
166 
167 //----------------------------------------------------------------------------------------------------
168 
170 {
171  // get conditions
172  edm::ESHandle<LHCInfo> hLHCInfo;
173  iSetup.get<LHCInfoRcd>().get(lhcInfoLabel_, hLHCInfo);
174 
175  // get input
176  edm::Handle<edm::HepMCProduct> hHepMCBeforeSmearing;
177  iEvent.getByToken(tokenHepMCBeforeSmearing_, hHepMCBeforeSmearing);
178  HepMC::GenEvent *hepMCEventBeforeSmearing = (HepMC::GenEvent *) hHepMCBeforeSmearing->GetEvent();
179 
180  edm::Handle<edm::HepMCProduct> hHepMCAfterSmearing;
181  iEvent.getByToken(tokenHepMCAfterSmearing_, hHepMCAfterSmearing);
182  HepMC::GenEvent *hepMCEventAfterSmearing = (HepMC::GenEvent *) hHepMCAfterSmearing->GetEvent();
183 
184  Handle<reco::ForwardProtonCollection> hRecoProtonsSingleRP;
185  iEvent.getByToken(tokenRecoProtonsSingleRP_, hRecoProtonsSingleRP);
186 
187  Handle<reco::ForwardProtonCollection> hRecoProtonsMultiRP;
188  iEvent.getByToken(tokenRecoProtonsMultiRP_, hRecoProtonsMultiRP);
189 
190  // extract vertex position
191  bool vertex_set = false;
192  FourVector vtx;
193  for (auto it = hepMCEventAfterSmearing->vertices_begin(); it != hepMCEventAfterSmearing->vertices_end(); ++it) {
194  if (vertex_set) {
195  LogError("CTPPSProtonReconstructionSimulationValidator") << "Multiple vertices found.";
196  return;
197  }
198 
199  vertex_set = true;
200  vtx = (*it)->position();
201  }
202 
203  // extract forward protons
204  bool proton_45_set = false;
205  bool proton_56_set = false;
206  FourVector mom_45, mom_56;
207 
208  for (auto it = hepMCEventBeforeSmearing->particles_begin(); it != hepMCEventBeforeSmearing->particles_end(); ++it) {
209  const auto &part = *it;
210 
211  // accept only stable non-beam protons
212  if (part->pdg_id() != 2212)
213  continue;
214 
215  if (part->status() != 1)
216  continue;
217 
218  if (part->is_beam())
219  continue;
220 
221  const auto &mom = part->momentum();
222 
223  if (mom.e() < 4500.)
224  continue;
225 
226  if (mom.z() > 0) {
227  // 45
228  if (proton_45_set) {
229  LogError("CTPPSProtonReconstructionSimulationValidator") << "Found multiple protons in sector 45.";
230  return;
231  }
232 
233  proton_45_set = true;
234  mom_45 = mom;
235  } else {
236  // 56
237  if (proton_56_set) {
238  LogError("CTPPSProtonReconstructionSimulationValidator") << "Found multiple protons in sector 56.";
239  return;
240  }
241 
242  proton_56_set = true;
243  mom_56 = mom;
244  }
245  }
246 
247  // do comparison
248  for (const auto& handle : { hRecoProtonsSingleRP, hRecoProtonsMultiRP } ) {
249  for (const auto& rec_pr : *handle) {
250  if (! rec_pr.validFit())
251  continue;
252 
253  unsigned int idx;
254 
255  bool mom_set = false;
256  FourVector mom;
257 
258  if (rec_pr.lhcSector() == reco::ForwardProton::LHCSector::sector45)
259  {
260  idx = 0;
261  mom_set = proton_45_set;
262  mom = mom_45;
263  } else {
264  idx = 1;
265  mom_set = proton_56_set;
266  mom = mom_56;
267  }
268 
269  if (! mom_set)
270  continue;
271 
272  unsigned int meth_idx = 1234;
273 
275  meth_idx = 0;
276 
277  CTPPSDetId rpId((*rec_pr.contributingLocalTracks().begin())->getRPId());
278  idx = 100*rpId.arm() + 10*rpId.station() + rpId.rp();
279  }
280 
282  meth_idx = 1;
283 
284  fillPlots(meth_idx, idx, rec_pr, vtx, mom, *hLHCInfo);
285  }
286  }
287 }
288 
289 //----------------------------------------------------------------------------------------------------
290 
291 void CTPPSProtonReconstructionSimulationValidator::fillPlots(unsigned int meth_idx, unsigned int idx,
292  const reco::ForwardProton &rec_pr, const HepMC::FourVector &vtx, const HepMC::FourVector &mom, const LHCInfo &lhcInfo)
293 {
294  const double p_nom = lhcInfo.energy();
295  const double xi_simu = (p_nom - mom.rho()) / p_nom;
296  const double th_x_simu = mom.x() / mom.rho();
297  const double th_y_simu = mom.y() / mom.rho();
298  const double vtx_y_simu = vtx.y();
299  const double th_simu = sqrt(th_x_simu*th_x_simu + th_y_simu*th_y_simu);
300  const double t_simu = - reco::ForwardProton::calculateT(p_nom, mom.rho(), th_simu);
301 
302  const double xi_reco = rec_pr.xi();
303  const double th_x_reco = rec_pr.thetaX();
304  const double th_y_reco = rec_pr.thetaY();
305  const double vtx_y_reco = rec_pr.vertex().y() * 10.; // conversion: cm --> mm
306  const double t_reco = - rec_pr.t();
307 
308  auto& plt = plots_[meth_idx][idx];
309 
310  plt.h_xi_reco_vs_xi_simu->Fill(xi_simu, xi_reco);
311  plt.h_de_xi->Fill(xi_reco - xi_simu);
312  plt.p_de_xi_vs_xi_simu->Fill(xi_simu, xi_reco - xi_simu);
313 
314  plt.h_de_th_x->Fill(th_x_reco - th_x_simu);
315  plt.p_de_th_x_vs_xi_simu->Fill(xi_simu, th_x_reco - th_x_simu);
316 
317  plt.h_de_th_y->Fill(th_y_reco - th_y_simu);
318  plt.p_de_th_y_vs_xi_simu->Fill(xi_simu, th_y_reco - th_y_simu);
319 
320  plt.h_de_vtx_y->Fill(vtx_y_reco - vtx_y_simu);
321  plt.p_de_vtx_y_vs_xi_simu->Fill(xi_simu, vtx_y_reco - vtx_y_simu);
322 
323  plt.h_de_t->Fill(t_reco - t_simu);
324  plt.p_de_t_vs_xi_simu->Fill(xi_simu, t_reco - t_simu);
325  plt.p_de_t_vs_t_simu->Fill(t_simu, t_reco - t_simu);
326 }
327 
328 //----------------------------------------------------------------------------------------------------
329 
331 {
332  auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
333 
334  for (const auto& mit : plots_) {
335  const char* method = (mit.first == 0) ? "single rp" : "multi rp";
336  TDirectory *d_method = f_out->mkdir(method);
337 
338  for (const auto& eit : mit.second) {
339  gDirectory = d_method->mkdir(Form("%i", eit.first));
340  eit.second.write();
341  }
342  }
343 }
344 
345 //----------------------------------------------------------------------------------------------------
346 
348 
const Point & vertex() const
fitted vertex position
Definition: ForwardProton.h:46
const double w
Definition: UKUtility.cc:23
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...
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
float t() const
four-momentum transfer squared, in GeV^2
edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsSingleRP_
float thetaX() const
vertical scattering angle, in rad
Definition: ForwardProton.h:78
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
T sqrt(T t)
Definition: SSEVec.h:18
edm::EDGetTokenT< reco::ForwardProtonCollection > tokenRecoProtonsMultiRP_
float xi() const
longitudinal fractional momentum loss
Definition: ForwardProton.h:76
#define N
Definition: blowfish.cc:9
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
static TGraphErrors profileToRMSGraph(TProfile *p, const char *name="")
part
Definition: HCALResponse.h:20
void analyze(const edm::Event &, const edm::EventSetup &) override
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
fixed size matrix
HLT enums.
T get() const
Definition: EventSetup.h:71
float const energy() const
Definition: LHCInfo.cc:192
float thetaY() const
horizontal scattering angle, in rad
Definition: ForwardProton.h:80
std::vector< ForwardProton > ForwardProtonCollection
Collection of ForwardProton objects.
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)
std::map< unsigned int, std::map< unsigned int, PlotGroup > > plots_