#include <RPCDeadChannelTest.h>
Public Member Functions | |
void | analyze (const edm::Event &, const edm::EventSetup &) |
Analyze. | |
void | beginJob (DQMStore *, std::string) |
BeginJob. | |
void | beginLuminosityBlock (edm::LuminosityBlock const &, edm::EventSetup const &) |
Begin Lumi block. | |
void | beginRun (const edm::Run &, const edm::EventSetup &) |
void | clientOperation (edm::EventSetup const &c) |
void | endJob () |
Endjob. | |
void | endLuminosityBlock (edm::LuminosityBlock const &, edm::EventSetup const &) |
End Lumi Block. | |
void | endRun (const edm::Run &, const edm::EventSetup &) |
void | getMonitorElements (std::vector< MonitorElement * > &, std::vector< RPCDetId > &) |
RPCDeadChannelTest (const edm::ParameterSet &ps) | |
Constructor. | |
virtual | ~RPCDeadChannelTest () |
Destructor. | |
Private Attributes | |
DQMStore * | dbe_ |
MonitorElement * | DEADDisk [10] |
MonitorElement * | DEADWheel [5] |
std::string | globalFolder_ |
std::vector< RPCDetId > | myDetIds_ |
std::vector< MonitorElement * > | myOccupancyMe_ |
int | numberOfDisks_ |
int | numberOfRings_ |
int | prescaleFactor_ |
bool | useRollInfo_ |
Definition at line 10 of file RPCDeadChannelTest.h.
RPCDeadChannelTest::RPCDeadChannelTest | ( | const edm::ParameterSet & | ps | ) |
Constructor.
Definition at line 16 of file RPCDeadChannelTest.cc.
References edm::ParameterSet::getUntrackedParameter(), numberOfDisks_, numberOfRings_, prescaleFactor_, and useRollInfo_.
{ edm::LogVerbatim ("rpcdeadchanneltest") << "[RPCDeadChannelTest]: Constructor"; useRollInfo_ = ps.getUntrackedParameter<bool>("UseRollInfo", false); prescaleFactor_ = ps.getUntrackedParameter<int>("DiagnosticPrescale", 1); numberOfDisks_ = ps.getUntrackedParameter<int>("NumberOfEndcapDisks", 3); numberOfRings_ = ps.getUntrackedParameter<int>("NumberOfEndcapRings", 2); }
RPCDeadChannelTest::~RPCDeadChannelTest | ( | ) | [virtual] |
void RPCDeadChannelTest::analyze | ( | const edm::Event & | iEvent, |
const edm::EventSetup & | c | ||
) | [virtual] |
void RPCDeadChannelTest::beginJob | ( | DQMStore * | dbe, |
std::string | workingFolder | ||
) | [virtual] |
BeginJob.
Implements RPCClient.
Definition at line 29 of file RPCDeadChannelTest.cc.
References dbe_, and globalFolder_.
{ edm::LogVerbatim ("rpcdeadchanneltest") << "[RPCDeadChannelTest]: Begin Job"; globalFolder_ = workingFolder; dbe_=dbe; }
void RPCDeadChannelTest::beginLuminosityBlock | ( | edm::LuminosityBlock const & | lumiSeg, |
edm::EventSetup const & | context | ||
) | [virtual] |
void RPCDeadChannelTest::beginRun | ( | const edm::Run & | r, |
const edm::EventSetup & | c | ||
) | [virtual] |
Implements RPCClient.
Definition at line 139 of file RPCDeadChannelTest.cc.
References DQMStore::book2D(), dbe_, DEADDisk, DEADWheel, DQMStore::get(), MonitorElement::getName(), globalFolder_, i, rpcdqm::utils::labelXAxisSector(), rpcdqm::utils::labelXAxisSegment(), rpcdqm::utils::labelYAxisRing(), rpcdqm::utils::labelYAxisRoll(), MessageLogger_cff::limit, numberOfDisks_, numberOfRings_, evf::evtn::offset(), DQMStore::removeElement(), DQMStore::setCurrentFolder(), useRollInfo_, x, and detailsBasic3DVector::y.
{ MonitorElement* me; dbe_->setCurrentFolder( globalFolder_); std::stringstream histoName; rpcdqm::utils rpcUtils; int limit = numberOfDisks_; if(numberOfDisks_ < 2) limit = 2; for (int i = -1 * limit; i<= limit;i++ ){//loop on wheels and disks if (i>-3 && i<3){//wheels histoName.str(""); histoName<<"DeadChannelFraction_Roll_vs_Sector_Wheel"<<i; me = 0; me = dbe_->get(globalFolder_ +"/"+ histoName.str()); if (0!=me ) { dbe_->removeElement(me->getName()); } DEADWheel[i+2] = dbe_->book2D(histoName.str().c_str(), histoName.str().c_str(), 12, 0.5, 12.5, 21, 0.5, 21.5); for (int x = 1; x<=12; x++) for(int y=1; y<=21; y++) DEADWheel[i+2]->setBinContent(x,y,-1); rpcUtils.labelXAxisSector( DEADWheel[i+2]); rpcUtils.labelYAxisRoll( DEADWheel[i+2], 0, i, useRollInfo_); }//end wheels if (i == 0 || i > numberOfDisks_ || i< (-1 * numberOfDisks_))continue; int offset = numberOfDisks_; if (i>0) offset --; //used to skip case equale to zero histoName.str(""); histoName<<"DeadChannelFraction_Ring_vs_Segment_Disk"<<i; me = 0; me = dbe_->get(globalFolder_ +"/"+ histoName.str()); if (0!=me ) { dbe_->removeElement(me->getName()); } DEADDisk[i+offset] = dbe_->book2D(histoName.str().c_str(), histoName.str().c_str(),36, 0.5, 36.5, 3*numberOfRings_, 0.5,3*numberOfRings_+ 0.5); rpcUtils.labelXAxisSegment(DEADDisk[i+offset]); rpcUtils.labelYAxisRing(DEADDisk[i+offset], numberOfRings_ ,useRollInfo_); }//end loop on wheels and disks }
void RPCDeadChannelTest::clientOperation | ( | edm::EventSetup const & | c | ) | [virtual] |
Implements RPCClient.
Definition at line 74 of file RPCDeadChannelTest.cc.
References DEADDisk, DEADWheel, rpcdqm::utils::detId2RollNr(), MonitorElement::getBinContent(), MonitorElement::getNbinsX(), MonitorElement::getQReport(), QReport::getQTresult(), i, myDetIds_, myOccupancyMe_, NULL, numberOfDisks_, numberOfRings_, RPCDetId::region(), RPCDetId::ring(), RPCDetId::roll(), RPCDetId::sector(), RPCGeomServ::segment(), MonitorElement::setBinContent(), RPCDetId::station(), and x.
{ edm::LogVerbatim ("rpcdeadchanneltest") <<"[RPCDeadChannelTest]:Client Operation"; MonitorElement * DEAD = NULL; //Loop on chambers for (unsigned int i = 0 ; i<myOccupancyMe_.size();i++){ RPCDetId & detId = myDetIds_[i]; MonitorElement * myMe = myOccupancyMe_[i]; if (! myMe ) continue; const QReport * theOccupancyQReport = myMe->getQReport("DeadChannel_0"); float deadFraction = 0.0 ; if(theOccupancyQReport) { float qtresult = theOccupancyQReport->getQTresult(); // std::vector<dqm::me_util::Channel> badChannels = theOccupancyQReport->getBadChannels(); deadFraction = 1.0 - qtresult; }else{ int xBins = myMe->getNbinsX(); float emptyBins = 0.0; for(int x = 1 ; x<= xBins ; x++){if(myMe->getBinContent(x) == 0 ) {emptyBins++;}} if (xBins != 0){ deadFraction = emptyBins/xBins;} } if (detId.region()==0) DEAD = DEADWheel[detId.ring() + 2] ; else{ if(-detId.station()+ numberOfDisks_ >= 0 ){ if(detId.region()<0){ DEAD = DEADDisk[-detId.station() + numberOfDisks_]; }else{ DEAD = DEADDisk[detId.station() + numberOfDisks_-1]; } } } if (DEAD){ int xBin,yBin; if(detId.region()==0){//Barrel xBin= detId.sector(); rpcdqm::utils rollNumber; yBin = rollNumber.detId2RollNr(detId); }else{//Endcap //get segment number RPCGeomServ RPCServ(detId); xBin = RPCServ.segment(); (numberOfRings_ == 3 ? yBin= detId.ring()*3-detId.roll()+1 : yBin= (detId.ring()-1)*3-detId.roll()+1); } DEAD->setBinContent(xBin,yBin, deadFraction ); } }//End loop on rolls in given chambers }
void RPCDeadChannelTest::endJob | ( | void | ) | [virtual] |
void RPCDeadChannelTest::endLuminosityBlock | ( | edm::LuminosityBlock const & | lumiSeg, |
edm::EventSetup const & | iSetup | ||
) | [virtual] |
void RPCDeadChannelTest::endRun | ( | const edm::Run & | r, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
Implements RPCClient.
Definition at line 35 of file RPCDeadChannelTest.cc.
{ edm::LogVerbatim ("rpcdeadchanneltest") << "[RPCDeadChannelTest]: End run"; }
void RPCDeadChannelTest::getMonitorElements | ( | std::vector< MonitorElement * > & | meVector, |
std::vector< RPCDetId > & | detIdVector | ||
) | [virtual] |
Implements RPCClient.
Definition at line 42 of file RPCDeadChannelTest.cc.
References i, myDetIds_, myOccupancyMe_, rpcdqm::OCCUPANCY, and packageDocSplitter::tagList.
{ //Get Occupancy ME for each roll for (unsigned int i = 0 ; i<meVector.size(); i++){ bool flag= false; DQMNet::TagList tagList; tagList = meVector[i]->getTags(); DQMNet::TagList::iterator tagItr = tagList.begin(); while (tagItr != tagList.end() && !flag ) { if((*tagItr) == rpcdqm::OCCUPANCY) flag= true; tagItr++; } if(flag){ myOccupancyMe_.push_back(meVector[i]); myDetIds_.push_back(detIdVector[i]); } } }
DQMStore* RPCDeadChannelTest::dbe_ [private] |
Definition at line 56 of file RPCDeadChannelTest.h.
Referenced by beginJob(), beginRun(), and ~RPCDeadChannelTest().
MonitorElement* RPCDeadChannelTest::DEADDisk[10] [private] |
Definition at line 62 of file RPCDeadChannelTest.h.
Referenced by beginRun(), and clientOperation().
MonitorElement* RPCDeadChannelTest::DEADWheel[5] [private] |
Definition at line 61 of file RPCDeadChannelTest.h.
Referenced by beginRun(), and clientOperation().
std::string RPCDeadChannelTest::globalFolder_ [private] |
Definition at line 52 of file RPCDeadChannelTest.h.
Referenced by beginJob(), and beginRun().
std::vector<RPCDetId> RPCDeadChannelTest::myDetIds_ [private] |
Definition at line 54 of file RPCDeadChannelTest.h.
Referenced by clientOperation(), and getMonitorElements().
std::vector<MonitorElement *> RPCDeadChannelTest::myOccupancyMe_ [private] |
Definition at line 53 of file RPCDeadChannelTest.h.
Referenced by clientOperation(), and getMonitorElements().
int RPCDeadChannelTest::numberOfDisks_ [private] |
Definition at line 59 of file RPCDeadChannelTest.h.
Referenced by beginRun(), clientOperation(), and RPCDeadChannelTest().
int RPCDeadChannelTest::numberOfRings_ [private] |
Definition at line 60 of file RPCDeadChannelTest.h.
Referenced by beginRun(), clientOperation(), and RPCDeadChannelTest().
int RPCDeadChannelTest::prescaleFactor_ [private] |
Definition at line 51 of file RPCDeadChannelTest.h.
Referenced by RPCDeadChannelTest().
bool RPCDeadChannelTest::useRollInfo_ [private] |
Definition at line 55 of file RPCDeadChannelTest.h.
Referenced by beginRun(), and RPCDeadChannelTest().