CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RPCBxTest.cc
Go to the documentation of this file.
3 
4 // Framework
6 
7 // //Geometry
9 # include "TH1F.h"
10 
11 
12 
14  edm::LogVerbatim ("rpcbxtest") << "[RPCBxTest]: Constructor";
15 
16  prescaleFactor_ = ps.getUntrackedParameter<int>("DiagnosticPrescale", 1);
17 
18  //Nome della dir per gli istogrammi nuovi . Cominciare sempre con RPC/RecHits/
19  globalFolder_ = ps.getUntrackedParameter<std::string>("RPCGlobalFolder", "RPC/RecHits/SummaryHistograms/");
20 
21  entriesCut_ = ps.getUntrackedParameter<int>("EntriesCut");
22  rmsCut_ = ps.getUntrackedParameter<double>("RMSCut");
23  distanceMean_ = ps.getUntrackedParameter<double>("DistanceFromZeroBx");
24 
25  numberOfDisks_ = ps.getUntrackedParameter<int>("NumberOfEndcapDisks", 3);
26  numberOfRings_ = ps.getUntrackedParameter<int>("NumberOfEndcapRings", 2);
27 }
28 
30 
32  edm::LogVerbatim ("rpcbxtest") << "[RPCBxTest]: Begin job ";
33  dbe_ = dbe;
34 }
35 
36 //Qui puoi definitre gli istogrammi nuovi che vuoi riempire
37 void RPCBxTest::beginRun(const edm::Run& r, const edm::EventSetup& c, const std::vector<MonitorElement *>& meVector, const std::vector<RPCDetId>& detIdVector){
38  edm::LogVerbatim ("rpcbxtest") << "[RPCBxTest]: Begin run";
39 
40  MonitorElement* me;
42 
43  std::stringstream histoName;
44 
45  histoName.str("");
46  histoName<<"BX_Mean_Distribution_Barrel";
47  BXMeanBarrel = dbe_->book1D(histoName.str().c_str(), histoName.str().c_str(), 11, -5.5, 5.5);
48  BXMeanBarrel->setAxisTitle("Bx",1);
49 
50  histoName.str("");
51  histoName<<"BX_Mean_Distribution_EndcapP";
52  BXMeanEndcapP = dbe_->book1D(histoName.str().c_str(), histoName.str().c_str(), 11, -5.5, 5.5);
53 
54  histoName.str("");
55  histoName<<"BX_Mean_Distribution_EndcapN";
56  BXMeanEndcapN = dbe_->book1D(histoName.str().c_str(), histoName.str().c_str(), 11, -5.5, 5.5);
57 
58 
59  histoName.str("");
60  histoName<<"BX_Entries_Distribution_Barrel";
61  BXEntriesBarrel = dbe_->book1D(histoName.str().c_str(), histoName.str().c_str(), 1000, -0.5, 999.5);
62 
63  histoName.str("");
64  histoName<<"BX_Entries_Distribution_EndcapP";
65  BXEntriesEndcapP = dbe_->book1D(histoName.str().c_str(), histoName.str().c_str(),1000, -0.5, 999.5);
66 
67  histoName.str("");
68  histoName<<"BX_Entries_Distribution_EndcapN";
69  BXEntriesEndcapN = dbe_->book1D(histoName.str().c_str(), histoName.str().c_str(), 1000, -0.5, 999.5);
70 
71 
72  histoName.str("");
73  histoName<<"BX_RMS_Distribution_Barrel";
74  BXRmsBarrel = dbe_->book1D(histoName.str().c_str(), histoName.str().c_str(), 21, -0.1, 4.1);
75 
76  histoName.str("");
77  histoName<<"BX_RMS_Distribution_EndcapP";
78  BXRmsEndcapP = dbe_->book1D(histoName.str().c_str(), histoName.str().c_str(), 21, -0.1, 4.1);
79 
80  histoName.str("");
81  histoName<<"BX_RMS_Distribution_EndcapN";
82  BXRmsEndcapN = dbe_->book1D(histoName.str().c_str(), histoName.str().c_str(), 21, -0.1, 4.1);
83 
84  rpcdqm::utils rpcUtils;
85 
86  int limit = numberOfDisks_;
87  if(numberOfDisks_ < 2) limit = 2;
88 
89  for (int w = -1 * limit; w<=limit;w++ ){//loop on wheels and disks
90  if (w>-3 && w<3){//wheels
91  histoName.str("");
92  histoName<<"BX_Mean_Distribution_Wheel"<<w;
93  me = 0;
94  me = dbe_->get(globalFolder_ + histoName.str()) ;
95  if ( 0!=me ) {
96  dbe_->removeElement(me->getName());
97  }
98  BXMeanWheel[w+2] = dbe_->book1D(histoName.str().c_str(), histoName.str().c_str(), 10, -0.5, 9.5);
99 
100  histoName.str("");
101  histoName<<"BX_RMS_Distribution_Wheel"<<w;
102  me = 0;
103  me = dbe_->get(globalFolder_ + histoName.str()) ;
104  if ( 0!=me){
105  dbe_->removeElement(me->getName());
106  }
107  BXRmsWheel[w+2] = dbe_->book1D(histoName.str().c_str(), histoName.str().c_str(), 50, -0.5, 4.5);
108  }//end loop on wheels
109 
110  if (w == 0 || w< (-1 * numberOfDisks_) || w > numberOfDisks_)continue;
111  //Endcap
112  int offset = numberOfDisks_;
113  if (w>0) offset --; //used to skip case equale to zero
114 
115  histoName.str("");
116  histoName<<"BX_Mean_Distribution_Disk"<<w;
117  me = 0;
118  me = dbe_->get(globalFolder_ + histoName.str()) ;
119  if ( 0!=me){
120  dbe_->removeElement(me->getName());
121  }
122  BXMeanDisk[w+offset] = dbe_->book1D(histoName.str().c_str(), histoName.str().c_str(),10, -0.5, 9.5);
123 
124 
125  histoName.str("");
126  histoName<<"BX_RMS_Distribution_Disk"<<w;
127  me = 0;
128  me = dbe_->get(globalFolder_ + histoName.str()) ;
129  if ( 0!=me){
130  dbe_->removeElement(me->getName());
131  }
132  BXRmsDisk[w+offset] = dbe_->book1D(histoName.str().c_str(), histoName.str().c_str(), 50, -0.5, 4.5);
133  }
134 
135 
136  //Qui prende gli isto del BX per roll definiti nel client
137  //crea due vettore ordinati myBXMe_(istogrammi) e myDetIds_(rpcDetId)
138  //Get ME for each roll
139  for (unsigned int i = 0 ; i<meVector.size(); i++){
140 
141  bool flag= false;
142 
143  DQMNet::TagList tagList;
144  tagList = meVector[i]->getTags();
145  DQMNet::TagList::iterator tagItr = tagList.begin();
146 
147  while (tagItr != tagList.end() && !flag ) {
148  if((*tagItr) == rpcdqm::OCCUPANCY)
149  flag= true;
150  tagItr++;
151  }
152 
153  if(flag){
154  myBXMe_.push_back(meVector[i]);
155  myDetIds_.push_back(detIdVector[i]);
156  }
157  }
158 
159 }
160 
162 }
163 
165 
167 
168 void RPCBxTest::endJob(void) {
169  edm::LogVerbatim ("rpcbxtest") << "[RPCBxTest]: end job ";
170 
171 }
172 
173 
174 //ME = monitor Element
175 
176 
177 //Link a public getters
178 // http://cmssw.cvs.cern.ch/cgi-bin/cmssw.cgi/CMSSW/DQMServices/Core/interface/MonitorElement.h?revision=1.39&view=markup
179 
180 //Link per ROOT
181 //http://root.cern.ch/root/html522/TH1.html
182 
183 //Loop sul vettore degli istogrammi (myBXMe_) e prendi le info
185 
186  MonitorElement * myMe;
187  RPCDetId detId;
188  TH1F * myTH1F;
189 
190  MonitorElement * ENTRIES =NULL;
191  MonitorElement * MEAN =NULL;
192  MonitorElement * MEANRing =NULL;
193  MonitorElement * RMS =NULL;
194  MonitorElement * RMSRing =NULL;
195 
196  for (unsigned int i = 0 ; i<myBXMe_.size();i++){
197 
198  myMe = myBXMe_[i];
199  detId = myDetIds_[i];
200 
201  //Prendi TH1F corrispondente al Monitor Element
202  myTH1F = myMe->getTH1F();
203  // //Spegni Overflow
204 // myTH1F->StatOverflows(false); // per accendere overflow mettere true
205 // //Ricalcola la media e l'RMS //commentare le 3 righe seguenti. Ricorda di ricompilare
206 // myTH1F->GetXaxis()->SetRangeUser(-9.5,9.5);
207 // Double_t stat[4];
208 // myTH1F->GetStats(stat);
209 
210 
211  float mean = myTH1F->GetMean();
212  float rms = myTH1F->GetRMS();
213  float entries = myTH1F->GetEntries();
214 
215  //Get Occupancy ME for roll
216  RPCGeomServ RPCname(detId);
217  // if(rms==0) cout<<RPCname.name()<<endl;
218 
219  if(detId.region()== 0){
220  ENTRIES = BXEntriesBarrel;
221  MEAN = BXMeanBarrel; //nome istogramma definito in beginRun
222  MEANRing = BXMeanWheel[detId.ring()+2];
223  RMS = BXRmsBarrel;
224  RMSRing = BXRmsWheel[detId.ring()+2];
225  }else if(detId.region()==1){
226  ENTRIES = BXEntriesEndcapP;
227  MEAN = BXMeanEndcapP;
228  MEANRing = BXMeanDisk[detId.station()+2];
229  RMS = BXRmsEndcapP;
230  RMSRing = BXRmsDisk[detId.station()+2];
231  }else if(detId.region()==-1){
232  ENTRIES = BXEntriesEndcapN;
233  MEAN = BXMeanEndcapN;
234  MEANRing = BXMeanDisk[3-detId.station()];
235  RMS = BXRmsEndcapN;
236  RMSRing = BXRmsDisk[3-detId.station()];
237  }
238 
239  ENTRIES->Fill(entries);
240 
241  if(entries >= entriesCut_){
242  RMSRing->Fill(rms);
243  RMS->Fill(rms);
244 
245  if(rms <= rmsCut_){
246 
247  //if(mean> distanceMean_ || mean<-distanceMean_ ) cout<<RPCname.name()<<endl;
248 
249  MEAN->Fill(mean);
250  MEANRing->Fill(mean);
251  }
252  }
253 
254  }
255 }
void analyze(const edm::Event &iEvent, const edm::EventSetup &c)
Analyze.
Definition: RPCBxTest.cc:164
T getUntrackedParameter(std::string const &, T const &) const
const std::string & getName(void) const
get name of ME
int i
Definition: DBlmapReader.cc:9
MonitorElement * BXRmsEndcapP
Definition: RPCBxTest.h:66
double distanceMean_
Definition: RPCBxTest.h:43
MonitorElement * BXRmsEndcapN
Definition: RPCBxTest.h:65
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:873
MonitorElement * BXRmsBarrel
Definition: RPCBxTest.h:67
MonitorElement * BXRmsWheel[5]
Definition: RPCBxTest.h:69
MonitorElement * BXMeanEndcapN
Definition: RPCBxTest.h:59
void beginJob(DQMStore *)
BeginJob.
Definition: RPCBxTest.cc:31
std::vector< RPCDetId > myDetIds_
Definition: RPCBxTest.h:53
#define NULL
Definition: scimark2.h:8
MonitorElement * BXEntriesBarrel
Definition: RPCBxTest.h:57
void Fill(long long x)
MonitorElement * BXMeanDisk[10]
Definition: RPCBxTest.h:63
int prescaleFactor_
Definition: RPCBxTest.h:46
int iEvent
Definition: GenABIO.cc:243
virtual ~RPCBxTest()
Destructor.
Definition: RPCBxTest.cc:29
virtual void endJob(void)
Definition: RPCBxTest.cc:168
void beginLuminosityBlock(edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &context)
Begin Lumi block.
Definition: RPCBxTest.cc:161
DQMNet::TagList getTags(void) const
int ring() const
Definition: RPCDetId.h:72
MonitorElement * BXRmsDisk[10]
Definition: RPCBxTest.h:68
void removeElement(const std::string &name)
Definition: DQMStore.cc:2773
virtual void endRun(const edm::Run &r, const edm::EventSetup &c)
Definition: RPCBxTest.cc:184
std::string globalFolder_
Definition: RPCBxTest.h:44
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:1624
DQMStore * dbe_
Definition: RPCBxTest.h:49
std::vector< uint32_t > TagList
Definition: DQMNet.h:83
int entriesCut_
Definition: RPCBxTest.h:48
RPCBxTest(const edm::ParameterSet &ps)
Constructor.
Definition: RPCBxTest.cc:13
TH1F * getTH1F(void) const
void endLuminosityBlock(edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &c)
End Lumi Block.
Definition: RPCBxTest.cc:166
MonitorElement * BXMeanEndcapP
Definition: RPCBxTest.h:60
int numberOfDisks_
Definition: RPCBxTest.h:45
void beginRun(const edm::Run &r, const edm::EventSetup &c, const std::vector< MonitorElement * > &, const std::vector< RPCDetId > &)
Definition: RPCBxTest.cc:37
MonitorElement * BXMeanWheel[5]
Definition: RPCBxTest.h:62
MonitorElement * BXEntriesEndcapN
Definition: RPCBxTest.h:55
MonitorElement * BXEntriesEndcapP
Definition: RPCBxTest.h:56
MonitorElement * BXMeanBarrel
Definition: RPCBxTest.h:61
T w() const
int numberOfRings_
Definition: RPCBxTest.h:45
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
std::vector< MonitorElement * > myBXMe_
Definition: RPCBxTest.h:52
double rmsCut_
Definition: RPCBxTest.h:47
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:585
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