CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CandidateBenchmark.cc
Go to the documentation of this file.
2 
4 
5 
6 #include <TROOT.h>
7 #include <TFile.h>
8 #include <TH1.h>
9 #include <TH2.h>
10 
11 
12 using namespace std;
13 
14 
16  pt_ = 0;
17  eta_ = 0;
18  phi_ = 0;
19  charge_ = 0;
20  pdgId_ = 0;
21 
22  histogramBooked_ = false;
23 
24 }
25 
27 
28 
30 
31  if (!histogramBooked_) {
32  PhaseSpace ptPS(100,0,100);
33  PhaseSpace phiPS(360, -3.1416, 3.1416);
34  PhaseSpace etaPS(100, -5,5);
35  switch(mode_) {
36  case DQMOFFLINE:
37  default:
38  ptPS = PhaseSpace(50, 0, 100);
39  phiPS.n = 50;
40  etaPS.n = 20;
41  break;
42  }
43 
44  pt_ = book1D("pt_", "pt_;p_{T} (GeV)", ptPS.n, ptPS.m, ptPS.M);
45 
46  eta_ = book1D("eta_", "eta_;#eta", etaPS.n, etaPS.m, etaPS.M);
47 
48  phi_ = book1D("phi_", "phi_;#phi", phiPS.n, phiPS.m, phiPS.M);
49 
50  charge_ = book1D("charge_", "charge_;charge", 3,-1.5,1.5);
51 
52  histogramBooked_ = true;
53  }
54 }
55 
56 void CandidateBenchmark::setup(const edm::ParameterSet& parameterSet) {
57  if (!histogramBooked_) {
58  edm::ParameterSet ptPS = parameterSet.getParameter<edm::ParameterSet>("PtHistoParameter");
59  edm::ParameterSet etaPS = parameterSet.getParameter<edm::ParameterSet>("EtaHistoParameter");
60  edm::ParameterSet phiPS = parameterSet.getParameter<edm::ParameterSet>("PhiHistoParameter");
61  edm::ParameterSet chPS = parameterSet.getParameter<edm::ParameterSet>("ChargeHistoParameter");
62 
63  if (ptPS.getParameter<bool>("switchOn")) {
64  pt_ = book1D("pt_", "pt_;p_{T} (GeV)", ptPS.getParameter<int32_t>("nBin"),
65  ptPS.getParameter<double>("xMin"),
66  ptPS.getParameter<double>("xMax"));
67  }
68 
69  if (etaPS.getParameter<bool>("switchOn")) {
70  eta_ = book1D("eta_", "eta_;#eta", etaPS.getParameter<int32_t>("nBin"),
71  etaPS.getParameter<double>("xMin"),
72  etaPS.getParameter<double>("xMax"));
73  }
74  if (phiPS.getParameter<bool>("switchOn")) {
75  phi_ = book1D("phi_", "phi_;#phi", phiPS.getParameter<int32_t>("nBin"),
76  phiPS.getParameter<double>("xMin"),
77  phiPS.getParameter<double>("xMax"));
78  }
79  if (chPS.getParameter<bool>("switchOn")) {
80  charge_ = book1D("charge_","charge_;charge",chPS.getParameter<int32_t>("nBin"),
81  chPS.getParameter<double>("xMin"),
82  chPS.getParameter<double>("xMax"));
83  }
84  histogramBooked_ = true;
85  }
86 }
87 
89 
90  if( !isInRange(cand.pt(), cand.eta(), cand.phi() ) ) return;
91 
92  if (histogramBooked_) {
93  if (pt_) pt_->Fill( cand.pt() );
94  if (eta_) eta_->Fill( cand.eta() );
95  if (phi_) phi_->Fill( cand.phi() );
96  if (charge_) charge_->Fill( cand.charge() );
97  }
98 }
T getParameter(std::string const &) const
virtual float eta() const =0
momentum pseudorapidity
void fillOne(const reco::Candidate &candidate)
fill histograms with a given particle
void setup()
book histograms
abstract base class
Definition: Benchmark.h:20
virtual float phi() const =0
momentum azimuthal angle
TH1F * pdgId_
COLIN add this histo.
virtual float pt() const =0
transverse momentum
virtual int charge() const =0
electric charge
Mode mode_
Definition: Benchmark.h:86
bool isInRange(float pt, float eta, float phi) const
Definition: Benchmark.h:58
TH1F * book1D(const char *histname, const char *title, int nbins, float xmin, float xmax)
book a 1D histogram, either with DQM or plain root.
Definition: Benchmark.cc:25