CMS 3D CMS Logo

List of all members | Public Member Functions | Public Attributes
CTPPSProtonReconstructionPlotter::MultiRPPlots Struct Reference

Public Member Functions

void fill (const reco::ForwardProton &p, unsigned int nTracks, bool n1f1)
 
 MultiRPPlots ()
 
void write () const
 

Public Attributes

std::unique_ptr< TH2D > h2_t_vs_xi
 
std::unique_ptr< TH2D > h2_th_x_vs_xi
 
std::unique_ptr< TH2D > h2_th_y_vs_xi
 
std::unique_ptr< TH2D > h2_timing_tracks_vs_prot_mult
 
std::unique_ptr< TH2D > h2_vtx_y_vs_xi
 
std::unique_ptr< TH2D > h2_y_vs_x_tt0_ClCo
 
std::unique_ptr< TH2D > h2_y_vs_x_tt1_ClCo
 
std::unique_ptr< TH2D > h2_y_vs_x_ttm_ClCo
 
std::unique_ptr< TH1D > h_chi_sq
 
std::unique_ptr< TH1D > h_chi_sq_norm
 
std::unique_ptr< TH1D > h_de_x_match_timing_vs_tracking
 
std::unique_ptr< TH1D > h_de_x_match_timing_vs_tracking_ClCo
 
std::unique_ptr< TH1D > h_de_x_rel_timing_vs_tracking
 
std::unique_ptr< TH1D > h_de_x_rel_timing_vs_tracking_ClCo
 
std::unique_ptr< TH1D > h_de_x_timing_vs_tracking
 
std::unique_ptr< TH1D > h_de_x_timing_vs_tracking_ClCo
 
std::unique_ptr< TH1D > h_log_chi_sq
 
std::unique_ptr< TH1D > h_multiplicity
 
std::unique_ptr< TH1D > h_n_contrib_timing_tracks
 
std::unique_ptr< TH1D > h_n_contrib_tracking_tracks
 
std::unique_ptr< TH1D > h_t
 
std::unique_ptr< TH1D > h_t_unif
 
std::unique_ptr< TH1D > h_t_xi_range1
 
std::unique_ptr< TH1D > h_t_xi_range2
 
std::unique_ptr< TH1D > h_t_xi_range3
 
std::unique_ptr< TH1D > h_th_x
 
std::unique_ptr< TH1D > h_th_y
 
std::unique_ptr< TH1D > h_time
 
std::unique_ptr< TH1D > h_vtx_y
 
std::unique_ptr< TH1D > h_xi
 
std::unique_ptr< TH1D > h_xi_n1f1
 
std::map< unsigned int, TH1D * > m_h_xi_nTracks
 
std::unique_ptr< TProfile > p_th_x_vs_xi
 
std::unique_ptr< TProfile > p_th_y_vs_xi
 
std::unique_ptr< TProfile > p_vtx_y_vs_xi
 

Detailed Description

Definition at line 150 of file CTPPSProtonReconstructionPlotter.cc.

Constructor & Destructor Documentation

CTPPSProtonReconstructionPlotter::MultiRPPlots::MultiRPPlots ( )
inline

Definition at line 170 of file CTPPSProtonReconstructionPlotter.cc.

References gen::n, and lumiQTWidget::t.

