CMS 3D CMS Logo

RPCClusterSizeTest.cc
Go to the documentation of this file.
3 
4 // Framework
6 
7 // //Geometry
9 
11  edm::LogVerbatim("rpceventsummary") << "[RPCClusterSizeTest]: Constructor";
12 
13  prescaleFactor_ = ps.getUntrackedParameter<int>("DiagnosticPrescale", 1);
14 
15  numberOfDisks_ = ps.getUntrackedParameter<int>("NumberOfEndcapDisks", 4);
16  numberOfRings_ = ps.getUntrackedParameter<int>("NumberOfEndcapRings", 2);
17  testMode_ = ps.getUntrackedParameter<bool>("testMode", false);
18  useRollInfo_ = ps.getUntrackedParameter<bool>("useRollInfo", false);
19 
20  resetMEArrays();
21 }
22 
24 
26  edm::LogVerbatim("rpceventsummary") << "[RPCClusterSizeTest]: Begin job ";
27 
28  globalFolder_ = workingFolder;
29 }
30 
31 void RPCClusterSizeTest::getMonitorElements(std::vector<MonitorElement*>& meVector,
32  std::vector<RPCDetId>& detIdVector,
33  std::string& clientHistoName) {
34  //Get ME for each roll
35  for (unsigned int i = 0; i < meVector.size(); i++) {
36  std::string meName = meVector[i]->getName();
37 
38  if (meName.find(clientHistoName) != std::string::npos) {
39  myClusterMe_.push_back(meVector[i]);
40  myDetIds_.push_back(detIdVector[i]);
41  }
42  }
43 }
44 
46  edm::LogVerbatim("rpceventsummary") << "[RPCClusterSizeTest]:Client Operation";
47 
48  //check some statements and prescale Factor
49  if (myClusterMe_.empty() || myDetIds_.empty())
50  return;
51 
52  MonitorElement* CLS = nullptr; // ClusterSize in 1 bin, Roll vs Sector
53  MonitorElement* CLSD = nullptr; // ClusterSize in 1 bin, Distribution
54  MonitorElement* MEAN = nullptr; // Mean ClusterSize, Roll vs Sector
55  MonitorElement* MEAND = nullptr; // Mean ClusterSize, Distribution
56 
57  std::stringstream meName;
58  RPCDetId detId;
59  MonitorElement* myMe;
60 
61  //Loop on chambers
62  for (unsigned int i = 0; i < myClusterMe_.size(); i++) {
63  myMe = myClusterMe_[i];
64  if (!myMe || myMe->getEntries() == 0)
65  continue;
66 
67  detId = myDetIds_[i];
68 
69  if (detId.region() == 0) {
70  CLS = CLSWheel[detId.ring() + 2];
71  MEAN = MEANWheel[detId.ring() + 2];
72  if (testMode_) {
73  CLSD = CLSDWheel[detId.ring() + 2];
74  MEAND = MEANDWheel[detId.ring() + 2];
75  }
76  } else {
77  if (((detId.station() * detId.region()) + numberOfDisks_) >= 0) {
78  if (detId.region() < 0) {
79  CLS = CLSDisk[(detId.station() * detId.region()) + numberOfDisks_];
80  MEAN = MEANDisk[(detId.station() * detId.region()) + numberOfDisks_];
81  if (testMode_) {
82  CLSD = CLSDDisk[(detId.station() * detId.region()) + numberOfDisks_];
83  MEAND = MEANDDisk[(detId.station() * detId.region()) + numberOfDisks_];
84  }
85  } else {
86  CLS = CLSDisk[(detId.station() * detId.region()) + numberOfDisks_ - 1];
87  MEAN = MEANDisk[(detId.station() * detId.region()) + numberOfDisks_ - 1];
88  if (testMode_) {
89  CLSD = CLSDDisk[(detId.station() * detId.region()) + numberOfDisks_ - 1];
90  MEAND = MEANDDisk[(detId.station() * detId.region()) + numberOfDisks_ - 1];
91  }
92  }
93  }
94  }
95 
96  int xBin, yBin;
97 
98  if (detId.region() == 0) { //Barrel
99 
100  rpcdqm::utils rollNumber;
101  yBin = rollNumber.detId2RollNr(detId);
102  xBin = detId.sector();
103  } else { //Endcap
104 
105  //get segment number
106  RPCGeomServ RPCServ(detId);
107  xBin = RPCServ.segment();
108  (numberOfRings_ == 3 ? yBin = detId.ring() * 3 - detId.roll() + 1
109  : yBin = (detId.ring() - 1) * 3 - detId.roll() + 1);
110  }
111 
112  // Normalization -> # of Entries in first Bin normalaized by total Entries
113 
114  float NormCLS = myMe->getBinContent(1) / myMe->getEntries();
115  float meanCLS = myMe->getMean();
116 
117  if (CLS)
118  CLS->setBinContent(xBin, yBin, NormCLS);
119  if (MEAN)
120  MEAN->setBinContent(xBin, yBin, meanCLS);
121 
122  if (testMode_) {
123  if (MEAND)
124  MEAND->Fill(meanCLS);
125  if (CLSD)
126  CLSD->Fill(NormCLS);
127  }
128 
129  } //End loop on chambers
130 }
131 
133  memset((void*)CLSWheel, 0, sizeof(MonitorElement*) * kWheels);
134  memset((void*)CLSDWheel, 0, sizeof(MonitorElement*) * kWheels);
135  memset((void*)MEANWheel, 0, sizeof(MonitorElement*) * kWheels);
136  memset((void*)MEANDWheel, 0, sizeof(MonitorElement*) * kWheels);
137 
138  memset((void*)CLSDisk, 0, sizeof(MonitorElement*) * kDisks);
139  memset((void*)CLSDDisk, 0, sizeof(MonitorElement*) * kDisks);
140  memset((void*)MEANDisk, 0, sizeof(MonitorElement*) * kDisks);
141  memset((void*)MEANDDisk, 0, sizeof(MonitorElement*) * kDisks);
142 }
143 
145  resetMEArrays();
146 
148 
149  std::stringstream histoName;
150 
151  rpcdqm::utils rpcUtils;
152 
153  // Loop over wheels
154  for (int w = -2; w <= 2; w++) {
155  histoName.str("");
156  histoName << "ClusterSizeIn1Bin_Roll_vs_Sector_Wheel"
157  << w; // ClusterSize in first bin norm. by Entries (2D Roll vs Sector)
158  CLSWheel[w + 2] = ibooker.book2D(histoName.str().c_str(), histoName.str().c_str(), 12, 0.5, 12.5, 21, 0.5, 21.5);
159  rpcUtils.labelXAxisSector(CLSWheel[w + 2]);
160  rpcUtils.labelYAxisRoll(CLSWheel[w + 2], 0, w, useRollInfo_);
161 
162  histoName.str("");
163  histoName << "ClusterSizeMean_Roll_vs_Sector_Wheel" << w; // Avarage ClusterSize (2D Roll vs Sector)
164  MEANWheel[w + 2] = ibooker.book2D(histoName.str().c_str(), histoName.str().c_str(), 12, 0.5, 12.5, 21, 0.5, 21.5);
165 
166  rpcUtils.labelXAxisSector(MEANWheel[w + 2]);
167  rpcUtils.labelYAxisRoll(MEANWheel[w + 2], 0, w, useRollInfo_);
168 
169  if (testMode_) {
170  histoName.str("");
171  histoName << "ClusterSizeIn1Bin_Distribution_Wheel" << w; // ClusterSize in first bin, distribution
172  CLSDWheel[w + 2] = ibooker.book1D(histoName.str().c_str(), histoName.str().c_str(), 20, 0.0, 1.0);
173 
174  histoName.str("");
175  histoName << "ClusterSizeMean_Distribution_Wheel" << w; // Avarage ClusterSize Distribution
176  MEANDWheel[w + 2] = ibooker.book1D(histoName.str().c_str(), histoName.str().c_str(), 100, 0.5, 10.5);
177  }
178  } //end loop on wheels
179 
180  for (int d = -numberOfDisks_; d <= numberOfDisks_; d++) {
181  if (d == 0)
182  continue;
183  //Endcap
184  int offset = numberOfDisks_;
185  if (d > 0)
186  offset--;
187 
188  histoName.str("");
189  histoName << "ClusterSizeIn1Bin_Ring_vs_Segment_Disk"
190  << d; // ClusterSize in first bin norm. by Entries (2D Roll vs Sector)
191  CLSDisk[d + offset] = ibooker.book2D(histoName.str().c_str(),
192  histoName.str().c_str(),
193  36,
194  0.5,
195  36.5,
196  3 * numberOfRings_,
197  0.5,
198  3 * numberOfRings_ + 0.5);
199  rpcUtils.labelXAxisSegment(CLSDisk[d + offset]);
200  rpcUtils.labelYAxisRing(CLSDisk[d + offset], numberOfRings_, useRollInfo_);
201 
202  if (testMode_) {
203  histoName.str("");
204  histoName << "ClusterSizeIn1Bin_Distribution_Disk" << d; // ClusterSize in first bin, distribution
205  CLSDDisk[d + offset] = ibooker.book1D(histoName.str().c_str(), histoName.str().c_str(), 20, 0.0, 1.0);
206 
207  histoName.str("");
208  histoName << "ClusterSizeMean_Distribution_Disk" << d; // Avarage ClusterSize Distribution
209  MEANDDisk[d + offset] = ibooker.book1D(histoName.str().c_str(), histoName.str().c_str(), 100, 0.5, 10.5);
210  }
211 
212  histoName.str("");
213  histoName << "ClusterSizeMean_Ring_vs_Segment_Disk" << d; // Avarage ClusterSize (2D Roll vs Sector)
214  MEANDisk[d + offset] = ibooker.book2D(histoName.str().c_str(),
215  histoName.str().c_str(),
216  36,
217  0.5,
218  36.5,
219  3 * numberOfRings_,
220  0.5,
221  3 * numberOfRings_ + 0.5);
222  rpcUtils.labelXAxisSegment(MEANDisk[d + offset]);
223  rpcUtils.labelYAxisRing(MEANDisk[d + offset], numberOfRings_, useRollInfo_);
224  }
225 }
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX)
Definition: DQMStore.cc:239
T getUntrackedParameter(std::string const &, T const &) const
const double w
Definition: UKUtility.cc:23
MonitorElement * MEANWheel[kWheels]
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:418
void myBooker(DQMStore::IBooker &) override
void labelXAxisSegment(MonitorElement *myMe)
Definition: utils.h:247
void labelYAxisRoll(MonitorElement *myMe, int region, int ring, bool useRollInfo)
Definition: utils.h:265
MonitorElement * MEANDWheel[kWheels]
MonitorElement * CLSDisk[kDisks]
void clientOperation() override
MonitorElement * CLSDWheel[kWheels]
void beginJob(std::string &) override
void Fill(long long x)
void labelYAxisRing(MonitorElement *myMe, int numberOfRings, bool useRollInfo)
Definition: utils.h:291
int roll() const
Definition: RPCDetId.h:92
void getMonitorElements(std::vector< MonitorElement * > &, std::vector< RPCDetId > &, std::string &) override
virtual double getEntries() const
get # of entries
int ring() const
Definition: RPCDetId.h:59
MonitorElement * MEANDisk[kDisks]
MonitorElement * MEANDDisk[kDisks]
std::vector< MonitorElement * > myClusterMe_
virtual double getBinContent(int binx) const
get content of bin (1-D)
void labelXAxisSector(MonitorElement *myMe)
Definition: utils.h:231
d
Definition: ztail.py:151
virtual int segment()
Definition: RPCGeomServ.cc:361
int detId2RollNr(const RPCDetId &_id)
Definition: utils.h:31
~RPCClusterSizeTest() override
Destructor.
virtual double getMean(int axis=1) const
get mean value of histogram along x, y or z axis (axis=1, 2, 3 respectively)
RPCClusterSizeTest(const edm::ParameterSet &ps)
Constructor.
virtual void setBinContent(int binx, double content)
set content of bin (1-D)
std::vector< RPCDetId > myDetIds_
int sector() const
Sector id: the group of chambers at same phi (and increasing r)
Definition: RPCDetId.h:81
MonitorElement * CLSDDisk[kDisks]
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
Definition: DQMStore.cc:266
MonitorElement * CLSWheel[kWheels]
int region() const
Region id: 0 for Barrel, +/-1 For +/- Endcap.
Definition: RPCDetId.h:53
int station() const
Definition: RPCDetId.h:78