CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RPCMultiplicityTest.cc
Go to the documentation of this file.
1 /*
2  * \author Anna Cimmino
3  */
6 
7 // Framework
9 //DQMServices
11 // Geometry
13 
14 #include <sstream>
15 
17  edm::LogVerbatim ("multiplicity") << "[RPCMultiplicityTest]: Constructor";
18  useRollInfo_ = ps.getUntrackedParameter<bool>("UseRollInfo", false);
19  prescaleFactor_ = ps.getUntrackedParameter<int>("DiagnosticPrescale", 1);
20  numberOfDisks_ = ps.getUntrackedParameter<int>("NumberOfEndcapDisks", 3);
21  numberOfRings_ = ps.getUntrackedParameter<int>("NumberOfEndcapRings", 2);
22  testMode_ = ps.getUntrackedParameter<bool>("testMode", false);
23 }
24 
26  dbe_ = 0;
27 }
28 
29 
31  edm::LogVerbatim ("multiplicity") << "[RPCMultiplicityTest]: Begin job";
32 
33  globalFolder_ = workingFolder;
34  dbe_=dbe;
35 }
36 
37 
39 
40  edm::LogVerbatim ("multiplicity") << "[RPCMultiplicityTest]: End run";
41 
42 }
43 
44 void RPCMultiplicityTest::getMonitorElements(std::vector<MonitorElement *> & meVector, std::vector<RPCDetId> & detIdVector){
45 
46  //Get NumberOfDigi ME for each roll
47  for (unsigned int i = 0 ; i<meVector.size(); i++){
48 
49  bool flag= false;
50 
51  DQMNet::TagList tagList;
52  tagList = meVector[i]->getTags();
53  DQMNet::TagList::iterator tagItr = tagList.begin();
54 
55  while (tagItr != tagList.end() && !flag ) {
56  if((*tagItr) == rpcdqm::MULTIPLICITY)
57  flag= true;
58 
59  tagItr++;
60  }
61 
62  if(flag){
63  myNumDigiMe_.push_back(meVector[i]);
64  myDetIds_.push_back(detIdVector[i]);
65  }
66  }
67 }
68 
69 
71 
73 
75 
77 
78  edm::LogVerbatim ("multiplicity") <<"[RPCMultiplicityTest]: Client Operation";
79 
80  //Loop on MEs
81  for (unsigned int i = 0 ; i<myNumDigiMe_.size();i++){
83  }//End loop on MEs
84 }
85 
87 
88  MonitorElement* me=NULL;
90 
91  std::stringstream histoName;
92 
93  rpcdqm::utils rpcUtils;
94 
95  for (int i = -2; i<=2;i++ ){//loop on wheels and disks
96 
97  histoName.str("");
98  histoName<<"NumberOfDigi_Mean_Roll_vs_Sector_Wheel"<<i;
99  me = 0;
100  me = dbe_->get(globalFolder_ +"/"+ histoName.str());
101  if ( 0!=me) {
102  dbe_->removeElement(me->getName());
103  }
104 
105  MULTWheel[i+2] = dbe_->book2D(histoName.str().c_str(), histoName.str().c_str(), 12, 0.5, 12.5, 21, 0.5, 21.5);
106 
107  rpcUtils.labelXAxisSector( MULTWheel[i+2]);
108  rpcUtils.labelYAxisRoll( MULTWheel[i+2], 0, i,useRollInfo_ );
109 
110  if(testMode_){
111  histoName.str("");
112  histoName<<"NumberOfDigi_Mean_Distribution_Wheel"<<i;
113  me = 0;
114  me = dbe_->get(globalFolder_ +"/"+ histoName.str());
115  if ( 0!=me) {
116  dbe_->removeElement(me->getName());
117  }
118 
119  MULTDWheel[i+2] = dbe_->book1D(histoName.str().c_str(), histoName.str().c_str(), 100, 0.5, 50.5);
120  }
121 
122  }//end wheels
123 
124  for(int d = -numberOfDisks_; d<=numberOfDisks_; d++ ){
125  if (d == 0 )continue;
126 
127  int offset = numberOfDisks_;
128  if (d>0) offset --; //used to skip case equale to zero
129 
130  histoName.str("");
131  histoName<<"NumberOfDigi_Mean_Ring_vs_Segment_Disk"<<d;
132  me = 0;
133  me = dbe_->get(globalFolder_ +"/"+ histoName.str());
134  if ( 0!=me) {
135  dbe_->removeElement(me->getName());
136  }
137  MULTDisk[d+offset] = dbe_->book2D(histoName.str().c_str(), histoName.str().c_str(),36, 0.5, 36.5, 3*numberOfRings_, 0.5,3*numberOfRings_+ 0.5);
138  rpcUtils.labelXAxisSegment(MULTDisk[d+offset]);
139  rpcUtils.labelYAxisRing(MULTDisk[d+offset], numberOfRings_,useRollInfo_ );
140 
141  if(testMode_){
142  histoName.str("");
143  histoName<<"NumberOfDigi_Mean_Distribution_Disk"<<d;
144  me = 0;
145  me = dbe_->get(globalFolder_ +"/"+ histoName.str());
146  if ( 0!=me) {
147  dbe_->removeElement(me->getName());
148  }
149 
150  MULTDDisk[d+offset] = dbe_->book1D(histoName.str().c_str(), histoName.str().c_str(), 100, 0.5, 50.5);
151  }
152  }//end loop on wheels and disks
153 
154 
155 }
156 
158 
160 
161  MonitorElement * MULT =NULL;
162  MonitorElement * MULTD = NULL;
163 
164  if (detId.region()==0) {
165  MULT = MULTWheel[detId.ring()+2];
166  if(testMode_) MULTD = MULTDWheel[detId.ring()+2];
167  }else{
168  if(-detId.station() + numberOfDisks_ >= 0 ){
169 
170  if(detId.region()<0){
171  MULT = MULTDisk[-detId.station() + numberOfDisks_];
172  if(testMode_) MULTD = MULTDDisk[-detId.station()+ numberOfDisks_];
173  }else{
174  MULT = MULTDisk[detId.station()+ numberOfDisks_ -1];
175  if(testMode_) MULTD = MULTDDisk[detId.station()+ numberOfDisks_-1];
176  }
177  }
178  }
179 
180 
181 
182  int xBin,yBin;
183  if(detId.region()==0){//Barrel
184  xBin= detId.sector();
185  rpcdqm::utils rollNumber;
186  yBin = rollNumber.detId2RollNr(detId);
187  }else{//Endcap
188  //get segment number
189  RPCGeomServ RPCServ(detId);
190  xBin = RPCServ.segment();
191  (numberOfRings_ == 3 ? yBin= detId.ring()*3-detId.roll()+1 : yBin= (detId.ring()-1)*3-detId.roll()+1);
192  }
193 
194  float mean = myMe->getMean();
195 
196  if(MULT) MULT->setBinContent(xBin,yBin, mean );
197  if(testMode_ && MULTD) MULTD->Fill(mean);
198 
199 
200 
201 }
T getUntrackedParameter(std::string const &, T const &) const
const std::string & getName(void) const
get name of ME
int i
Definition: DBlmapReader.cc:9
void setBinContent(int binx, double content)
set content of bin (1-D)
MonitorElement * MULTDDisk[10]
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:872
tuple yBin
Definition: cuy.py:891
MonitorElement * MULTDisk[10]
void labelXAxisSegment(MonitorElement *myMe)
Definition: utils.h:251
void labelYAxisRoll(MonitorElement *myMe, int region, int ring, bool useRollInfo)
Definition: utils.h:269
#define NULL
Definition: scimark2.h:8
RPCMultiplicityTest(const edm::ParameterSet &ps)
Constructor.
double getMean(int axis=1) const
get mean value of histogram along x, y or z axis (axis=1, 2, 3 respectively)
void beginJob(DQMStore *, std::string)
BeginJob.
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:243
void labelYAxisRing(MonitorElement *myMe, int numberOfRings, bool useRollInfo)
Definition: utils.h:296
void clientOperation(edm::EventSetup const &c)
int roll() const
Definition: RPCDetId.h:120
void beginLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &)
Begin Lumi block.
int ring() const
Definition: RPCDetId.h:72
void removeElement(const std::string &name)
Definition: DQMStore.cc:2772
void fillGlobalME(RPCDetId &detId, MonitorElement *myMe)
MonitorElement * MULTDWheel[5]
MonitorElement * MULTWheel[5]
virtual ~RPCMultiplicityTest()
Destructor.
void beginRun(const edm::Run &, const edm::EventSetup &)
unsigned int offset(bool)
MonitorElement * get(const std::string &path) const
get ME from full pathname (e.g. &quot;my/long/dir/my_histo&quot;)
Definition: DQMStore.cc:1623
void labelXAxisSector(MonitorElement *myMe)
Definition: utils.h:236
virtual int segment()
Definition: RPCGeomServ.cc:467
int detId2RollNr(const RPCDetId &_id)
Definition: utils.h:18
void endRun(const edm::Run &, const edm::EventSetup &)
std::vector< uint32_t > TagList
Definition: DQMNet.h:83
int sector() const
Sector id: the group of chambers at same phi (and increasing r)
Definition: RPCDetId.h:102
std::vector< RPCDetId > myDetIds_
void getMonitorElements(std::vector< MonitorElement * > &, std::vector< RPCDetId > &)
MonitorElement * book2D(const char *name, const char *title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
Book 2D histogram.
Definition: DQMStore.cc:1000
void analyze(const edm::Event &, const edm::EventSetup &)
Analyze.
void endLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &)
End Lumi Block.
std::vector< MonitorElement * > myNumDigiMe_
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:584
Definition: Run.h:41
int region() const
Region id: 0 for Barrel, +/-1 For +/- Endcap.
Definition: RPCDetId.h:63
int station() const
Definition: RPCDetId.h:96