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 
45  struct RPPlots {
46  std::unique_ptr<TH2D> h2_y_vs_x;
47  std::unique_ptr<TProfile> p_y_vs_x;
48  std::unique_ptr<TH1D> h_x;
49 
51  : h2_y_vs_x(new TH2D("", "", 300, -10., +70., 300, -30, +30.)),
52  p_y_vs_x(new TProfile("", "", 300, -10., +70.)),
53  h_x(new TH1D("", "", 600, -10., +70.)) {}
54 
55  void fill(double x, double y) {
56  h2_y_vs_x->Fill(x, y);
57  p_y_vs_x->Fill(x, y);
58  h_x->Fill(x);
59  }
60 
61  void write() const {
62  h2_y_vs_x->Write("h2_y_vs_x");
63  p_y_vs_x->Write("p_y_vs_x");
64  h_x->Write("h_x");
65  }
66  };
67 
68  std::map<unsigned int, RPPlots> rpPlots;
69 
70  struct ArmPlots {
71  std::unique_ptr<TProfile2D> p2_de_x_vs_x_y, p2_de_y_vs_x_y;
72 
74  : p2_de_x_vs_x_y(new TProfile2D("", ";x;y", 40, 0., 40., 40, -20., +20.)),
75  p2_de_y_vs_x_y(new TProfile2D("", ";x;y", 40, 0., 40., 40, -20., +20.)) {}
76 
77  void fill(double x_N, double y_N, double x_F, double y_F) {
78  p2_de_x_vs_x_y->Fill(x_N, y_N, x_F - x_N);
79  p2_de_y_vs_x_y->Fill(x_N, y_N, y_F - y_N);
80  }
81 
82  void write() const {
83  p2_de_x_vs_x_y->Write("p2_de_x_vs_x_y");
84  p2_de_y_vs_x_y->Write("p2_de_y_vs_x_y");
85  }
86  };
87 
88  std::map<unsigned int, ArmPlots> armPlots;
89 };
90 
91 //----------------------------------------------------------------------------------------------------
92 
94  : tracksToken_(consumes<CTPPSLocalTrackLiteCollection>(iConfig.getParameter<edm::InputTag>("tagTracks"))),
95  outputFile_(iConfig.getParameter<std::string>("outputFile")) {}
96 
97 //----------------------------------------------------------------------------------------------------
98 
100  // get input
102  iEvent.getByToken(tracksToken_, tracks);
103 
104  // process tracks
105  for (const auto& trk : *tracks) {
106  CTPPSDetId rpId(trk.rpId());
107  unsigned int rpDecId = rpId.arm() * 100 + rpId.station() * 10 + rpId.rp();
108  rpPlots[rpDecId].fill(trk.x(), trk.y());
109  }
110 
111  for (const auto& t1 : *tracks) {
112  CTPPSDetId rpId1(t1.rpId());
113 
114  for (const auto& t2 : *tracks) {
115  CTPPSDetId rpId2(t2.rpId());
116 
117  if (rpId1.arm() != rpId2.arm())
118  continue;
119 
120  if (rpId1.station() == 0 && rpId2.station() == 2)
121  armPlots[rpId1.arm()].fill(t1.x(), t1.y(), t2.x(), t2.y());
122  }
123  }
124 }
125 
126 //----------------------------------------------------------------------------------------------------
127 
129  auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
130 
131  for (const auto& it : rpPlots) {
132  gDirectory = f_out->mkdir(Form("RP %u", it.first));
133  it.second.write();
134  }
135 
136  for (const auto& it : armPlots) {
137  gDirectory = f_out->mkdir(Form("arm %u", it.first));
138  it.second.write();
139  }
140 }
141 
142 //----------------------------------------------------------------------------------------------------
143 
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
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
edm::EDGetTokenT< CTPPSLocalTrackLiteCollection > tracksToken_
uint32_t arm() const
Definition: CTPPSDetId.h:51
std::map< unsigned int, RPPlots > rpPlots
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)