CMS 3D CMS Logo

CTPPSHepMCDistributionPlotter.cc
Go to the documentation of this file.
1 /****************************************************************************
2  * Authors:
3  * Jan Kašpar
4  ****************************************************************************/
5 
10 
13 
15 
17 
20 
21 #include "TFile.h"
22 #include "TH1D.h"
23 
24 #include <string>
25 
26 //----------------------------------------------------------------------------------------------------
27 
29 public:
31 
32 private:
33  void analyze(const edm::Event &, const edm::EventSetup &) override;
34  void endJob() override;
35 
39 
40  std::unique_ptr<TH1D> h_vtx_x_, h_vtx_y_, h_vtx_z_, h_vtx_t_;
41  std::unique_ptr<TH1D> h_xi_, h_th_x_, h_th_y_;
42 };
43 
44 //----------------------------------------------------------------------------------------------------
45 
46 using namespace std;
47 using namespace edm;
48 using namespace HepMC;
49 
50 //----------------------------------------------------------------------------------------------------
51 
53  : tokenHepMC_(consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("tagHepMC"))),
54  lhcInfoESToken_(esConsumes(ESInputTag("", iConfig.getParameter<std::string>("lhcInfoLabel")))),
55  outputFile_(iConfig.getParameter<string>("outputFile")),
56 
57  h_vtx_x_(new TH1D("h_vtx_x", ";vtx_x (mm)", 100, 0., 0.)),
58  h_vtx_y_(new TH1D("h_vtx_y", ";vtx_y (mm)", 100, 0., 0.)),
59  h_vtx_z_(new TH1D("h_vtx_z", ";vtx_z (mm)", 100, 0., 0.)),
60  h_vtx_t_(new TH1D("h_vtx_t", ";vtx_t (mm)", 100, 0., 0.)),
61 
62  h_xi_(new TH1D("h_xi", ";#xi", 100, 0., 0.30)),
63  h_th_x_(new TH1D("h_th_x", ";#theta^{*}_{x}", 100, -300E-6, +300E-6)),
64  h_th_y_(new TH1D("h_th_y", ";#theta^{*}_{y}", 100, -300E-6, +300E-6)) {}
65 
66 //----------------------------------------------------------------------------------------------------
67 
69  // get conditions
70  const auto &lhcInfo = iSetup.getData(lhcInfoESToken_);
71 
72  // get input
74  iEvent.getByToken(tokenHepMC_, hHepMC);
75  HepMC::GenEvent *hepMCEvent = (HepMC::GenEvent *)hHepMC->GetEvent();
76 
77  // plot vertices
78  for (HepMC::GenEvent::vertex_iterator vit = hepMCEvent->vertices_begin(); vit != hepMCEvent->vertices_end(); ++vit) {
79  const auto pos = (*vit)->position();
80  h_vtx_x_->Fill(pos.x());
81  h_vtx_y_->Fill(pos.y());
82  h_vtx_z_->Fill(pos.z());
83  h_vtx_t_->Fill(pos.t());
84  }
85 
86  // extract protons
87  for (auto it = hepMCEvent->particles_begin(); it != hepMCEvent->particles_end(); ++it) {
88  const auto &part = *it;
89 
90  // accept only stable non-beam protons
91  if (part->pdg_id() != 2212)
92  continue;
93 
94  if (part->status() != 1)
95  continue;
96 
97  if (part->is_beam())
98  continue;
99 
100  const auto &mom = part->momentum();
101  const double p_nom = lhcInfo.energy();
102 
103  if (mom.rho() / p_nom < 0.7)
104  continue;
105 
106  const double xi_simu = (p_nom - mom.e()) / p_nom;
107  const double th_x_simu = mom.x() / mom.rho();
108  const double th_y_simu = mom.y() / mom.rho();
109 
110  h_xi_->Fill(xi_simu);
111  h_th_x_->Fill(th_x_simu);
112  h_th_y_->Fill(th_y_simu);
113  }
114 }
115 
116 //----------------------------------------------------------------------------------------------------
117 
119  auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
120 
121  h_vtx_x_->Write();
122  h_vtx_y_->Write();
123  h_vtx_z_->Write();
124  h_vtx_t_->Write();
125 
126  h_xi_->Write();
127  h_th_x_->Write();
128  h_th_y_->Write();
129 }
130 
131 //----------------------------------------------------------------------------------------------------
132 
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
edm::EDGetTokenT< edm::HepMCProduct > tokenHepMC_
edm::ESGetToken< LHCInfo, LHCInfoRcd > lhcInfoESToken_
CTPPSHepMCDistributionPlotter(const edm::ParameterSet &)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
bool getData(T &iHolder) const
Definition: EventSetup.h:122
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:37
part
Definition: HCALResponse.h:20
void analyze(const edm::Event &, const edm::EventSetup &) override
HLT enums.
VertexRefVector::iterator vertex_iterator
iterator over a vector of references to Vertex objects in the same collection
Definition: VertexFwd.h:19