171  : h_multiplicity(new TH1D("", ";reconstructed protons per event", 11, -0.5, 10.5)),
172  h_xi(new TH1D("", ";#xi", 100, 0., 0.3)),
173  h_th_x(new TH1D("", ";#theta_{x} (rad)", 250, -500E-6, +500E-6)),
174  h_th_y(new TH1D("", ";#theta_{y} (rad)", 500, -1000E-6, +1000E-6)),
175  h_vtx_y(new TH1D("", ";vtx_{y} (cm)", 100, -100E-3, +100E-3)),
176  h_chi_sq(new TH1D("", ";#chi^{2}", 100, 0., 10.)),
177  h_log_chi_sq(new TH1D("", ";log_{10} #chi^{2}", 100, -20., 5.)),
178  h_chi_sq_norm(new TH1D("", ";#chi^{2}/ndf", 100, 0., 5.)),
179  h_time(new TH1D("", ";time", 100, -2., +2.)),
180  h_n_contrib_tracking_tracks(new TH1D("", ";n of contrib. tracking tracks per reco proton", 4, -0.5, +3.5)),
181  h_n_contrib_timing_tracks(new TH1D("", ";n of contrib. timing tracks per reco proton", 4, -0.5, +3.5)),
182  h2_th_x_vs_xi(new TH2D("", ";#xi;#theta_{x} (rad)", 100, 0., 0.3, 100, -500E-6, +500E-6)),
183  h2_th_y_vs_xi(new TH2D("", ";#xi;#theta_{y} (rad)", 100, 0., 0.3, 200, -1000E-6, +1000E-6)),
184  h2_vtx_y_vs_xi(new TH2D("", ";#xi;vtx_{y} (cm)", 100, 0., 0.3, 100, -100E-3, +100E-3)),
185  p_th_x_vs_xi(new TProfile("", ";#xi;#theta_{x} (rad)", 100, 0., 0.3)),
186  p_th_y_vs_xi(new TProfile("", ";#xi;#theta_{y} (rad)", 100, 0., 0.3)),
187  p_vtx_y_vs_xi(new TProfile("", ";#xi;vtx_{y} (cm)", 100, 0., 0.3)),
189  new TH2D("", ";reco protons per event;timing tracks per event", 11, -0.5, 10.5, 11, -0.5, 10.5)),
190  h_xi_n1f1(new TH1D("", ";#xi", 100, 0., 0.3)),
191 
192  h_de_x_timing_vs_tracking(new TH1D("", ";#Delta x (mm)", 200, -1., +1.)),
193  h_de_x_rel_timing_vs_tracking(new TH1D("", ";#Delta x / #sigma(x)", 200, -20., +20.)),
194  h_de_x_match_timing_vs_tracking(new TH1D("", ";match between tracking and timing tracks", 2, -0.5, +1.5)),
195 
196  h_de_x_timing_vs_tracking_ClCo(new TH1D("", ";#Delta x (mm)", 200, -1., +1.)),
197  h_de_x_rel_timing_vs_tracking_ClCo(new TH1D("", ";#Delta x / #sigma(x)", 200, -20., +20.)),
199  new TH1D("", ";match between tracking and timing tracks", 2, -0.5, +1.5)),
200 
201  h2_y_vs_x_tt0_ClCo(new TH2D("", ";x (mm);y (mm)", 100, -5., 25., 100, -15., +15.)),
202  h2_y_vs_x_tt1_ClCo(new TH2D("", ";x (mm);y (mm)", 100, -5., 25., 100, -15., +15.)),
203  h2_y_vs_x_ttm_ClCo(new TH2D("", ";x (mm);y (mm)", 100, -5., 25., 100, -15., +15.)) {
204  std::vector<double> v_t_bin_edges;
205  for (double t = 0; t <= 5.;) {
206  v_t_bin_edges.push_back(t);
207  const double de_t = 0.05 + 0.09 * t + 0.02 * t * t;
208  t += de_t;
209  }
210  h_t_unif.reset(new TH1D("", ";|t| (GeV^2)", 100, 0., 5.));
211  h_t.reset(new TH1D("", ";|t| (GeV^2)", v_t_bin_edges.size() - 1, v_t_bin_edges.data()));
212  h_t_xi_range1.reset(new TH1D("", ";|t| (GeV^2)", v_t_bin_edges.size() - 1, v_t_bin_edges.data()));
213  h_t_xi_range2.reset(new TH1D("", ";|t| (GeV^2)", v_t_bin_edges.size() - 1, v_t_bin_edges.data()));
214  h_t_xi_range3.reset(new TH1D("", ";|t| (GeV^2)", v_t_bin_edges.size() - 1, v_t_bin_edges.data()));
215  h2_t_vs_xi.reset(
216  new TH2D("", ";#xi;|t| (GeV^2)", 100, 0., 0.3, v_t_bin_edges.size() - 1, v_t_bin_edges.data()));
217 
218  for (unsigned int n = 2; n <= 10; ++n)
219  m_h_xi_nTracks[n] = new TH1D(*h_xi);
220  }

