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