CMS 3D CMS Logo

CTPPSTrackDistributionPlotter.cc
Go to the documentation of this file.
1 /****************************************************************************
2  *
3  * This is a part of CTPPS validation software
4  * Authors:
5  * Jan Kašpar
6  * Laurent Forthomme
7  *
8  ****************************************************************************/
9 
12 
15 
17 
21 
22 #include "TFile.h"
23 #include "TH2D.h"
24 #include "TProfile.h"
25 #include "TProfile2D.h"
26 #include "TGraph.h"
27 
28 #include <map>
29 
30 //----------------------------------------------------------------------------------------------------
31 
33 public:
35 
37 
38 private:
39  void analyze(const edm::Event&, const edm::EventSetup&) override;
40  void endJob() override;
41 
43 
45 
47 
48  unsigned int events_total_;
49  std::map<unsigned int, unsigned int> events_per_arm_;
50 
51  struct RPPlots {
53 
54  std::unique_ptr<TH2D> h2_y_vs_x;
55  std::unique_ptr<TProfile> p_y_vs_x;
56  std::unique_ptr<TH1D> h_x;
57  std::unique_ptr<TH1D> h_y;
58  std::unique_ptr<TH1D> h_time;
59 
60  RPPlots() : initialized(false) {}
61 
62  void init(bool pixel, double pitch) {
63  const double bin_size_x = (pixel) ? pitch * cos(18.4 / 180. * M_PI) : 100E-3;
64 
65  h2_y_vs_x.reset(new TH2D("", "", 300, -10., +70., 600, -30., +30.));
66  p_y_vs_x.reset(new TProfile("", "", 300, -10., +70.));
67 
68  int n_mi = ceil(10. / bin_size_x);
69  int n_pl = ceil(70. / bin_size_x);
70 
71  h_x.reset(new TH1D("", "", n_mi + n_pl, -n_mi * bin_size_x, +n_pl * bin_size_x));
72 
73  h_y.reset(new TH1D("", "", 300, -15., +15.));
74 
75  h_time.reset(new TH1D("", ";time", 500, -50., +50.));
76 
77  initialized = true;
78  }
79 
80  void fill(double x, double y, double time) {
81  h2_y_vs_x->Fill(x, y);
82  p_y_vs_x->Fill(x, y);
83  h_x->Fill(x);
84  h_y->Fill(y);
85  h_time->Fill(time);
86  }
87 
88  void write() const {
89  h2_y_vs_x->Write("h2_y_vs_x");
90  p_y_vs_x->Write("p_y_vs_x");
91  h_x->Write("h_x");
92  h_y->Write("h_y");
93  h_time->Write("h_time");
94  }
95  };
96 
97  std::map<unsigned int, RPPlots> rpPlots;
98 
99  struct ArmPlots {
100  unsigned int rpId_N, rpId_F;
101 
102  std::unique_ptr<TH1D> h_de_x, h_de_y;
103  std::unique_ptr<TProfile> p_de_x_vs_x, p_de_y_vs_x;
104  std::unique_ptr<TProfile2D> p2_de_x_vs_x_y, p2_de_y_vs_x_y;
105  std::unique_ptr<TH2D> h2_de_x_vs_x, h2_de_y_vs_x;
106  std::unique_ptr<TH2D> h2_de_y_vs_de_x;
107 
108  struct Stat {
109  unsigned int sN = 0, s1 = 0;
110  };
111 
112  std::map<unsigned int, std::map<unsigned int, Stat>> m_stat;
113 
115  : h_de_x(new TH1D("", ";x^{F} - x^{N}", 100, -1., +1.)),
116  h_de_y(new TH1D("", ";y^{F} - y^{N}", 100, -1., +1.)),
117  p_de_x_vs_x(new TProfile("", ";x^{N};x^{F} - x^{N}", 40, 0., 40.)),
118  p_de_y_vs_x(new TProfile("", ";x^{N};y^{F} - y^{N}", 40, 0., 40.)),
119  p2_de_x_vs_x_y(new TProfile2D("", ";x;y", 40, 0., 40., 40, -20., +20.)),
120  p2_de_y_vs_x_y(new TProfile2D("", ";x;y", 40, 0., 40., 40, -20., +20.)),
121  h2_de_x_vs_x(new TH2D("", ";x^{N};x^{F} - x^{N}", 80, 0., 40., 100, -1., +1.)),
122  h2_de_y_vs_x(new TH2D("", ";x^{N};y^{F} - y^{N}", 80, 0., 40., 100, -1., +1.)),
123  h2_de_y_vs_de_x(new TH2D("", ";x^{F} - x^{N};y^{F} - y^{N}", 100, -1., +1., 100, -1., +1.)) {}
124 
125  void fill(double x_N, double y_N, double x_F, double y_F) {
126  h_de_x->Fill(x_F - x_N);
127  h_de_y->Fill(y_F - y_N);
128 
129  p_de_x_vs_x->Fill(x_N, x_F - x_N);
130  p_de_y_vs_x->Fill(x_N, y_F - y_N);
131 
132  p2_de_x_vs_x_y->Fill(x_N, y_N, x_F - x_N);
133  p2_de_y_vs_x_y->Fill(x_N, y_N, y_F - y_N);
134 
135  h2_de_x_vs_x->Fill(x_N, x_F - x_N);
136  h2_de_y_vs_x->Fill(x_N, y_F - y_N);
137 
138  h2_de_y_vs_de_x->Fill(x_F - x_N, y_F - y_N);
139  }
140 
141  void write() const {
142  h_de_x->Write("h_de_x");
143  h_de_y->Write("h_de_y");
144 
145  p_de_x_vs_x->Write("p_de_x_vs_x");
146  p_de_y_vs_x->Write("p_de_y_vs_x");
147 
148  p2_de_x_vs_x_y->Write("p2_de_x_vs_x_y");
149  p2_de_y_vs_x_y->Write("p2_de_y_vs_x_y");
150 
151  h2_de_x_vs_x->Write("h2_de_x_vs_x");
152  h2_de_y_vs_x->Write("h2_de_y_vs_x");
153 
154  h2_de_y_vs_de_x->Write("h2_de_y_vs_de_x");
155 
156  for (const auto& rp : m_stat) {
157  TGraph* g = new TGraph();
158 
159  char buf[100];
160  sprintf(buf, "g_mean_track_mult_run_%u", rp.first);
161  g->SetName(buf);
162 
163  for (const auto& lsp : rp.second) {
164  const int idx = g->GetN();
165  const double m = (lsp.second.s1 > 0) ? double(lsp.second.sN) / lsp.second.s1 : 0.;
166  g->SetPoint(idx, lsp.first, m);
167  }
168 
169  g->Write();
170  }
171  }
172  };
173 
174  std::map<unsigned int, ArmPlots> armPlots;
175 };
176 
177 //----------------------------------------------------------------------------------------------------
178 
180  : tracksToken_(consumes<CTPPSLocalTrackLiteCollection>(iConfig.getParameter<edm::InputTag>("tagTracks"))),
181  x_pitch_pixels_(iConfig.getUntrackedParameter<double>("x_pitch_pixels", 150E-3)),
182  outputFile_(iConfig.getParameter<std::string>("outputFile")),
183  events_total_(0) {
184  armPlots[0].rpId_N = iConfig.getParameter<unsigned int>("rpId_45_N");
185  armPlots[0].rpId_F = iConfig.getParameter<unsigned int>("rpId_45_F");
186 
187  armPlots[1].rpId_N = iConfig.getParameter<unsigned int>("rpId_56_N");
188  armPlots[1].rpId_F = iConfig.getParameter<unsigned int>("rpId_56_F");
189 }
190 
191 //----------------------------------------------------------------------------------------------------
192 
194  // get input
196  iEvent.getByToken(tracksToken_, tracks);
197 
198  // process tracks
199  std::map<unsigned int, unsigned int> m_mult;
200 
201  for (const auto& trk : *tracks) {
202  CTPPSDetId rpId(trk.getRPId());
203  unsigned int rpDecId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
204  bool rpPixel = (rpId.subdetId() == CTPPSDetId::sdTrackingPixel);
205 
206  auto& pl = rpPlots[rpDecId];
207  if (!pl.initialized)
208  pl.init(rpPixel, x_pitch_pixels_);
209 
210  pl.fill(trk.getX(), trk.getY(), trk.getTime());
211 
212  m_mult[rpId.arm()]++;
213  }
214 
215  for (unsigned int arm = 0; arm < 2; ++arm) {
216  auto& st = armPlots[arm].m_stat[iEvent.id().run()][iEvent.id().luminosityBlock()];
217  st.s1++;
218  st.sN += m_mult[arm];
219  }
220 
221  for (const auto& t1 : *tracks) {
222  const CTPPSDetId rpId1(t1.getRPId());
223 
224  for (const auto& t2 : *tracks) {
225  const CTPPSDetId rpId2(t2.getRPId());
226 
227  if (rpId1.arm() != rpId2.arm())
228  continue;
229 
230  auto& ap = armPlots[rpId1.arm()];
231 
232  const unsigned int rpDecId1 = rpId1.arm() * 100 + rpId1.station() * 10 + rpId1.rp();
233  const unsigned int rpDecId2 = rpId2.arm() * 100 + rpId2.station() * 10 + rpId2.rp();
234 
235  if (rpDecId1 == ap.rpId_N && rpDecId2 == ap.rpId_F)
236  ap.fill(t1.getX(), t1.getY(), t2.getX(), t2.getY());
237  }
238  }
239 
240  // update counters
241  events_total_++;
242 
243  if (m_mult[0] > 0)
244  events_per_arm_[0]++;
245  if (m_mult[1] > 0)
246  events_per_arm_[1]++;
247 }
248 
249 //----------------------------------------------------------------------------------------------------
250 
252  edm::LogInfo("CTPPSTrackDistributionPlotter")
253  << " events processed: " << events_total_ << " (" << std::scientific << double(events_total_) << ")\n"
254  << " events with tracks in sector 45: " << events_per_arm_[0] << " (" << double(events_per_arm_[0]) << ")\n"
255  << " events with tracks in sector 56: " << events_per_arm_[1] << " (" << double(events_per_arm_[1]) << ")";
256 
257  auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
258 
259  for (const auto& it : rpPlots) {
260  gDirectory = f_out->mkdir(Form("RP %u", it.first));
261  it.second.write();
262  }
263 
264  for (const auto& it : armPlots) {
265  gDirectory = f_out->mkdir(Form("arm %u", it.first));
266  it.second.write();
267  }
268 }
269 
270 //----------------------------------------------------------------------------------------------------
271 
RunNumber_t run() const
Definition: EventID.h:39
T getParameter(std::string const &) const
std::map< unsigned int, unsigned int > events_per_arm_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
void analyze(const edm::Event &, const edm::EventSetup &) override
std::map< unsigned int, ArmPlots > armPlots
LuminosityBlockNumber_t luminosityBlock() const
Definition: EventID.h:40
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
std::map< unsigned int, std::map< unsigned int, Stat > > m_stat
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
edm::EDGetTokenT< CTPPSLocalTrackLiteCollection > tracksToken_
std::map< unsigned int, RPPlots > rpPlots
#define M_PI
std::vector< CTPPSLocalTrackLite > CTPPSLocalTrackLiteCollection
Collection of CTPPSLocalTrackLite objects.
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
edm::EventID id() const
Definition: EventBase.h:59
HLT enums.
CTPPSTrackDistributionPlotter(const edm::ParameterSet &)
void fill(double x_N, double y_N, double x_F, double y_F)
void fill(double x, double y, double time)