#include <HLTriggerOffline/Tau/interface/HLTMuonRate.h>
Definition at line 30 of file HLTMuonRate.h.
HLTMuonRate::HLTMuonRate | ( | const edm::ParameterSet & | pset, | |
int | Index | |||
) |
Definition at line 33 of file HLTMuonRate.cc.
References DQMStore::cd(), dbe_, e, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), i, InputLabel, NULL, DQMStore::setCurrentFolder(), DQMStore::setVerbose(), theCrossSection, theHLTCollectionLabels, theHLTReferenceThreshold, theL1CollectionLabel, theL1ReferenceThreshold, theLuminosity, theNbins, theNSigmas, theNumberOfEvents, theNumberOfL1Events, theNumberOfObjects, thePtMax, thePtMin, and theRootFileName.
00034 { 00035 InputLabel = pset.getUntrackedParameter<InputTag>("InputLabel"); 00036 Parameters TriggerLists=pset.getParameter<Parameters>("TriggerCollection"); 00037 int i=0; 00038 for(Parameters::iterator itTrigger = TriggerLists.begin(); itTrigger != TriggerLists.end(); ++itTrigger) { 00039 if ( i == Index ) { 00040 theL1CollectionLabel = itTrigger->getParameter<InputTag>("L1CollectionLabel"); 00041 theHLTCollectionLabels = itTrigger->getParameter<std::vector<InputTag> >("HLTCollectionLabels"); 00042 theL1ReferenceThreshold = itTrigger->getParameter<double>("L1ReferenceThreshold"); 00043 theHLTReferenceThreshold = itTrigger->getParameter<double>("HLTReferenceThreshold"); 00044 theNumberOfObjects = itTrigger->getParameter<unsigned int>("NumberOfObjects"); 00045 break; 00046 } 00047 ++i; 00048 } 00049 theNSigmas = pset.getUntrackedParameter<std::vector<double> >("NSigmas90"); 00050 theCrossSection = pset.getParameter<double>("CrossSection"); 00051 // Convert it already into /nb/s) 00052 theLuminosity = pset.getUntrackedParameter<double>("Luminosity",1.e+32)*1.e-33; 00053 thePtMin = pset.getUntrackedParameter<double>("PtMin",0.); 00054 thePtMax = pset.getUntrackedParameter<double>("PtMax",40.); 00055 theNbins = pset.getUntrackedParameter<unsigned int>("Nbins",40); 00056 theNumberOfEvents = 0; 00057 theNumberOfL1Events = 0; 00058 dbe_ = 0 ; 00059 if (pset.getUntrackedParameter < bool > ("DQMStore", false)) { 00060 dbe_ = Service < DQMStore > ().operator->(); 00061 dbe_->setVerbose(0); 00062 } 00063 bool disable =pset.getUntrackedParameter < bool > ("disableROOToutput", false); 00064 if (disable) theRootFileName=""; 00065 else theRootFileName = pset.getUntrackedParameter<std::string>("RootFileName"); 00066 if (dbe_ != NULL) { 00067 dbe_->cd(); 00068 dbe_->setCurrentFolder("HLT/Muon"); 00069 dbe_->setCurrentFolder("HLT/Muon/RateEfficiencies"); 00070 dbe_->setCurrentFolder("HLT/Muon/Distributions"); 00071 } 00072 }
HLTMuonRate::~HLTMuonRate | ( | ) | [virtual] |
void HLTMuonRate::analyze | ( | const edm::Event & | event | ) |
Definition at line 75 of file HLTMuonRate.cc.
References d, geometryDiff::epsilon, err0, MonitorElement::Fill(), edm::Ref< C, T, F >::get(), getAngle(), hEtaNor, hHLTeff, hHLTpt, hL1DR, hL1eff, hL1Eta, hL1Phi, hL1pt, hL2DR, hL3DR, hPhiNor, hPtNor, hSteps, i, InputLabel, j, k, LogTrace, NumberOfEvents, funct::sqrt(), theHLTCollectionLabels, theHLTReferenceThreshold, theL1CollectionLabel, theL1ReferenceThreshold, theNbins, theNSigmas, theNumberOfEvents, theNumberOfObjects, thePtMax, thePtMin, and this_event_weight.
00075 { 00076 this_event_weight=1; 00077 NumberOfEvents->Fill(++theNumberOfEvents); 00078 LogTrace("HLTMuonVal")<<"In analyze for L1 trigger "<<theL1CollectionLabel<<" Event:"<<theNumberOfEvents; 00079 bool refmuon_found = false; 00080 double ptuse = -1; 00081 Handle<vector<XYZTLorentzVectorD> > refVector; 00082 event.getByLabel(InputLabel,refVector); 00083 for(unsigned int i = 0; i < refVector->size(); i++){ 00084 double pt1 = refVector->at(i).pt(); 00085 hEtaNor->Fill(refVector->at(i).eta()); 00086 hPhiNor->Fill(refVector->at(i).phi()); 00087 if (pt1>ptuse) { 00088 refmuon_found = true; 00089 ptuse = pt1; 00090 } 00091 } 00092 if (ptuse > 0 ) hPtNor->Fill(ptuse,this_event_weight); 00093 // Get the L1 collection 00094 Handle<TriggerFilterObjectWithRefs> l1candsHandle; 00095 event.getByLabel(theL1CollectionLabel, l1candsHandle); 00096 // Get the HLT collections 00097 std::vector<Handle<TriggerFilterObjectWithRefs> > hltcandsHandle(theHLTCollectionLabels.size()); 00098 for (unsigned int i=0; i<theHLTCollectionLabels.size(); i++) { 00099 event.getByLabel(theHLTCollectionLabels[i], hltcandsHandle[i]); 00100 if (hltcandsHandle[i].failedToGet()) break; 00101 } 00102 vector<L1MuonParticleRef> l1cands; 00103 l1candsHandle->getObjects(81,l1cands); 00104 unsigned hltsize=theHLTCollectionLabels.size(); 00105 vector<vector<RecoChargedCandidateRef> > hltcands(hltsize); 00106 unsigned int modules_in_this_event = 0; 00107 for (unsigned int i=0; i<theHLTCollectionLabels.size(); i++) { 00108 hltcandsHandle[i]->getObjects(93, hltcands[i]); 00109 modules_in_this_event++; 00110 } 00111 // Fix L1 thresholds to obtain HLT plots 00112 unsigned int nL1FoundRef = 0; 00113 double epsilon = 0.001; 00114 for (unsigned int k=0; k<l1cands.size(); k++) { 00115 double ptLUT = l1cands.at(k)->pt(); 00116 // Add "epsilon" to avoid rounding errors when ptLUT==L1Threshold 00117 if (ptLUT+epsilon>theL1ReferenceThreshold) { 00118 nL1FoundRef++; 00119 hL1pt->Fill(ptLUT); 00120 pair<double,double> angularInfo=getAngle(l1cands.at(k)->eta(),l1cands.at(k)->phi(), refVector ); 00121 LogTrace("HLTMuonVal")<<"Filling L1 histos...."; 00122 if ( angularInfo.first < 999.){ 00123 hL1Eta->Fill(angularInfo.first); 00124 hL1Phi->Fill(angularInfo.second); 00125 float dphi=angularInfo.second-l1cands.at(k)->phi(); 00126 float deta=angularInfo.first-l1cands.at(k)->eta(); 00127 hL1DR->Fill(sqrt(dphi*dphi+deta*deta)); 00128 LogTrace("HLTMuonVal")<<"Filling done"; 00129 } 00130 } 00131 } 00132 if (nL1FoundRef>=theNumberOfObjects){ 00133 hSteps->Fill(1.); 00134 } 00135 if (ptuse>0){ 00136 int last_module = modules_in_this_event - 1; 00137 for ( int i=0; i<=last_module; i++) { 00138 double ptcut = theHLTReferenceThreshold; 00139 unsigned nFound = 0; 00140 for (unsigned int k=0; k<hltcands[i].size(); k++) { 00141 RecoChargedCandidateRef candref = RecoChargedCandidateRef(hltcands[i][k]); 00142 reco::TrackRef tk = candref->get<reco::TrackRef>(); 00143 double pt = tk->pt(); 00144 if (pt>ptcut) nFound++; 00145 } 00146 if (nFound>=theNumberOfObjects){ 00147 hSteps->Fill(2+i); 00148 } 00149 } 00150 } 00151 for (unsigned int j=0; j<theNbins; j++) { 00152 double ptcut = thePtMin + j*(thePtMax-thePtMin)/theNbins; 00153 // L1 filling 00154 unsigned int nFound = 0; 00155 for (unsigned int k=0; k<l1cands.size(); k++) { 00156 L1MuonParticleRef candref = L1MuonParticleRef(l1cands[k]); 00157 double pt = candref->pt(); 00158 if (pt>ptcut) nFound++; 00159 } 00160 if (nFound>=theNumberOfObjects) hL1eff->Fill(ptcut,this_event_weight); 00161 // Stop here if L1 reference cuts were not satisfied 00162 if (nL1FoundRef<theNumberOfObjects) continue; 00163 // HLT filling 00164 for (unsigned int i=0; i<modules_in_this_event; i++) { 00165 unsigned nFound = 0; 00166 for (unsigned int k=0; k<hltcands[i].size(); k++) { 00167 RecoChargedCandidateRef candref = RecoChargedCandidateRef(hltcands[i][k]); 00168 reco::TrackRef tk = candref->get<reco::TrackRef>(); 00169 double pt = tk->pt(); 00170 if ( ptcut == thePtMin ) { 00171 hHLTpt[i]->Fill(pt); 00172 pair<double,double> angularInfo=getAngle(candref->eta(),candref->phi(), refVector ); 00173 if ( angularInfo.first < 999.){ 00174 float dphi=angularInfo.second-candref->phi(); 00175 float deta=angularInfo.first-candref->eta(); 00176 float d=sqrt(dphi*dphi+deta*deta); 00177 if (i==0)hL2DR->Fill(d); 00178 if ((modules_in_this_event==2 && i==1)||i==2 )hL3DR->Fill(d); 00179 LogTrace("HLTMuonVal")<<"Filling done"; 00180 } 00181 } 00182 double err0 = tk->error(0); 00183 double abspar0 = fabs(tk->parameter(0)); 00184 // convert to 90% efficiency threshold 00185 if (abspar0>0) pt += theNSigmas[i]*err0/abspar0*pt; 00186 if (pt>ptcut) nFound++; 00187 } 00188 if (nFound>=theNumberOfObjects) { 00189 hHLTeff[i]->Fill(ptcut,this_event_weight); 00190 } else { 00191 break; 00192 } 00193 } 00194 } 00195 }
void HLTMuonRate::BookHistograms | ( | ) |
Definition at line 217 of file HLTMuonRate.cc.
References DQMStore::bookInt(), BookIt(), DQMStore::cd(), dbe_, edm::InputTag::encode(), encode(), h, hEtaNor, hHLTeff, hHLTpt, hHLTrate, hL1DR, hL1eff, hL1Eta, hL1Phi, hL1pt, hL1rate, hL2DR, hL3DR, hPhiNor, hPtNor, hSteps, i, NumberOfEvents, NumberOfL1Events, MonitorElement::setAxisTitle(), DQMStore::setCurrentFolder(), theHLTCollectionLabels, theL1CollectionLabel, theLuminosity, theNbins, thePtMax, and thePtMin.
00217 { 00218 char chname[256]; 00219 char chtitle[256]; 00220 char * mylabel ; 00221 char * mydirlabel ; 00222 char str[100],str2[100],str3[100]; 00223 vector<TH1F*> h; 00224 if (dbe_) { 00225 dbe_->cd(); 00226 dbe_->setCurrentFolder("HLT/Muon"); 00227 snprintf(str2, 99, "%s",theL1CollectionLabel.encode().c_str() ); 00228 mydirlabel = strtok(str2,"L1"); 00229 snprintf(str3,99, "HLT/Muon/RateEfficiencies/%s",mydirlabel); 00230 dbe_->setCurrentFolder(str3); 00231 NumberOfEvents=dbe_->bookInt("NumberOfEvents"); 00232 NumberOfL1Events=dbe_->bookInt("NumberOfL1Events"); 00233 snprintf(str, 99, "%s",theL1CollectionLabel.encode().c_str() ); 00234 mylabel = strtok(str,":"); 00235 snprintf(chname, 255, "L1DR_%s", mylabel); 00236 snprintf(chtitle, 255, "L1 association #DR, label=%s", mylabel); 00237 hL1DR= BookIt(chname, chtitle,theNbins, 0., 0.5); 00238 snprintf(chname, 255, "L2DR_%s", mylabel); 00239 snprintf(chtitle, 255, "L2 association #DR, label=%s", mylabel); 00240 hL2DR= BookIt(chname, chtitle,theNbins, 0., 0.5); 00241 snprintf(chname, 255, "L3DR_%s", mylabel); 00242 snprintf(chtitle, 255, "L3 association #DR, label=%s", mylabel); 00243 hL3DR= BookIt(chname, chtitle,theNbins, 0., 0.5); 00244 snprintf(chname, 255, "HLTSteps_%s", mylabel); 00245 snprintf(chtitle, 255, "Events passing the HLT filters, label=%s", mylabel); 00246 hSteps= BookIt(chname, chtitle,5, 0.5, 5.5); 00247 snprintf(chname, 255, "eff_%s", mylabel); 00248 snprintf(chtitle, 255, "Efficiency (%%) vs L1 Pt threshold (GeV/c), label=%s", mylabel); 00249 hL1eff = BookIt(chname, chtitle, theNbins, thePtMin, thePtMax); 00250 snprintf(chname, 255, "rate_%s", mylabel); 00251 snprintf(chtitle, 255, "Rate (Hz) vs L1 Pt threshold (GeV/c), label=%s, L=%.2E (cm^{-2} s^{-1})", mylabel, theLuminosity*1.e33); 00252 hL1rate = BookIt(chname, chtitle, theNbins, thePtMin, thePtMax); 00253 dbe_->cd(); 00254 snprintf(str3,99, "HLT/Muon/Distributions/%s",mydirlabel); 00255 dbe_->setCurrentFolder(str3); 00256 snprintf(chname, 255, "pt_%s", mylabel); 00257 snprintf(chtitle, 255, "L1 Pt distribution label=%s", mylabel); 00258 hL1pt = BookIt(chname, chtitle, theNbins, thePtMin, thePtMax); 00259 snprintf(chtitle, 255, "L1 max ref pt distribution label=%s", mylabel); 00260 snprintf(chname, 255, "ptNorm_%s", mylabel); 00261 hPtNor = BookIt(chname, chtitle, theNbins, thePtMin, thePtMax); 00262 snprintf(chname, 255, "etaNorm_%s",mylabel); 00263 snprintf(chtitle, 255, "Norm #eta "); 00264 hEtaNor = BookIt(chname, chtitle, 50, -2.5, 2.5); 00265 snprintf(chname, 255, "phiNorm_%s",mylabel); 00266 snprintf(chtitle, 255, "Norm #phi"); 00267 hPhiNor = BookIt(chname, chtitle, 50, -3.3, 3.3); 00268 snprintf(chname, 255, "eta_%s", mylabel); 00269 snprintf(chtitle, 255, "L1 Eta distribution label=%s", mylabel); 00270 hL1Eta = BookIt(chname, chtitle, 50, -2.5, 2.5); 00271 snprintf(chname, 255, "phi_%s", mylabel); 00272 snprintf(chtitle, 255, "L1 Phi distribution label=%s", mylabel); 00273 hL1Phi = BookIt(chname, chtitle, 50, -3.3, 3.3); 00274 for (unsigned int i=0; i<theHLTCollectionLabels.size(); i++) { 00275 dbe_->cd(); 00276 snprintf(str3,99, "HLT/Muon/RateEfficiencies/%s",mydirlabel); 00277 dbe_->setCurrentFolder(str3); 00278 snprintf(str, 99, "%s",theHLTCollectionLabels[i].encode().c_str() ); 00279 mylabel = strtok(str,":"); 00280 snprintf(chname, 255, "eff_%s",mylabel ); 00281 snprintf(chtitle, 255, "Efficiency (%%) vs HLT Pt threshold (GeV/c), label=%s", mylabel); 00282 hHLTeff.push_back(BookIt(chname, chtitle, theNbins, thePtMin, thePtMax)); 00283 snprintf(chname, 255, "rate_%s", mylabel); 00284 snprintf(chtitle, 255, "Rate (Hz) vs HLT Pt threshold (GeV/c), label=%s, L=%.2E (cm^{-2} s^{-1})", mylabel, theLuminosity*1.e33); 00285 hHLTrate.push_back(BookIt(chname, chtitle, theNbins, thePtMin, thePtMax)); 00286 dbe_->cd(); 00287 snprintf(str3,99, "HLT/Muon/Distributions/%s",mydirlabel); 00288 dbe_->setCurrentFolder(str3); 00289 snprintf(chname, 255, "pt_%s",mylabel ); 00290 snprintf(chtitle, 255, "Pt distribution label=%s", mylabel); 00291 hHLTpt.push_back(BookIt(chname, chtitle, theNbins, thePtMin, thePtMax)); 00292 } 00293 hSteps->setAxisTitle("Trigger Filtering Step"); 00294 hSteps->setAxisTitle("Events passing Trigger Step",2); 00295 hL1eff->setAxisTitle("90% Muon Pt threshold (GeV/c)"); 00296 hL1rate->setAxisTitle("90% Muon Pt threshold (GeV/c)"); 00297 hL1rate->setAxisTitle("Rate (Hz)",2); 00298 hL1pt->setAxisTitle("Muon Pt (GeV/c)"); 00299 for (unsigned int i=0; i<theHLTCollectionLabels.size(); i++) { 00300 hHLTeff[i]->setAxisTitle("90% Muon Pt threshold (GeV/c)"); 00301 hHLTrate[i]->setAxisTitle("Rate (Hz)",2); 00302 hHLTrate[i]->setAxisTitle("90% Muon Pt threshold (GeV/c)",1); 00303 hHLTpt[i]->setAxisTitle("HLT Muon Pt (GeV/c)",1); 00304 } 00305 } 00306 }
MonitorElement* HLTMuonRate::BookIt | ( | char | name[], | |
char | title[], | |||
int | bins, | |||
float | min, | |||
float | max | |||
) |
Referenced by BookHistograms().
pair< double, double > HLTMuonRate::getAngle | ( | double | eta, | |
double | phi, | |||
Handle< vector< XYZTLorentzVectorD > > & | refVector | |||
) | [private] |
Definition at line 197 of file HLTMuonRate.cc.
References angle(), deltaR(), i, LogTrace, and funct::sqrt().
Referenced by analyze().
00198 { 00199 LogTrace("HLTMuonVal")<< "in getAngle"; 00200 double candDeltaR = 0.4; 00201 pair<double,double> angle(999.,999.); 00202 LogTrace("HLTMuonVal")<< " candidate eta="<<eta<<" and phi="<<phi; 00203 for (unsigned int i = 0; i < refVector->size(); i++ ) { 00204 double Deta=eta - refVector->at(i).eta(); 00205 double Dphi=phi - refVector->at(i).phi(); 00206 double deltaR = sqrt(Deta*Deta+Dphi*Dphi); 00207 if ( deltaR < candDeltaR ) { 00208 candDeltaR=deltaR; 00209 angle.first = refVector->at(i).eta(); 00210 angle.second = refVector->at(i).phi(); 00211 } 00212 } 00213 LogTrace("HLTMuonVal")<< "Best deltaR="<<candDeltaR; 00214 return angle; 00215 }
void HLTMuonRate::WriteHistograms | ( | ) |
Definition at line 316 of file HLTMuonRate.cc.
References dbe_, DQMStore::save(), and theRootFileName.
00316 { 00317 if (theRootFileName.size() != 0 && dbe_) dbe_->save(theRootFileName); 00318 return; 00319 }
DQMStore* HLTMuonRate::dbe_ [private] |
Definition at line 50 of file HLTMuonRate.h.
Referenced by BookHistograms(), HLTMuonRate(), and WriteHistograms().
MonitorElement* HLTMuonRate::hEtaNor [private] |
std::vector<MonitorElement*> HLTMuonRate::hHLTeff [private] |
std::vector<MonitorElement*> HLTMuonRate::hHLTeta [private] |
Definition at line 64 of file HLTMuonRate.h.
std::vector<MonitorElement*> HLTMuonRate::hHLTphi [private] |
Definition at line 65 of file HLTMuonRate.h.
std::vector<MonitorElement*> HLTMuonRate::hHLTpt [private] |
std::vector<MonitorElement*> HLTMuonRate::hHLTrate [private] |
MonitorElement* HLTMuonRate::hL1DR [private] |
MonitorElement* HLTMuonRate::hL1eff [private] |
MonitorElement* HLTMuonRate::hL1Eta [private] |
MonitorElement* HLTMuonRate::hL1Phi [private] |
MonitorElement* HLTMuonRate::hL1pt [private] |
MonitorElement* HLTMuonRate::hL1rate [private] |
MonitorElement * HLTMuonRate::hL2DR [private] |
MonitorElement * HLTMuonRate::hL3DR [private] |
MonitorElement* HLTMuonRate::hPhiNor [private] |
MonitorElement* HLTMuonRate::hPtNor [private] |
MonitorElement* HLTMuonRate::hSteps [private] |
edm::InputTag HLTMuonRate::InputLabel [private] |
MonitorElement* HLTMuonRate::NumberOfEvents [private] |
MonitorElement * HLTMuonRate::NumberOfL1Events [private] |
double HLTMuonRate::theCrossSection [private] |
std::vector<edm::InputTag> HLTMuonRate::theHLTCollectionLabels [private] |
Definition at line 40 of file HLTMuonRate.h.
Referenced by analyze(), BookHistograms(), and HLTMuonRate().
double HLTMuonRate::theHLTReferenceThreshold [private] |
Definition at line 39 of file HLTMuonRate.h.
Referenced by analyze(), BookHistograms(), and HLTMuonRate().
double HLTMuonRate::theL1ReferenceThreshold [private] |
double HLTMuonRate::theLuminosity [private] |
unsigned int HLTMuonRate::theNbins [private] |
Definition at line 48 of file HLTMuonRate.h.
Referenced by analyze(), BookHistograms(), and HLTMuonRate().
std::vector<double> HLTMuonRate::theNSigmas [private] |
int HLTMuonRate::theNumberOfEvents [private] |
int HLTMuonRate::theNumberOfL1Events [private] |
unsigned int HLTMuonRate::theNumberOfObjects [private] |
double HLTMuonRate::thePtMax [private] |
Definition at line 47 of file HLTMuonRate.h.
Referenced by analyze(), BookHistograms(), and HLTMuonRate().
double HLTMuonRate::thePtMin [private] |
Definition at line 46 of file HLTMuonRate.h.
Referenced by analyze(), BookHistograms(), and HLTMuonRate().
std::string HLTMuonRate::theRootFileName [private] |
int HLTMuonRate::this_event_weight [private] |