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_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  lhcInfoLabel_(iConfig.getParameter<std::string>("lhcInfoLabel")),
55  outputFile_(iConfig.getParameter<string>("outputFile")),
56  h_xi_(new TH1D("h_xi", ";#xi", 100, 0., 0.30)),
57  h_th_x_(new TH1D("h_th_x", ";#theta^{*}_{x}", 100, -300E-6, +300E-6)),
58  h_th_y_(new TH1D("h_th_y", ";#theta^{*}_{y}", 100, -300E-6, +300E-6))
59 {}
60 
61 //----------------------------------------------------------------------------------------------------
62 
64 {
65  // get conditions
66  edm::ESHandle<LHCInfo> hLHCInfo;
67  iSetup.get<LHCInfoRcd>().get(lhcInfoLabel_, hLHCInfo);
68 
69  // get input
71  iEvent.getByToken(tokenHepMC_, hHepMC);
72  HepMC::GenEvent *hepMCEvent = (HepMC::GenEvent *) hHepMC->GetEvent();
73 
74  // extract protons
75  for (auto it = hepMCEvent->particles_begin(); it != hepMCEvent->particles_end(); ++it)
76  {
77  const auto &part = *it;
78 
79  // accept only stable non-beam protons
80  if (part->pdg_id() != 2212)
81  continue;
82 
83  if (part->status() != 1)
84  continue;
85 
86  if (part->is_beam())
87  continue;
88 
89  const auto &mom = part->momentum();
90  const double p_nom = hLHCInfo->energy();
91 
92  if (mom.rho() / p_nom < 0.7)
93  continue;
94 
95  const double xi_simu = (p_nom - mom.e()) / p_nom;
96  const double th_x_simu = mom.x() / mom.rho();
97  const double th_y_simu = mom.y() / mom.rho();
98 
99  h_xi_->Fill(xi_simu);
100  h_th_x_->Fill(th_x_simu);
101  h_th_y_->Fill(th_y_simu);
102  }
103 }
104 
105 //----------------------------------------------------------------------------------------------------
106 
108 {
109  auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
110 
111  h_xi_->Write();
112  h_th_x_->Write();
113  h_th_y_->Write();
114 }
115 
116 //----------------------------------------------------------------------------------------------------
117 
119 
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