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 
10 
13 
16 
18 
22 
23 #include "TFile.h"
24 #include "TH2D.h"
25 #include "TProfile.h"
26 #include "TProfile2D.h"
27 
28 #include <map>
29 
30 //----------------------------------------------------------------------------------------------------
31 
33 {
34  public:
36 
38 
39  private:
40  void analyze( const edm::Event&, const edm::EventSetup& ) override;
41  void endJob() override;
42 
44 
46 
47  struct RPPlots
48  {
49  std::unique_ptr<TH2D> h2_y_vs_x;
50  std::unique_ptr<TProfile> p_y_vs_x;
51  std::unique_ptr<TH1D> h_x;
52 
53  RPPlots() :
54  h2_y_vs_x(new TH2D("", "", 300, -10., +70., 300, -30, +30.)),
55  p_y_vs_x(new TProfile("", "", 300, -10., +70.)),
56  h_x(new TH1D("", "", 600, -10., +70.))
57  {}
58 
59  void fill(double x, double y)
60  {
61  h2_y_vs_x->Fill(x, y);
62  p_y_vs_x->Fill(x, y);
63  h_x->Fill(x);
64  }
65 
66  void write() const
67  {
68  h2_y_vs_x->Write("h2_y_vs_x");
69  p_y_vs_x->Write("p_y_vs_x");
70  h_x->Write("h_x");
71  }
72  };
73 
74  std::map<unsigned int, RPPlots> rpPlots;
75 
76 
77  struct ArmPlots
78  {
79  std::unique_ptr<TProfile2D> p2_de_x_vs_x_y, p2_de_y_vs_x_y;
80 
82  p2_de_x_vs_x_y(new TProfile2D("", ";x;y", 40, 0., 40., 40, -20., +20.)),
83  p2_de_y_vs_x_y(new TProfile2D("", ";x;y", 40, 0., 40., 40, -20., +20.))
84  {}
85 
86  void fill(double x_N, double y_N, double x_F, double y_F)
87  {
88  p2_de_x_vs_x_y->Fill(x_N, y_N, x_F - x_N);
89  p2_de_y_vs_x_y->Fill(x_N, y_N, y_F - y_N);
90  }
91 
92  void write() const
93  {
94  p2_de_x_vs_x_y->Write("p2_de_x_vs_x_y");
95  p2_de_y_vs_x_y->Write("p2_de_y_vs_x_y");
96  }
97  };
98 
99  std::map<unsigned int, ArmPlots> armPlots;
100 };
101 
102 //----------------------------------------------------------------------------------------------------
103 
105  tracksToken_( consumes<CTPPSLocalTrackLiteCollection>( iConfig.getParameter<edm::InputTag>( "tagTracks" ) ) ),
106  outputFile_( iConfig.getParameter<std::string>("outputFile") )
107 {
108 }
109 
110 //----------------------------------------------------------------------------------------------------
111 
113 {
114  // get input
116  iEvent.getByToken( tracksToken_, tracks );
117 
118  // process tracks
119  for (const auto& trk : *tracks) {
120  CTPPSDetId rpId(trk.getRPId());
121  unsigned int rpDecId = rpId.arm()*100 + rpId.station()*10 + rpId.rp();
122  rpPlots[rpDecId].fill(trk.getX(), trk.getY());
123  }
124 
125  for (const auto& t1 : *tracks) {
126  CTPPSDetId rpId1(t1.getRPId());
127 
128  for (const auto& t2 : *tracks) {
129  CTPPSDetId rpId2(t2.getRPId());
130 
131  if (rpId1.arm() != rpId2.arm())
132  continue;
133 
134  if (rpId1.station() == 0 && rpId2.station() == 2)
135  armPlots[rpId1.arm()].fill(t1.getX(), t1.getY(), t2.getX(), t2.getY());
136  }
137  }
138 }
139 
140 //----------------------------------------------------------------------------------------------------
141 
143 {
144  auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
145 
146  for (const auto& it : rpPlots) {
147  gDirectory = f_out->mkdir(Form("RP %u", it.first));
148  it.second.write();
149  }
150 
151  for (const auto& it : armPlots) {
152  gDirectory = f_out->mkdir(Form("arm %u", it.first));
153  it.second.write();
154  }
155 }
156 
157 //----------------------------------------------------------------------------------------------------
158 
160 
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
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)