00001
00002
00003
00004
00005
00006
00007 #include <memory>
00008
00009
00010
00011 #include "FWCore/Framework/interface/ESHandle.h"
00012
00013 #include "FWCore/Framework/interface/Frameworkfwd.h"
00014 #include "FWCore/Framework/interface/EDAnalyzer.h"
00015
00016 #include "FWCore/Framework/interface/Event.h"
00017 #include "FWCore/Framework/interface/MakerMacros.h"
00018
00019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00020
00021 #include "DataFormats/HLTReco/interface/TriggerEvent.h"
00022 #include "DataFormats/L1Trigger/interface/L1JetParticle.h"
00023 #include "DataFormats/L1Trigger/interface/L1JetParticleFwd.h"
00024
00025 #include "CondFormats/L1TObjects/interface/L1GtTriggerMenu.h"
00026 #include "CondFormats/DataRecord/interface/L1GtTriggerMenuRcd.h"
00027 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutSetupFwd.h"
00028 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutSetup.h"
00029 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutRecord.h"
00030
00031 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00032 #include "DataFormats/TrackReco/interface/Track.h"
00033
00034 #include "CondFormats/L1TObjects/interface/L1GtPrescaleFactors.h"
00035 #include "CondFormats/DataRecord/interface/L1GtPrescaleFactorsAlgoTrigRcd.h"
00036 #include "CondFormats/DataRecord/interface/L1GtPrescaleFactorsTechTrigRcd.h"
00037
00038 #include "DQMServices/Core/interface/DQMStore.h"
00039 #include "DQMServices/Core/interface/MonitorElement.h"
00040 #include "FWCore/ServiceRegistry/interface/Service.h"
00041
00042 #include "DataFormats/JetReco/interface/GenJetCollection.h"
00043
00044 #include "DataFormats/HcalIsolatedTrack/interface/IsolatedPixelTrackCandidate.h"
00045 #include "DataFormats/HcalIsolatedTrack/interface/IsolatedPixelTrackCandidateFwd.h"
00046
00047 #include "TH1F.h"
00048 #include "TH2F.h"
00049
00050 #include <fstream>
00051
00052 class ValHcalIsoTrackHLT : public edm::EDAnalyzer {
00053 public:
00054 explicit ValHcalIsoTrackHLT(const edm::ParameterSet&);
00055 ~ValHcalIsoTrackHLT();
00056 double getDist(double,double,double,double);
00057
00058 private:
00059
00060 int evtBuf;
00061
00062 DQMStore* dbe_;
00063
00064 virtual void beginJob() ;
00065 virtual void analyze(const edm::Event&, const edm::EventSetup&);
00066 virtual void endJob() ;
00067
00068 bool produceRates_;
00069 double sampleXsec_;
00070 double lumi_;
00071 std::string outTxtFileName_;
00072 std::string folderName_;
00073 bool saveToRootFile_;
00074 std::string outRootFileName_;
00075 edm::InputTag hltEventTag_;
00076 edm::InputTag hltFilterTag_;
00077 std::vector<edm::InputTag> l1extraJetTag_;
00078 std::vector<std::string> l1seedNames_;
00079 edm::InputTag gtDigiLabel_;
00080 bool useReco_;
00081 edm::InputTag recoTracksLabel_;
00082 bool checkL2_;
00083 edm::InputTag l2colLabel_;
00084
00085 bool checkL1eff_;
00086 edm::InputTag genJetsLabel_;
00087
00088 bool produceRatePdep_;
00089
00090 bool doL1Prescaling_;
00091
00092 MonitorElement* hl3Pt;
00093 MonitorElement* hL3L2trackMatch;
00094 MonitorElement* hL3L2pTrat;
00095 MonitorElement* hl3Pt0005;
00096 MonitorElement* hl3Pt0510;
00097 MonitorElement* hl3Pt1015;
00098 MonitorElement* hl3Pt1520;
00099 MonitorElement* hl3P0005;
00100 MonitorElement* hl3P0510;
00101 MonitorElement* hl3P1015;
00102 MonitorElement* hl3P1520;
00103 MonitorElement* hl3eta;
00104 MonitorElement* hl3phi;
00105
00106 MonitorElement* hl3pVsEta;
00107
00108 MonitorElement* hOffL3TrackMatch;
00109 MonitorElement* hOffL3TrackPtRat;
00110
00111 MonitorElement* hl2eta;
00112 MonitorElement* hl2phi;
00113 MonitorElement* hl2pT;
00114 MonitorElement* hl2pVsEta;
00115 MonitorElement* hisopT;
00116 MonitorElement* hisopTvsEta;
00117
00118 MonitorElement* haccepts;
00119 MonitorElement* hOffPvsEta;
00120 MonitorElement* hpTgenLead;
00121 MonitorElement* hpTgenLeadL1;
00122 MonitorElement* hpTgenNext;
00123 MonitorElement* hpTgenNextL1;
00124 MonitorElement* hLeadTurnOn;
00125 MonitorElement* hNextToLeadTurnOn;
00126
00127 MonitorElement* hRateVsThr;
00128
00129 MonitorElement* hPhiToGJ;
00130 MonitorElement* hDistToGJ;
00131
00132 std::vector<int> l1counter;
00133
00134 std::ofstream txtout;
00135
00136 int nTotal;
00137 int nL1accepts;
00138 int nHLTL2accepts;
00139 int nHLTL3accepts;
00140 int nHLTL3acceptsPure;
00141
00142 int nl3_0005;
00143 int nl3_0510;
00144 int nl3_1015;
00145 int nl3_1520;
00146
00147 int purnl3_0005;
00148 int purnl3_0510;
00149 int purnl3_1015;
00150 int purnl3_1520;
00151
00152 double hltPThr_;
00153 };
00154
00155 double ValHcalIsoTrackHLT::getDist(double eta1, double phi1, double eta2, double phi2)
00156 {
00157 double dphi = fabs(phi1 - phi2);
00158 if(dphi>acos(-1)) dphi = 2*acos(-1)-dphi;
00159 double dr = sqrt(dphi*dphi + pow(eta1-eta2,2));
00160 return dr;
00161 }
00162
00163 ValHcalIsoTrackHLT::ValHcalIsoTrackHLT(const edm::ParameterSet& iConfig)
00164
00165 {
00166 produceRates_=iConfig.getParameter<bool>("produceRates");
00167 sampleXsec_=iConfig.getParameter<double>("sampleCrossSection");
00168 lumi_=iConfig.getParameter<double>("luminosity");
00169 outTxtFileName_=iConfig.getParameter<std::string>("outputTxtFileName");
00170
00171 folderName_ = iConfig.getParameter<std::string>("folderName");
00172 saveToRootFile_=iConfig.getParameter<bool>("SaveToRootFile");
00173 outRootFileName_=iConfig.getParameter<std::string>("outputRootFileName");
00174
00175 hltEventTag_=iConfig.getParameter<edm::InputTag>("hltTriggerEventLabel");
00176 hltFilterTag_=iConfig.getParameter<edm::InputTag>("hltL3FilterLabel");
00177
00178 l1extraJetTag_=iConfig.getParameter<std::vector<edm::InputTag> >("hltL1extraJetLabel");
00179 gtDigiLabel_=iConfig.getParameter<edm::InputTag>("gtDigiLabel");
00180 checkL1eff_=iConfig.getParameter<bool>("CheckL1TurnOn");
00181 genJetsLabel_=iConfig.getParameter<edm::InputTag>("genJetsLabel");
00182 l1seedNames_=iConfig.getParameter<std::vector<std::string> >("l1seedNames");
00183 useReco_=iConfig.getParameter<bool>("useReco");
00184 recoTracksLabel_=iConfig.getParameter<edm::InputTag>("recoTracksLabel");
00185 checkL2_=iConfig.getParameter<bool>("DebugL2");
00186 l2colLabel_=iConfig.getParameter<edm::InputTag>("L2producerLabel");
00187
00188 produceRatePdep_=iConfig.getParameter<bool>("produceRatePdep");
00189
00190 hltPThr_=iConfig.getUntrackedParameter<double>("l3momThreshold",10);
00191
00192 doL1Prescaling_=iConfig.getParameter<bool>("doL1Prescaling");
00193
00194 for (unsigned int i=0; i<l1seedNames_.size(); i++)
00195 {
00196 l1counter.push_back(0);
00197 }
00198
00199 if (produceRates_) txtout.open(outTxtFileName_.c_str());
00200
00201 nTotal=0;
00202 nL1accepts=0;
00203 nHLTL2accepts=0;
00204 nHLTL3accepts=0;
00205 nHLTL3acceptsPure=0;
00206
00207 nl3_0005=0;
00208 nl3_0510=0;
00209 nl3_1015=0;
00210 nl3_1520=0;
00211
00212 purnl3_0005=0;
00213 purnl3_0510=0;
00214 purnl3_1015=0;
00215 purnl3_1520=0;
00216 }
00217
00218
00219 ValHcalIsoTrackHLT::~ValHcalIsoTrackHLT()
00220 {
00221 if (produceRates_)
00222 {
00223 double sampleRate=(lumi_)*(sampleXsec_*1E-36);
00224 double l1Rate=nL1accepts*pow(nTotal,-1)*sampleRate;
00225 double hltRate=nHLTL3accepts*pow(nL1accepts,-1)*l1Rate;
00226 double hltRatePure=nHLTL3acceptsPure*pow(nL1accepts,-1)*l1Rate;
00227
00228 double l1rateError=l1Rate/sqrt(nL1accepts);
00229 double hltRateError=hltRate/sqrt(nHLTL3accepts);
00230 double hltRatePureError=hltRatePure/sqrt(nHLTL3acceptsPure);
00231
00232 double rate_0005=nl3_0005*pow(nL1accepts,-1)*l1Rate;
00233 double rate_0510=nl3_0510*pow(nL1accepts,-1)*l1Rate;
00234 double rate_1015=nl3_1015*pow(nL1accepts,-1)*l1Rate;
00235 double rate_1520=nl3_1520*pow(nL1accepts,-1)*l1Rate;
00236
00237 double prate_0005=purnl3_0005*pow(nL1accepts,-1)*l1Rate;
00238 double prate_0510=purnl3_0510*pow(nL1accepts,-1)*l1Rate;
00239 double prate_1015=purnl3_1015*pow(nL1accepts,-1)*l1Rate;
00240 double prate_1520=purnl3_1520*pow(nL1accepts,-1)*l1Rate;
00241
00242 txtout<<std::setw(50)<<std::left<<"sample xsec(pb)"<<sampleXsec_<<std::endl;
00243 txtout<<std::setw(50)<<std::left<<"lumi(cm^-2*s^-1)"<<lumi_<<std::endl;
00244 txtout<<std::setw(50)<<std::left<<"Events processed/rate(Hz)"<<nTotal<<"/"<<sampleRate<<std::endl;
00245 txtout<<std::setw(50)<<std::left<<"L1 accepts/(rate+-error (Hz))"<<nL1accepts<<"/("<<l1Rate<<"+-"<<l1rateError<<")"<<std::endl;
00246 txtout<<std::setw(50)<<std::left<<"HLTL3accepts/(rate+-error (Hz))"<<nHLTL3accepts<<"/("<<hltRate<<"+-"<<hltRateError<<")"<<std::endl;
00247 txtout<<std::setw(50)<<std::left<<"L3 acc. |eta|<0.5 / rate"<<nl3_0005<<" / "<<rate_0005<<std::endl;
00248 txtout<<std::setw(50)<<std::left<<"L3 acc. |eta|>0.5 && |eta|<1.0 / rate"<<nl3_0510<<" / "<<rate_0510<<std::endl;
00249 txtout<<std::setw(50)<<std::left<<"L3 acc. |eta|>1.0 && |eta|<1.5 / rate"<<nl3_1015<<" / "<<rate_1015<<std::endl;
00250 txtout<<std::setw(50)<<std::left<<"L3 acc. |eta|>1.5 && |eta|<2.0 / rate"<<nl3_1520<<" / "<<rate_1520<<std::endl;
00251 txtout<<"\n"<<std::endl;
00252 txtout<<std::setw(50)<<std::left<<"HLTL3acceptsPure/(rate+-error (Hz))"<<nHLTL3acceptsPure<<"/("<<hltRatePure<<"+-"<<hltRatePureError<<")"<<std::endl;
00253 txtout<<std::setw(50)<<std::left<<"pure L3 acc. |eta|<0.5 / rate"<<purnl3_0005<<" / "<<prate_0005<<std::endl;
00254 txtout<<std::setw(50)<<std::left<<"pure L3 acc. |eta|>0.5 && |eta|<1.0 / rate"<<purnl3_0510<<" / "<<prate_0510<<std::endl;
00255 txtout<<std::setw(50)<<std::left<<"pure L3 acc. |eta|>1.0 && |eta|<1.5 / rate"<<purnl3_1015<<" / "<<prate_1015<<std::endl;
00256 txtout<<std::setw(50)<<std::left<<"pure L3 acc. |eta|>1.5 && |eta|<2.0 / rate"<<purnl3_1520<<" / "<<prate_1520<<std::endl;
00257 }
00258 }
00259
00260 void ValHcalIsoTrackHLT::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00261 {
00262 nTotal++;
00263 haccepts->Fill(0.0001,1);
00264
00265 double phiGJLead=-10000;
00266 double etaGJLead=-10000;
00267
00268 bool l1pass=false;
00269
00270 edm::ESHandle<L1GtTriggerMenu> menuRcd;
00271 iSetup.get<L1GtTriggerMenuRcd>().get(menuRcd) ;
00272 const L1GtTriggerMenu* menu = menuRcd.product();
00273 const AlgorithmMap& bitMap = menu->gtAlgorithmMap();
00274
00275 edm::ESHandle<L1GtPrescaleFactors> l1GtPfAlgo;
00276 iSetup.get<L1GtPrescaleFactorsAlgoTrigRcd>().get(l1GtPfAlgo);
00277 const L1GtPrescaleFactors* preFac = l1GtPfAlgo.product();
00278 const std::vector< std::vector< int > > prescaleSet=preFac->gtPrescaleFactors();
00279
00280 edm::Handle< L1GlobalTriggerReadoutRecord > gtRecord;
00281 iEvent.getByLabel(gtDigiLabel_ ,gtRecord);
00282 const DecisionWord dWord = gtRecord->decisionWord();
00283
00284 for (unsigned int i=0; i<l1seedNames_.size(); i++)
00285 {
00286 int l1seedBitNumber=-10;
00287 for (CItAlgo itAlgo = bitMap.begin(); itAlgo != bitMap.end(); itAlgo++)
00288 {
00289 if (itAlgo->first==l1seedNames_[i]) l1seedBitNumber = (itAlgo->second).algoBitNumber();
00290 }
00291 int prescale=prescaleSet[0][l1seedBitNumber];
00292
00293 if (menu->gtAlgorithmResult( l1seedNames_[i], dWord))
00294 {
00295 l1counter[i]++;
00296 if ((doL1Prescaling_&&l1counter[i]%prescale==0)||(!doL1Prescaling_))
00297 {
00298 l1pass=true;
00299 break;
00300 }
00301 }
00302 }
00303 if (checkL1eff_)
00304 {
00305 edm::Handle<reco::GenJetCollection> gjcol;
00306 iEvent.getByLabel(genJetsLabel_,gjcol);
00307
00308 reco::GenJetCollection::const_iterator gjit=gjcol->begin();
00309 if (gjit!=gjcol->end())
00310 {
00311 etaGJLead=gjit->eta();
00312 phiGJLead=gjit->phi();
00313 hpTgenLead->Fill(gjit->pt(),1);
00314 if (l1pass) hpTgenLeadL1->Fill(gjit->pt(),1);
00315 gjit++;
00316 }
00317 if (gjit!=gjcol->end())
00318 {
00319 hpTgenNext->Fill(gjit->pt(),1);
00320 if (l1pass) hpTgenNextL1->Fill(gjit->pt(),1);
00321 }
00322 }
00323
00324 if (!l1pass) return;
00325 else haccepts->Fill(1+0.0001,1);
00326
00327 nL1accepts++;
00328
00329 edm::Handle<reco::IsolatedPixelTrackCandidateCollection> l2col;
00330 if (checkL2_)
00331 {
00332 iEvent.getByLabel(l2colLabel_,l2col);
00333 }
00334
00335 edm::Handle<trigger::TriggerEvent> trEv;
00336 iEvent.getByLabel(hltEventTag_,trEv);
00337
00338 edm::Handle<reco::TrackCollection> recoTr;
00339
00340 if (useReco_)
00341 {
00342 iEvent.getByLabel(recoTracksLabel_,recoTr);
00343 }
00344
00345 const trigger::TriggerObjectCollection& TOCol(trEv->getObjects());
00346
00347 trigger::Keys KEYS;
00348 const trigger::size_type nFilt(trEv->sizeFilters());
00349
00350 int nFired=0;
00351
00352 bool passl3=false;
00353
00354 for (trigger::size_type iFilt=0; iFilt!=nFilt; iFilt++)
00355 {
00356 trigger::Keys KEYS1=trEv->filterKeys(iFilt);
00357 if (KEYS1.size()>0) nFired++;
00358 if (trEv->filterTag(iFilt)==hltFilterTag_) KEYS=trEv->filterKeys(iFilt);
00359 }
00360
00361 trigger::size_type nReg=KEYS.size();
00362
00363 if (nFired==2&&nReg>0)
00364 {
00365 nHLTL3acceptsPure++;
00366 for (trigger::size_type iReg=0; iReg<nReg; iReg++)
00367 {
00368 const trigger::TriggerObject& TObj(TOCol[KEYS[iReg]]);
00369 if (TObj.pt()*cosh(TObj.eta())<10) continue;
00370 if (fabs(TObj.eta())<0.5) purnl3_0005++;
00371 if (fabs(TObj.eta())>0.5&&fabs(TObj.eta())<1.0) purnl3_0510++;
00372 if (fabs(TObj.eta())>1.0&&fabs(TObj.eta())<1.5) purnl3_1015++;
00373 if (fabs(TObj.eta())<2.0&&fabs(TObj.eta())>1.5) purnl3_1520++;
00374 }
00375 }
00376
00377 for (trigger::size_type iReg=0; iReg<nReg; iReg++)
00378 {
00379 const trigger::TriggerObject& TObj(TOCol[KEYS[iReg]]);
00380
00381 if (produceRatePdep_){
00382 for (int i=0; i<50; i++)
00383 {
00384 double pthr=5+i;
00385 if (TObj.pt()*cosh(TObj.eta())>pthr&&fabs(TObj.eta())<1.0) hRateVsThr->Fill(pthr+0.001,1);
00386 }
00387 }
00388
00389 if (TObj.pt()*cosh(TObj.eta())<hltPThr_) continue;
00390
00391 passl3=true;
00392
00393 double dphiGJ=fabs(TObj.phi()-phiGJLead);
00394 if (dphiGJ>acos(-1)) dphiGJ=2*acos(-1)-dphiGJ;
00395 double dR=sqrt(dphiGJ*dphiGJ+pow(TObj.eta()-etaGJLead,2));
00396 hPhiToGJ->Fill(dphiGJ,1);
00397 hDistToGJ->Fill(dR,1);
00398
00399 hl3eta->Fill(TObj.eta(),1);
00400 hl3phi->Fill(TObj.phi(),1);
00401 if (fabs(TObj.eta())<0.5)
00402 {
00403 hl3P0005->Fill(cosh(TObj.eta())*TObj.pt(),1);
00404 hl3Pt0005->Fill(TObj.pt(),1);
00405 nl3_0005++;
00406 }
00407 if (fabs(TObj.eta())>0.5&&fabs(TObj.eta())<1.0)
00408 {
00409 nl3_0510++;
00410 hl3P0510->Fill(cosh(TObj.eta())*TObj.pt(),1);
00411 hl3Pt0510->Fill(TObj.pt(),1);
00412 }
00413 if (fabs(TObj.eta())>1.0&&fabs(TObj.eta())<1.5)
00414 {
00415 nl3_1015++;
00416 hl3P1015->Fill(cosh(TObj.eta())*TObj.pt(),1);
00417 hl3Pt1015->Fill(TObj.pt(),1);
00418 }
00419 if (fabs(TObj.eta())<2.0&&fabs(TObj.eta())>1.5)
00420 {
00421 nl3_1520++;
00422 hl3P1520->Fill(cosh(TObj.eta())*TObj.pt(),1);
00423 hl3Pt1520->Fill(TObj.pt(),1);
00424 }
00425 if (l2col->size()==0) continue;
00426 double l2l3d=100;
00427 reco::IsolatedPixelTrackCandidateCollection::const_iterator l2match;
00428 for (reco::IsolatedPixelTrackCandidateCollection::const_iterator eptit=l2col->begin(); eptit!=l2col->end(); eptit++)
00429 {
00430 hl2pT->Fill(eptit->track()->pt(),1);
00431 hl2eta->Fill(eptit->track()->eta(),1);
00432 hl2phi->Fill(eptit->track()->phi(),1);
00433 hl2pVsEta->Fill(eptit->track()->eta(),eptit->track()->p(),1);
00434 hisopT->Fill(eptit->maxPtPxl(),1);
00435 hisopTvsEta->Fill(eptit->track()->eta(),eptit->maxPtPxl(),1);
00436
00437 double R=getDist(eptit->eta(), eptit->phi(), TObj.eta(), TObj.phi());
00438 if (R<l2l3d)
00439 {
00440 l2match=eptit;
00441 l2l3d=R;
00442 }
00443 }
00444
00445 hL3L2trackMatch->Fill(l2l3d,1);
00446 hL3L2pTrat->Fill(l2match->pt()/TObj.pt(),1);
00447 hl3Pt->Fill(TObj.pt(),1);
00448
00449 if (useReco_&&recoTr->size()>0)
00450 {
00451 double minRecoL3dist=100;
00452 reco::TrackCollection::const_iterator mrtr;
00453 for (reco::TrackCollection::const_iterator rtrit=recoTr->begin(); rtrit!=recoTr->end(); rtrit++)
00454 {
00455 double R=getDist(rtrit->eta(),rtrit->phi(),TObj.eta(),TObj.phi());
00456 if (R<minRecoL3dist)
00457 {
00458 mrtr=rtrit;
00459 minRecoL3dist=R;
00460 }
00461 }
00462 hOffL3TrackMatch->Fill(minRecoL3dist,1);
00463 hOffL3TrackPtRat->Fill(TObj.pt()/mrtr->pt(),1);
00464 hOffPvsEta->Fill(mrtr->eta(),mrtr->p(),1);
00465 }
00466 }
00467
00468 if (passl3) nHLTL3accepts++;
00469 if (passl3) haccepts->Fill(2+0.0001,1);
00470
00471 if (!l1pass||!passl3) return;
00472
00473 }
00474
00475 void ValHcalIsoTrackHLT::beginJob()
00476 {
00477 dbe_ = edm::Service<DQMStore>().operator->();
00478 dbe_->setCurrentFolder(folderName_);
00479
00480 hRateVsThr=dbe_->book1D("hRatesVsThr","hRateVsThr",100,0,100);
00481
00482 hPhiToGJ=dbe_->book1D("hPhiToGJ","hPhiToGJ",100,0,4);
00483
00484 hDistToGJ=dbe_->book1D("hDistToGJ","hDistToGJ",100,0,10);
00485
00486 hL3L2trackMatch=dbe_->book1D("hL3L2trackMatch","R from L3 object to L2 object ",1000,0,1);
00487 hL3L2trackMatch->setAxisTitle("R(eta,phi)",1);
00488
00489 hL3L2pTrat=dbe_->book1D("hL3L2pTrat","ratio of L2 to L3 measurement",1000,0,10);
00490 hL3L2pTrat->setAxisTitle("pT_L2/pT_L3",1);
00491
00492 hl3Pt=dbe_->book1D("hl3Pt","pT of L3 objects",1000,0,100);
00493 hl3Pt->setAxisTitle("pT(GeV)",1);
00494
00495 hl3Pt0005=dbe_->book1D("hl3Pt0005","hl3Pt0005",1000,0,100);
00496 hl3Pt0005->setAxisTitle("pT(GeV)",1);
00497
00498 hl3Pt0510=dbe_->book1D("hl3Pt0510","hl3Pt0510",1000,0,100);
00499 hl3Pt0510->setAxisTitle("pT(GeV)",1);
00500
00501 hl3Pt1015=dbe_->book1D("hl3Pt1015","hl3Pt1015",1000,0,100);
00502 hl3Pt1015->setAxisTitle("pT(GeV)",1);
00503
00504 hl3Pt1520=dbe_->book1D("hl3Pt1520","hl3Pt1520",1000,0,100);
00505 hl3Pt1520->setAxisTitle("pT(GeV)",1);
00506
00507 hl3P0005=dbe_->book1D("hl3P0005","hl3P0005",1000,0,100);
00508 hl3P0005->setAxisTitle("P(GeV)",1);
00509
00510 hl3P0510=dbe_->book1D("hl3P0510","hl3P0510",1000,0,100);
00511 hl3P0510->setAxisTitle("P(GeV)",1);
00512
00513 hl3P1015=dbe_->book1D("hl3P1015","hl3P1015",1000,0,100);
00514 hl3P1015->setAxisTitle("P(GeV)",1);
00515
00516 hl3P1520=dbe_->book1D("hl3P1520","hl3P1520",1000,0,100);
00517 hl3P1520->setAxisTitle("P(GeV)",1);
00518
00519 hl3eta=dbe_->book1D("hl3eta","eta of L3 objects",50,-2.5,2.5);
00520 hl3eta->setAxisTitle("eta",1);
00521
00522 hl3phi=dbe_->book1D("hl3phi","phi of L3 objects",70,-3.5,3.5);
00523 hl3phi->setAxisTitle("phi(rad)",1);
00524
00525 hl2pT=dbe_->book1D("hl2pT","pT of L2 objects",1000,0,1000);
00526 hl2pT->setAxisTitle("pT(GeV)",1);
00527
00528 hl2eta=dbe_->book1D("hl2eta","eta of L2 objects",50,-2.5,2.5);
00529 hl2eta->setAxisTitle("eta",1);
00530
00531 hl2phi=dbe_->book1D("hl2phi","phi of L2 objects",70,-3.5,3.5);
00532 hl2phi->setAxisTitle("phi(rad)",1);
00533
00534 hisopT=dbe_->book1D("hisopT","isolation pT",100,0,5.5);
00535 hisopT->setAxisTitle("iso pT (GeV)",1);
00536
00537 hisopTvsEta=dbe_->book2D("hisopTvsEta","isolation pT vs Eta",8,-2,2,100,0,5.5);
00538 hisopTvsEta->setAxisTitle("eta",1);
00539 hisopTvsEta->setAxisTitle("iso pT (GeV)",2);
00540
00541 hl2pVsEta=dbe_->book2D("hl2pVsEta","Distribution of l2 track energy vs eta",25,-2.5,2.5,100,0,100);
00542 hl2pVsEta->setAxisTitle("eta",1);
00543 hl2pVsEta->setAxisTitle("E(GeV)",2);
00544
00545 haccepts=dbe_->book1D("haccepts","Number of accepts at each level",3,0,3);
00546 haccepts->setAxisTitle("selection level",1);
00547
00548 hOffL3TrackMatch=dbe_->book1D("hOffL3TrackMatch","Distance from L3 object to offline track",200,0,0.5);
00549 hOffL3TrackMatch->setAxisTitle("R(eta,phi)",1);
00550
00551 hOffL3TrackPtRat=dbe_->book1D("hOffL3TrackPtRat","Ratio of pT: L3/offline",100,0,10);
00552 hOffL3TrackPtRat->setAxisTitle("ratio L3/offline",1);
00553
00554 hOffPvsEta=dbe_->book2D("hOffPvsEta","Distribution of offline track energy vs eta",25,-2.5,2.5,100,0,100);
00555 hOffPvsEta->setAxisTitle("eta",1);
00556 hOffPvsEta->setAxisTitle("E(GeV)",2);
00557
00558 hpTgenLead=dbe_->book1D("hpTgenLead","hpTgenLead",100,0,100);
00559
00560 hpTgenLeadL1=dbe_->book1D("hpTgenLeadL1","hpTgenLeadL1",100,0,100);
00561
00562 hpTgenNext=dbe_->book1D("hpTgenNext","hpTgenNext",100,0,100);
00563
00564 hpTgenNextL1=dbe_->book1D("hpTgenNextL1","hpTgenNextL1",100,0,100);
00565 }
00566
00567 void ValHcalIsoTrackHLT::endJob() {
00568
00569 if(dbe_)
00570 {
00571 if (saveToRootFile_) dbe_->save(outRootFileName_);
00572 }
00573 }
00574
00575 DEFINE_FWK_MODULE(ValHcalIsoTrackHLT);