CMS 3D CMS Logo

RPCDeadChannelTest.cc
Go to the documentation of this file.
1 /* * \author Anna Cimmino*/
4 // Framework
6 // Geometry
10 
11 #include <sstream>
12 
14 
15  edm::LogVerbatim ("rpcdeadchanneltest") << "[RPCDeadChannelTest]: Constructor";
16 
17  useRollInfo_ = ps.getUntrackedParameter<bool>("UseRollInfo", false);
18 
19  prescaleFactor_ = ps.getUntrackedParameter<int>("DiagnosticPrescale", 1);
20  numberOfDisks_ = ps.getUntrackedParameter<int>("NumberOfEndcapDisks", 4);
21  numberOfRings_ = ps.getUntrackedParameter<int>("NumberOfEndcapRings", 2);
22 }
23 
25 
27  edm::LogVerbatim ("rpcdeadchanneltest") << "[RPCDeadChannelTest]: Begin Job";
28  globalFolder_ = workingFolder;
29 
30 }
31 
32 
33 void RPCDeadChannelTest::getMonitorElements(std::vector<MonitorElement *> & meVector, std::vector<RPCDetId> & detIdVector, std::string & clientHistoName){
34 
35  for (unsigned int i = 0 ; i<meVector.size(); i++){
36 
37  std::string meName = meVector[i]->getName();
38 
39  if(meName.find(clientHistoName) != std::string::npos){
40  myOccupancyMe_.push_back(meVector[i]);
41  myDetIds_.push_back(detIdVector[i]);
42  }
43  }
44 }
45 
46 
48 
49  edm::LogVerbatim ("rpcdeadchanneltest") <<"[RPCDeadChannelTest]:Client Operation";
50 
51 
52  MonitorElement * DEAD = nullptr;
53 
54  //Loop on chambers
55  for (unsigned int i = 0 ; i<myOccupancyMe_.size();i++){
56 
57  RPCDetId & detId = myDetIds_[i];
59 
60  if (! myMe ) continue;
61 
62  const QReport * theOccupancyQReport = myMe->getQReport("DeadChannel_0");
63 
64  float deadFraction = 0.0 ;
65 
66  if(theOccupancyQReport) {
67 
68  float qtresult = theOccupancyQReport->getQTresult();
69  // std::vector<dqm::me_util::Channel> badChannels = theOccupancyQReport->getBadChannels();
70  deadFraction = 1.0 - qtresult;
71 
72  }else{
73  int xBins = myMe->getNbinsX();
74  float emptyBins = 0.0;
75  for(int x = 1 ; x<= xBins ; x++){if(myMe->getBinContent(x) == 0 ) {emptyBins++;}}
76  if (xBins != 0){ deadFraction = emptyBins/xBins;}
77  }
78 
79  if (detId.region()==0) DEAD = DEADWheel[detId.ring() + 2] ;
80  else{
81  if(-detId.station()+ numberOfDisks_ >= 0 ){
82 
83  if(detId.region()<0){
84  DEAD = DEADDisk[-detId.station() + numberOfDisks_];
85  }else{
86  DEAD = DEADDisk[detId.station() + numberOfDisks_-1];
87  }
88  }
89  }
90 
91  if (DEAD){
92  int xBin,yBin;
93  if(detId.region()==0){//Barrel
94  xBin= detId.sector();
95  rpcdqm::utils rollNumber;
96  yBin = rollNumber.detId2RollNr(detId);
97  }else{//Endcap
98  //get segment number
99  RPCGeomServ RPCServ(detId);
100  xBin = RPCServ.segment();
101  (numberOfRings_ == 3 ? yBin= detId.ring()*3-detId.roll()+1 : yBin= (detId.ring()-1)*3-detId.roll()+1);
102  }
103  DEAD->setBinContent(xBin,yBin, deadFraction );
104 
105  }
106 
107  }//End loop on rolls in given chambers
108 
109 }
110 
112 
114 
115  std::stringstream histoName;
116 
117  rpcdqm::utils rpcUtils;
118 
119  int limit = numberOfDisks_;
120  if(numberOfDisks_ < 2) limit = 2;
121 
122  for (int i = -1 * limit; i<= limit;i++ ){//loop on wheels and disks
123  if (i>-3 && i<3){//wheels
124  histoName.str("");
125  histoName<<"DeadChannelFraction_Roll_vs_Sector_Wheel"<<i;
126  DEADWheel[i+2] = ibooker.book2D(histoName.str().c_str(), histoName.str().c_str(), 12, 0.5, 12.5, 21, 0.5, 21.5);
127 
128  for (int x = 1; x<=12; x++){
129  for(int y=1; y<=21; y++){
130  DEADWheel[i+2]->setBinContent(x,y,-1);
131  }
132  }
133 
134  rpcUtils.labelXAxisSector( DEADWheel[i+2]);
135  rpcUtils.labelYAxisRoll( DEADWheel[i+2], 0, i, useRollInfo_);
136  }//end wheels
137 
138  if (i == 0 || i > numberOfDisks_ || i< (-1 * numberOfDisks_)){continue;}
139 
140  int offset = numberOfDisks_;
141  if (i>0) {offset --;} //used to skip case equale to zero
142 
143  histoName.str("");
144  histoName<<"DeadChannelFraction_Ring_vs_Segment_Disk"<<i;
145  DEADDisk[i+offset] = ibooker.book2D(histoName.str().c_str(), histoName.str().c_str(),36, 0.5, 36.5, 3*numberOfRings_, 0.5,3*numberOfRings_+ 0.5);
146 
147  rpcUtils.labelXAxisSegment(DEADDisk[i+offset]);
148  rpcUtils.labelYAxisRing(DEADDisk[i+offset], numberOfRings_ ,useRollInfo_);
149 
150  }//end loop on wheels and disks
151 
152 
153 }
154 
155 
void clientOperation() override
T getUntrackedParameter(std::string const &, T const &) const
const QReport * getQReport(const std::string &qtname) const
get QReport corresponding to <qtname> (null pointer if QReport does not exist)
void setBinContent(int binx, double content)
set content of bin (1-D)
~RPCDeadChannelTest() override
Destructor.
void getMonitorElements(std::vector< MonitorElement * > &, std::vector< RPCDetId > &, std::string &) override
MonitorElement * DEADWheel[5]
void beginJob(std::string &) override
void labelXAxisSegment(MonitorElement *myMe)
Definition: utils.h:264
void labelYAxisRoll(MonitorElement *myMe, int region, int ring, bool useRollInfo)
Definition: utils.h:283
void labelYAxisRing(MonitorElement *myMe, int numberOfRings, bool useRollInfo)
Definition: utils.h:311
int roll() const
Definition: RPCDetId.h:120
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
int ring() const
Definition: RPCDetId.h:72
void myBooker(DQMStore::IBooker &) override
void labelXAxisSector(MonitorElement *myMe)
Definition: utils.h:249
std::vector< MonitorElement * > myOccupancyMe_
virtual int segment()
Definition: RPCGeomServ.cc:469
int detId2RollNr(const RPCDetId &_id)
Definition: utils.h:31
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:109
std::vector< RPCDetId > myDetIds_
int sector() const
Sector id: the group of chambers at same phi (and increasing r)
Definition: RPCDetId.h:102
double getBinContent(int binx) const
get content of bin (1-D)
float getQTresult() const
get test result i.e. prob value
Definition: QReport.h:20
RPCDeadChannelTest(const edm::ParameterSet &ps)
Constructor.
int getNbinsX() const
get # of bins in X-axis
MonitorElement * DEADDisk[10]
yBin
Definition: cuy.py:893
int region() const
Region id: 0 for Barrel, +/-1 For +/- Endcap.
Definition: RPCDetId.h:63
int station() const
Definition: RPCDetId.h:96