00001 #include "DQM/DataScouting/plugins/RazorVarAnalyzer.h"
00002
00003 #include "DataFormats/JetReco/interface/CaloJet.h"
00004 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00005
00006 #include "DataFormats/EgammaCandidates/interface/Electron.h"
00007 #include "DataFormats/EgammaCandidates/interface/ElectronFwd.h"
00008
00009 #include "DataFormats/MuonReco/interface/Muon.h"
00010 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
00011
00012 #include <cmath>
00013
00014
00015
00016 RazorVarAnalyzer::RazorVarAnalyzer( const edm::ParameterSet & conf ):
00017 ScoutingAnalyzerBase(conf),
00018 m_eleCollectionTag(conf.getUntrackedParameter<edm::InputTag>("eleCollectionName",edm::InputTag("hltPixelMatchElectronsActivity"))),
00019 m_jetCollectionTag(conf.getUntrackedParameter<edm::InputTag>("jetCollectionName",edm::InputTag("hltCaloJetIDPassed"))),
00020 m_muCollectionTag(conf.getUntrackedParameter<edm::InputTag>("muCollectionName",edm::InputTag("hltL3MuonCandidates"))),
00021 m_razorVarCollectionTag(conf.getUntrackedParameter<edm::InputTag>("razorVarCollectionName")){
00022 }
00023
00024
00025
00026 RazorVarAnalyzer::~RazorVarAnalyzer(){}
00027
00028
00029
00030 void RazorVarAnalyzer::analyze( const edm::Event & iEvent, const edm::EventSetup & c ){
00031
00032
00033 edm::Handle<reco::CaloJetCollection> calojets_handle;
00034 iEvent.getByLabel(m_jetCollectionTag,calojets_handle);
00035
00036 unsigned int njets = 0;
00037 for(reco::CaloJetCollection::const_iterator it = calojets_handle->begin(); it != calojets_handle->end(); ++it){
00038 if(it->pt() >= 30. && fabs(it->eta()) <= 3.0){
00039 njets++;
00040 }
00041 }
00042
00043
00044 edm::Handle<std::vector<reco::RecoChargedCandidate> > muon_handle;
00045 iEvent.getByLabel(m_muCollectionTag,muon_handle);
00046
00047 unsigned int nmu_loose = 0;
00048 unsigned int nmu_tight = 0;
00049 if(muon_handle.isValid()){
00050 for(std::vector<reco::RecoChargedCandidate>::const_iterator it = muon_handle->begin(); it != muon_handle->end(); ++it){
00051 if(it->pt() >= 15 && fabs(it->eta()) <= 2.1) nmu_tight++;
00052 if(it->pt() >= 10 && fabs(it->eta()) <= 2.4) nmu_loose++;
00053 }
00054 }
00055
00056
00057 edm::Handle<reco::ElectronCollection> ele_handle;
00058 iEvent.getByLabel(m_eleCollectionTag,ele_handle);
00059
00060 unsigned int nele_loose = 0;
00061 unsigned int nele_tight = 0;
00062 if(ele_handle.isValid()){
00063 for(reco::ElectronCollection::const_iterator it = ele_handle->begin(); it != ele_handle->end(); ++it){
00064 if(it->pt() >= 20 && fabs(it->eta()) <= 2.5) nele_tight++;
00065 if(it->pt() >= 10 && fabs(it->eta()) <= 2.5) nele_loose++;
00066 }
00067 }
00068
00069
00070 unsigned int box_num = 5;
00071 if(nmu_tight > 0 && nele_tight > 0){
00072 box_num = 0;
00073 }else if(nmu_tight > 0 && nmu_loose > 1){
00074 box_num = 1;
00075 }else if(nele_tight > 0 && nele_loose > 1){
00076 box_num = 2;
00077 }else if(nmu_tight > 0){
00078 box_num = 3;
00079 }else if(nele_tight > 0){
00080 box_num = 4;
00081 }
00082
00083 edm::Handle<std::vector<double> > razorvar_handle;
00084 iEvent.getByLabel(m_razorVarCollectionTag,razorvar_handle);
00085
00086 if(razorvar_handle->size() > 1){
00087 const double MR = razorvar_handle->at(0);
00088 const double R = razorvar_handle->at(1);
00089 m_rsqMRFullyInc->Fill(MR,R*R);
00090 if(njets >= 4) m_rsqMRInc4J->Fill(MR,R*R);
00091 if(njets >= 6) m_rsqMRInc6J->Fill(MR,R*R);
00092 if(njets >= 8) m_rsqMRInc8J->Fill(MR,R*R);
00093 if(njets >= 10) m_rsqMRInc10J->Fill(MR,R*R);
00094 if(njets >= 12) m_rsqMRInc12J->Fill(MR,R*R);
00095 if(njets >= 14) m_rsqMRInc14J->Fill(MR,R*R);
00096
00097
00098 if(box_num == 0) m_rsqMREleMu->Fill(MR,R*R);
00099 if(box_num == 1) m_rsqMRMuMu->Fill(MR,R*R);
00100 if(box_num == 2) m_rsqMREleEle->Fill(MR,R*R);
00101 if(box_num == 3) m_rsqMRMu->Fill(MR,R*R);
00102 if(box_num == 4) m_rsqMREle->Fill(MR,R*R);
00103 if(box_num == 5) m_rsqMRHad->Fill(MR,R*R);
00104
00105
00106
00107 if( box_num == 3 && njets >= 4) m_rsqMRMuMJ->Fill(MR,R*R);
00108
00109 else if( box_num == 4 && njets >= 5) m_rsqMREleMJ->Fill(MR,R*R);
00110
00111 else if( box_num == 5 && njets >= 6) m_rsqMRHadMJ->Fill(MR,R*R);
00112 }
00113
00114 }
00115
00116 void RazorVarAnalyzer::endRun( edm::Run const &, edm::EventSetup const & ){
00117 }
00118
00119
00120
00121 void RazorVarAnalyzer::bookMEs(){
00122
00123
00124 m_rsqMRFullyInc = bookH2withSumw2("rsqMRFullyInc",
00125 "M_{R} vs R^{2} (All Events)",
00126 400,0.,4000.,
00127 50,0.,1.,
00128 "M_{R} [GeV]",
00129 "R^{2}");
00130 m_rsqMRInc4J = bookH2withSumw2("rsqMRInc4J",
00131 "M_{R} vs R^{2} (>= 4j)",
00132 400,0.,4000.,
00133 50,0.,1.,
00134 "M_{R} [GeV]",
00135 "R^{2}");
00136 m_rsqMRInc6J = bookH2withSumw2("rsqMRInc6J",
00137 "M_{R} vs R^{2} (>= 6j)",
00138 400,0.,4000.,
00139 50,0.,1.,
00140 "M_{R} [GeV]",
00141 "R^{2}");
00142 m_rsqMRInc8J = bookH2withSumw2("rsqMRInc8J",
00143 "M_{R} vs R^{2} (>= 8j)",
00144 400,0.,4000.,
00145 50,0.,1.,
00146 "M_{R} [GeV]",
00147 "R^{2}");
00148 m_rsqMRInc10J = bookH2withSumw2("rsqMRInc10J",
00149 "M_{R} vs R^{2} (>= 10j)",
00150 400,0.,4000.,
00151 50,0.,1.,
00152 "M_{R} [GeV]",
00153 "R^{2}");
00154 m_rsqMRInc12J = bookH2withSumw2("rsqMRInc12J",
00155 "M_{R} vs R^{2} (>= 12j)",
00156 400,0.,4000.,
00157 50,0.,1.,
00158 "M_{R} [GeV]",
00159 "R^{2}");
00160 m_rsqMRInc14J = bookH2withSumw2("rsqMRInc14J",
00161 "M_{R} vs R^{2} (>= 14j)",
00162 400,0.,4000.,
00163 50,0.,1.,
00164 "M_{R} [GeV]",
00165 "R^{2}");
00166
00167
00168 m_rsqMREleMu = bookH2withSumw2("rsqMREleMu",
00169 "M_{R} vs R^{2} (EleMu box)",
00170 400,0.,4000.,
00171 50,0.,1.,
00172 "M_{R} [GeV]",
00173 "R^{2}");
00174 m_rsqMRMuMu = bookH2withSumw2("rsqMRMuMu",
00175 "M_{R} vs R^{2} (MuMu box)",
00176 400,0.,4000.,
00177 50,0.,1.,
00178 "M_{R} [GeV]",
00179 "R^{2}");
00180 m_rsqMREleEle = bookH2withSumw2("rsqMREleEle",
00181 "M_{R} vs R^{2} (EleEle box)",
00182 400,0.,4000.,
00183 50,0.,1.,
00184 "M_{R} [GeV]",
00185 "R^{2}");
00186 m_rsqMRMu = bookH2withSumw2("rsqMRMu",
00187 "M_{R} vs R^{2} (Mu box)",
00188 400,0.,4000.,
00189 50,0.,1.,
00190 "M_{R} [GeV]",
00191 "R^{2}");
00192 m_rsqMREle = bookH2withSumw2("rsqMREle",
00193 "M_{R} vs R^{2} (Ele box)",
00194 400,0.,4000.,
00195 50,0.,1.,
00196 "M_{R} [GeV]",
00197 "R^{2}");
00198 m_rsqMRHad = bookH2withSumw2("rsqMRHad",
00199 "M_{R} vs R^{2} (Had box)",
00200 400,0.,4000.,
00201 50,0.,1.,
00202 "M_{R} [GeV]",
00203 "R^{2}");
00204
00205
00206 m_rsqMRMuMJ = bookH2withSumw2("rsqMRMuMJ",
00207 "M_{R} vs R^{2} (Mu box >= 4j)",
00208 400,0.,4000.,
00209 50,0.,1.,
00210 "M_{R} [GeV]",
00211 "R^{2}");
00212 m_rsqMREleMJ = bookH2withSumw2("rsqMREleMJ",
00213 "M_{R} vs R^{2} (Ele box >= 5j)",
00214 400,0.,4000.,
00215 50,0.,1.,
00216 "M_{R} [GeV]",
00217 "R^{2}");
00218 m_rsqMRHadMJ = bookH2withSumw2("rsqMRHadMJ",
00219 "M_{R} vs R^{2} (Had box >= 6j)",
00220 400,0.,4000.,
00221 50,0.,1.,
00222 "M_{R} [GeV]",
00223 "R^{2}");
00224
00225
00226 }
00227