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 {
30  public:
32 
33  private:
34  void analyze(const edm::Event&, const edm::EventSetup&) override;
35  void endJob() override;
36 
40 
41  std::unique_ptr<TH1D> h_vtx_x_, h_vtx_y_, h_vtx_z_;
42  std::unique_ptr<TH1D> h_xi_, h_th_x_, h_th_y_;
43 };
44 
45 //----------------------------------------------------------------------------------------------------
46 
47 using namespace std;
48 using namespace edm;
49 using namespace HepMC;
50 
51 //----------------------------------------------------------------------------------------------------
52 
54  tokenHepMC_( consumes<edm::HepMCProduct>(iConfig.getParameter<edm::InputTag>("tagHepMC")) ),
55  lhcInfoLabel_(iConfig.getParameter<std::string>("lhcInfoLabel")),
56  outputFile_(iConfig.getParameter<string>("outputFile")),
57 
58  h_vtx_x_(new TH1D("h_vtx_x", ";vtx_x (mm)", 100, 0., 0.)),
59  h_vtx_y_(new TH1D("h_vtx_y", ";vtx_y (mm)", 100, 0., 0.)),
60  h_vtx_z_(new TH1D("h_vtx_z", ";vtx_y (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 //----------------------------------------------------------------------------------------------------
68 
70 {
71  // get conditions
72  edm::ESHandle<LHCInfo> hLHCInfo;
73  iSetup.get<LHCInfoRcd>().get(lhcInfoLabel_, hLHCInfo);
74 
75  // get input
77  iEvent.getByToken(tokenHepMC_, hHepMC);
78  HepMC::GenEvent *hepMCEvent = (HepMC::GenEvent *) hHepMC->GetEvent();
79 
80  // plot vertices
81  for (HepMC::GenEvent::vertex_iterator vit = hepMCEvent->vertices_begin(); vit != hepMCEvent->vertices_end(); ++vit)
82  {
83  const auto pos = (*vit)->position();
84  h_vtx_x_->Fill(pos.x());
85  h_vtx_y_->Fill(pos.y());
86  h_vtx_z_->Fill(pos.z());
87  }
88 
89  // extract protons
90  for (auto it = hepMCEvent->particles_begin(); it != hepMCEvent->particles_end(); ++it)
91  {
92  const auto &part = *it;
93 
94  // accept only stable non-beam protons
95  if (part->pdg_id() != 2212)
96  continue;
97 
98  if (part->status() != 1)
99  continue;
100 
101  if (part->is_beam())
102  continue;
103 
104  const auto &mom = part->momentum();
105  const double p_nom = hLHCInfo->energy();
106 
107  if (mom.rho() / p_nom < 0.7)
108  continue;
109 
110  const double xi_simu = (p_nom - mom.e()) / p_nom;
111  const double th_x_simu = mom.x() / mom.rho();
112  const double th_y_simu = mom.y() / mom.rho();
113 
114  h_xi_->Fill(xi_simu);
115  h_th_x_->Fill(th_x_simu);
116  h_th_y_->Fill(th_y_simu);
117  }
118 }
119 
120 //----------------------------------------------------------------------------------------------------
121 
123 {
124  auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
125 
126  h_vtx_x_->Write();
127  h_vtx_y_->Write();
128  h_vtx_z_->Write();
129 
130  h_xi_->Write();
131  h_th_x_->Write();
132  h_th_y_->Write();
133 }
134 
135 //----------------------------------------------------------------------------------------------------
136 
138 
edm::EDGetTokenT< edm::HepMCProduct > tokenHepMC_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
CTPPSHepMCDistributionPlotter(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
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
part
Definition: HCALResponse.h:20
void analyze(const edm::Event &, const edm::EventSetup &) override
HLT enums.
T get() const
Definition: EventSetup.h:71
float const energy() const
Definition: LHCInfo.cc:192
VertexRefVector::iterator vertex_iterator
iterator over a vector of references to Vertex objects in the same collection
Definition: VertexFwd.h:19