CMS 3D CMS Logo

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