CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MatchCandidateBenchmark.cc
Go to the documentation of this file.
2 
4 
5 
6 // #include "DQMServices/Core/interface/MonitorElement.h"
7 // #include "DQMServices/Core/interface/DQMStore.h"
8 
9 #include <TROOT.h>
10 #include <TFile.h>
11 #include <TH1.h>
12 #include <TH2.h>
13 
14 
15 using namespace std;
16 
19  delta_et_VS_et_ = 0;
20  delta_eta_VS_et_ = 0;
21  delta_phi_VS_et_ = 0;
22 
23  histogramBooked_ = false;
24 }
25 
27 
28 
30  if (!histogramBooked_) {
31  PhaseSpace ptPS;
32  PhaseSpace dptOvptPS;
33  PhaseSpace dptPS;
34  PhaseSpace detaPS;
35  PhaseSpace dphiPS;
36  switch(mode_) {
37  case VALIDATION:
38  ptPS = PhaseSpace(100,0,1000);
39  dptOvptPS = PhaseSpace( 200, -1, 1);
40  dphiPS = PhaseSpace( 200, -1, 1);
41  detaPS = PhaseSpace( 200, -1, 1);
42  dptPS = PhaseSpace( 100, -100, 100);
43  break;
44  case DQMOFFLINE:
45  default:
46  ptPS = PhaseSpace(50,0,100);
47  dptOvptPS = PhaseSpace( 50, -1, 1);
48  dphiPS = PhaseSpace( 50, -1, 1);
49  detaPS = PhaseSpace( 50, -1, 1);
50  dptPS = PhaseSpace( 50, -50, 50);
51  break;
52  }
53  float ptBins[11] = {0, 1, 2, 5, 10, 20, 50, 100, 200, 400, 1000};
54 
55  delta_et_Over_et_VS_et_ = book2D("delta_et_Over_et_VS_et_",
56  ";E_{T, true} (GeV);#DeltaE_{T}/E_{T}",
57  10, ptBins,
58  dptOvptPS.n, dptOvptPS.m, dptOvptPS.M );
59 
60 
61  delta_et_VS_et_ = book2D("delta_et_VS_et_",
62  ";E_{T, true} (GeV);#DeltaE_{T}",
63  10, ptBins,
64  dptPS.n, dptPS.m, dptPS.M );
65 
66  delta_eta_VS_et_ = book2D("delta_eta_VS_et_",
67  ";#E_{T, true} (GeV);#Delta#eta",
68  10, ptBins,
69  detaPS.n, detaPS.m, detaPS.M );
70 
71  delta_phi_VS_et_ = book2D("delta_phi_VS_et_",
72  ";E_{T, true} (GeV);#Delta#phi",
73  10, ptBins,
74  dphiPS.n, dphiPS.m, dphiPS.M );
75  histogramBooked_ = true;
76  }
77 }
79 
80  if (histogramBooked_) return;
81 
82  edm::ParameterSet dptPS = parameterSet.getParameter<edm::ParameterSet>("DeltaPtHistoParameter");
83  edm::ParameterSet dptOvptPS = parameterSet.getParameter<edm::ParameterSet>("DeltaPtOvPtHistoParameter");
84  edm::ParameterSet detaPS = parameterSet.getParameter<edm::ParameterSet>("DeltaEtaHistoParameter");
85  edm::ParameterSet dphiPS = parameterSet.getParameter<edm::ParameterSet>("DeltaPhiHistoParameter");
86 
87  std::vector<double> ptBinsPS = parameterSet.getParameter< std::vector<double> >( "VariablePtBins" );
88  float* ptBins = new float[ptBinsPS.size()];
89  for (size_t i = 0; i < ptBinsPS.size(); i++) {
90  ptBins[i] = ptBinsPS[i];
91  }
92 
93  if (dptOvptPS.getParameter<bool>("switchOn")) {
94  delta_et_Over_et_VS_et_ = book2D("delta_et_Over_et_VS_et_",
95  ";E_{T, true} (GeV);#DeltaE_{T}/E_{T}",
96  ptBinsPS.size()-1, ptBins,
97  dptOvptPS.getParameter<int32_t>("nBin"),
98  dptOvptPS.getParameter<double>("xMin"),
99  dptOvptPS.getParameter<double>("xMax"));
100  }
101 
102  if (dptPS.getParameter<bool>("switchOn")) {
103  delta_et_VS_et_ = book2D("delta_et_VS_et_",
104  ";E_{T, true} (GeV);#DeltaE_{T}",
105  ptBinsPS.size()-1, ptBins,
106  dptPS.getParameter<int32_t>("nBin"),
107  dptPS.getParameter<double>("xMin"),
108  dptPS.getParameter<double>("xMax"));
109  }
110 
111  if (detaPS.getParameter<bool>("switchOn")) {
112  delta_eta_VS_et_ = book2D("delta_eta_VS_et_",
113  ";#E_{T, true} (GeV);#Delta#eta",
114  ptBinsPS.size()-1, ptBins,
115  detaPS.getParameter<int32_t>("nBin"),
116  detaPS.getParameter<double>("xMin"),
117  detaPS.getParameter<double>("xMax"));
118  }
119 
120  if (dphiPS.getParameter<bool>("switchOn")) {
121  delta_phi_VS_et_ = book2D("delta_phi_VS_et_",
122  ";E_{T, true} (GeV);#Delta#phi",
123  ptBinsPS.size()-1, ptBins,
124  dphiPS.getParameter<int32_t>("nBin"),
125  dphiPS.getParameter<double>("xMin"),
126  dphiPS.getParameter<double>("xMax"));
127  }
128  histogramBooked_ = true;
129 }
130 
132  const reco::Candidate& matchedCand) {
133 
134  if( !isInRange(cand.pt(), cand.eta(), cand.phi() ) ) return;
135 
136  if (histogramBooked_) {
137  if (delta_et_Over_et_VS_et_) delta_et_Over_et_VS_et_->Fill( matchedCand.pt(), (cand.pt() - matchedCand.pt())/matchedCand.pt() );
138  if (delta_et_VS_et_) delta_et_VS_et_->Fill( matchedCand.pt(), cand.pt() - matchedCand.pt() );
139  if (delta_eta_VS_et_) delta_eta_VS_et_->Fill( matchedCand.pt(), cand.eta() - matchedCand.eta() );
140  if (delta_phi_VS_et_) delta_phi_VS_et_->Fill( matchedCand.pt(), cand.phi() - matchedCand.phi() );
141  }
142 
143 }
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
void fillOne(const reco::Candidate &candidate, const reco::Candidate &matchedCandidate)
fill histograms with a given particle
virtual double pt() const =0
transverse momentum
void setup()
book histograms
abstract base class
Definition: Benchmark.h:20
TH2F * book2D(const char *histname, const char *title, int nbinsx, float xmin, float xmax, int nbinsy, float ymin, float ymax)
book a 2D histogram, either with DQM or plain root.
Definition: Benchmark.cc:43
int mode
Definition: AMPTWrapper.h:139
Mode mode_
Definition: Benchmark.h:86
bool isInRange(float pt, float eta, float phi) const
Definition: Benchmark.h:58
virtual double phi() const =0
momentum azimuthal angle
virtual double eta() const =0
momentum pseudorapidity