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 
27 #include <map>
28 
29 //----------------------------------------------------------------------------------------------------
30 
32 public:
34 
36 
37 private:
38  void analyze(const edm::Event&, const edm::EventSetup&) override;
39  void endJob() override;
40 
42 
44 
46 
47  struct RPPlots {
49 
50  std::unique_ptr<TH2D> h2_y_vs_x;
51  std::unique_ptr<TProfile> p_y_vs_x;
52  std::unique_ptr<TH1D> h_x;
53  std::unique_ptr<TH1D> h_y;
54 
55  RPPlots() : initialized(false) {}
56 
57  void init(bool pixel, double pitch) {
58  const double bin_size_x = (pixel) ? pitch * cos(18.4 / 180. * M_PI) : 100E-3;
59 
60  h2_y_vs_x.reset(new TH2D("", "", 300, -10., +70., 300, -30., +30.));
61  p_y_vs_x.reset(new TProfile("", "", 300, -10., +70.));
62 
63  int n_mi = ceil(10. / bin_size_x);
64  int n_pl = ceil(70. / bin_size_x);
65 
66  h_x.reset(new TH1D("", "", n_mi + n_pl, -n_mi * bin_size_x, +n_pl * bin_size_x));
67 
68  h_y.reset(new TH1D("", "", 300, -15., +15.));
69 
70  initialized = true;
71  }
72 
73  void fill(double x, double y) {
74  h2_y_vs_x->Fill(x, y);
75  p_y_vs_x->Fill(x, y);
76  h_x->Fill(x);
77  h_y->Fill(y);
78  }
79 
80  void write() const {
81  h2_y_vs_x->Write("h2_y_vs_x");
82  p_y_vs_x->Write("p_y_vs_x");
83  h_x->Write("h_x");
84  h_y->Write("h_y");
85  }
86  };
87 
88  std::map<unsigned int, RPPlots> rpPlots;
89 
90  struct ArmPlots {
91  std::unique_ptr<TH1D> h_de_x, h_de_y;
92  std::unique_ptr<TProfile> p_de_x_vs_x, p_de_y_vs_x;
93  std::unique_ptr<TProfile2D> p2_de_x_vs_x_y, p2_de_y_vs_x_y;
94 
96  : h_de_x(new TH1D("", ";x^{F} - x^{N}", 100, -1., +1.)),
97  h_de_y(new TH1D("", ";y^{F} - y^{N}", 100, -1., +1.)),
98  p_de_x_vs_x(new TProfile("", ";x^{N};x^{F} - x^{N}", 40, 0., 40.)),
99  p_de_y_vs_x(new TProfile("", ";x^{N};y^{F} - y^{N}", 40, 0., 40.)),
100  p2_de_x_vs_x_y(new TProfile2D("", ";x;y", 40, 0., 40., 40, -20., +20.)),
101  p2_de_y_vs_x_y(new TProfile2D("", ";x;y", 40, 0., 40., 40, -20., +20.)) {}
102 
103  void fill(double x_N, double y_N, double x_F, double y_F) {
104  h_de_x->Fill(x_F - x_N);
105  h_de_y->Fill(y_F - y_N);
106 
107  p_de_x_vs_x->Fill(x_N, x_F - x_N);
108  p_de_y_vs_x->Fill(x_N, y_F - y_N);
109 
110  p2_de_x_vs_x_y->Fill(x_N, y_N, x_F - x_N);
111  p2_de_y_vs_x_y->Fill(x_N, y_N, y_F - y_N);
112  }
113 
114  void write() const {
115  h_de_x->Write("h_de_x");
116  h_de_y->Write("h_de_y");
117 
118  p_de_x_vs_x->Write("p_de_x_vs_x");
119  p_de_y_vs_x->Write("p_de_y_vs_x");
120 
121  p2_de_x_vs_x_y->Write("p2_de_x_vs_x_y");
122  p2_de_y_vs_x_y->Write("p2_de_y_vs_x_y");
123  }
124  };
125 
126  std::map<unsigned int, ArmPlots> armPlots;
127 };
128 
129 //----------------------------------------------------------------------------------------------------
130 
132  : tracksToken_(consumes<CTPPSLocalTrackLiteCollection>(iConfig.getParameter<edm::InputTag>("tagTracks"))),
133  x_pitch_pixels_(iConfig.getUntrackedParameter<double>("x_pitch_pixels", 150E-3)),
134  outputFile_(iConfig.getParameter<std::string>("outputFile")) {}
135 
136 //----------------------------------------------------------------------------------------------------
137 
139  // get input
141  iEvent.getByToken(tracksToken_, tracks);
142 
143  // process tracks
144  for (const auto& trk : *tracks) {
145  CTPPSDetId rpId(trk.getRPId());
146  unsigned int rpDecId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
147  bool rpPixel = (rpId.subdetId() == CTPPSDetId::sdTrackingPixel);
148 
149  auto& pl = rpPlots[rpDecId];
150  if (!pl.initialized)
151  pl.init(rpPixel, x_pitch_pixels_);
152 
153  pl.fill(trk.getX(), trk.getY());
154  }
155 
156  for (const auto& t1 : *tracks) {
157  CTPPSDetId rpId1(t1.getRPId());
158 
159  for (const auto& t2 : *tracks) {
160  CTPPSDetId rpId2(t2.getRPId());
161 
162  if (rpId1.arm() != rpId2.arm())
163  continue;
164 
165  if (rpId1.station() == 0 && rpId2.station() == 2)
166  armPlots[rpId1.arm()].fill(t1.getX(), t1.getY(), t2.getX(), t2.getY());
167  }
168  }
169 }
170 
171 //----------------------------------------------------------------------------------------------------
172 
174  auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
175 
176  for (const auto& it : rpPlots) {
177  gDirectory = f_out->mkdir(Form("RP %u", it.first));
178  it.second.write();
179  }
180 
181  for (const auto& it : armPlots) {
182  gDirectory = f_out->mkdir(Form("arm %u", it.first));
183  it.second.write();
184  }
185 }
186 
187 //----------------------------------------------------------------------------------------------------
188 
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
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
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
HLT enums.
CTPPSTrackDistributionPlotter(const edm::ParameterSet &)
void fill(double x_N, double y_N, double x_F, double y_F)