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)
 
 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_vtx_y_vs_xi
 
std::unique_ptr< TH1D > h_chi_sq
 
std::unique_ptr< TH1D > h_chi_sq_norm
 
std::unique_ptr< TH1D > h_n_timing_RPs
 
std::unique_ptr< TH1D > h_n_tracking_RPs
 
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_vtx_y
 
std::unique_ptr< TH1D > h_xi
 
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 99 of file CTPPSProtonReconstructionPlotter.cc.

Constructor & Destructor Documentation

CTPPSProtonReconstructionPlotter::MultiRPPlots::MultiRPPlots ( )
inline

Definition at line 106 of file CTPPSProtonReconstructionPlotter.cc.

References OrderedSet::t.

107  : h_xi(new TH1D("", ";#xi", 100, 0., 0.3)),
108  h_th_x(new TH1D("", ";#theta_{x} (rad)", 100, -500E-6, +500E-6)),
109  h_th_y(new TH1D("", ";#theta_{y} (rad)", 100, -500E-6, +500E-6)),
110  h_vtx_y(new TH1D("", ";vtx_{y} (cm)", 100, -2., +2.)),
111  h_chi_sq(new TH1D("", ";#chi^{2}", 100, 0., 0.)),
112  h_chi_sq_norm(new TH1D("", ";#chi^{2}/ndf", 100, 0., 5.)),
113  h_n_tracking_RPs(new TH1D("", ";n of tracking RPs", 4, -0.5, +3.5)),
114  h_n_timing_RPs(new TH1D("", ";n of timing RPs", 4, -0.5, +3.5)),
115  h2_th_x_vs_xi(new TH2D("", ";#xi;#theta_{x} (rad)", 100, 0., 0.3, 100, -500E-6, +500E-6)),
116  h2_th_y_vs_xi(new TH2D("", ";#xi;#theta_{y} (rad)", 100, 0., 0.3, 100, -500E-6, +500E-6)),
117  h2_vtx_y_vs_xi(new TH2D("", ";#xi;vtx_{y} (cm)", 100, 0., 0.3, 100, -500E-3, +500E-3)),
118  p_th_x_vs_xi(new TProfile("", ";#xi;#theta_{x} (rad)", 100, 0., 0.3)),
119  p_th_y_vs_xi(new TProfile("", ";#xi;#theta_{y} (rad)", 100, 0., 0.3)),
120  p_vtx_y_vs_xi(new TProfile("", ";#xi;vtx_{y} (cm)", 100, 0., 0.3)) {
121  std::vector<double> v_t_bin_edges;
122  for (double t = 0; t <= 5.;) {
123  v_t_bin_edges.push_back(t);
124  const double de_t = 0.05 + 0.09 * t + 0.02 * t * t;
125  t += de_t;
126  }
127  h_t_unif.reset(new TH1D("", ";|t| (GeV^2)", 100, 0., 5.));
128  h_t.reset(new TH1D("", ";|t| (GeV^2)", v_t_bin_edges.size() - 1, v_t_bin_edges.data()));
129  h_t_xi_range1.reset(new TH1D("", ";|t| (GeV^2)", v_t_bin_edges.size() - 1, v_t_bin_edges.data()));
130  h_t_xi_range2.reset(new TH1D("", ";|t| (GeV^2)", v_t_bin_edges.size() - 1, v_t_bin_edges.data()));
131  h_t_xi_range3.reset(new TH1D("", ";|t| (GeV^2)", v_t_bin_edges.size() - 1, v_t_bin_edges.data()));
132  h2_t_vs_xi.reset(
133  new TH2D("", ";#xi;|t| (GeV^2)", 100, 0., 0.3, v_t_bin_edges.size() - 1, v_t_bin_edges.data()));
134  }

Member Function Documentation

void CTPPSProtonReconstructionPlotter::MultiRPPlots::fill ( const reco::ForwardProton p)
inline

Definition at line 136 of file CTPPSProtonReconstructionPlotter.cc.

References reco::ForwardProton::chi2(), reco::ForwardProton::contributingLocalTracks(), 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::validFit(), reco::ForwardProton::vertex(), and reco::ForwardProton::xi().

136  {
137  if (!p.validFit())
138  return;
139 
140  unsigned int n_tracking_RPs = 0, n_timing_RPs = 0;
141  for (const auto &tr : p.contributingLocalTracks()) {
142  CTPPSDetId detId(tr->rpId());
143  if (detId.subdetId() == CTPPSDetId::sdTrackingStrip || detId.subdetId() == CTPPSDetId::sdTrackingPixel)
144  n_tracking_RPs++;
145  if (detId.subdetId() == CTPPSDetId::sdTimingDiamond || detId.subdetId() == CTPPSDetId::sdTimingFastSilicon)
146  n_timing_RPs++;
147  }
148 
149  const double th_x = p.thetaX();
150  const double th_y = p.thetaY();
151  const double mt = -p.t();
152 
153  h_chi_sq->Fill(p.chi2());
154  if (p.ndof() > 0)
155  h_chi_sq_norm->Fill(p.normalizedChi2());
156 
157  h_n_tracking_RPs->Fill(n_tracking_RPs);
158  h_n_timing_RPs->Fill(n_timing_RPs);
159 
160  h_xi->Fill(p.xi());
161 
162  h_th_x->Fill(th_x);
163  h_th_y->Fill(th_y);
164 
165  h_vtx_y->Fill(p.vertex().y());
166 
167  h_t_unif->Fill(mt);
168  h_t->Fill(mt);
169  if (p.xi() > 0.04 && p.xi() < 0.07)
170  h_t_xi_range1->Fill(mt);
171  if (p.xi() > 0.07 && p.xi() < 0.10)
172  h_t_xi_range2->Fill(mt);
173  if (p.xi() > 0.10 && p.xi() < 0.13)
174  h_t_xi_range3->Fill(mt);
175 
176  h2_th_x_vs_xi->Fill(p.xi(), th_x);
177  h2_th_y_vs_xi->Fill(p.xi(), th_y);
178  h2_vtx_y_vs_xi->Fill(p.xi(), p.vertex().y());
179  h2_t_vs_xi->Fill(p.xi(), mt);
180 
181  p_th_x_vs_xi->Fill(p.xi(), th_x);
182  p_th_y_vs_xi->Fill(p.xi(), th_y);
183  p_vtx_y_vs_xi->Fill(p.xi(), p.vertex().y());
184  }
const Point & vertex() const
fitted vertex position
Definition: ForwardProton.h:51
float t() const
four-momentum transfer squared, in GeV^2
float chi2() const
chi-squared of the fit
Definition: ForwardProton.h:72
float thetaX() const
vertical scattering angle, in rad
Definition: ForwardProton.h:81
float normalizedChi2() const
chi-squared divided by ndof (or chi-squared * 1e6 if ndof is zero)
Definition: ForwardProton.h:76
float xi() const
longitudinal fractional momentum loss
Definition: ForwardProton.h:79
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:83
unsigned int ndof() const
number of degrees of freedom for the track fit
Definition: ForwardProton.h:74
void CTPPSProtonReconstructionPlotter::MultiRPPlots::write ( ) const
inline

