CMS 3D CMS Logo

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", 4);
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
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 
85  rpcdqm::utils rpcUtils;
86 
87  int limit = numberOfDisks_;
88  if(numberOfDisks_ < 2) limit = 2;
89 
90  for (int w = -1 * limit; w<=limit;w++ ){//loop on wheels and disks
91  if (w>-3 && w<3){//wheels
92  histoName.str("");
93  histoName<<"BX_Mean_Distribution_Wheel"<<w;
94  me = nullptr;
95  me = dbe_->get(globalFolder_ + histoName.str()) ;
96  if ( nullptr!=me ) {
97  dbe_->removeElement(me->getName());
98  }
99  BXMeanWheel[w+2] = dbe_->book1D(histoName.str().c_str(), histoName.str().c_str(), 10, -0.5, 9.5);
100 
101  histoName.str("");
102  histoName<<"BX_RMS_Distribution_Wheel"<<w;
103  me = nullptr;
104  me = dbe_->get(globalFolder_ + histoName.str()) ;
105  if ( nullptr!=me){
106  dbe_->removeElement(me->getName());
107  }
108  BXRmsWheel[w+2] = dbe_->book1D(histoName.str().c_str(), histoName.str().c_str(), 50, -0.5, 4.5);
109  }//end loop on wheels
110 
111  if (w == 0 || w< (-1 * numberOfDisks_) || w > numberOfDisks_)continue;
112  //Endcap
113  int offset = numberOfDisks_;
114  if (w>0) offset --; //used to skip case equale to zero
115 
116  histoName.str("");
117  histoName<<"BX_Mean_Distribution_Disk"<<w;
118  me = nullptr;
119  me = dbe_->get(globalFolder_ + histoName.str()) ;
120  if ( nullptr!=me){
121  dbe_->removeElement(me->getName());
122  }
123  BXMeanDisk[w+offset] = dbe_->book1D(histoName.str().c_str(), histoName.str().c_str(),10, -0.5, 9.5);
124 
125 
126  histoName.str("");
127  histoName<<"BX_RMS_Distribution_Disk"<<w;
128  me = nullptr;
129  me = dbe_->get(globalFolder_ + histoName.str()) ;
130  if ( nullptr!=me){
131  dbe_->removeElement(me->getName());
132  }
133  BXRmsDisk[w+offset] = dbe_->book1D(histoName.str().c_str(), histoName.str().c_str(), 50, -0.5, 4.5);
134  }
135 
136 }
137 
138 void RPCBxTest::getMonitorElements(std::vector<MonitorElement *> & meVector, std::vector<RPCDetId> & detIdVector){
139 
140  //Qui prende gli isto del BX per roll definiti nel client
141  //crea due vettore ordinati myBXMe_(istogrammi) e myDetIds_(rpcDetId)
142  //Get ME for each roll
143  for (unsigned int i = 0 ; i<meVector.size(); i++){
144 
145  bool flag= false;
146 
147  DQMNet::TagList tagList;
148  tagList = meVector[i]->getTags();
149  DQMNet::TagList::iterator tagItr = tagList.begin();
150 
151  while (tagItr != tagList.end() && !flag ) {
152  if((*tagItr) == rpcdqm::BX){ flag= true;}
153  tagItr++;
154  }
155 
156  if(flag){
157  myBXMe_.push_back(meVector[i]);
158  myDetIds_.push_back(detIdVector[i]);
159  }
160  }
161 
162 }
163 
164 
166 
167 
168 void RPCBxTest::endJob(void) {
169  edm::LogVerbatim ("rpcbxtest") << "[RPCBxTest]: end job ";
170 }
171 
172 
173 //Loop sul vettore degli istogrammi (myBXMe_) e prendi le info
175 
176  MonitorElement * myMe;
177  RPCDetId detId;
178  TH1F * myTH1F;
179 
180  MonitorElement * ENTRIES =nullptr;
181  MonitorElement * MEAN =nullptr;
182  MonitorElement * MEANRing =nullptr;
183  MonitorElement * RMS =nullptr;
184  MonitorElement * RMSRing =nullptr;
185 
186  for (unsigned int i = 0 ; i<myBXMe_.size();i++){
187 
188  myMe = myBXMe_[i];
189  detId = myDetIds_[i];
190 
191  //Prendi TH1F corrispondente al Monitor Element
192  myTH1F = myMe->getTH1F();
193  // //Spegni Overflow
194 // myTH1F->StatOverflows(false); // per accendere overflow mettere true
195 // //Ricalcola la media e l'RMS //commentare le 3 righe seguenti. Ricorda di ricompilare
196 // myTH1F->GetXaxis()->SetRangeUser(-9.5,9.5);
197 // Double_t stat[4];
198 // myTH1F->GetStats(stat);
199 
200 
201  float mean = myTH1F->GetMean();
202  float rms = myTH1F->GetRMS();
203  float entries = myTH1F->GetEntries();
204 
205  //Get Occupancy ME for roll
206  RPCGeomServ RPCname(detId);
207  // if(rms==0) cout<<RPCname.name()<<endl;
208 
209  if(detId.region()== 0){
210  ENTRIES = BXEntriesBarrel;
211  MEAN = BXMeanBarrel; //nome istogramma definito in beginRun
212  MEANRing = BXMeanWheel[detId.ring()+2];
213  RMS = BXRmsBarrel;
214  RMSRing = BXRmsWheel[detId.ring()+2];
215  }else if(detId.region()==1){
216  ENTRIES = BXEntriesEndcapP;
217  MEAN = BXMeanEndcapP;
218  MEANRing = BXMeanDisk[detId.station()+2];
219  RMS = BXRmsEndcapP;
220  RMSRing = BXRmsDisk[detId.station()+2];
221  }else if(detId.region()==-1){
222  ENTRIES = BXEntriesEndcapN;
223  MEAN = BXMeanEndcapN;
224  MEANRing = BXMeanDisk[3-detId.station()];
225  RMS = BXRmsEndcapN;
226  RMSRing = BXRmsDisk[3-detId.station()];
227  }
228 
229  ENTRIES->Fill(entries);
230 
231  if(entries >= entriesCut_){
232  RMSRing->Fill(rms);
233  RMS->Fill(rms);
234 
235  if(rms <= rmsCut_){
236 
237  //if(mean> distanceMean_ || mean<-distanceMean_ ) cout<<RPCname.name()<<endl;
238 
239  MEAN->Fill(mean);
240  MEANRing->Fill(mean);
241  }
242  }
243 
244  }
245 }
void analyze(const edm::Event &iEvent, const edm::EventSetup &c)
Analyze.
Definition: RPCBxTest.cc:165
T getUntrackedParameter(std::string const &, T const &) const
MonitorElement * BXRmsEndcapP
Definition: RPCBxTest.h:67
double distanceMean_
Definition: RPCBxTest.h:44
const double w
Definition: UKUtility.cc:23
MonitorElement * BXRmsEndcapN
Definition: RPCBxTest.h:66
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:837
TH1F * getTH1F() const
MonitorElement * BXRmsBarrel
Definition: RPCBxTest.h:68
MonitorElement * BXRmsWheel[5]
Definition: RPCBxTest.h:70
MonitorElement * BXMeanEndcapN
Definition: RPCBxTest.h:60
void beginJob(DQMStore *)
BeginJob.
Definition: RPCBxTest.cc:31
std::vector< RPCDetId > myDetIds_
Definition: RPCBxTest.h:54
const std::string & getName() const
get name of ME
void beginRun(const edm::Run &r, const edm::EventSetup &c)
Definition: RPCBxTest.cc:37
MonitorElement * BXEntriesBarrel
Definition: RPCBxTest.h:58
void Fill(long long x)
MonitorElement * BXMeanDisk[10]
Definition: RPCBxTest.h:64
int prescaleFactor_
Definition: RPCBxTest.h:47
int iEvent
Definition: GenABIO.cc:230
virtual void endJob(void)
Definition: RPCBxTest.cc:168
int ring() const
Definition: RPCDetId.h:72
MonitorElement * BXRmsDisk[10]
Definition: RPCBxTest.h:69
void removeElement(const std::string &name)
Definition: DQMStore.cc:3183
virtual void endRun(const edm::Run &r, const edm::EventSetup &c)
Definition: RPCBxTest.cc:174
std::string globalFolder_
Definition: RPCBxTest.h:45
void getMonitorElements(std::vector< MonitorElement * > &, std::vector< RPCDetId > &)
Definition: RPCBxTest.cc:138
MonitorElement * get(const std::string &path) const
get ME from full pathname (e.g. "my/long/dir/my_histo")
Definition: DQMStore.cc:1610
DQMStore * dbe_
Definition: RPCBxTest.h:50
int entriesCut_
Definition: RPCBxTest.h:49
RPCBxTest(const edm::ParameterSet &ps)
Constructor.
Definition: RPCBxTest.cc:13
MonitorElement * BXMeanEndcapP
Definition: RPCBxTest.h:61
int numberOfDisks_
Definition: RPCBxTest.h:46
MonitorElement * BXMeanWheel[5]
Definition: RPCBxTest.h:63
MonitorElement * BXEntriesEndcapN
Definition: RPCBxTest.h:56
MonitorElement * BXEntriesEndcapP
Definition: RPCBxTest.h:57
MonitorElement * BXMeanBarrel
Definition: RPCBxTest.h:62
std::vector< uint32_t > TagList
Definition: DQMNet.h:85
int numberOfRings_
Definition: RPCBxTest.h:46
~RPCBxTest() override
Destructor.
Definition: RPCBxTest.cc:29
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:53
double rmsCut_
Definition: RPCBxTest.h:48
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:545
Definition: Run.h:44
int region() const
Region id: 0 for Barrel, +/-1 For +/- Endcap.
Definition: RPCDetId.h:63
int station() const
Definition: RPCDetId.h:96