00001
00002
00003
00004
00005
00006
00007
00008 #include "DataFormats/Common/interface/Handle.h"
00009
00010 #include "DQMServices/Core/interface/DQMStore.h"
00011 #include "DQMServices/Core/interface/MonitorElement.h"
00012
00013 #include "FWCore/Framework/interface/Frameworkfwd.h"
00014 #include "FWCore/Framework/interface/Event.h"
00015 #include "FWCore/Framework/interface/MakerMacros.h"
00016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00017 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00018 #include "FWCore/ServiceRegistry/interface/Service.h"
00019
00020
00021
00022
00023 #include "DataFormats/CaloTowers/interface/CaloTowerCollection.h"
00024 #include "DataFormats/CaloTowers/interface/CaloTowerDetId.h"
00025 #include "DataFormats/EcalDetId/interface/EcalSubdetector.h"
00026 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
00027
00028 #include "DataFormats/JetReco/interface/CaloJet.h"
00029 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00030 #include "DataFormats/JetReco/interface/PFJet.h"
00031 #include "DataFormats/JetReco/interface/PFJetCollection.h"
00032 #include "DataFormats/JetReco/interface/JPTJet.h"
00033 #include "DataFormats/JetReco/interface/JPTJetCollection.h"
00034
00035 #include "DataFormats/METReco/interface/MET.h"
00036 #include "DataFormats/METReco/interface/METCollection.h"
00037 #include "DataFormats/METReco/interface/CaloMETCollection.h"
00038 #include "DataFormats/METReco/interface/CaloMET.h"
00039 #include "DataFormats/METReco/interface/PFMET.h"
00040 #include "DataFormats/METReco/interface/PFMETCollection.h"
00041 #include "DataFormats/METReco/interface/PFMETFwd.h"
00042
00043
00044
00045
00046 #include "DataFormats/Math/interface/LorentzVector.h"
00047 #include "DQMOffline/JetMET/interface/SUSYDQMAnalyzer.h"
00048 #include "DQMOffline/JetMET/interface/SusyDQM/alpha_T.h"
00049 #include "DQMOffline/JetMET/interface/SusyDQM/HT.h"
00050
00051
00052
00053
00054 #include "TH1.h"
00055 #include "TVector2.h"
00056 #include "TLorentzVector.h"
00057
00058
00059
00060
00061 #include <cmath>
00062 #include <vector>
00063 #include <fstream>
00064 #include <sstream>
00065 #include <string>
00066
00067 using namespace edm;
00068 using namespace reco;
00069 using namespace math;
00070 using namespace std;
00071
00072 SUSYDQMAnalyzer::SUSYDQMAnalyzer( const edm::ParameterSet& pSet)
00073 {
00074 iConfig = pSet;
00075
00076 SUSYFolder = iConfig.getParameter<std::string>("folderName");
00077 dqm = edm::Service<DQMStore>().operator->();
00078 }
00079
00080 const char* SUSYDQMAnalyzer::messageLoggerCatregory = "SUSYDQM";
00081
00082 void SUSYDQMAnalyzer::beginJob(void){
00083
00084 thePFMETCollectionLabel = iConfig.getParameter<edm::InputTag>("PFMETCollectionLabel");
00085 theCaloMETCollectionLabel = iConfig.getParameter<edm::InputTag>("CaloMETCollectionLabel");
00086 theTCMETCollectionLabel = iConfig.getParameter<edm::InputTag>("TCMETCollectionLabel");
00087
00088 theCaloJetCollectionLabel = iConfig.getParameter<edm::InputTag>("CaloJetCollectionLabel");
00089 theJPTJetCollectionLabel = iConfig.getParameter<edm::InputTag>("JPTJetCollectionLabel");
00090 thePFJetCollectionLabel = iConfig.getParameter<edm::InputTag>("PFJetCollectionLabel");
00091
00092 _ptThreshold = iConfig.getParameter<double>("ptThreshold");
00093 _maxNJets = iConfig.getParameter<double>("maxNJets");
00094 _maxAbsEta = iConfig.getParameter<double>("maxAbsEta");
00095 }
00096
00097 void SUSYDQMAnalyzer::beginRun(const edm::Run& iRun, const edm::EventSetup& iSetup){
00098 if( dqm ) {
00099
00100
00101
00102 std::string dir=SUSYFolder;
00103 dir+="HT";
00104 dqm->setCurrentFolder(dir);
00105 hCaloHT = dqm->book1D("Calo_HT", "", 500, 0., 2000);
00106 hPFHT = dqm->book1D("PF_HT" , "", 500, 0., 2000);
00107 hJPTHT = dqm->book1D("JPT_HT" , "", 500, 0., 2000);
00108
00109
00110
00111
00112 dir=SUSYFolder;
00113 dir+="MET";
00114 dqm->setCurrentFolder(dir);
00115 hCaloMET = dqm->book1D("Calo_MET", "", 500, 0., 1000);
00116 hPFMET = dqm->book1D("PF_MET" , "", 500, 0., 1000);
00117 hTCMET = dqm->book1D("TC_MET" , "", 500, 0., 1000);
00118
00119
00120
00121
00122 dir=SUSYFolder;
00123 dir+="MHT";
00124 dqm->setCurrentFolder(dir);
00125 hCaloMHT = dqm->book1D("Calo_MHT", "", 500, 0., 1000);
00126 hPFMHT = dqm->book1D("PF_MHT" , "", 500, 0., 1000);
00127 hJPTMHT = dqm->book1D("JPT_MHT" , "", 500, 0., 1000);
00128
00129
00130
00131
00132 dir=SUSYFolder;
00133 dir+="Alpha_T";
00134 dqm->setCurrentFolder(dir);
00135 hCaloAlpha_T = dqm->book1D("Calo_AlphaT", "", 100, 0., 1.);
00136 hJPTAlpha_T = dqm->book1D("PF_AlphaT" , "", 100, 0., 1.);
00137 hPFAlpha_T = dqm->book1D("JPT_AlphaT" , "", 100, 0., 1.);
00138 }
00139 }
00140
00141 void SUSYDQMAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00142 {
00143
00144
00145
00146
00147
00148
00149 edm::Handle<reco::CaloJetCollection> CaloJetcoll;
00150
00151 iEvent.getByLabel(theCaloJetCollectionLabel, CaloJetcoll);
00152
00153 if(!CaloJetcoll.isValid()) return;
00154
00155 std::vector<math::XYZTLorentzVector> Ps;
00156 for (reco::CaloJetCollection::const_iterator jet = CaloJetcoll->begin(); jet!=CaloJetcoll->end(); ++jet){
00157 if ((jet->pt()>_ptThreshold) && (abs(jet->eta()) < _maxAbsEta)){
00158 if(Ps.size()>_maxNJets) {
00159 edm::LogInfo(messageLoggerCatregory)<<"NMax Jets exceded..";
00160 break;
00161 }
00162 Ps.push_back(jet->p4());
00163 }
00164 }
00165
00166 hCaloAlpha_T->Fill( alpha_T()(Ps));
00167
00168 HT< reco::CaloJetCollection > CaloHT(CaloJetcoll, _ptThreshold, _maxAbsEta);
00169
00170 hCaloHT->Fill(CaloHT.ScalarSum);
00171 hCaloMHT->Fill(CaloHT.v.Mod());
00172
00173
00174
00175
00176 edm::Handle<reco::PFJetCollection> PFjetcoll;
00177
00178 iEvent.getByLabel(thePFJetCollectionLabel, PFjetcoll);
00179
00180 if(!PFjetcoll.isValid()) return;
00181
00182 Ps.clear();
00183 for (reco::PFJetCollection::const_iterator jet = PFjetcoll->begin(); jet!=PFjetcoll->end(); ++jet){
00184 if ((jet->pt()>_ptThreshold) && (abs(jet->eta()) < _maxAbsEta)){
00185 if(Ps.size()>_maxNJets) {
00186 edm::LogInfo(messageLoggerCatregory)<<"NMax Jets exceded..";
00187 break;
00188 }
00189 Ps.push_back(jet->p4());
00190 }
00191 }
00192 hPFAlpha_T->Fill( alpha_T()(Ps));
00193
00194 HT<reco::PFJetCollection> PFHT(PFjetcoll, _ptThreshold, _maxAbsEta);
00195
00196 hPFHT->Fill(PFHT.ScalarSum);
00197 hPFMHT->Fill(PFHT.v.Mod());
00198
00199
00200
00201
00202 edm::Handle<reco::JPTJetCollection> JPTjetcoll;
00203
00204 iEvent.getByLabel(theJPTJetCollectionLabel, JPTjetcoll);
00205
00206 if(!JPTjetcoll.isValid()) return;
00207
00208 Ps.clear();
00209 for (reco::JPTJetCollection::const_iterator jet = JPTjetcoll->begin(); jet!=JPTjetcoll->end(); ++jet){
00210 if ((jet->pt()>_ptThreshold) && (abs(jet->eta())<_maxAbsEta)){
00211 if(Ps.size()>_maxNJets) {
00212 edm::LogInfo(messageLoggerCatregory)<<"NMax Jets exceded..";
00213 break;
00214 }
00215 Ps.push_back(jet->p4());
00216 }
00217 }
00218 hJPTAlpha_T->Fill( alpha_T()(Ps));
00219
00220 HT<reco::JPTJetCollection> JPTHT(JPTjetcoll, _ptThreshold, _maxAbsEta);
00221
00222 hJPTHT->Fill(JPTHT.ScalarSum);
00223 hJPTMHT->Fill(JPTHT.v.Mod());
00224
00225
00226
00227
00228
00229
00230
00231 edm::Handle<reco::CaloMETCollection> calometcoll;
00232 iEvent.getByLabel(theCaloMETCollectionLabel, calometcoll);
00233
00234 if(!calometcoll.isValid()) return;
00235
00236 const CaloMETCollection *calometcol = calometcoll.product();
00237 const CaloMET *calomet;
00238 calomet = &(calometcol->front());
00239
00240 hCaloMET->Fill(calomet->pt());
00241
00242
00243
00244
00245 edm::Handle<reco::PFMETCollection> pfmetcoll;
00246 iEvent.getByLabel(thePFMETCollectionLabel, pfmetcoll);
00247
00248 if(!pfmetcoll.isValid()) return;
00249
00250 const PFMETCollection *pfmetcol = pfmetcoll.product();
00251 const PFMET *pfmet;
00252 pfmet = &(pfmetcol->front());
00253
00254 hPFMET->Fill(pfmet->pt());
00255
00256
00257
00258
00259 edm::Handle<reco::METCollection> tcmetcoll;
00260 iEvent.getByLabel(theTCMETCollectionLabel, tcmetcoll);
00261
00262 if(!tcmetcoll.isValid()) return;
00263
00264 const METCollection *tcmetcol = tcmetcoll.product();
00265 const MET *tcmet;
00266 tcmet = &(tcmetcol->front());
00267
00268 hTCMET->Fill(tcmet->pt());
00269
00270 }
00271
00272 void SUSYDQMAnalyzer::endRun(const edm::Run&, const edm::EventSetup&){
00273 }
00274
00275 SUSYDQMAnalyzer::~SUSYDQMAnalyzer(){
00276 }