Member Function Documentation

void CTPPSProtonReconstructionPlotter::MultiRPPlots::fill ( const reco::ForwardProton p,
unsigned int  nTracks,
bool  n1f1 
)
inline

Definition at line 222 of file CTPPSProtonReconstructionPlotter.cc.

References reco::ForwardProton::chi2(), reco::ForwardProton::contributingLocalTracks(), CTPPSLocalTrackLite::getRPId(), TtSemiLepEvtBuilder_cfi::mt, reco::ForwardProton::ndof(), reco::ForwardProton::normalizedChi2(), CTPPSDetId::sdTimingDiamond, CTPPSDetId::sdTimingFastSilicon, CTPPSDetId::sdTrackingPixel, CTPPSDetId::sdTrackingStrip, reco::ForwardProton::t(), reco::ForwardProton::thetaX(), reco::ForwardProton::thetaY(), reco::ForwardProton::time(), reco::ForwardProton::validFit(), reco::ForwardProton::vertex(), and reco::ForwardProton::xi().

222  {
223  if (!p.validFit())
224  return;
225 
226  unsigned int n_contrib_tracking_tracks = 0, n_contrib_timing_tracks = 0;
227  for (const auto &tr : p.contributingLocalTracks()) {
228  CTPPSDetId detId(tr->getRPId());
229  if (detId.subdetId() == CTPPSDetId::sdTrackingStrip || detId.subdetId() == CTPPSDetId::sdTrackingPixel)
230  n_contrib_tracking_tracks++;
231  if (detId.subdetId() == CTPPSDetId::sdTimingDiamond || detId.subdetId() == CTPPSDetId::sdTimingFastSilicon)
232  n_contrib_timing_tracks++;
233  }
234 
235  const double th_x = p.thetaX();
236  const double th_y = p.thetaY();
237  const double mt = -p.t();
238 
239  h_chi_sq->Fill(p.chi2());
240  h_log_chi_sq->Fill(log10(p.chi2()));
241  if (p.ndof() > 0)
242  h_chi_sq_norm->Fill(p.normalizedChi2());
243 
244  h_n_contrib_tracking_tracks->Fill(n_contrib_tracking_tracks);
245  h_n_contrib_timing_tracks->Fill(n_contrib_timing_tracks);
246 
247  h_xi->Fill(p.xi());
248 
249  h_th_x->Fill(th_x);
250  h_th_y->Fill(th_y);
251 
252  h_vtx_y->Fill(p.vertex().y());
253 
254  h_t_unif->Fill(mt);
255  h_t->Fill(mt);
256  if (p.xi() > 0.04 && p.xi() < 0.07)
257  h_t_xi_range1->Fill(mt);
258  if (p.xi() > 0.07 && p.xi() < 0.10)
259  h_t_xi_range2->Fill(mt);
260  if (p.xi() > 0.10 && p.xi() < 0.13)
261  h_t_xi_range3->Fill(mt);
262 
263  h_time->Fill(p.time());
264 
265  h2_th_x_vs_xi->Fill(p.xi(), th_x);
266  h2_th_y_vs_xi->Fill(p.xi(), th_y);
267  h2_vtx_y_vs_xi->Fill(p.xi(), p.vertex().y());
268  h2_t_vs_xi->Fill(p.xi(), mt);
269 
270  p_th_x_vs_xi->Fill(p.xi(), th_x);
271  p_th_y_vs_xi->Fill(p.xi(), th_y);
272  p_vtx_y_vs_xi->Fill(p.xi(), p.vertex().y());
273 
274  auto it = m_h_xi_nTracks.find(nTracks);
275  if (it != m_h_xi_nTracks.end())
276  it->second->Fill(p.xi());
277 
278  if (n1f1)
279  h_xi_n1f1->Fill(p.xi());
280  }
const Point & vertex() const
fitted vertex position
Definition: ForwardProton.h:46
const unsigned int nTracks(const reco::Vertex &sv)
float t() const
four-momentum transfer squared, in GeV^2
float chi2() const
chi-squared of the fit
Definition: ForwardProton.h:67
float thetaX() const
vertical scattering angle, in rad
Definition: ForwardProton.h:78
float normalizedChi2() const
chi-squared divided by ndof (or chi-squared * 1e6 if ndof is zero)
Definition: ForwardProton.h:71
float xi() const
longitudinal fractional momentum loss
Definition: ForwardProton.h:76
const CTPPSLocalTrackLiteRefVector & contributingLocalTracks() const
list of RP tracks that contributed to this global track
bool validFit() const
flag for the fit validity
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
float thetaY() const
horizontal scattering angle, in rad
Definition: ForwardProton.h:80
float time() const
time of proton arrival at forward stations
unsigned int ndof() const
number of degrees of freedom for the track fit
Definition: ForwardProton.h:69
void CTPPSProtonReconstructionPlotter::MultiRPPlots::write ( ) const
inline

