#include <DTLocalTriggerEfficiencyTest.h>
Public Member Functions | |
DTLocalTriggerEfficiencyTest (const edm::ParameterSet &ps) | |
Constructor. | |
virtual | ~DTLocalTriggerEfficiencyTest () |
Destructor. | |
Protected Member Functions | |
void | beginRun (const edm::Run &r, const edm::EventSetup &c) |
BeginRun. | |
void | bookChambHistos (DTChamberId chambId, std::string htype) |
Book the new MEs (for each chamber) | |
void | makeEfficiencyME (TH1D *numerator, TH1D *denominator, MonitorElement *result) |
Compute efficiency plots. | |
void | makeEfficiencyME2D (TH2F *numerator, TH2F *denominator, MonitorElement *result) |
Compute 2D efficiency plots. | |
void | runClientDiagnostic () |
DQM Client Diagnostic. | |
Private Attributes | |
std::map< uint32_t, std::map < std::string, MonitorElement * > > | chambME |
DTTrigGeomUtils * | trigGeomUtils |
* DQM Test Client
Definition at line 21 of file DTLocalTriggerEfficiencyTest.h.
DTLocalTriggerEfficiencyTest::DTLocalTriggerEfficiencyTest | ( | const edm::ParameterSet & | ps | ) |
Constructor.
Definition at line 37 of file DTLocalTriggerEfficiencyTest.cc.
References DTLocalTriggerBaseTest::baseFolderDCC, DTLocalTriggerBaseTest::baseFolderDDU, and DTLocalTriggerBaseTest::setConfig().
: trigGeomUtils(0) { setConfig(ps,"DTLocalTriggerEfficiency"); baseFolderDCC = "DT/03-LocalTrigger-DCC/"; baseFolderDDU = "DT/04-LocalTrigger-DDU/"; }
DTLocalTriggerEfficiencyTest::~DTLocalTriggerEfficiencyTest | ( | ) | [virtual] |
Destructor.
Definition at line 46 of file DTLocalTriggerEfficiencyTest.cc.
References trigGeomUtils.
{ if (trigGeomUtils) { delete trigGeomUtils; } }
void DTLocalTriggerEfficiencyTest::beginRun | ( | const edm::Run & | r, |
const edm::EventSetup & | c | ||
) | [protected, virtual] |
BeginRun.
Reimplemented from DTLocalTriggerBaseTest.
Definition at line 55 of file DTLocalTriggerEfficiencyTest.cc.
References bookChambHistos(), DTLocalTriggerBaseTest::bookSectorHistos(), DTLocalTriggerBaseTest::bookWheelHistos(), edm::ParameterSet::getUntrackedParameter(), DTLocalTriggerBaseTest::hwSource, DTLocalTriggerBaseTest::hwSources, DTLocalTriggerBaseTest::muonGeom, DTLocalTriggerBaseTest::parameters, trigGeomUtils, DTLocalTriggerBaseTest::trigSource, and DTLocalTriggerBaseTest::trigSources.
{ DTLocalTriggerBaseTest::beginRun(r,c); trigGeomUtils = new DTTrigGeomUtils(muonGeom); vector<string>::const_iterator iTr = trigSources.begin(); vector<string>::const_iterator trEnd = trigSources.end(); vector<string>::const_iterator iHw = hwSources.begin(); vector<string>::const_iterator hwEnd = hwSources.end(); //Booking if(parameters.getUntrackedParameter<bool>("staticBooking", true)){ for (; iTr != trEnd; ++iTr){ trigSource = (*iTr); for (; iHw != hwEnd; ++iHw){ hwSource = (*iHw); // Loop over the TriggerUnits for (int wh=-2; wh<=2; ++wh){ for (int sect=1; sect<=12; ++sect){ for (int stat=1; stat<=4; ++stat){ DTChamberId chId(wh,stat,sect); bookChambHistos(chId,"TrigEffPosvsAnglePhi"); bookChambHistos(chId,"TrigEffPosvsAngleHHHLPhi"); bookChambHistos(chId,"TrigEffPosPhi"); bookChambHistos(chId,"TrigEffPosHHHLPhi"); bookChambHistos(chId,"TrigEffAnglePhi"); bookChambHistos(chId,"TrigEffAngleHHHLPhi"); if (stat<=3) { bookChambHistos(chId,"TrigEffPosvsAngleTheta"); bookChambHistos(chId,"TrigEffPosvsAngleHTheta"); bookChambHistos(chId,"TrigEffPosTheta"); bookChambHistos(chId,"TrigEffPosHTheta"); bookChambHistos(chId,"TrigEffAngleTheta"); bookChambHistos(chId,"TrigEffAngleHTheta"); } } bookSectorHistos(wh,sect,"TrigEffPhi"); bookSectorHistos(wh,sect,"TrigEffTheta"); } bookWheelHistos(wh,"TrigEffPhi"); bookWheelHistos(wh,"TrigEffHHHLPhi"); bookWheelHistos(wh,"TrigEffTheta"); bookWheelHistos(wh,"TrigEffHTheta"); } } } } }
void DTLocalTriggerEfficiencyTest::bookChambHistos | ( | DTChamberId | chambId, |
std::string | htype | ||
) | [protected] |
Book the new MEs (for each chamber)
Referenced by beginRun(), and runClientDiagnostic().
void DTLocalTriggerEfficiencyTest::makeEfficiencyME | ( | TH1D * | numerator, |
TH1D * | denominator, | ||
MonitorElement * | result | ||
) | [protected] |
Compute efficiency plots.
Definition at line 245 of file DTLocalTriggerEfficiencyTest.cc.
References newFWLiteAna::bin, error, MonitorElement::getTH1F(), ExpressReco_HICollisions_FallBack::nbins, and mathSSE::sqrt().
Referenced by runClientDiagnostic().
{ TH1F* efficiency = result->getTH1F(); efficiency->Divide(numerator,denominator,1,1,""); int nbins = efficiency->GetNbinsX(); for (int bin=1; bin<=nbins; ++bin){ float error = 0; float bineff = efficiency->GetBinContent(bin); if (denominator->GetBinContent(bin)){ error = sqrt(bineff*(1-bineff)/denominator->GetBinContent(bin)); } else { error = 1; efficiency->SetBinContent(bin,1.); } efficiency->SetBinError(bin,error); } }
void DTLocalTriggerEfficiencyTest::makeEfficiencyME2D | ( | TH2F * | numerator, |
TH2F * | denominator, | ||
MonitorElement * | result | ||
) | [protected] |
Compute 2D efficiency plots.
Definition at line 269 of file DTLocalTriggerEfficiencyTest.cc.
References error, MonitorElement::getTH2F(), and mathSSE::sqrt().
Referenced by runClientDiagnostic().
{ TH2F* efficiency = result->getTH2F(); efficiency->Divide(numerator,denominator,1,1,""); int nbinsx = efficiency->GetNbinsX(); int nbinsy = efficiency->GetNbinsY(); for (int binx=1; binx<=nbinsx; ++binx){ for (int biny=1; biny<=nbinsy; ++biny){ float error = 0; float bineff = efficiency->GetBinContent(binx,biny); if (denominator->GetBinContent(binx,biny)){ error = sqrt(bineff*(1-bineff)/denominator->GetBinContent(binx,biny)); } else { error = 1; efficiency->SetBinContent(binx,biny,0.); } efficiency->SetBinError(binx,biny,error); } } }
void DTLocalTriggerEfficiencyTest::runClientDiagnostic | ( | ) | [protected, virtual] |
DQM Client Diagnostic.
Implements DTLocalTriggerBaseTest.
Definition at line 107 of file DTLocalTriggerEfficiencyTest.cc.
References bookChambHistos(), DTLocalTriggerBaseTest::bookSectorHistos(), DTLocalTriggerBaseTest::bookWheelHistos(), chambME, DTLocalTriggerBaseTest::dbe, end, spr::find(), DTLocalTriggerBaseTest::fullName(), DQMStore::get(), DTLocalTriggerBaseTest::getMEName(), DTLocalTriggerBaseTest::hwSource, DTLocalTriggerBaseTest::hwSources, makeEfficiencyME(), makeEfficiencyME2D(), DetId::rawId(), DTLocalTriggerBaseTest::secME, MonitorElement::setBinContent(), mathSSE::sqrt(), reco::tau::disc::TrackAngle(), DTLocalTriggerBaseTest::trigSource, DTLocalTriggerBaseTest::trigSources, and DTLocalTriggerBaseTest::whME.
{ // Loop over Trig & Hw sources for (vector<string>::const_iterator iTr = trigSources.begin(); iTr != trigSources.end(); ++iTr){ trigSource = (*iTr); for (vector<string>::const_iterator iHw = hwSources.begin(); iHw != hwSources.end(); ++iHw){ hwSource = (*iHw); // Loop over the TriggerUnits for (int stat=1; stat<=4; ++stat){ for (int wh=-2; wh<=2; ++wh){ for (int sect=1; sect<=12; ++sect){ DTChamberId chId(wh,stat,sect); int sector_id = (wh+3)+(sect-1)*5; uint32_t indexCh = chId.rawId(); // Perform Efficiency analysis (Phi+Segments 2D) TH2F * TrackPosvsAngle = getHisto<TH2F>(dbe->get(getMEName("TrackPosvsAngle","Segment", chId))); TH2F * TrackPosvsAngleandTrig = getHisto<TH2F>(dbe->get(getMEName("TrackPosvsAngleandTrig","Segment", chId))); TH2F * TrackPosvsAngleandTrigHHHL = getHisto<TH2F>(dbe->get(getMEName("TrackPosvsAngleandTrigHHHL","Segment", chId))); if (TrackPosvsAngle && TrackPosvsAngleandTrig && TrackPosvsAngleandTrigHHHL && TrackPosvsAngle->GetEntries()>1) { if( chambME[indexCh].find(fullName("TrigEffAnglePhi")) == chambME[indexCh].end()){ bookChambHistos(chId,"TrigEffPosvsAnglePhi"); bookChambHistos(chId,"TrigEffPosvsAngleHHHLPhi"); bookChambHistos(chId,"TrigEffPosPhi"); bookChambHistos(chId,"TrigEffPosHHHLPhi"); bookChambHistos(chId,"TrigEffAnglePhi"); bookChambHistos(chId,"TrigEffAngleHHHLPhi"); } if( secME[sector_id].find(fullName("TrigEffPhi")) == secME[sector_id].end() ){ bookSectorHistos(wh,sect,"TrigEffPhi"); } if( whME[wh].find(fullName("TrigEffPhi")) == whME[wh].end() ){ bookWheelHistos(wh,"TrigEffPhi"); bookWheelHistos(wh,"TrigEffHHHLPhi"); } std::map<std::string,MonitorElement*> *innerME = &(secME[sector_id]); TH1D* TrackPos = TrackPosvsAngle->ProjectionY(); TH1D* TrackAngle = TrackPosvsAngle->ProjectionX(); TH1D* TrackPosandTrig = TrackPosvsAngleandTrig->ProjectionY(); TH1D* TrackAngleandTrig = TrackPosvsAngleandTrig->ProjectionX(); TH1D* TrackPosandTrigHHHL = TrackPosvsAngleandTrigHHHL->ProjectionY(); TH1D* TrackAngleandTrigHHHL = TrackPosvsAngleandTrigHHHL->ProjectionX(); float binEff = float(TrackPosandTrig->GetEntries())/TrackPos->GetEntries(); float binEffHHHL = float(TrackPosandTrigHHHL->GetEntries())/TrackPos->GetEntries(); float binErr = sqrt(binEff*(1-binEff)/TrackPos->GetEntries()); float binErrHHHL = sqrt(binEffHHHL*(1-binEffHHHL)/TrackPos->GetEntries()); MonitorElement* globalEff = innerME->find(fullName("TrigEffPhi"))->second; globalEff->setBinContent(stat,binEff); globalEff->setBinError(stat,binErr); innerME = &(whME[wh]); globalEff = innerME->find(fullName("TrigEffPhi"))->second; globalEff->setBinContent(sect,stat,binEff); globalEff->setBinError(sect,stat,binErr); globalEff = innerME->find(fullName("TrigEffHHHLPhi"))->second; globalEff->setBinContent(sect,stat,binEffHHHL); globalEff->setBinError(sect,stat,binErrHHHL); innerME = &(chambME[indexCh]); makeEfficiencyME(TrackPosandTrig,TrackPos,innerME->find(fullName("TrigEffPosPhi"))->second); makeEfficiencyME(TrackPosandTrigHHHL,TrackPos,innerME->find(fullName("TrigEffPosHHHLPhi"))->second); makeEfficiencyME(TrackAngleandTrig,TrackAngle,innerME->find(fullName("TrigEffAnglePhi"))->second); makeEfficiencyME(TrackAngleandTrigHHHL,TrackAngle,innerME->find(fullName("TrigEffAngleHHHLPhi"))->second); makeEfficiencyME2D(TrackPosvsAngleandTrig,TrackPosvsAngle,innerME->find(fullName("TrigEffPosvsAnglePhi"))->second); makeEfficiencyME2D(TrackPosvsAngleandTrigHHHL,TrackPosvsAngle,innerME->find(fullName("TrigEffPosvsAngleHHHLPhi"))->second); } // Perform Efficiency analysis (Theta+Segments) CB FIXME -> no DCC theta qual info TH2F * TrackThetaPosvsAngle = getHisto<TH2F>(dbe->get(getMEName("TrackThetaPosvsAngle","Segment", chId))); TH2F * TrackThetaPosvsAngleandTrig = getHisto<TH2F>(dbe->get(getMEName("TrackThetaPosvsAngleandTrig","Segment", chId))); TH2F * TrackThetaPosvsAngleandTrigH = getHisto<TH2F>(dbe->get(getMEName("TrackThetaPosvsAngleandTrigH","Segment", chId))); if (TrackThetaPosvsAngle && TrackThetaPosvsAngleandTrig && TrackThetaPosvsAngleandTrigH && TrackThetaPosvsAngle->GetEntries()>1) { if( chambME[indexCh].find(fullName("TrigEffAngleTheta")) == chambME[indexCh].end()){ bookChambHistos(chId,"TrigEffPosvsAngleTheta"); bookChambHistos(chId,"TrigEffPosvsAngleHTheta"); bookChambHistos(chId,"TrigEffPosTheta"); bookChambHistos(chId,"TrigEffPosHTheta"); bookChambHistos(chId,"TrigEffAngleTheta"); bookChambHistos(chId,"TrigEffAngleHTheta"); } if( secME[sector_id].find(fullName("TrigEffTheta")) == secME[sector_id].end() ){ bookSectorHistos(wh,sect,"TrigEffTheta"); } if( whME[wh].find(fullName("TrigEffTheta")) == whME[wh].end() ){ bookWheelHistos(wh,"TrigEffTheta"); bookWheelHistos(wh,"TrigEffHTheta"); } std::map<std::string,MonitorElement*> *innerME = &(secME[sector_id]); TH1D* TrackThetaPos = TrackThetaPosvsAngle->ProjectionY(); TH1D* TrackThetaAngle = TrackThetaPosvsAngle->ProjectionX(); TH1D* TrackThetaPosandTrig = TrackThetaPosvsAngleandTrig->ProjectionY(); TH1D* TrackThetaAngleandTrig = TrackThetaPosvsAngleandTrig->ProjectionX(); TH1D* TrackThetaPosandTrigH = TrackThetaPosvsAngleandTrigH->ProjectionY(); TH1D* TrackThetaAngleandTrigH = TrackThetaPosvsAngleandTrigH->ProjectionX(); float binEff = float(TrackThetaPosandTrig->GetEntries())/TrackThetaPos->GetEntries(); float binErr = sqrt(binEff*(1-binEff)/TrackThetaPos->GetEntries()); float binEffH = float(TrackThetaPosandTrigH->GetEntries())/TrackThetaPos->GetEntries(); float binErrH = sqrt(binEffH*(1-binEffH)/TrackThetaPos->GetEntries()); MonitorElement* globalEff = innerME->find(fullName("TrigEffTheta"))->second; globalEff->setBinContent(stat,binEff); globalEff->setBinError(stat,binErr); innerME = &(whME[wh]); globalEff = innerME->find(fullName("TrigEffTheta"))->second; globalEff->setBinContent(sect,stat,binEff); globalEff->setBinError(sect,stat,binErr); globalEff = innerME->find(fullName("TrigEffHTheta"))->second; globalEff->setBinContent(sect,stat,binEffH); globalEff->setBinError(sect,stat,binErrH); innerME = &(chambME[indexCh]); makeEfficiencyME(TrackThetaPosandTrig,TrackThetaPos,innerME->find(fullName("TrigEffPosTheta"))->second); makeEfficiencyME(TrackThetaPosandTrigH,TrackThetaPos,innerME->find(fullName("TrigEffPosHTheta"))->second); makeEfficiencyME(TrackThetaAngleandTrig,TrackThetaAngle,innerME->find(fullName("TrigEffAngleTheta"))->second); makeEfficiencyME(TrackThetaAngleandTrigH,TrackThetaAngle,innerME->find(fullName("TrigEffAngleHTheta"))->second); makeEfficiencyME2D(TrackThetaPosvsAngleandTrig,TrackThetaPosvsAngle,innerME->find(fullName("TrigEffPosvsAngleTheta"))->second); makeEfficiencyME2D(TrackThetaPosvsAngleandTrigH,TrackThetaPosvsAngle,innerME->find(fullName("TrigEffPosvsAngleHTheta"))->second); } } } } } } }
std::map<uint32_t,std::map<std::string,MonitorElement*> > DTLocalTriggerEfficiencyTest::chambME [private] |
Definition at line 52 of file DTLocalTriggerEfficiencyTest.h.
Referenced by DTTriggerEfficiencyTest::bookChambHistos(), and runClientDiagnostic().
Definition at line 53 of file DTLocalTriggerEfficiencyTest.h.
Referenced by beginRun(), DTTriggerEfficiencyTest::bookChambHistos(), and ~DTLocalTriggerEfficiencyTest().