CMS 3D CMS Logo

DijetRatio.h
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: DijetRatio
4 // Class: DijetRatio
5 //
13 //
14 // Original Author: Manoj Jha
15 // Created: Thu Apr 12 15:04:37 CDT 2007
16 // Kalanand Mishra (November 22, 2009):
17 // Modified and cleaned up to work in 3.3.X
18 //
19 //
20 
21 // system include files
22 #ifndef DIJETRATIO_HH
23 #define DIJETRATIO_HH
24 
25 #include <memory>
26 #include <string>
27 #include <iostream>
28 #include <map>
29 #include <algorithm>
30 
31 // user include files
36 
41 
48 #include "CLHEP/Vector/LorentzVector.h"
49 
50 #include "TFile.h"
51 #include "TH1.h"
52 #include "TH2.h"
53 
54 const int histoSize = 5;
55 
56 //Histo Initializations
57 inline void hInit(TH1F* hJet[], const char* name) {
58  int const binSize = 35;
59  float massBin[binSize + 1] = {100, 113, 132, 153, 176, 201, 229, 259, 292, 327, 366, 400,
60  453, 501, 553, 609, 669, 733, 802, 875, 954, 1038, 1127, 1222,
61  1323, 1431, 1546, 1667, 1796, 1934, 2079, 2233, 2396, 2569, 2752, 3000};
62 
63  // (jetEta1 > 0 && jetEta1 < 0.7), (jetEta2 > 0 && jetEta2 < 0.7 )
64  std::string tit = std::string(name) + "_Eta_innerEtaCut_outerEtaCut";
65  hJet[0] = new TH1F(tit.c_str(), "DiJet Mass", binSize, massBin);
66 
67  // (jetEta1 > 0.7 && jetEta1 < 1.3), (jetEta2 > 0.7 && jetEta2 < 1.3 )
68  tit = std::string(name) + "_Eta_0_innerEtaCut";
69  hJet[1] = new TH1F(tit.c_str(), "DiJet Mass", binSize, massBin);
70 
71  tit = std::string(name) + "_LeadJetEta";
72  hJet[2] = new TH1F(tit.c_str(), "1^{st} Leading Jet #eta", 120, -6., 6.);
73  tit = std::string(name) + "_SecondJetEta";
74  hJet[3] = new TH1F(tit.c_str(), "2^{nd} Leading Jet #eta", 120, -6., 6.);
75  tit = std::string(name) + "_numEvents";
76  hJet[4] = new TH1F(tit.c_str(), "No. of events", 10, 0., 10.);
77 
78  return;
79 }
80 
81 template <class R>
82 void histoFill(TH1F* jetHisto[], edm::Handle<R> jetsRec, double eta1, double eta2) {
83  //For no. of events
84  jetHisto[4]->Fill(1.);
85 
86  if ((*jetsRec).size() >= 2) {
87  double px1 = (*jetsRec)[0].px();
88  double py1 = (*jetsRec)[0].py();
89  double pz1 = (*jetsRec)[0].pz();
90  double e1 = (*jetsRec)[0].energy();
91  double jetEta1 = (*jetsRec)[0].eta();
92  jetHisto[2]->Fill(jetEta1);
93 
94  double px2 = (*jetsRec)[1].px();
95  double py2 = (*jetsRec)[1].py();
96  double pz2 = (*jetsRec)[1].pz();
97  double e2 = (*jetsRec)[1].energy();
98  double jetEta2 = (*jetsRec)[1].eta();
99  jetHisto[3]->Fill(jetEta2);
100 
101  CLHEP::HepLorentzVector v1(px1, py1, pz1, e1);
102  CLHEP::HepLorentzVector v2(px2, py2, pz2, e2);
103  CLHEP::HepLorentzVector v(0., 0., 0., 0.);
104  v = v1 + v2;
105  float mass = v.m();
106 
107  if (fabs(jetEta1) > 0.0 && fabs(jetEta1) < eta1)
108  if (fabs(jetEta2) > 0.0 && fabs(jetEta2) < eta1)
109  jetHisto[0]->Fill(mass);
110 
111  if (fabs(jetEta1) > eta1 && fabs(jetEta1) < eta2)
112  if (fabs(jetEta2) > eta1 && fabs(jetEta2) < eta2)
113  jetHisto[1]->Fill(mass);
114  }
115 } //histoFill
116 
117 //
118 // class decleration
119 //
120 template <class Jet>
121 class DijetRatio : public edm::EDAnalyzer {
122 public:
123  explicit DijetRatio(const edm::ParameterSet&);
124  ~DijetRatio() override;
125 
126  typedef std::vector<Jet> JetCollection;
127  void beginJob() override;
128  void analyze(const edm::Event&, const edm::EventSetup&) override;
129  void endJob() override;
130 
131  // ----------member data ---------------------------
133 
134  // names of modules, producing object collections
137 
138  // eta limit
139  double m_eta3; // eta limit for numerator
140  double m_eta4; // eta limit for denominator
141 
142  //Jet Kinematics for leading Jet
143  static const int hisotNumber = 10;
144 
147 
148  // Root file for saving histo
149  TFile* hOutputFile;
150 };
151 
152 #endif
std::string m_Mid5CaloJetsSrc
Definition: DijetRatio.h:136
static const int hisotNumber
Definition: DijetRatio.h:143
std::string fOutputFileName
Definition: DijetRatio.h:132
std::string m_Mid5CorRecJetsSrc
Definition: DijetRatio.h:135
~DijetRatio() override
Definition: DijetRatio.cc:38
DijetRatio(const edm::ParameterSet &)
Definition: DijetRatio.cc:24
void endJob() override
Definition: DijetRatio.cc:79
double m_eta3
Definition: DijetRatio.h:139
TFile * hOutputFile
Definition: DijetRatio.h:149
const int histoSize
Definition: DijetRatio.h:54
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: DijetRatio.cc:50
void beginJob() override
Definition: DijetRatio.cc:66
double m_eta4
Definition: DijetRatio.h:140
void histoFill(TH1F *jetHisto[], edm::Handle< R > jetsRec, double eta1, double eta2)
Definition: DijetRatio.h:82
TH1F * hCalo[hisotNumber]
Definition: DijetRatio.h:145
TH1F * hCor[hisotNumber]
Definition: DijetRatio.h:146
std::vector< Jet > JetCollection
Definition: DijetRatio.h:126
void hInit(TH1F *hJet[], const char *name)
Definition: DijetRatio.h:57