CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RPCMonitorSync.cc
Go to the documentation of this file.
1 
12 
13 
14 
15 // RPC Digi
18 
21 
22 #include <vector>
23 #include <ctime>
24 #include <sstream>
25 
26 
27 // Detector overview historam booking routines
28 
30  std::stringstream htit;
31  htit<<"Synchronization offset - "<<title;
32  dbe->setCurrentFolder("Global");
33  MonitorElement *hist = dbe->book2D(name,htit.str(),18,1,7,24,1,13);
34  return hist;
35 }
36 
38  std::stringstream htit;
39  htit<<"Synchronization offset - "<<title;
40  dbe->setCurrentFolder("Global");
41  MonitorElement *hist = dbe->book2D(name,htit.str(),12,1,4,36,1,7);
42  return hist;
43 }
44 
46  std::stringstream htit;
47  htit<<"Synchronization width - "<<title;
48  dbe->setCurrentFolder("Global");
49  MonitorElement *hist = dbe->book2D(name,htit.str(),18,1,7,24,1,13);
50  return hist;
51 }
52 
54  std::stringstream htit;
55  htit<<"Synchronization width - "<<title;
56  dbe->setCurrentFolder("Global");
57  MonitorElement *hist = dbe->book2D(name,htit.str(),12,1,4,36,1,7);
58  return hist;
59 }
60 
61 
62 // Class constructor and destructor
63 
65 {
66  counter=0;
67 
68  nameInLog = pset.getUntrackedParameter<std::string>("moduleLogName", "RPC_DQM");
69 
70  saveRootFile = pset.getUntrackedParameter<bool>("SyncDQMSaveRootFile", false);
71  saveRootFileEventsInterval = pset.getUntrackedParameter<int>("SyncEventsInterval", 10000);
72  RootFileName = pset.getUntrackedParameter<std::string>("RootFileNameSync", "RPCMonitorSync.root");
75  dbe->showDirStructure();
76 }
77 
79 
80 }
81 
82 
83 // Histogram for individual DetUnit - booking routine
84 
85 std::map<std::string, MonitorElement*> RPCMonitorSync::bookDetUnitME(RPCDetId & detId) {
86 
87  std::map<std::string, MonitorElement*> meMap;
88  std::string regionName;
89  std::string ringType;
90 
91  if(detId.region() == 0) {
92  regionName="Barrel";
93  ringType="Wheel";
94  } else {
95  ringType="Disk";
96  if(detId.region() == -1) regionName="Encap-";
97  if(detId.region() == 1) regionName="Encap+";
98  }
99 
100  char folder[120];
101  sprintf(folder,"%s/%s_%d/station_%d/sector_%d",regionName.c_str(),ringType.c_str(),
102  detId.ring(),detId.station(),detId.sector());
103 
104  dbe->setCurrentFolder(folder);
105 
106  // Name components common to current RPDDetId
107  char detUnitLabel[128];
108  char layerLabel[128];
109  sprintf(detUnitLabel ,"%d",detId());
110  sprintf(layerLabel ,"layer%d_subsector%d_roll%d",detId.layer(),detId.subsector(),detId.roll());
111 
112  char meId [128];
113  char meTitle [128];
114 
115  sprintf(meId,"Sync_%s",detUnitLabel);
116  sprintf(meTitle,"Sychronization_for_%s",layerLabel);
117  meMap[meId] = dbe->book1D(meId, meTitle, 9, -4.5, 4.5);
118 
119  return meMap;
120 }
121 
122 
123 // Fill histograms and save data
124 
126 {
127  int station_map[8]={1,2,3,4,5,0,6,0}; // map RPC layer numbers to histogram bins
128 
129  MonitorElement *hOffsetMBp1 = barrelOffsetHist("hOffsetMBp1","Barrell Wheel +1");
130  MonitorElement *hOffsetMBp2 = barrelOffsetHist("hOffsetMBp2","Barrell Wheel +2");
131  MonitorElement *hOffsetMB0 = barrelOffsetHist("hOffsetMB0" ,"Barrell Wheel 0");
132  MonitorElement *hOffsetMBm1 = barrelOffsetHist("hOffsetMBm1","Barrell Wheel -1");
133  MonitorElement *hOffsetMBm2 = barrelOffsetHist("hOffsetMBm2","Barrell Wheel -2");
134  MonitorElement *hOffsetMEp1 = endcapOffsetHist("hOffsetMEp1","Endcap Disk +1");
135  MonitorElement *hOffsetMEp2 = endcapOffsetHist("hOffsetMEp2","Endcap Disk +2");
136  MonitorElement *hOffsetMEp3 = endcapOffsetHist("hOffsetMEp3","Endcap Disk +3");
137  MonitorElement *hOffsetMEm1 = endcapOffsetHist("hOffsetMEm1","Endcap Disk -1");
138  MonitorElement *hOffsetMEm2 = endcapOffsetHist("hOffsetMEm2","Endcap Disk -2");
139  MonitorElement *hOffsetMEm3 = endcapOffsetHist("hOffsetMEm3","Endcap Disk -3");
140 
141  MonitorElement *hWidthMBp1 = barrelWidthHist("hWidthMBp1","Barrell Wheel +1");
142  MonitorElement *hWidthMBp2 = barrelWidthHist("hWidthMBp2","Barrell Wheel +2");
143  MonitorElement *hWidthMB0 = barrelWidthHist("hWidthMB0" ,"Barrell Wheel 0");
144  MonitorElement *hWidthMBm1 = barrelWidthHist("hWidthMBm1","Barrell Wheel -1");
145  MonitorElement *hWidthMBm2 = barrelWidthHist("hWidthMBm2","Barrell Wheel -2");
146  MonitorElement *hWidthMEp1 = endcapWidthHist("hWidthMEp1","Endcap Disk +1");
147  MonitorElement *hWidthMEp2 = endcapWidthHist("hWidthMEp2","Endcap Disk +2");
148  MonitorElement *hWidthMEp3 = endcapWidthHist("hWidthMEp3","Endcap Disk +3");
149  MonitorElement *hWidthMEm1 = endcapWidthHist("hWidthMEm1","Endcap Disk -1");
150  MonitorElement *hWidthMEm2 = endcapWidthHist("hWidthMEm2","Endcap Disk -2");
151  MonitorElement *hWidthMEm3 = endcapWidthHist("hWidthMEm3","Endcap Disk -3");
152 
153  std::map <uint32_t,timing>::const_iterator ci;
154  float offset,width;
155  float xf=0,yf=0;
156 
157  for(ci=synchroMap.begin();ci!=synchroMap.end();ci++){
158  uint32_t id = ci->first;
159  RPCDetId detId(id);
160 
161  offset = ci->second.offset();
162  width = ci->second.width();
163 
164  // RPCDetId *tempDetId=new RPCDetId(ci->first);
165 
166 
167  if( detId.region()==0 ) {
168  // Fill barrel histogram
169  xf=station_map[detId.station()*2+detId.layer()-3]+((float)(detId.roll()-0.5)/3.);
170  yf=detId.sector()+((float)(detId.subsector()-0.5)/2.);
171  if ((detId.sector()==4) && (detId.station()==4))
172  yf=detId.sector()+((float)(detId.subsector()-0.5)/4.);
173 
174  if (detId.ring()==1) hOffsetMBp1->Fill(xf,yf,offset);
175  if (detId.ring()==2) hOffsetMBp2->Fill(xf,yf,offset);
176  if (detId.ring()==-1) hOffsetMBm1->Fill(xf,yf,offset);
177  if (detId.ring()==-2) hOffsetMBm2->Fill(xf,yf,offset);
178  if (detId.ring()==0) hOffsetMB0->Fill(xf,yf,offset);
179  if (detId.ring()==1) hWidthMBp1->Fill(xf,yf,width);
180  if (detId.ring()==2) hWidthMBp2->Fill(xf,yf,width);
181  if (detId.ring()==-1) hWidthMBm1->Fill(xf,yf,width);
182  if (detId.ring()==-2) hWidthMBm2->Fill(xf,yf,width);
183  if (detId.ring()==0) hWidthMB0->Fill(xf,yf,width);
184  } else {
185  // Fill endcap histogram
186  xf=detId.ring() +((float)(detId.roll()-0.5)/4.);
187  yf=detId.sector()+((float)(detId.subsector()-0.5)/6.);
188 
189  if (detId.region()==1) {
190  if (detId.station()==1) hOffsetMEp1->Fill(xf,yf,offset);
191  if (detId.station()==2) hOffsetMEp2->Fill(xf,yf,offset);
192  if (detId.station()==3) hOffsetMEp3->Fill(xf,yf,offset);
193  if (detId.station()==1) hWidthMEp1->Fill(xf,yf,width);
194  if (detId.station()==2) hWidthMEp2->Fill(xf,yf,width);
195  if (detId.station()==3) hWidthMEp3->Fill(xf,yf,width);
196  }
197  if (detId.region()==-1) {
198  if (detId.station()==1) hOffsetMEm1->Fill(xf,yf,offset);
199  if (detId.station()==2) hOffsetMEm2->Fill(xf,yf,offset);
200  if (detId.station()==3) hOffsetMEm3->Fill(xf,yf,offset);
201  if (detId.station()==1) hWidthMEm1->Fill(xf,yf,width);
202  if (detId.station()==2) hWidthMEm2->Fill(xf,yf,width);
203  if (detId.station()==3) hWidthMEm3->Fill(xf,yf,width);
204  }
205  } // end histogram filling
206 
207  } // end loop over synchroMap
208 
209  if(saveRootFile)
211 
212 }
213 
214 
215 // Loop over Digis and read synchronization information
216 
218 
219  timing aTiming;
220  aTiming.inTime = 0;
221  for (int i=0;i<4;i++) {
222  aTiming.early_all[i] = 0;
223  aTiming.late_all[i] = 0;
224  }
225 
226  char detUnitLabel[128];
227  char layerLabel[128];
228  char meId [128];
229 
231 // LogInfo (nameInLog) << " RPC Digis: " << rpcDigis->size() << endl;
232  iEvent.getByType(rpcDigis);
234  for(rpcDigiCI = rpcDigis->begin();rpcDigiCI!=rpcDigis->end();rpcDigiCI++){
235 
236  RPCDetId detId=(*rpcDigiCI).first;
237  uint32_t id=detId();
238 
239  sprintf(detUnitLabel ,"%d",detId());
240  sprintf(layerLabel ,"layer%d_subsector%d_roll%d",detId.layer(),detId.subsector(),detId.roll());
241 
242  std::map<uint32_t, std::map<std::string,MonitorElement*> >::iterator meItr = meCollection.find(id);
243  if (meItr == meCollection.end() || (meCollection.size()==0)) {
244  meCollection[id]=bookDetUnitME(detId);
245  }
246  std::map<std::string, MonitorElement*> meMap=meCollection[id];
247 
248  if(synchroMap.find(id)==synchroMap.end()) {
249  synchroMap[id] = aTiming;
250  }
251  const RPCDigiCollection::Range& range = (*rpcDigiCI).second;
252  int aBX=4;
253  for (RPCDigiCollection::const_iterator digiIt = range.first; digiIt!=range.second; ++digiIt){
254  if( digiIt->bx()<aBX) aBX=digiIt->bx();
255  }
256  if (abs(aBX)<4) {
257  if(aBX<0) synchroMap[id].early_all[-aBX]++;
258  if(aBX==0) synchroMap[id].inTime++;
259  if(aBX>0) synchroMap[id].late_all[aBX]++;
260  }
261  sprintf(meId,"Sync_%s",detUnitLabel);
262  meMap[meId]->Fill(aBX);
263  }
265 }
266 
267 
268 // Process event
269 
271 
272  edm::LogInfo (nameInLog) << "Beginning analyzing event " << counter++;
273 
274  // time_t aTime = iEvent.time().value();
275 
276 
277  readRPCDAQStrips(iEvent);
278 
279  return;
280 }
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
int late_all[4]
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:717
MonitorElement * endcapOffsetHist(std::string name, std::string title)
void readRPCDAQStrips(const edm::Event &iEvent)
virtual void endJob(void)
void save(const std::string &filename, const std::string &path="", const std::string &pattern="", const std::string &rewrite="", SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, const std::string &fileupdate="RECREATE")
Definition: DQMStore.cc:2113
std::map< uint32_t, std::map< std::string, MonitorElement * > > meCollection
#define abs(x)
Definition: mlp_lapack.h:159
DQMStore * dbe
back-end interface
bool getByType(Handle< PROD > &result) const
Definition: Event.h:398
virtual void analyze(const edm::Event &, const edm::EventSetup &)
void Fill(long long x)
int iEvent
Definition: GenABIO.cc:243
int roll() const
Definition: RPCDetId.h:124
MonitorElement * endcapWidthHist(std::string name, std::string title)
int ring() const
Definition: RPCDetId.h:76
std::string nameInLog
RPCMonitorSync(const edm::ParameterSet &)
unsigned int offset(bool)
int early_all[4]
int inTime
int layer() const
Definition: RPCDetId.h:112
int saveRootFileEventsInterval
std::vector< DigiType >::const_iterator const_iterator
std::map< uint32_t, timing > synchroMap
int sector() const
Sector id: the group of chambers at same phi (and increasing r)
Definition: RPCDetId.h:106
std::string RootFileName
int subsector() const
SubSector id : some sectors are divided along the phi direction in subsectors (from 1 to 4 in Barrel...
Definition: RPCDetId.h:118
std::pair< const_iterator, const_iterator > Range
MonitorElement * barrelWidthHist(std::string name, std::string title)
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:845
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:429
MonitorElement * barrelOffsetHist(std::string name, std::string title)
Log messages.
std::map< std::string, MonitorElement * > bookDetUnitME(RPCDetId &detId)
Booking of MonitoringElemnt for one RPCDetId (= roll)
int region() const
Region id: 0 for Barrel, +/-1 For +/- Endcap.
Definition: RPCDetId.h:67
int station() const
Definition: RPCDetId.h:100