CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Types | Private Member Functions | Private Attributes
DijetMass< Jet > Class Template Reference

#include <DijetMass.h>

Inheritance diagram for DijetMass< Jet >:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

 DijetMass (const edm::ParameterSet &)
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Types

typedef std::vector< Jet > JetCollection
 

Private Member Functions

void analyze (const edm::Event &, const edm::EventSetup &)
 
void beginJob ()
 
void endJob ()
 

Private Attributes

std::string AKCorJets
 
std::string AKJets
 
TH1F etaAKcor
 
TH1F etaAKunc
 
TH1F etaICcor
 
TH1F etaICunc
 
TH1F etaKTcor
 
TH1F etaKTunc
 
double EtaMax
 
TH1F etaSCcor
 
TH1F etaSCunc
 
int evtCount
 
std::string histogramFile
 
std::string ICCorJets
 
std::string ICJets
 
std::string KTCorJets
 
std::string KTJets
 
TH1F m2jAKcor
 
TH1F m2jAKunc
 
TH1F m2jICcor
 
TH1F m2jICunc
 
TH1F m2jKTcor
 
TH1F m2jKTunc
 
TH1F m2jSCcor
 
TH1F m2jSCunc
 
TFile * m_file
 
int numJets
 
TH1F phiAKcor
 
TH1F phiAKunc
 
TH1F phiICcor
 
TH1F phiICunc
 
TH1F phiKTcor
 
TH1F phiKTunc
 
TH1F phiSCcor
 
TH1F phiSCunc
 
TH1F ptAKcor
 
TH1F ptAKunc
 
double PtHistMax
 
TH1F ptICcor
 
TH1F ptICunc
 
TH1F ptKTcor
 
TH1F ptKTunc
 
TH1F ptSCcor
 
TH1F ptSCunc
 
std::string SCCorJets
 
std::string SCJets
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

template<class Jet>
class DijetMass< Jet >

Definition at line 24 of file DijetMass.h.

Member Typedef Documentation

template<class Jet >
typedef std::vector<Jet> DijetMass< Jet >::JetCollection
private

Definition at line 29 of file DijetMass.h.

Constructor & Destructor Documentation

template<class Jet >
DijetMass< Jet >::DijetMass ( const edm::ParameterSet cfg)

Definition at line 30 of file DijetMass.cc.

References HLT_25ns14e33_v3_cff::EtaMax, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), and AlCaHLTBitMon_QueryRunRegistry::string.

