File indexing completed on 2024-04-06 12:09:32
0001
0002
0003
0004
0005
0006
0007
0008 #include "DataFormats/Common/interface/Handle.h"
0009
0010 #include "DQMServices/Core/interface/DQMStore.h"
0011
0012 #include "FWCore/Framework/interface/Frameworkfwd.h"
0013 #include "FWCore/Framework/interface/Event.h"
0014 #include "FWCore/Framework/interface/MakerMacros.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0017 #include "FWCore/ServiceRegistry/interface/Service.h"
0018
0019
0020
0021
0022 #include "DataFormats/CaloTowers/interface/CaloTowerCollection.h"
0023 #include "DataFormats/CaloTowers/interface/CaloTowerDetId.h"
0024 #include "DataFormats/EcalDetId/interface/EcalSubdetector.h"
0025 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0026
0027 #include "DataFormats/JetReco/interface/CaloJet.h"
0028 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
0029 #include "DataFormats/JetReco/interface/PFJet.h"
0030 #include "DataFormats/JetReco/interface/PFJetCollection.h"
0031
0032
0033
0034 #include "DataFormats/METReco/interface/MET.h"
0035 #include "DataFormats/METReco/interface/METCollection.h"
0036 #include "DataFormats/METReco/interface/CaloMETCollection.h"
0037 #include "DataFormats/METReco/interface/CaloMET.h"
0038 #include "DataFormats/METReco/interface/PFMET.h"
0039 #include "DataFormats/METReco/interface/PFMETCollection.h"
0040 #include "DataFormats/METReco/interface/PFMETFwd.h"
0041
0042
0043
0044
0045 #include "DataFormats/Math/interface/LorentzVector.h"
0046 #include "DQMOffline/JetMET/interface/SUSYDQMAnalyzer.h"
0047 #include "DQMOffline/JetMET/interface/SusyDQM/alpha_T.h"
0048 #include "DQMOffline/JetMET/interface/SusyDQM/HT.h"
0049
0050
0051
0052
0053 #include "TH1.h"
0054 #include "TVector2.h"
0055 #include "TLorentzVector.h"
0056
0057
0058
0059
0060 #include <cmath>
0061 #include <vector>
0062 #include <fstream>
0063 #include <sstream>
0064 #include <string>
0065
0066 using namespace edm;
0067 using namespace reco;
0068 using namespace math;
0069 using namespace std;
0070
0071 SUSYDQMAnalyzer::SUSYDQMAnalyzer(const edm::ParameterSet& pSet) {
0072 iConfig = pSet;
0073
0074 SUSYFolder = iConfig.getParameter<std::string>("folderName");
0075
0076 thePFMETCollectionToken =
0077 consumes<reco::PFMETCollection>(iConfig.getParameter<edm::InputTag>("PFMETCollectionLabel"));
0078 theCaloMETCollectionToken =
0079 consumes<reco::CaloMETCollection>(iConfig.getParameter<edm::InputTag>("CaloMETCollectionLabel"));
0080
0081
0082
0083
0084 theCaloJetCollectionToken =
0085 consumes<reco::CaloJetCollection>(iConfig.getParameter<edm::InputTag>("CaloJetCollectionLabel"));
0086
0087 thePFJetCollectionToken =
0088 consumes<std::vector<reco::PFJet> >(iConfig.getParameter<edm::InputTag>("PFJetCollectionLabel"));
0089
0090 _ptThreshold = iConfig.getParameter<double>("ptThreshold");
0091 _maxNJets = iConfig.getParameter<double>("maxNJets");
0092 _maxAbsEta = iConfig.getParameter<double>("maxAbsEta");
0093 }
0094
0095 const char* SUSYDQMAnalyzer::messageLoggerCatregory = "SUSYDQM";
0096
0097 void SUSYDQMAnalyzer::bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& iRun, edm::EventSetup const&) {
0098
0099
0100
0101 std::string dir = SUSYFolder;
0102 dir += "HT";
0103 ibooker.setCurrentFolder(dir);
0104 hCaloHT = ibooker.book1D("Calo_HT", "", 500, 0., 2000);
0105 hPFHT = ibooker.book1D("PF_HT", "", 500, 0., 2000);
0106
0107
0108
0109
0110 dir = SUSYFolder;
0111 dir += "MET";
0112 ibooker.setCurrentFolder(dir);
0113 hCaloMET = ibooker.book1D("Calo_MET", "", 500, 0., 1000);
0114 hPFMET = ibooker.book1D("PF_MET", "", 500, 0., 1000);
0115
0116
0117
0118
0119
0120 dir = SUSYFolder;
0121 dir += "MHT";
0122 ibooker.setCurrentFolder(dir);
0123 hCaloMHT = ibooker.book1D("Calo_MHT", "", 500, 0., 1000);
0124 hPFMHT = ibooker.book1D("PF_MHT", "", 500, 0., 1000);
0125
0126
0127
0128
0129
0130 dir = SUSYFolder;
0131 dir += "Alpha_T";
0132 ibooker.setCurrentFolder(dir);
0133 hCaloAlpha_T = ibooker.book1D("Calo_AlphaT", "", 100, 0., 1.);
0134
0135 hPFAlpha_T = ibooker.book1D("PF_AlphaT", "", 100, 0., 1.);
0136
0137 }
0138
0139 void SUSYDQMAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0140
0141
0142
0143
0144
0145
0146 edm::Handle<reco::CaloJetCollection> CaloJetcoll;
0147
0148 iEvent.getByToken(theCaloJetCollectionToken, CaloJetcoll);
0149
0150 if (!CaloJetcoll.isValid())
0151 return;
0152
0153 std::vector<math::XYZTLorentzVector> Ps;
0154 for (reco::CaloJetCollection::const_iterator jet = CaloJetcoll->begin(); jet != CaloJetcoll->end(); ++jet) {
0155 if ((jet->pt() > _ptThreshold) && (abs(jet->eta()) < _maxAbsEta)) {
0156 if (Ps.size() > _maxNJets) {
0157 edm::LogInfo(messageLoggerCatregory) << "NMax Jets exceded..";
0158 break;
0159 }
0160 Ps.push_back(jet->p4());
0161 }
0162 }
0163
0164 hCaloAlpha_T->Fill(alpha_T()(Ps));
0165
0166 HT<reco::CaloJetCollection> CaloHT(CaloJetcoll, _ptThreshold, _maxAbsEta);
0167
0168 hCaloHT->Fill(CaloHT.ScalarSum);
0169 hCaloMHT->Fill(CaloHT.v.Mod());
0170
0171
0172
0173
0174 edm::Handle<reco::PFJetCollection> PFjetcoll;
0175
0176 iEvent.getByToken(thePFJetCollectionToken, PFjetcoll);
0177
0178 if (!PFjetcoll.isValid())
0179 return;
0180
0181 Ps.clear();
0182 for (reco::PFJetCollection::const_iterator jet = PFjetcoll->begin(); jet != PFjetcoll->end(); ++jet) {
0183 if ((jet->pt() > _ptThreshold) && (abs(jet->eta()) < _maxAbsEta)) {
0184 if (Ps.size() > _maxNJets) {
0185 edm::LogInfo(messageLoggerCatregory) << "NMax Jets exceded..";
0186 break;
0187 }
0188 Ps.push_back(jet->p4());
0189 }
0190 }
0191 hPFAlpha_T->Fill(alpha_T()(Ps));
0192
0193 HT<reco::PFJetCollection> PFHT(PFjetcoll, _ptThreshold, _maxAbsEta);
0194
0195 hPFHT->Fill(PFHT.ScalarSum);
0196 hPFMHT->Fill(PFHT.v.Mod());
0197
0198
0199
0200
0201
0202
0203
0204
0205
0206
0207
0208
0209
0210
0211
0212
0213
0214
0215
0216
0217
0218
0219
0220
0221
0222
0223
0224
0225
0226
0227
0228
0229
0230 edm::Handle<reco::CaloMETCollection> calometcoll;
0231 iEvent.getByToken(theCaloMETCollectionToken, calometcoll);
0232
0233 if (!calometcoll.isValid())
0234 return;
0235
0236 const CaloMETCollection* calometcol = calometcoll.product();
0237 const CaloMET* calomet;
0238 calomet = &(calometcol->front());
0239
0240 hCaloMET->Fill(calomet->pt());
0241
0242
0243
0244
0245 edm::Handle<reco::PFMETCollection> pfmetcoll;
0246 iEvent.getByToken(thePFMETCollectionToken, pfmetcoll);
0247
0248 if (!pfmetcoll.isValid())
0249 return;
0250
0251 const PFMETCollection* pfmetcol = pfmetcoll.product();
0252 const PFMET* pfmet;
0253 pfmet = &(pfmetcol->front());
0254
0255 hPFMET->Fill(pfmet->pt());
0256
0257
0258
0259
0260
0261
0262
0263
0264
0265
0266
0267
0268
0269
0270 }
0271
0272 SUSYDQMAnalyzer::~SUSYDQMAnalyzer() {}