File indexing completed on 2024-04-06 12:25:26
0001
0002
0003
0004
0005
0006
0007
0008 #include "RecoJets/JetAnalyzers/interface/DijetMass.h"
0009
0010 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
0011 #include "DataFormats/JetReco/interface/PFJetCollection.h"
0012 #include "DataFormats/JetReco/interface/GenJetCollection.h"
0013 #include "DataFormats/JetReco/interface/CaloJet.h"
0014 #include "DataFormats/JetReco/interface/PFJet.h"
0015 #include "DataFormats/JetReco/interface/GenJet.h"
0016 #include "DataFormats/Math/interface/deltaR.h"
0017 #include "FWCore/Framework/interface/Event.h"
0018 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0019 #include <TROOT.h>
0020 #include <TSystem.h>
0021 #include <TFile.h>
0022 #include <TCanvas.h>
0023 #include <cmath>
0024 using namespace edm;
0025 using namespace reco;
0026 using namespace std;
0027
0028 template <class Jet>
0029 DijetMass<Jet>::DijetMass(const edm::ParameterSet& cfg) {
0030 PtHistMax = cfg.getUntrackedParameter<double>("PtHistMax", 3000.0);
0031 EtaMax = cfg.getUntrackedParameter<double>("EtaMax", 1.3);
0032 histogramFile = cfg.getUntrackedParameter<std::string>("HistoFileName", "DijetMassHistos.root");
0033
0034 AKJets = cfg.getParameter<std::string>("AKJets");
0035 AKCorJets = cfg.getParameter<std::string>("AKCorrectedJets");
0036 ICJets = cfg.getParameter<std::string>("ICJets");
0037 ICCorJets = cfg.getParameter<std::string>("ICCorrectedJets");
0038 SCJets = cfg.getParameter<std::string>("SCJets");
0039 SCCorJets = cfg.getParameter<std::string>("SCCorrectedJets");
0040 KTJets = cfg.getParameter<std::string>("KTJets");
0041 KTCorJets = cfg.getParameter<std::string>("KTCorrectedJets");
0042 }
0043
0044 template <class Jet>
0045 void DijetMass<Jet>::beginJob() {
0046 cout << "DijetMass: Maximum bin edge for Pt Hists = " << PtHistMax << endl;
0047 numJets = 2;
0048
0049
0050 evtCount = 0;
0051
0052
0053 m_file = new TFile(histogramFile.c_str(), "RECREATE");
0054
0055
0056
0057
0058 ptAKunc = TH1F("ptAKunc", "p_{T} of leading Jets (AK)", 50, 0.0, PtHistMax);
0059 etaAKunc = TH1F("etaAKunc", "#eta of leading Jets (AK)", 23, -1.0, 1.0);
0060 phiAKunc = TH1F("phiAKunc", "#phi of leading Jets (AK)", 72, -M_PI, M_PI);
0061 m2jAKunc = TH1F("m2jAKunc", "Dijet Mass of leading Jets (AK)", 100, 0.0, 2 * PtHistMax);
0062
0063
0064 ptAKcor = TH1F("ptAKcor", "p_{T} of leading Corrected Jets (AK)", 50, 0.0, PtHistMax);
0065 etaAKcor = TH1F("etaAKcor", "#eta of leading Corrected Jets (AK)", 23, -1.0, 1.0);
0066 phiAKcor = TH1F("phiAKcor", "#phi of leading Corrected Jets (AK)", 72, -M_PI, M_PI);
0067 m2jAKcor = TH1F("m2jAKcor", "Dijet Mass of leading Corrected Jets (AK)", 100, 0.0, 2 * PtHistMax);
0068
0069
0070 ptICunc = TH1F("ptICunc", "p_{T} of leading Jets (IC)", 50, 0.0, PtHistMax);
0071 etaICunc = TH1F("etaICunc", "#eta of leading Jets (IC)", 23, -1.0, 1.0);
0072 phiICunc = TH1F("phiICunc", "#phi of leading Jets (IC)", 72, -M_PI, M_PI);
0073 m2jICunc = TH1F("m2jICunc", "Dijet Mass of leading Jets (IC)", 100, 0.0, 2 * PtHistMax);
0074
0075
0076 ptICcor = TH1F("ptICcor", "p_{T} of leading Corrected Jets (IC)", 50, 0.0, PtHistMax);
0077 etaICcor = TH1F("etaICcor", "#eta of leading Corrected Jets (IC)", 23, -1.0, 1.0);
0078 phiICcor = TH1F("phiICcor", "#phi of leading Corrected Jets (IC)", 72, -M_PI, M_PI);
0079 m2jICcor = TH1F("m2jICcor", "Dijet Mass of leading Corrected Jets (IC)", 100, 0.0, 2 * PtHistMax);
0080
0081
0082 ptKTunc = TH1F("ptKTunc", "p_{T} of leading Jets (KT)", 50, 0.0, PtHistMax);
0083 etaKTunc = TH1F("etaKTunc", "#eta of leading Jets (KT)", 23, -1.0, 1.0);
0084 phiKTunc = TH1F("phiKTunc", "#phi of leading Jets (KT)", 72, -M_PI, M_PI);
0085 m2jKTunc = TH1F("m2jKTunc", "Dijet Mass of leading Jets (KT)", 100, 0.0, 2 * PtHistMax);
0086
0087
0088 ptKTcor = TH1F("ptKTcor", "p_{T} of leading Corrected Jets (KT)", 50, 0.0, PtHistMax);
0089 etaKTcor = TH1F("etaKTcor", "#eta of leading Corrected Jets (KT)", 23, -1.0, 1.0);
0090 phiKTcor = TH1F("phiKTcor", "#phi of leading Corrected Jets (KT)", 72, -M_PI, M_PI);
0091 m2jKTcor = TH1F("m2jKTcor", "Dijet Mass of leading Corrected Jets (KT)", 100, 0.0, 2 * PtHistMax);
0092
0093
0094 ptSCunc = TH1F("ptSCunc", "p_{T} of leading Jets (SC)", 50, 0.0, PtHistMax);
0095 etaSCunc = TH1F("etaSCunc", "#eta of leading Jets (SC)", 23, -1.0, 1.0);
0096 phiSCunc = TH1F("phiSCunc", "#phi of leading Jets (SC)", 72, -M_PI, M_PI);
0097 m2jSCunc = TH1F("m2jSCunc", "Dijet Mass of leading Jets (SC)", 100, 0.0, 2 * PtHistMax);
0098
0099
0100 ptSCcor = TH1F("ptSCcor", "p_{T} of leading Corrected Jets (SC)", 50, 0.0, PtHistMax);
0101 etaSCcor = TH1F("etaSCcor", "#eta of leading Corrected Jets (SC)", 23, -1.0, 1.0);
0102 phiSCcor = TH1F("phiSCcor", "#phi of leading Corrected Jets (SC)", 72, -M_PI, M_PI);
0103 m2jSCcor = TH1F("m2jSCcor", "Dijet Mass of leading Corrected Jets (SC)", 100, 0.0, 2 * PtHistMax);
0104 }
0105
0106 template <class Jet>
0107 void DijetMass<Jet>::analyze(const Event& evt, const EventSetup& es) {
0108 evtCount++;
0109 math::XYZTLorentzVector p4jet[2];
0110 int jetInd;
0111 Handle<JetCollection> Jets;
0112
0113
0114 typename JetCollection::const_iterator i_jet;
0115
0116
0117 evt.getByLabel(AKJets, Jets);
0118 jetInd = 0;
0119 for (i_jet = Jets->begin(); i_jet != Jets->end() && jetInd < 2; ++i_jet) {
0120 p4jet[jetInd] = i_jet->p4();
0121 jetInd++;
0122 }
0123 if (jetInd == 2 && abs(p4jet[0].eta()) < EtaMax && abs(p4jet[1].eta()) < EtaMax) {
0124 m2jAKunc.Fill((p4jet[0] + p4jet[1]).mass());
0125 ptAKunc.Fill(p4jet[0].Pt());
0126 ptAKunc.Fill(p4jet[1].Pt());
0127 etaAKunc.Fill(p4jet[0].eta());
0128 etaAKunc.Fill(p4jet[1].eta());
0129 phiAKunc.Fill(p4jet[0].phi());
0130 phiAKunc.Fill(p4jet[1].phi());
0131 }
0132
0133
0134 evt.getByLabel(AKCorJets, Jets);
0135 jetInd = 0;
0136 for (i_jet = Jets->begin(); i_jet != Jets->end() && jetInd < 2; ++i_jet) {
0137 p4jet[jetInd] = i_jet->p4();
0138 jetInd++;
0139 }
0140 if (jetInd == 2 && abs(p4jet[0].eta()) < EtaMax && abs(p4jet[1].eta()) < EtaMax) {
0141 m2jAKcor.Fill((p4jet[0] + p4jet[1]).mass());
0142 ptAKcor.Fill(p4jet[0].Pt());
0143 ptAKcor.Fill(p4jet[1].Pt());
0144 etaAKcor.Fill(p4jet[0].eta());
0145 etaAKcor.Fill(p4jet[1].eta());
0146 phiAKcor.Fill(p4jet[0].phi());
0147 phiAKcor.Fill(p4jet[1].phi());
0148 }
0149
0150
0151 evt.getByLabel(ICJets, Jets);
0152 jetInd = 0;
0153 for (i_jet = Jets->begin(); i_jet != Jets->end() && jetInd < 2; ++i_jet) {
0154 p4jet[jetInd] = i_jet->p4();
0155 jetInd++;
0156 }
0157 if (jetInd == 2 && abs(p4jet[0].eta()) < EtaMax && abs(p4jet[1].eta()) < EtaMax) {
0158 m2jICunc.Fill((p4jet[0] + p4jet[1]).mass());
0159 ptICunc.Fill(p4jet[0].Pt());
0160 ptICunc.Fill(p4jet[1].Pt());
0161 etaICunc.Fill(p4jet[0].eta());
0162 etaICunc.Fill(p4jet[1].eta());
0163 phiICunc.Fill(p4jet[0].phi());
0164 phiICunc.Fill(p4jet[1].phi());
0165 }
0166
0167
0168 evt.getByLabel(ICCorJets, Jets);
0169 jetInd = 0;
0170 for (i_jet = Jets->begin(); i_jet != Jets->end() && jetInd < 2; ++i_jet) {
0171 p4jet[jetInd] = i_jet->p4();
0172 jetInd++;
0173 }
0174 if (jetInd == 2 && abs(p4jet[0].eta()) < EtaMax && abs(p4jet[1].eta()) < EtaMax) {
0175 m2jICcor.Fill((p4jet[0] + p4jet[1]).mass());
0176 ptICcor.Fill(p4jet[0].Pt());
0177 ptICcor.Fill(p4jet[1].Pt());
0178 etaICcor.Fill(p4jet[0].eta());
0179 etaICcor.Fill(p4jet[1].eta());
0180 phiICcor.Fill(p4jet[0].phi());
0181 phiICcor.Fill(p4jet[1].phi());
0182 }
0183
0184
0185 evt.getByLabel(KTJets, Jets);
0186 jetInd = 0;
0187 for (i_jet = Jets->begin(); i_jet != Jets->end() && jetInd < 2; ++i_jet) {
0188 p4jet[jetInd] = i_jet->p4();
0189 jetInd++;
0190 }
0191 if (jetInd == 2 && abs(p4jet[0].eta()) < EtaMax && abs(p4jet[1].eta()) < EtaMax) {
0192 m2jKTunc.Fill((p4jet[0] + p4jet[1]).mass());
0193 ptKTunc.Fill(p4jet[0].Pt());
0194 ptKTunc.Fill(p4jet[1].Pt());
0195 etaKTunc.Fill(p4jet[0].eta());
0196 etaKTunc.Fill(p4jet[1].eta());
0197 phiKTunc.Fill(p4jet[0].phi());
0198 phiKTunc.Fill(p4jet[1].phi());
0199 }
0200
0201
0202 evt.getByLabel(KTCorJets, Jets);
0203 jetInd = 0;
0204 for (i_jet = Jets->begin(); i_jet != Jets->end() && jetInd < 2; ++i_jet) {
0205 p4jet[jetInd] = i_jet->p4();
0206 jetInd++;
0207 }
0208 if (jetInd == 2 && abs(p4jet[0].eta()) < EtaMax && abs(p4jet[1].eta()) < EtaMax) {
0209 m2jKTcor.Fill((p4jet[0] + p4jet[1]).mass());
0210 ptKTcor.Fill(p4jet[0].Pt());
0211 ptKTcor.Fill(p4jet[1].Pt());
0212 etaKTcor.Fill(p4jet[0].eta());
0213 etaKTcor.Fill(p4jet[1].eta());
0214 phiKTcor.Fill(p4jet[0].phi());
0215 phiKTcor.Fill(p4jet[1].phi());
0216 }
0217
0218
0219 evt.getByLabel(SCJets, Jets);
0220 jetInd = 0;
0221 for (i_jet = Jets->begin(); i_jet != Jets->end() && jetInd < 2; ++i_jet) {
0222 p4jet[jetInd] = i_jet->p4();
0223 jetInd++;
0224 }
0225 if (jetInd == 2 && abs(p4jet[0].eta()) < EtaMax && abs(p4jet[1].eta()) < EtaMax) {
0226 m2jSCunc.Fill((p4jet[0] + p4jet[1]).mass());
0227 ptSCunc.Fill(p4jet[0].Pt());
0228 ptSCunc.Fill(p4jet[1].Pt());
0229 etaSCunc.Fill(p4jet[0].eta());
0230 etaSCunc.Fill(p4jet[1].eta());
0231 phiSCunc.Fill(p4jet[0].phi());
0232 phiSCunc.Fill(p4jet[1].phi());
0233 }
0234
0235
0236 evt.getByLabel(SCCorJets, Jets);
0237 jetInd = 0;
0238 for (i_jet = Jets->begin(); i_jet != Jets->end() && jetInd < 2; ++i_jet) {
0239 p4jet[jetInd] = i_jet->p4();
0240 jetInd++;
0241 }
0242 if (jetInd == 2 && abs(p4jet[0].eta()) < EtaMax && abs(p4jet[1].eta()) < EtaMax) {
0243 m2jSCcor.Fill((p4jet[0] + p4jet[1]).mass());
0244 ptSCcor.Fill(p4jet[0].Pt());
0245 ptSCcor.Fill(p4jet[1].Pt());
0246 etaSCcor.Fill(p4jet[0].eta());
0247 etaSCcor.Fill(p4jet[1].eta());
0248 phiSCcor.Fill(p4jet[0].phi());
0249 phiSCcor.Fill(p4jet[1].phi());
0250 }
0251 }
0252
0253 template <class Jet>
0254 void DijetMass<Jet>::endJob() {
0255
0256 m_file->cd();
0257
0258 ptAKunc.Write();
0259 etaAKunc.Write();
0260 phiAKunc.Write();
0261 m2jAKunc.Write();
0262
0263 ptAKcor.Write();
0264 etaAKcor.Write();
0265 phiAKcor.Write();
0266 m2jAKcor.Write();
0267
0268 ptICunc.Write();
0269 etaICunc.Write();
0270 phiICunc.Write();
0271 m2jICunc.Write();
0272
0273 ptICcor.Write();
0274 etaICcor.Write();
0275 phiICcor.Write();
0276 m2jICcor.Write();
0277
0278 ptKTunc.Write();
0279 etaKTunc.Write();
0280 phiKTunc.Write();
0281 m2jKTunc.Write();
0282
0283 ptKTcor.Write();
0284 etaKTcor.Write();
0285 phiKTcor.Write();
0286 m2jKTcor.Write();
0287
0288 ptSCunc.Write();
0289 etaSCunc.Write();
0290 phiSCunc.Write();
0291 m2jSCunc.Write();
0292
0293 ptSCcor.Write();
0294 etaSCcor.Write();
0295 phiSCcor.Write();
0296 m2jSCcor.Write();
0297
0298 m_file->Close();
0299 }
0300 #include "FWCore/Framework/interface/MakerMacros.h"
0301
0302 typedef DijetMass<CaloJet> DijetMassCaloJets;
0303 DEFINE_FWK_MODULE(DijetMassCaloJets);
0304
0305 typedef DijetMass<GenJet> DijetMassGenJets;
0306 DEFINE_FWK_MODULE(DijetMassGenJets);
0307
0308 typedef DijetMass<PFJet> DijetMassPFJets;
0309 DEFINE_FWK_MODULE(DijetMassPFJets);