30  {
31  PtHistMax = cfg.getUntrackedParameter<double>( "PtHistMax", 3000.0);
32  EtaMax = cfg.getUntrackedParameter<double> ("EtaMax", 1.3);
33  histogramFile = cfg.getUntrackedParameter<std::string>("HistoFileName",
34  "DijetMassHistos.root");
35 
36  AKJets = cfg.getParameter<std::string> ("AKJets");
37  AKCorJets = cfg.getParameter<std::string> ("AKCorrectedJets");
38  ICJets = cfg.getParameter<std::string> ("ICJets");
39  ICCorJets = cfg.getParameter<std::string> ("ICCorrectedJets");
40  SCJets = cfg.getParameter<std::string> ("SCJets");
41  SCCorJets = cfg.getParameter<std::string> ("SCCorrectedJets");
42  KTJets = cfg.getParameter<std::string> ("KTJets");
43  KTCorJets = cfg.getParameter<std::string> ("KTCorrectedJets");
44 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::string AKCorJets
Definition: DijetMass.h:40
std::string KTJets
Definition: DijetMass.h:45
std::string SCJets
Definition: DijetMass.h:43
std::string KTCorJets
Definition: DijetMass.h:46
std::string AKJets
Definition: DijetMass.h:39
double PtHistMax
Definition: DijetMass.h:36
double EtaMax
Definition: DijetMass.h:37
std::string ICJets
Definition: DijetMass.h:41
std::string histogramFile
Definition: DijetMass.h:38
std::string SCCorJets
Definition: DijetMass.h:44
std::string ICCorJets
Definition: DijetMass.h:42

Member Function Documentation

template<class Jet >
void DijetMass< Jet >::analyze ( const edm::Event evt,
const edm::EventSetup es 
)
privatevirtual

Implements edm::EDAnalyzer.

Definition at line 115 of file DijetMass.cc.

References funct::abs(), eta(), HLT_25ns14e33_v3_cff::EtaMax, edm::Event::getByLabel(), METSkim_cff::Jets, phi, and reco::tau::disc::Pt().

115  {
116 
117  evtCount++;
118  math::XYZTLorentzVector p4jet[2];
119  int jetInd;
121 
122 
123  //Fill Simple Histos
124  typename JetCollection::const_iterator i_jet;
125 
126  //AK unc
127  evt.getByLabel( AKJets, Jets );
128  jetInd = 0;
129  for( i_jet = Jets->begin(); i_jet != Jets->end() && jetInd<2; ++ i_jet ) {
130  p4jet[jetInd] = i_jet->p4();
131  jetInd++;
132  }
133  if(jetInd==2&&abs(p4jet[0].eta())<EtaMax &&abs(p4jet[1].eta())<EtaMax ){
134  m2jAKunc.Fill( (p4jet[0]+p4jet[1]).mass() );
135  ptAKunc.Fill( p4jet[0].Pt() ); ptAKunc.Fill( p4jet[1].Pt() );
136  etaAKunc.Fill( p4jet[0].eta() ); etaAKunc.Fill( p4jet[1].eta() );
137  phiAKunc.Fill( p4jet[0].phi() ); phiAKunc.Fill( p4jet[1].phi() );
138  }
139 
140 
141  //AK corrected
142  evt.getByLabel( AKCorJets, Jets );
143  jetInd = 0;
144  for( i_jet = Jets->begin(); i_jet != Jets->end() && jetInd<2; ++ i_jet ) {
145  p4jet[jetInd] = i_jet->p4();
146  jetInd++;
147  }
148  if(jetInd==2&&abs(p4jet[0].eta())<EtaMax &&abs(p4jet[1].eta())<EtaMax ){
149  m2jAKcor.Fill( (p4jet[0]+p4jet[1]).mass() );
150  ptAKcor.Fill( p4jet[0].Pt() ); ptAKcor.Fill( p4jet[1].Pt() );
151  etaAKcor.Fill( p4jet[0].eta() ); etaAKcor.Fill( p4jet[1].eta() );
152  phiAKcor.Fill( p4jet[0].phi() ); phiAKcor.Fill( p4jet[1].phi() );
153  }
154 
155 
156  //IC unc
157  evt.getByLabel( ICJets, Jets );
158  jetInd = 0;
159  for( i_jet = Jets->begin(); i_jet != Jets->end() && jetInd<2; ++ i_jet ) {
160  p4jet[jetInd] = i_jet->p4();
161  jetInd++;
162  }
163  if(jetInd==2&&abs(p4jet[0].eta())<EtaMax &&abs(p4jet[1].eta())<EtaMax ){
164  m2jICunc.Fill( (p4jet[0]+p4jet[1]).mass() );
165  ptICunc.Fill( p4jet[0].Pt() ); ptICunc.Fill( p4jet[1].Pt() );
166  etaICunc.Fill( p4jet[0].eta() ); etaICunc.Fill( p4jet[1].eta() );
167  phiICunc.Fill( p4jet[0].phi() ); phiICunc.Fill( p4jet[1].phi() );
168  }
169 
170 
171  //IC corrected
172  evt.getByLabel( ICCorJets, Jets );
173  jetInd = 0;
174  for( i_jet = Jets->begin(); i_jet != Jets->end() && jetInd<2; ++ i_jet ) {
175  p4jet[jetInd] = i_jet->p4();
176  jetInd++;
177  }
178  if(jetInd==2&&abs(p4jet[0].eta())<EtaMax &&abs(p4jet[1].eta())<EtaMax ){
179  m2jICcor.Fill( (p4jet[0]+p4jet[1]).mass() );
180  ptICcor.Fill( p4jet[0].Pt() ); ptICcor.Fill( p4jet[1].Pt() );
181  etaICcor.Fill( p4jet[0].eta() ); etaICcor.Fill( p4jet[1].eta() );
182  phiICcor.Fill( p4jet[0].phi() ); phiICcor.Fill( p4jet[1].phi() );
183  }
184 
185  //KT unc
186  evt.getByLabel( KTJets, Jets );
187  jetInd = 0;
188  for( i_jet = Jets->begin(); i_jet != Jets->end() && jetInd<2; ++ i_jet ) {
189  p4jet[jetInd] = i_jet->p4();
190  jetInd++;
191  }
192  if(jetInd==2&&abs(p4jet[0].eta())<EtaMax &&abs(p4jet[1].eta())<EtaMax ){
193  m2jKTunc.Fill( (p4jet[0]+p4jet[1]).mass() );
194  ptKTunc.Fill( p4jet[0].Pt() ); ptKTunc.Fill( p4jet[1].Pt() );
195  etaKTunc.Fill( p4jet[0].eta() ); etaKTunc.Fill( p4jet[1].eta() );
196  phiKTunc.Fill( p4jet[0].phi() ); phiKTunc.Fill( p4jet[1].phi() );
197  }
198 
199 
200  //KT corrected
201  evt.getByLabel( KTCorJets, Jets );
202  jetInd = 0;
203  for( i_jet = Jets->begin(); i_jet != Jets->end() && jetInd<2; ++i_jet ) {
204  p4jet[jetInd] = i_jet->p4();
205  jetInd++;
206  }
207  if(jetInd==2&&abs(p4jet[0].eta())<EtaMax &&abs(p4jet[1].eta())<EtaMax ){
208  m2jKTcor.Fill( (p4jet[0]+p4jet[1]).mass() );
209  ptKTcor.Fill( p4jet[0].Pt() ); ptKTcor.Fill( p4jet[1].Pt() );
210  etaKTcor.Fill( p4jet[0].eta() ); etaKTcor.Fill( p4jet[1].eta() );
211  phiKTcor.Fill( p4jet[0].phi() ); phiKTcor.Fill( p4jet[1].phi() );
212  }
213 
214 
215  //SC unc
216  evt.getByLabel( SCJets, Jets );
217  jetInd = 0;
218  for( i_jet = Jets->begin(); i_jet != Jets->end() && jetInd<2; ++i_jet ) {
219  p4jet[jetInd] = i_jet->p4();
220  jetInd++;
221  }
222  if(jetInd==2&&abs(p4jet[0].eta())<EtaMax &&abs(p4jet[1].eta())<EtaMax ){
223  m2jSCunc.Fill( (p4jet[0]+p4jet[1]).mass() );
224  ptSCunc.Fill( p4jet[0].Pt() ); ptSCunc.Fill( p4jet[1].Pt() );
225  etaSCunc.Fill( p4jet[0].eta() ); etaSCunc.Fill( p4jet[1].eta() );
226  phiSCunc.Fill( p4jet[0].phi() ); phiSCunc.Fill( p4jet[1].phi() );
227  }
228 
229 
230  //SC corrected
231  evt.getByLabel( SCCorJets, Jets );
232  jetInd = 0;
233  for( i_jet = Jets->begin(); i_jet != Jets->end() && jetInd<2; ++i_jet ) {
234  p4jet[jetInd] = i_jet->p4();
235  jetInd++;
236  }
237  if(jetInd==2&&abs(p4jet[0].eta())<EtaMax &&abs(p4jet[1].eta())<EtaMax ){
238  m2jSCcor.Fill( (p4jet[0]+p4jet[1]).mass() );
239  ptSCcor.Fill( p4jet[0].Pt() ); ptSCcor.Fill( p4jet[1].Pt() );
240  etaSCcor.Fill( p4jet[0].eta() ); etaSCcor.Fill( p4jet[1].eta() );
241  phiSCcor.Fill( p4jet[0].phi() ); phiSCcor.Fill( p4jet[1].phi() );
242  }
243 
244 }
TH1F phiSCcor
Definition: DijetMass.h:57
TH1F phiICunc
Definition: DijetMass.h:52
TH1F etaICcor
Definition: DijetMass.h:53
TH1F ptKTcor
Definition: DijetMass.h:55
TH1F m2jAKcor
Definition: DijetMass.h:51
TH1F ptKTunc
Definition: DijetMass.h:54
std::string AKCorJets
Definition: DijetMass.h:40
TH1F ptSCunc
Definition: DijetMass.h:56
TH1F phiICcor
Definition: DijetMass.h:53
std::string KTJets
Definition: DijetMass.h:45
TH1F ptICcor
Definition: DijetMass.h:53
TH1F phiSCunc
Definition: DijetMass.h:56
TH1F etaSCunc
Definition: DijetMass.h:56
std::string SCJets
Definition: DijetMass.h:43
TH1F m2jKTunc
Definition: DijetMass.h:54
T eta() const
TH1F etaICunc
Definition: DijetMass.h:52
TH1F ptAKunc
Definition: DijetMass.h:50
TH1F etaKTcor
Definition: DijetMass.h:55
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
TH1F m2jKTcor
Definition: DijetMass.h:55
std::string KTCorJets
Definition: DijetMass.h:46
TH1F m2jSCunc
Definition: DijetMass.h:56
TH1F m2jSCcor
Definition: DijetMass.h:57
TH1F phiKTcor
Definition: DijetMass.h:55
TH1F etaAKunc
Definition: DijetMass.h:50
tuple Jets
Definition: METSkim_cff.py:17
int evtCount
Definition: DijetMass.h:63
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
TH1F ptICunc
Definition: DijetMass.h:52
std::string AKJets
Definition: DijetMass.h:39
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:405
TH1F phiAKunc
Definition: DijetMass.h:50
TH1F etaSCcor
Definition: DijetMass.h:57
double EtaMax
Definition: DijetMass.h:37
TH1F etaAKcor
Definition: DijetMass.h:51
TH1F ptSCcor
Definition: DijetMass.h:57
TH1F phiKTunc
Definition: DijetMass.h:54
std::string ICJets
Definition: DijetMass.h:41
TH1F m2jICunc
Definition: DijetMass.h:52
TH1F phiAKcor
Definition: DijetMass.h:51
TH1F etaKTunc
Definition: DijetMass.h:54
std::string SCCorJets
Definition: DijetMass.h:44
TH1F m2jICcor
Definition: DijetMass.h:53
std::string ICCorJets
Definition: DijetMass.h:42
TH1F m2jAKunc
Definition: DijetMass.h:50
TH1F ptAKcor
Definition: DijetMass.h:51
Definition: DDAxes.h:10
template<class Jet >
void DijetMass< Jet >::beginJob ( void  )
privatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 47 of file DijetMass.cc.

References gather_cfg::cout, and M_PI.

47  {
48  cout << "DijetMass: Maximum bin edge for Pt Hists = " << PtHistMax << endl;
49  numJets=2;
50 
51  //Initialize some stuff
52  evtCount = 0;
53 
54  // Open the histogram file and book some associated histograms
55  m_file=new TFile( histogramFile.c_str(),"RECREATE" );
56 
57  //Simple histos
58 
59  //AK unc
60  ptAKunc = TH1F("ptAKunc","p_{T} of leading Jets (AK)",50,0.0,PtHistMax);
61  etaAKunc = TH1F("etaAKunc","#eta of leading Jets (AK)",23,-1.0,1.0);
62  phiAKunc = TH1F("phiAKunc","#phi of leading Jets (AK)",72,-M_PI, M_PI);
63  m2jAKunc = TH1F("m2jAKunc","Dijet Mass of leading Jets (AK)",100,0.0,2*PtHistMax);
64 
65 
66  //AK cor
67  ptAKcor = TH1F("ptAKcor","p_{T} of leading Corrected Jets (AK)",50,0.0,PtHistMax);
68  etaAKcor = TH1F("etaAKcor","#eta of leading Corrected Jets (AK)",23,-1.0,1.0);
69  phiAKcor = TH1F("phiAKcor","#phi of leading Corrected Jets (AK)",72,-M_PI, M_PI);
70  m2jAKcor = TH1F("m2jAKcor","Dijet Mass of leading Corrected Jets (AK)",100,0.0,2*PtHistMax);
71 
72  //IC unc
73  ptICunc = TH1F("ptICunc","p_{T} of leading Jets (IC)",50,0.0,PtHistMax);
74  etaICunc = TH1F("etaICunc","#eta of leading Jets (IC)",23,-1.0,1.0);
75  phiICunc = TH1F("phiICunc","#phi of leading Jets (IC)",72,-M_PI, M_PI);
76  m2jICunc = TH1F("m2jICunc","Dijet Mass of leading Jets (IC)",100,0.0,2*PtHistMax);
77 
78 
79  //IC cor
80  ptICcor = TH1F("ptICcor","p_{T} of leading Corrected Jets (IC)",50,0.0,PtHistMax);
81  etaICcor = TH1F("etaICcor","#eta of leading Corrected Jets (IC)",23,-1.0,1.0);
82  phiICcor = TH1F("phiICcor","#phi of leading Corrected Jets (IC)",72,-M_PI, M_PI);
83  m2jICcor = TH1F("m2jICcor","Dijet Mass of leading Corrected Jets (IC)",100,0.0,2*PtHistMax);
84 
85  //KT unc
86  ptKTunc = TH1F("ptKTunc","p_{T} of leading Jets (KT)",50,0.0,PtHistMax);
87  etaKTunc = TH1F("etaKTunc","#eta of leading Jets (KT)",23,-1.0,1.0);
88  phiKTunc = TH1F("phiKTunc","#phi of leading Jets (KT)",72,-M_PI, M_PI);
89  m2jKTunc = TH1F("m2jKTunc","Dijet Mass of leading Jets (KT)",100,0.0,2*PtHistMax);
90 
91 
92  //KT cor
93  ptKTcor = TH1F("ptKTcor","p_{T} of leading Corrected Jets (KT)",50,0.0,PtHistMax);
94  etaKTcor = TH1F("etaKTcor","#eta of leading Corrected Jets (KT)",23,-1.0,1.0);
95  phiKTcor = TH1F("phiKTcor","#phi of leading Corrected Jets (KT)",72,-M_PI, M_PI);
96  m2jKTcor = TH1F("m2jKTcor","Dijet Mass of leading Corrected Jets (KT)",100,0.0,2*PtHistMax);
97 
98  //SC unc
99  ptSCunc = TH1F("ptSCunc","p_{T} of leading Jets (SC)",50,0.0,PtHistMax);
100  etaSCunc = TH1F("etaSCunc","#eta of leading Jets (SC)",23,-1.0,1.0);
101  phiSCunc = TH1F("phiSCunc","#phi of leading Jets (SC)",72,-M_PI, M_PI);
102  m2jSCunc = TH1F("m2jSCunc","Dijet Mass of leading Jets (SC)",100,0.0,2*PtHistMax);
103 
104  //SC cor
105  ptSCcor = TH1F("ptSCcor","p_{T} of leading Corrected Jets (SC)",50,0.0,PtHistMax);
106  etaSCcor = TH1F("etaSCcor","#eta of leading Corrected Jets (SC)",23,-1.0,1.0);
107  phiSCcor = TH1F("phiSCcor","#phi of leading Corrected Jets (SC)",72,-M_PI, M_PI);
108  m2jSCcor = TH1F("m2jSCcor","Dijet Mass of leading Corrected Jets (SC)",100,0.0,2*PtHistMax);
109 
110 }
int numJets
Definition: DijetMass.h:64
TH1F phiSCcor
Definition: DijetMass.h:57
TH1F phiICunc
Definition: DijetMass.h:52
TH1F etaICcor
Definition: DijetMass.h:53
TH1F ptKTcor
Definition: DijetMass.h:55
TH1F m2jAKcor
Definition: DijetMass.h:51
TH1F ptKTunc
Definition: DijetMass.h:54
TH1F ptSCunc
Definition: DijetMass.h:56
TH1F phiICcor
Definition: DijetMass.h:53
TH1F ptICcor
Definition: DijetMass.h:53
TH1F phiSCunc
Definition: DijetMass.h:56
TH1F etaSCunc
Definition: DijetMass.h:56
TH1F m2jKTunc
Definition: DijetMass.h:54
TH1F etaICunc
Definition: DijetMass.h:52
TH1F ptAKunc
Definition: DijetMass.h:50
TH1F etaKTcor
Definition: DijetMass.h:55
TH1F m2jKTcor
Definition: DijetMass.h:55
TH1F m2jSCunc
Definition: DijetMass.h:56
TH1F m2jSCcor
Definition: DijetMass.h:57
TH1F phiKTcor
Definition: DijetMass.h:55
TH1F etaAKunc
Definition: DijetMass.h:50
int evtCount
Definition: DijetMass.h:63
TFile * m_file
Definition: DijetMass.h:60
TH1F ptICunc
Definition: DijetMass.h:52
#define M_PI
double PtHistMax
Definition: DijetMass.h:36
TH1F phiAKunc
Definition: DijetMass.h:50
TH1F etaSCcor
Definition: DijetMass.h:57
TH1F etaAKcor
Definition: DijetMass.h:51
TH1F ptSCcor
Definition: DijetMass.h:57
TH1F phiKTunc
Definition: DijetMass.h:54
TH1F m2jICunc
Definition: DijetMass.h:52
TH1F phiAKcor
Definition: DijetMass.h:51
std::string histogramFile
Definition: DijetMass.h:38
TH1F etaKTunc
Definition: DijetMass.h:54
tuple cout
Definition: gather_cfg.py:121
TH1F m2jICcor
Definition: DijetMass.h:53
TH1F m2jAKunc
Definition: DijetMass.h:50
TH1F ptAKcor
Definition: DijetMass.h:51
template<class Jet >
void DijetMass< Jet >::endJob ( void  )
privatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 247 of file DijetMass.cc.

247  {
248 
249  //Write out the histogram file.
250  m_file->cd();
251 
252  ptAKunc.Write();
253  etaAKunc.Write();
254  phiAKunc.Write();
255  m2jAKunc.Write();
256 
257  ptAKcor.Write();
258  etaAKcor.Write();
259  phiAKcor.Write();
260  m2jAKcor.Write();
261 
262  ptICunc.Write();
263  etaICunc.Write();
264  phiICunc.Write();
265  m2jICunc.Write();
266 
267 
268  ptICcor.Write();
269  etaICcor.Write();
270  phiICcor.Write();
271  m2jICcor.Write();
272 
273  ptKTunc.Write();
274  etaKTunc.Write();
275  phiKTunc.Write();
276  m2jKTunc.Write();
277 
278 
279  ptKTcor.Write();
280  etaKTcor.Write();
281  phiKTcor.Write();
282  m2jKTcor.Write();
283 
284  ptSCunc.Write();
285  etaSCunc.Write();
286  phiSCunc.Write();
287  m2jSCunc.Write();
288 
289  ptSCcor.Write();
290  etaSCcor.Write();
291  phiSCcor.Write();
292  m2jSCcor.Write();
293 
294  m_file->Close();
295 }
TH1F phiSCcor
Definition: DijetMass.h:57
TH1F phiICunc
Definition: DijetMass.h:52
TH1F etaICcor
Definition: DijetMass.h:53
TH1F ptKTcor
Definition: DijetMass.h:55
TH1F m2jAKcor
Definition: DijetMass.h:51
TH1F ptKTunc
Definition: DijetMass.h:54
TH1F ptSCunc
Definition: DijetMass.h:56
TH1F phiICcor
Definition: DijetMass.h:53
TH1F ptICcor
Definition: DijetMass.h:53
TH1F phiSCunc
Definition: DijetMass.h:56
TH1F etaSCunc
Definition: DijetMass.h:56
TH1F m2jKTunc
Definition: DijetMass.h:54
TH1F etaICunc
Definition: DijetMass.h:52
TH1F ptAKunc
Definition: DijetMass.h:50
TH1F etaKTcor
Definition: DijetMass.h:55
TH1F m2jKTcor
Definition: DijetMass.h:55
TH1F m2jSCunc
Definition: DijetMass.h:56
TH1F m2jSCcor
Definition: DijetMass.h:57
TH1F phiKTcor
Definition: DijetMass.h:55
TH1F etaAKunc
Definition: DijetMass.h:50
TFile * m_file
Definition: DijetMass.h:60
TH1F ptICunc
Definition: DijetMass.h:52
TH1F phiAKunc
Definition: DijetMass.h:50
TH1F etaSCcor
Definition: DijetMass.h:57
TH1F etaAKcor
Definition: DijetMass.h:51
TH1F ptSCcor
Definition: DijetMass.h:57
TH1F phiKTunc
Definition: DijetMass.h:54
TH1F m2jICunc
Definition: DijetMass.h:52
TH1F phiAKcor
Definition: DijetMass.h:51
TH1F etaKTunc
Definition: DijetMass.h:54
TH1F m2jICcor
Definition: DijetMass.h:53
TH1F m2jAKunc
Definition: DijetMass.h:50
TH1F ptAKcor
Definition: DijetMass.h:51

Member Data Documentation

template<class Jet >
std::string DijetMass< Jet >::AKCorJets
private

Definition at line 40 of file DijetMass.h.

template<class Jet >
std::string DijetMass< Jet >::AKJets
private

Definition at line 39 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::etaAKcor
private

Definition at line 51 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::etaAKunc
private

Definition at line 50 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::etaICcor
private

Definition at line 53 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::etaICunc
private

Definition at line 52 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::etaKTcor
private

Definition at line 55 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::etaKTunc
private

Definition at line 54 of file DijetMass.h.

template<class Jet >
double DijetMass< Jet >::EtaMax
private

Definition at line 37 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::etaSCcor
private

Definition at line 57 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::etaSCunc
private

Definition at line 56 of file DijetMass.h.

template<class Jet >
int DijetMass< Jet >::evtCount
private

Definition at line 63 of file DijetMass.h.

template<class Jet >
std::string DijetMass< Jet >::histogramFile
private

Definition at line 38 of file DijetMass.h.

template<class Jet >
std::string DijetMass< Jet >::ICCorJets
private

Definition at line 42 of file DijetMass.h.

template<class Jet >
std::string DijetMass< Jet >::ICJets
private

Definition at line 41 of file DijetMass.h.

template<class Jet >
std::string DijetMass< Jet >::KTCorJets
private

Definition at line 46 of file DijetMass.h.

template<class Jet >
std::string DijetMass< Jet >::KTJets
private

Definition at line 45 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::m2jAKcor
private

Definition at line 51 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::m2jAKunc
private

Definition at line 50 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::m2jICcor
private

Definition at line 53 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::m2jICunc
private

Definition at line 52 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::m2jKTcor
private

Definition at line 55 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::m2jKTunc
private

Definition at line 54 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::m2jSCcor
private

Definition at line 57 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::m2jSCunc
private

Definition at line 56 of file DijetMass.h.

template<class Jet >
TFile* DijetMass< Jet >::m_file
private

Definition at line 60 of file DijetMass.h.

template<class Jet >
int DijetMass< Jet >::numJets
private

Definition at line 64 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::phiAKcor
private

Definition at line 51 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::phiAKunc
private

Definition at line 50 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::phiICcor
private

Definition at line 53 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::phiICunc
private

Definition at line 52 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::phiKTcor
private

Definition at line 55 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::phiKTunc
private

Definition at line 54 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::phiSCcor
private

Definition at line 57 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::phiSCunc
private

Definition at line 56 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::ptAKcor
private

Definition at line 51 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::ptAKunc
private

Definition at line 50 of file DijetMass.h.

template<class Jet >
double DijetMass< Jet >::PtHistMax
private

Definition at line 36 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::ptICcor
private

Definition at line 53 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::ptICunc
private

Definition at line 52 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::ptKTcor
private

Definition at line 55 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::ptKTunc
private

Definition at line 54 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::ptSCcor
private

Definition at line 57 of file DijetMass.h.

template<class Jet >
TH1F DijetMass< Jet >::ptSCunc
private

Definition at line 56 of file DijetMass.h.

template<class Jet >
std::string DijetMass< Jet >::SCCorJets
private

Definition at line 44 of file DijetMass.h.

template<class Jet >
std::string DijetMass< Jet >::SCJets
private

Definition at line 43 of file DijetMass.h.