Definition at line 282 of file CTPPSProtonReconstructionPlotter.cc.

References AlCaHLTBitMon_ParallelJobs::p, and CTPPSProtonReconstructionPlotter::profileToRMSGraph().

282  {
283  h_multiplicity->Write("h_multiplicity");
284 
285  h_chi_sq->Write("h_chi_sq");
286  h_log_chi_sq->Write("h_log_chi_sq");
287  h_chi_sq_norm->Write("h_chi_sq_norm");
288 
289  h_n_contrib_tracking_tracks->Write("h_n_contrib_tracking_tracks");
290  h_n_contrib_timing_tracks->Write("h_n_contrib_timing_tracks");
291 
292  h2_timing_tracks_vs_prot_mult->Write("h2_timing_tracks_vs_prot_mult");
293 
294  h_xi->Write("h_xi");
295 
296  h_th_x->Write("h_th_x");
297  h2_th_x_vs_xi->Write("h2_th_x_vs_xi");
298  p_th_x_vs_xi->Write("p_th_x_vs_xi");
299  auto g_th_x_RMS_vs_xi = std::make_unique<TGraphErrors>();
300  profileToRMSGraph(p_th_x_vs_xi.get(), g_th_x_RMS_vs_xi.get());
301  g_th_x_RMS_vs_xi->Write("g_th_x_RMS_vs_xi");
302 
303  h_th_y->Write("h_th_y");
304  h2_th_y_vs_xi->Write("h2_th_y_vs_xi");
305  p_th_y_vs_xi->Write("p_th_y_vs_xi");
306  auto g_th_y_RMS_vs_xi = std::make_unique<TGraphErrors>();
307  profileToRMSGraph(p_th_y_vs_xi.get(), g_th_y_RMS_vs_xi.get());
308  g_th_y_RMS_vs_xi->Write("g_th_y_RMS_vs_xi");
309 
310  h_vtx_y->Write("h_vtx_y");
311  h2_vtx_y_vs_xi->Write("h2_vtx_y_vs_xi");
312  p_vtx_y_vs_xi->Write("p_vtx_y_vs_xi");
313  auto g_vtx_y_RMS_vs_xi = std::make_unique<TGraphErrors>();
314  profileToRMSGraph(p_vtx_y_vs_xi.get(), g_vtx_y_RMS_vs_xi.get());
315  g_vtx_y_RMS_vs_xi->Write("g_vtx_y_RMS_vs_xi");
316 
317  h_t->Scale(1., "width");
318 
319  h_t_unif->Write("h_t_unif");
320  h_t->Write("h_t");
321  h_t_xi_range1->Write("h_t_xi_range1");
322  h_t_xi_range2->Write("h_t_xi_range2");
323  h_t_xi_range3->Write("h_t_xi_range3");
324 
325  h2_t_vs_xi->Write("h2_t_vs_xi");
326 
327  h_time->Write("h_time");
328 
329  TDirectory *d_top = gDirectory;
330 
331  gDirectory = d_top->mkdir("h_xi_nTracks");
332  for (const auto p : m_h_xi_nTracks) {
333  char buf[100];
334  sprintf(buf, "h_xi_nTracks_%u", p.first);
335  p.second->Write(buf);
336  }
337 
338  gDirectory = d_top;
339 
340  h_xi_n1f1->Write("h_xi_n1f1");
341 
342  h_de_x_timing_vs_tracking->Write("h_de_x_timing_vs_tracking");
343  h_de_x_rel_timing_vs_tracking->Write("h_de_x_rel_timing_vs_tracking");
344  h_de_x_match_timing_vs_tracking->Write("h_de_x_match_timing_vs_tracking");
345 
346  h_de_x_timing_vs_tracking_ClCo->Write("h_de_x_timing_vs_tracking_ClCo");
347  h_de_x_rel_timing_vs_tracking_ClCo->Write("h_de_x_rel_timing_vs_tracking_ClCo");
348  h_de_x_match_timing_vs_tracking_ClCo->Write("h_de_x_match_timing_vs_tracking_ClCo");
349 
350  h2_y_vs_x_tt0_ClCo->Write("h2_y_vs_x_tt0_ClCo");
351  h2_y_vs_x_tt1_ClCo->Write("h2_y_vs_x_tt1_ClCo");
352  h2_y_vs_x_ttm_ClCo->Write("h2_y_vs_x_ttm_ClCo");
353  }
static void profileToRMSGraph(TProfile *p, TGraphErrors *g)

