CMS 3D CMS Logo

CTPPSAcceptancePlotter.cc
Go to the documentation of this file.
1 /****************************************************************************
2  * Authors:
3  * Jan Kašpar
4  ****************************************************************************/
5 
8 
11 
13 
17 
18 #include "TFile.h"
19 #include "TH1D.h"
20 #include "TH2D.h"
21 
22 #include <map>
23 #include <set>
24 #include <string>
25 
26 //----------------------------------------------------------------------------------------------------
27 
29 {
30  public:
31  explicit CTPPSAcceptancePlotter( const edm::ParameterSet& );
32 
33  private:
34  void beginJob() override {}
35  void analyze(const edm::Event&, const edm::EventSetup&) override;
36  void endJob() override;
37 
40 
42 
44 
46  {
47  std::unique_ptr<TH1D> h_xi_all, h_xi_acc;
49  h_xi_all(new TH1D("", ";#xi", 100, 0., 0.25)),
50  h_xi_acc(new TH1D("", ";#xi", 100, 0., 0.25))
51  {}
52 
53  void fill(double xi, bool acc)
54  {
55  h_xi_all->Fill(xi);
56  if (acc)
57  h_xi_acc->Fill(xi);
58  }
59 
60  void write() const
61  {
62  h_xi_all->Write("h_xi_all");
63  h_xi_acc->Write("h_xi_acc");
64 
65  auto h_xi_rat = std::make_unique<TH1D>(*h_xi_acc);
66  h_xi_rat->Divide(h_xi_all.get());
67  h_xi_rat->Write("h_xi_rat");
68  }
69  };
70 
71  std::vector<std::set<unsigned int>> singleArmConfigurations;
72  std::map<std::set<unsigned int>, SingleArmPlots> singleArmPlots;
73 
75  {
76  std::unique_ptr<TH1D> h_m_all, h_m_acc;
77  std::unique_ptr<TH2D> h2_xi_45_vs_xi_56_all, h2_xi_45_vs_xi_56_acc;
78  std::unique_ptr<TH2D> h2_y_vs_m_all, h2_y_vs_m_acc;
79 
81  h_m_all(new TH1D("", ";m (GeV)", 100, 0., 2500.)),
82  h_m_acc(new TH1D("", ";m (GeV)", 100, 0., 2500.)),
83  h2_xi_45_vs_xi_56_all(new TH2D("", ";xi_56;xi_45", 25, 0., 0.25, 25, 0., 0.25)),
84  h2_xi_45_vs_xi_56_acc(new TH2D("", ";xi_56;xi_45", 25, 0., 0.25, 25, 0., 0.25)),
85  h2_y_vs_m_all(new TH2D("", ";m (GeV);y", 25, 0., 2500., 25, -1.5, +1.5)),
86  h2_y_vs_m_acc(new TH2D("", ";m (GeV);y", 25, 0., 2500., 25, -1.5, +1.5))
87  {}
88 
89  void fill(double xi_45, double xi_56, bool acc)
90  {
91  const double p_nom = 6500.;
92  const double m = 2. * p_nom * sqrt(xi_45 * xi_56);
93  const double y = log(xi_45 / xi_56) / 2.;
94 
95  h_m_all->Fill(m);
96  h2_xi_45_vs_xi_56_all->Fill(xi_56, xi_45);
97  h2_y_vs_m_all->Fill(m, y);
98 
99  if (acc) {
100  h_m_acc->Fill(m);
101  h2_xi_45_vs_xi_56_acc->Fill(xi_56, xi_45);
102  h2_y_vs_m_acc->Fill(m, y);
103  }
104  }
105 
106  void write() const
107  {
108  h_m_all->Write("h_m_all");
109  h_m_acc->Write("h_m_acc");
110 
111  auto h_m_rat = std::make_unique<TH1D>(*h_m_acc);
112  h_m_rat->Divide(h_m_all.get());
113  h_m_rat->Write("h_m_rat");
114 
115  h2_xi_45_vs_xi_56_all->Write("h2_xi_45_vs_xi_56_all");
116  h2_xi_45_vs_xi_56_acc->Write("h2_xi_45_vs_xi_56_acc");
117 
118  auto h2_xi_45_vs_xi_56_rat = std::make_unique<TH2D>(*h2_xi_45_vs_xi_56_acc);
119  h2_xi_45_vs_xi_56_rat->Divide(h2_xi_45_vs_xi_56_all.get());
120  h2_xi_45_vs_xi_56_rat->Write("h2_xi_45_vs_xi_56_rat");
121 
122  h2_y_vs_m_all->Write("h2_y_vs_m_all");
123  h2_y_vs_m_acc->Write("h2_y_vs_m_acc");
124 
125  auto h2_y_vs_m_rat = std::make_unique<TH2D>(*h2_y_vs_m_acc);
126  h2_y_vs_m_rat->Divide(h2_y_vs_m_all.get());
127  h2_y_vs_m_rat->Write("h2_y_vs_m_rat");
128  }
129  };
130 
131  std::vector<std::set<unsigned int>> doubleArmConfigurations;
132  std::map<std::set<unsigned int>, DoubleArmPlots> doubleArmPlots;
133 };
134 
135 //----------------------------------------------------------------------------------------------------
136 
137 using namespace std;
138 using namespace edm;
139 using namespace HepMC;
140 
141 //----------------------------------------------------------------------------------------------------
142 
144  tokenHepMC_ (consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("tagHepMC"))),
145  tokenTracks_(consumes<CTPPSLocalTrackLiteCollection>(iConfig.getParameter<edm::InputTag>("tagTracks"))),
146  rpId_45_N_ (iConfig.getParameter<unsigned int>("rpId_45_N")),
147  rpId_45_F_ (iConfig.getParameter<unsigned int>("rpId_45_F")),
148  rpId_56_N_ (iConfig.getParameter<unsigned int>("rpId_56_N")),
149  rpId_56_F_ (iConfig.getParameter<unsigned int>("rpId_56_F")),
150  outputFile_(iConfig.getParameter<string>("outputFile"))
151 {
153  { rpId_45_N_ },
154  { rpId_45_F_ },
155  { rpId_56_N_ },
156  { rpId_56_F_ },
157  { rpId_45_N_, rpId_45_F_ },
158  { rpId_56_N_, rpId_56_F_ },
159  };
160 
162  { rpId_45_N_, rpId_56_N_ },
163  { rpId_45_F_, rpId_56_F_ },
164  { rpId_45_N_, rpId_45_F_, rpId_56_N_, rpId_56_F_ },
165  };
166 }
167 
168 //----------------------------------------------------------------------------------------------------
169 
171 {
172  // get input
174  iEvent.getByToken(tokenHepMC_, hHepMC);
175  HepMC::GenEvent *hepMCEvent = (HepMC::GenEvent *) hHepMC->GetEvent();
176 
178  iEvent.getByToken(tokenTracks_, hTracks);
179 
180  // extract protons
181  bool proton_45_set = false;
182  bool proton_56_set = false;
183  FourVector mom_45, mom_56;
184 
185  for (auto it = hepMCEvent->particles_begin(); it != hepMCEvent->particles_end(); ++it)
186  {
187  const auto &part = *it;
188 
189  // accept only stable non-beam protons
190  if (part->pdg_id() != 2212)
191  continue;
192 
193  if (part->status() != 1)
194  continue;
195 
196  if (part->is_beam())
197  continue;
198 
199  const auto &mom = part->momentum();
200 
201  if (mom.e() < 4500.)
202  continue;
203 
204  if (mom.z() > 0)
205  {
206  // 45
207  if (proton_45_set)
208  {
209  LogError("CTPPSAcceptancePlotter") << "Multiple protons found in sector 45.";
210  return;
211  }
212 
213  proton_45_set = true;
214  mom_45 = mom;
215  } else {
216  // 56
217  if (proton_56_set)
218  {
219  LogError("CTPPSAcceptancePlotter") << "Multiple protons found in sector 56.";
220  return;
221  }
222 
223  proton_56_set = true;
224  mom_56 = mom;
225  }
226  }
227 
228  // stop if protons missing
229  if (!proton_45_set || !proton_56_set)
230  return;
231 
232  // calculate xi's
233  const double p_nom = 6500.;
234  const double xi_45 = (p_nom - mom_45.e()) / p_nom;
235  const double xi_56 = (p_nom - mom_56.e()) / p_nom;
236 
237  // process tracks
238  map<unsigned int, bool> trackPresent;
239  for (const auto& trk : *hTracks)
240  {
241  CTPPSDetId rpId(trk.getRPId());
242  unsigned int rpDecId = rpId.arm()*100 + rpId.station()*10 + rpId.rp();
243  trackPresent[rpDecId] = true;
244  }
245 
246  // update plots
247  for (const auto rpIds : singleArmConfigurations)
248  {
249  bool acc = true;
250  signed int arm = -1;
251  for (const auto rpId : rpIds)
252  {
253  acc &= trackPresent[rpId];
254  arm = rpId / 100;
255  }
256 
257  if (arm < 0)
258  continue;
259 
260  const double xi = (arm == 0) ? xi_45 : xi_56;
261 
262  singleArmPlots[rpIds].fill(xi, acc);
263  }
264 
265  for (const auto rpIds : doubleArmConfigurations)
266  {
267  bool acc = true;
268  for (const auto rpId : rpIds)
269  acc &= trackPresent[rpId];
270 
271  doubleArmPlots[rpIds].fill(xi_45, xi_56, acc);
272  }
273 }
274 
275 
276 //----------------------------------------------------------------------------------------------------
277 
279 {
280  auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
281 
282  for (const auto &p : singleArmPlots) {
283  string dirName;
284  for (const auto &rpId : p.first) {
285  if (!dirName.empty())
286  dirName += ",";
287  dirName += Form("%u", rpId);
288  }
289 
290  gDirectory = f_out->mkdir(dirName.c_str());
291  p.second.write();
292  }
293 
294  for (const auto &p : doubleArmPlots) {
295  string dirName;
296  for (const auto &rpId : p.first) {
297  if (!dirName.empty())
298  dirName += ",";
299  dirName += Form("%u", rpId);
300  }
301 
302  gDirectory = f_out->mkdir(dirName.c_str());
303  p.second.write();
304  }
305 }
306 
307 //----------------------------------------------------------------------------------------------------
308 
310 
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
CTPPSAcceptancePlotter(const edm::ParameterSet &)
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
T sqrt(T t)
Definition: SSEVec.h:18
std::map< std::set< unsigned int >, SingleArmPlots > singleArmPlots
edm::EDGetTokenT< edm::HepMCProduct > tokenHepMC_
std::vector< std::set< unsigned int > > singleArmConfigurations
std::map< std::set< unsigned int >, DoubleArmPlots > doubleArmPlots
edm::EDGetTokenT< CTPPSLocalTrackLiteCollection > tokenTracks_
void fill(double xi_45, double xi_56, bool acc)
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
std::vector< CTPPSLocalTrackLite > CTPPSLocalTrackLiteCollection
Collection of CTPPSLocalTrackLite objects.
part
Definition: HCALResponse.h:20
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
HLT enums.
void analyze(const edm::Event &, const edm::EventSetup &) override
std::vector< std::set< unsigned int > > doubleArmConfigurations