CMS 3D CMS Logo

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::one::EDAnalyzer<> edm::one::EDAnalyzerBase edm::EDConsumerBase

Public Member Functions

 DijetMass (const edm::ParameterSet &)
 
- Public Member Functions inherited from edm::one::EDAnalyzer<>
 EDAnalyzer ()=default
 
 EDAnalyzer (const EDAnalyzer &)=delete
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
const EDAnalyzeroperator= (const EDAnalyzer &)=delete
 
bool wantsGlobalLuminosityBlocks () const noexcept final
 
bool wantsGlobalRuns () const noexcept final
 
bool wantsInputProcessBlocks () const noexcept final
 
bool wantsProcessBlocks () const noexcept final
 
- Public Member Functions inherited from edm::one::EDAnalyzerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const noexcept
 
bool wantsStreamRuns () const noexcept
 
 ~EDAnalyzerBase () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESResolverIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
std::vector< ESResolverIndex > const & esGetTokenIndicesVector (edm::Transition iTrans) const
 
std::vector< ESRecordIndex > const & esGetTokenRecordIndicesVector (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::array< std::vector< ModuleDescription const *> *, NumBranchTypes > &modulesAll, std::vector< ModuleProcessName > &modulesInPreviousProcesses, ProductRegistry const &preg, std::map< std::string, ModuleDescription const *> const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
void selectInputProcessBlocks (ProductRegistry const &productRegistry, ProcessBlockHelperBase const &processBlockHelperBase)
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProductResolverIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Types

typedef std::vector< Jet > JetCollection
 

Private Member Functions

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

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::one::EDAnalyzerBase
typedef EDAnalyzerBase ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::one::EDAnalyzerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
template<BranchType B = InEvent>
EDConsumerBaseAdaptor< Bconsumes (edm::InputTag tag) noexcept
 
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 ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<Transition Tr = Transition::Event>
constexpr auto esConsumes ()
 
template<Transition Tr = Transition::Event>
auto esConsumes (ESInputTag tag)
 
template<Transition Tr = Transition::Event>
ESGetTokenGeneric esConsumes (eventsetup::EventSetupRecordKey const &iRecord, eventsetup::DataKey const &iKey)
 Used with EventSetupRecord::doGet. More...
 
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)
 
void resetItemsToGetFrom (BranchType iType)
 

Detailed Description

template<class Jet>
class DijetMass< Jet >

Definition at line 23 of file DijetMass.h.

Member Typedef Documentation

◆ JetCollection

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

Definition at line 28 of file DijetMass.h.

Constructor & Destructor Documentation

◆ DijetMass()

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

Definition at line 29 of file DijetMass.cc.

References looper::cfg, MonitorTrackInnerTrackMuons_cff::EtaMax, and AlCaHLTBitMon_QueryRunRegistry::string.

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

Member Function Documentation

◆ analyze()

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

Implements edm::one::EDAnalyzerBase.

Definition at line 107 of file DijetMass.cc.

References funct::abs(), PVValHelper::eta, MonitorTrackInnerTrackMuons_cff::EtaMax, edm::Event::getByLabel(), METSkim_cff::Jets, and EgHLTOffHistBins_cfi::mass.

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

◆ beginJob()

template<class Jet >
void DijetMass< Jet >::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 45 of file DijetMass.cc.

References gather_cfg::cout, and M_PI.

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

◆ endJob()

template<class Jet >
void DijetMass< Jet >::endJob ( void  )
overrideprivatevirtual

Reimplemented from edm::one::EDAnalyzerBase.

Definition at line 254 of file DijetMass.cc.

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

Member Data Documentation

◆ AKCorJets

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

Definition at line 39 of file DijetMass.h.

◆ AKJets

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

Definition at line 38 of file DijetMass.h.

◆ etaAKcor

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

Definition at line 49 of file DijetMass.h.

◆ etaAKunc

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

Definition at line 48 of file DijetMass.h.

◆ etaICcor

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

Definition at line 51 of file DijetMass.h.

◆ etaICunc

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

Definition at line 50 of file DijetMass.h.

◆ etaKTcor

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

Definition at line 53 of file DijetMass.h.

◆ etaKTunc

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

Definition at line 52 of file DijetMass.h.

◆ EtaMax

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

Definition at line 36 of file DijetMass.h.

◆ etaSCcor

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

Definition at line 55 of file DijetMass.h.

◆ etaSCunc

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

Definition at line 54 of file DijetMass.h.

◆ evtCount

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

Definition at line 61 of file DijetMass.h.

◆ histogramFile

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

Definition at line 37 of file DijetMass.h.

◆ ICCorJets

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

Definition at line 41 of file DijetMass.h.

◆ ICJets

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

Definition at line 40 of file DijetMass.h.

◆ KTCorJets

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

Definition at line 45 of file DijetMass.h.

◆ KTJets

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

Definition at line 44 of file DijetMass.h.

◆ m2jAKcor

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

Definition at line 49 of file DijetMass.h.

◆ m2jAKunc

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

Definition at line 48 of file DijetMass.h.

◆ m2jICcor

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

Definition at line 51 of file DijetMass.h.

◆ m2jICunc

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

Definition at line 50 of file DijetMass.h.

◆ m2jKTcor

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

Definition at line 53 of file DijetMass.h.

◆ m2jKTunc

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

Definition at line 52 of file DijetMass.h.

◆ m2jSCcor

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

Definition at line 55 of file DijetMass.h.

◆ m2jSCunc

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

Definition at line 54 of file DijetMass.h.

◆ m_file

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

Definition at line 58 of file DijetMass.h.

◆ numJets

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

Definition at line 62 of file DijetMass.h.

◆ phiAKcor

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

Definition at line 49 of file DijetMass.h.

◆ phiAKunc

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

Definition at line 48 of file DijetMass.h.

◆ phiICcor

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

Definition at line 51 of file DijetMass.h.

◆ phiICunc

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

Definition at line 50 of file DijetMass.h.

◆ phiKTcor

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

Definition at line 53 of file DijetMass.h.

◆ phiKTunc

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

Definition at line 52 of file DijetMass.h.

◆ phiSCcor

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

Definition at line 55 of file DijetMass.h.

◆ phiSCunc

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

Definition at line 54 of file DijetMass.h.

◆ ptAKcor

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

Definition at line 49 of file DijetMass.h.

◆ ptAKunc

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

Definition at line 48 of file DijetMass.h.

◆ PtHistMax

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

Definition at line 35 of file DijetMass.h.

◆ ptICcor

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

Definition at line 51 of file DijetMass.h.

◆ ptICunc

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

Definition at line 50 of file DijetMass.h.

◆ ptKTcor

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

Definition at line 53 of file DijetMass.h.

◆ ptKTunc

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

Definition at line 52 of file DijetMass.h.

◆ ptSCcor

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

Definition at line 55 of file DijetMass.h.

◆ ptSCunc

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

Definition at line 54 of file DijetMass.h.

◆ SCCorJets

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

Definition at line 43 of file DijetMass.h.

◆ SCJets

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

Definition at line 42 of file DijetMass.h.