Member Data Documentation

std::unique_ptr<TH2D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h2_t_vs_xi

Definition at line 156 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH2D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h2_th_x_vs_xi

Definition at line 156 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH2D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h2_th_y_vs_xi

Definition at line 156 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH2D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h2_timing_tracks_vs_prot_mult

Definition at line 159 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH2D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h2_vtx_y_vs_xi

Definition at line 156 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH2D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h2_y_vs_x_tt0_ClCo

Definition at line 168 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH2D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h2_y_vs_x_tt1_ClCo

Definition at line 168 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH2D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h2_y_vs_x_ttm_ClCo

Definition at line 168 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_chi_sq

Definition at line 152 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_chi_sq_norm

Definition at line 152 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_de_x_match_timing_vs_tracking

Definition at line 164 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_de_x_match_timing_vs_tracking_ClCo

Definition at line 165 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_de_x_rel_timing_vs_tracking

Definition at line 164 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_de_x_rel_timing_vs_tracking_ClCo

Definition at line 165 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_de_x_timing_vs_tracking

Definition at line 164 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_de_x_timing_vs_tracking_ClCo

Definition at line 165 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_log_chi_sq

Definition at line 152 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_multiplicity

Definition at line 151 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_n_contrib_timing_tracks

Definition at line 155 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_n_contrib_tracking_tracks

Definition at line 155 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_t

Definition at line 152 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_t_unif

Definition at line 152 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_t_xi_range1

Definition at line 153 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_t_xi_range2

Definition at line 153 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_t_xi_range3

Definition at line 153 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_th_x

Definition at line 152 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_th_y

Definition at line 152 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_time

Definition at line 154 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_vtx_y

Definition at line 152 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_xi

Definition at line 152 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TH1D> CTPPSProtonReconstructionPlotter::MultiRPPlots::h_xi_n1f1

Definition at line 162 of file CTPPSProtonReconstructionPlotter.cc.

std::map<unsigned int, TH1D *> CTPPSProtonReconstructionPlotter::MultiRPPlots::m_h_xi_nTracks

Definition at line 161 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TProfile> CTPPSProtonReconstructionPlotter::MultiRPPlots::p_th_x_vs_xi

Definition at line 157 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TProfile> CTPPSProtonReconstructionPlotter::MultiRPPlots::p_th_y_vs_xi

Definition at line 157 of file CTPPSProtonReconstructionPlotter.cc.

std::unique_ptr<TProfile> CTPPSProtonReconstructionPlotter::MultiRPPlots::p_vtx_y_vs_xi

Definition at line 157 of file CTPPSProtonReconstructionPlotter.cc.