Definition at line 186 of file CTPPSProtonReconstructionPlotter.cc.

References CTPPSProtonReconstructionPlotter::profileToRMSGraph().

186  {
187  h_chi_sq->Write("h_chi_sq");
188  h_chi_sq_norm->Write("h_chi_sq_norm");
189 
190  h_n_tracking_RPs->Write("h_n_tracking_RPs");
191  h_n_timing_RPs->Write("h_n_timing_RPs");
192 
193  h_xi->Write("h_xi");
194 
195  h_th_x->Write("h_th_x");
196  h2_th_x_vs_xi->Write("h2_th_x_vs_xi");
197  p_th_x_vs_xi->Write("p_th_x_vs_xi");
198  auto g_th_x_RMS_vs_xi = std::make_unique<TGraphErrors>();
199  profileToRMSGraph(p_th_x_vs_xi.get(), g_th_x_RMS_vs_xi.get());
200  g_th_x_RMS_vs_xi->Write("g_th_x_RMS_vs_xi");
201 
202  h_th_y->Write("h_th_y");
203  h2_th_y_vs_xi->Write("h2_th_y_vs_xi");
204  p_th_y_vs_xi->Write("p_th_y_vs_xi");
205  auto g_th_y_RMS_vs_xi = std::make_unique<TGraphErrors>();
206  profileToRMSGraph(p_th_y_vs_xi.get(), g_th_y_RMS_vs_xi.get());
207  g_th_y_RMS_vs_xi->Write("g_th_y_RMS_vs_xi");
208 
209  h_vtx_y->Write("h_vtx_y");
210  h2_vtx_y_vs_xi->Write("h2_vtx_y_vs_xi");
211  p_vtx_y_vs_xi->Write("p_vtx_y_vs_xi");
212  auto g_vtx_y_RMS_vs_xi = std::make_unique<TGraphErrors>();
213  profileToRMSGraph(p_vtx_y_vs_xi.get(), g_vtx_y_RMS_vs_xi.get());
214  g_vtx_y_RMS_vs_xi->Write("g_vtx_y_RMS_vs_xi");
215 
216  h_t->Scale(1., "width");
217 
218  h_t_unif->Write("h_t_unif");
219  h_t->Write("h_t");
220  h_t_xi_range1->Write("h_t_xi_range1");
221  h_t_xi_range2->Write("h_t_xi_range2");
222  h_t_xi_range3->Write("h_t_xi_range3");
223 
224  h2_t_vs_xi->Write("h2_t_vs_xi");
225  }
static void profileToRMSGraph(TProfile *p, TGraphErrors *g)

Member Data Documentation

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

Definition at line 103 of file CTPPSProtonReconstructionPlotter.cc.

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

Definition at line 103 of file CTPPSProtonReconstructionPlotter.cc.

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

Definition at line 103 of file CTPPSProtonReconstructionPlotter.cc.

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

Definition at line 103 of file CTPPSProtonReconstructionPlotter.cc.

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

Definition at line 100 of file CTPPSProtonReconstructionPlotter.cc.

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

Definition at line 100 of file CTPPSProtonReconstructionPlotter.cc.

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

Definition at line 102 of file CTPPSProtonReconstructionPlotter.cc.

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

Definition at line 102 of file CTPPSProtonReconstructionPlotter.cc.

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

Definition at line 100 of file CTPPSProtonReconstructionPlotter.cc.

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

Definition at line 100 of file CTPPSProtonReconstructionPlotter.cc.

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

Definition at line 101 of file CTPPSProtonReconstructionPlotter.cc.

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

Definition at line 101 of file CTPPSProtonReconstructionPlotter.cc.

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

Definition at line 101 of file CTPPSProtonReconstructionPlotter.cc.

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

Definition at line 100 of file CTPPSProtonReconstructionPlotter.cc.

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

Definition at line 100 of file CTPPSProtonReconstructionPlotter.cc.

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

Definition at line 100 of file CTPPSProtonReconstructionPlotter.cc.

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

Definition at line 100 of file CTPPSProtonReconstructionPlotter.cc.

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

Definition at line 104 of file CTPPSProtonReconstructionPlotter.cc.

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

Definition at line 104 of file CTPPSProtonReconstructionPlotter.cc.

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

Definition at line 104 of file CTPPSProtonReconstructionPlotter.cc.