CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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,
60  229, 259, 292, 327, 366, 400,
61  453, 501, 553, 609, 669, 733,
62  802, 875, 954, 1038, 1127, 1222,
63  1323, 1431, 1546, 1667, 1796, 1934,
64  2079, 2233, 2396, 2569, 2752,3000};
65 
66 
67  // (jetEta1 > 0 && jetEta1 < 0.7), (jetEta2 > 0 && jetEta2 < 0.7 )
68  std::string tit = std::string(name) + "_Eta_innerEtaCut_outerEtaCut";
69  hJet[0] = new TH1F(tit.c_str(), "DiJet Mass", binSize, massBin);
70 
71 
72  // (jetEta1 > 0.7 && jetEta1 < 1.3), (jetEta2 > 0.7 && jetEta2 < 1.3 )
73  tit = std::string(name) + "_Eta_0_innerEtaCut";
74  hJet[1] = new TH1F(tit.c_str(), "DiJet Mass", binSize, massBin);
75 
76  tit = std::string(name) + "_LeadJetEta";
77  hJet[2] = new TH1F(tit.c_str(), "1^{st} Leading Jet #eta", 120, -6., 6.);
78  tit = std::string(name) + "_SecondJetEta";
79  hJet[3] = new TH1F(tit.c_str(), "2^{nd} Leading Jet #eta", 120, -6., 6.);
80  tit = std::string(name) + "_numEvents";
81  hJet[4] = new TH1F(tit.c_str(), "No. of events", 10, 0.,10.);
82 
83  return ;
84 }
85 
86 
87 
88 template <class R>
89 void histoFill(TH1F* jetHisto[], edm::Handle<R> jetsRec, double eta1, double eta2)
90 {
91  //For no. of events
92  jetHisto[4]->Fill(1.);
93 
94  if ((*jetsRec).size() >=2){
95  double px1 = (*jetsRec)[0].px();
96  double py1 = (*jetsRec)[0].py();
97  double pz1 = (*jetsRec)[0].pz();
98  double e1 = (*jetsRec)[0].energy();
99  double jetEta1 = (*jetsRec)[0].eta();
100  jetHisto[2]->Fill(jetEta1);
101 
102  double px2 = (*jetsRec)[1].px();
103  double py2 = (*jetsRec)[1].py();
104  double pz2 = (*jetsRec)[1].pz();
105  double e2 = (*jetsRec)[1].energy();
106  double jetEta2 = (*jetsRec)[1].eta();
107  jetHisto[3]->Fill(jetEta2);
108 
109  CLHEP::HepLorentzVector v1(px1,py1,pz1,e1);
110  CLHEP::HepLorentzVector v2(px2,py2,pz2,e2);
111  CLHEP::HepLorentzVector v(0.,0.,0.,0.);
112  v = v1 + v2;
113  float mass = v.m();
114 
115  if ( fabs(jetEta1) > 0.0 && fabs(jetEta1) < eta1)
116  if ( fabs(jetEta2) > 0.0 && fabs(jetEta2) < eta1)
117  jetHisto[0]->Fill(mass);
118 
119  if ( fabs(jetEta1) > eta1 && fabs(jetEta1) < eta2)
120  if ( fabs(jetEta2) > eta1 && fabs(jetEta2) < eta2)
121  jetHisto[1]->Fill(mass);
122 
123  }
124 }//histoFill
125 
126 
127 //
128 // class decleration
129 //
130 template<class Jet>
131 class DijetRatio : public edm::EDAnalyzer {
132 public:
133  explicit DijetRatio(const edm::ParameterSet&);
134  ~DijetRatio();
135 
136 
137  typedef std::vector<Jet> JetCollection;
138  virtual void beginJob() ;
139  virtual void analyze(const edm::Event&, const edm::EventSetup&);
140  virtual void endJob() ;
141 
142 
143  // ----------member data ---------------------------
145 
146  // names of modules, producing object collections
149 
150  // eta limit
151  double m_eta3; // eta limit for numerator
152  double m_eta4; // eta limit for denominator
153 
154  //Jet Kinematics for leading Jet
155  static const int hisotNumber = 10;
156 
159 
160  // Root file for saving histo
161  TFile* hOutputFile ;
162 };
163 
164 #endif
165 
std::string m_Mid5CaloJetsSrc
Definition: DijetRatio.h:148
virtual void analyze(const edm::Event &, const edm::EventSetup &)
Definition: DijetRatio.cc:60
static const int hisotNumber
Definition: DijetRatio.h:155
std::string fOutputFileName
Definition: DijetRatio.h:144
std::string m_Mid5CorRecJetsSrc
Definition: DijetRatio.h:147
DijetRatio(const edm::ParameterSet &)
Definition: DijetRatio.cc:24
double m_eta3
Definition: DijetRatio.h:151
virtual void endJob()
Definition: DijetRatio.cc:93
virtual void beginJob()
Definition: DijetRatio.cc:79
Float e1
Definition: deltaR.h:22
TFile * hOutputFile
Definition: DijetRatio.h:161
const int histoSize
Definition: DijetRatio.h:54
Float e2
Definition: deltaR.h:23
return(e1-e2)*(e1-e2)+dp *dp
double m_eta4
Definition: DijetRatio.h:152
void histoFill(TH1F *jetHisto[], edm::Handle< R > jetsRec, double eta1, double eta2)
Definition: DijetRatio.h:89
TH1F * hCalo[hisotNumber]
Definition: DijetRatio.h:157
TH1F * hCor[hisotNumber]
Definition: DijetRatio.h:158
std::vector< Jet > JetCollection
Definition: DijetRatio.h:137
void hInit(TH1F *hJet[], const char *name)
Definition: DijetRatio.h:57