CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiPixelRawDataErrorSource.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiPixelMonitorRawData
4 // Class: SiPixelRawDataErrorSource
5 //
21 //
22 // Original Author: Andrew York
23 //
25 // Framework
28 // DQM Framework
31 // Geometry
36 // DataFormats
42 //
43 #include <string>
44 #include <stdlib.h>
45 
46 using namespace std;
47 using namespace edm;
48 
50  conf_(iConfig),
51  src_( consumes<DetSetVector<SiPixelRawDataError> >( conf_.getParameter<edm::InputTag>( "src" ) ) ),
52  saveFile( conf_.getUntrackedParameter<bool>("saveFile",false) ),
53  isPIB( conf_.getUntrackedParameter<bool>("isPIB",false) ),
54  slowDown( conf_.getUntrackedParameter<bool>("slowDown",false) ),
55  reducedSet( conf_.getUntrackedParameter<bool>("reducedSet",false) ),
56  modOn( conf_.getUntrackedParameter<bool>("modOn",true) ),
57  ladOn( conf_.getUntrackedParameter<bool>("ladOn",false) ),
58  bladeOn( conf_.getUntrackedParameter<bool>("bladeOn",false) )
59 {
61  LogInfo ("PixelDQM") << "SiPixelRawDataErrorSource::SiPixelRawDataErrorSource: Got DQM BackEnd interface"<<endl;
62 }
63 
64 
66 {
67  // do anything here that needs to be done at desctruction time
68  // (e.g. close files, deallocate resources etc.)
69  LogInfo ("PixelDQM") << "SiPixelRawDataErrorSource::~SiPixelRawDataErrorSource: Destructor"<<endl;
70 }
71 
72 
74  firstRun = true;
75 }
76 
78 
79  LogInfo ("PixelDQM") << " SiPixelRawDataErrorSource::beginRun - Initialisation ... " << std::endl;
80  LogInfo ("PixelDQM") << "Mod/Lad/Blade " << modOn << "/" << ladOn << "/" << bladeOn << std::endl;
81 
82  if(firstRun){
83  eventNo = 0;
84  // Build map
85  buildStructure(iSetup);
86  // Book Monitoring Elements
87  bookMEs();
88 
89  firstRun = false;
90  }
91 }
92 
93 
95 
96  if(saveFile) {
97  LogInfo ("PixelDQM") << " SiPixelRawDataErrorSource::endJob - Saving Root File " << std::endl;
99  theDMBE->save( outputFile.c_str() );
100  }
101 
102 }
103 
104 //------------------------------------------------------------------
105 // Method called for every event
106 //------------------------------------------------------------------
108 {
109  eventNo++;
110  //std::cout<<"Event number: "<<eventNo<<std::endl;
111  // get input data
113  iEvent.getByToken( src_, input );
114  if (!input.isValid()) return;
115 
116  // Get DQM interface
118 
119  //float iOrbitSec = iEvent.orbitNumber()/11223.;
120  //int bx = iEvent.bunchCrossing();
121  //long long tbx = (long long)iEvent.orbitNumber() * 3564 + bx;
122  int lumiSection = (int)iEvent.luminosityBlock();
123 
124  int nEventBPIXModuleErrors = 0; int nEventFPIXModuleErrors = 0; int nEventBPIXFEDErrors = 0; int nEventFPIXFEDErrors = 0;
125  int nErrors = 0;
126 
127  std::map<uint32_t,SiPixelRawDataErrorModule*>::iterator struct_iter;
128  std::map<uint32_t,SiPixelRawDataErrorModule*>::iterator struct_iter2;
129  for (struct_iter = thePixelStructure.begin() ; struct_iter != thePixelStructure.end() ; struct_iter++) {
130 
131  int numberOfModuleErrors = (*struct_iter).second->fill(*input, modOn, ladOn, bladeOn);
132  if(DetId( (*struct_iter).first ).subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel)) nEventBPIXModuleErrors = nEventBPIXModuleErrors + numberOfModuleErrors;
133  if(DetId( (*struct_iter).first ).subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap)) nEventFPIXModuleErrors = nEventFPIXModuleErrors + numberOfModuleErrors;
134  //cout<<"NErrors: "<<nEventBPIXModuleErrors<<" "<<nEventFPIXModuleErrors<<endl;
135  nErrors = nErrors + numberOfModuleErrors;
136  //if(nErrors>0) cout<<"MODULES: nErrors: "<<nErrors<<endl;
137  }
138  for (struct_iter2 = theFEDStructure.begin() ; struct_iter2 != theFEDStructure.end() ; struct_iter2++) {
139 
140  int numberOfFEDErrors = (*struct_iter2).second->fillFED(*input);
141  if((*struct_iter2).first <= 31) nEventBPIXFEDErrors = nEventBPIXFEDErrors + numberOfFEDErrors; // (*struct_iter2).first >= 0, since (*struct_iter2).first is unsigned
142  if((*struct_iter2).first >= 32 && (*struct_iter2).first <= 39) nEventFPIXFEDErrors = nEventFPIXFEDErrors + numberOfFEDErrors;
143  //cout<<"NFEDErrors: "<<nEventBPIXFEDErrors<<" "<<nEventFPIXFEDErrors<<endl;
144  nErrors = nErrors + numberOfFEDErrors;
145  //if(nErrors>0) cout<<"FEDS: nErrors: "<<nErrors<<endl;
146  }
147  MonitorElement* me = theDMBE->get("Pixel/AdditionalPixelErrors/byLumiErrors");
148  if(me){
149  me->setBinContent(0,eventNo);
150  //cout<<"NErrors: "<<nEventBPIXModuleErrors<<" "<<nEventFPIXModuleErrors<<" "<<nEventBPIXFEDErrors<<" "<<nEventFPIXFEDErrors<<endl;
151  //if(nEventBPIXModuleErrors+nEventBPIXFEDErrors>0) me->setBinContent(1,me->getBinContent(1)+nEventBPIXModuleErrors+nEventBPIXFEDErrors);
152  //if(nEventFPIXModuleErrors+nEventFPIXFEDErrors>0) me->setBinContent(2,me->getBinContent(2)+nEventFPIXModuleErrors+nEventFPIXFEDErrors);
153  if(nEventBPIXModuleErrors+nEventBPIXFEDErrors>0) me->Fill(0,1.);
154  if(nEventFPIXModuleErrors+nEventFPIXFEDErrors>0) me->Fill(1,1.);
155  //cout<<"histo: "<<me->getBinContent(0)<<" "<<me->getBinContent(1)<<" "<<me->getBinContent(2)<<endl;
156  }
157 
158  // Rate of errors per lumi section:
159  MonitorElement* me1 = theDMBE->get("Pixel/AdditionalPixelErrors/errorRate");
160  if(me1) me1->Fill(lumiSection, nErrors);
161 
162  // slow down...
163  if(slowDown) usleep(100000);
164 
165 }
166 
167 //------------------------------------------------------------------
168 // Build data structure
169 //------------------------------------------------------------------
171 
172  LogInfo ("PixelDQM") <<" SiPixelRawDataErrorSource::buildStructure" ;
174  iSetup.get<TrackerDigiGeometryRecord>().get( pDD );
175 
176  LogVerbatim ("PixelDQM") << " *** Geometry node for TrackerGeom is "<<&(*pDD)<<std::endl;
177  LogVerbatim ("PixelDQM") << " *** I have " << pDD->dets().size() <<" detectors"<<std::endl;
178  LogVerbatim ("PixelDQM") << " *** I have " << pDD->detTypes().size() <<" types"<<std::endl;
179 
180  for(TrackerGeometry::DetContainer::const_iterator it = pDD->dets().begin(); it != pDD->dets().end(); it++){
181 
182  if( ((*it)->subDetector()==GeomDetEnumerators::PixelBarrel) || ((*it)->subDetector()==GeomDetEnumerators::PixelEndcap) ){
183  DetId detId = (*it)->geographicalId();
184  const GeomDetUnit * geoUnit = pDD->idToDetUnit( detId );
185  const PixelGeomDetUnit * pixDet = dynamic_cast<const PixelGeomDetUnit*>(geoUnit);
186  int nrows = (pixDet->specificTopology()).nrows();
187  int ncols = (pixDet->specificTopology()).ncolumns();
188 
189  if(detId.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel)) {
190  if(isPIB) continue;
191  LogDebug ("PixelDQM") << " ---> Adding Barrel Module " << detId.rawId() << endl;
192  uint32_t id = detId();
193  SiPixelRawDataErrorModule* theModule = new SiPixelRawDataErrorModule(id, ncols, nrows);
194  thePixelStructure.insert(pair<uint32_t,SiPixelRawDataErrorModule*> (id,theModule));
195 
196  } else if(detId.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap)) {
197  LogDebug ("PixelDQM") << " ---> Adding Endcap Module " << detId.rawId() << endl;
198  uint32_t id = detId();
199  SiPixelRawDataErrorModule* theModule = new SiPixelRawDataErrorModule(id, ncols, nrows);
200 
202  int disk = PixelEndcapName(DetId(id)).diskName();
203  int blade = PixelEndcapName(DetId(id)).bladeName();
204  int panel = PixelEndcapName(DetId(id)).pannelName();
206 
207  char sside[80]; sprintf(sside, "HalfCylinder_%i",side);
208  char sdisk[80]; sprintf(sdisk, "Disk_%i",disk);
209  char sblade[80]; sprintf(sblade, "Blade_%02i",blade);
210  char spanel[80]; sprintf(spanel, "Panel_%i",panel);
211  char smodule[80];sprintf(smodule,"Module_%i",module);
212  std::string side_str = sside;
213  std::string disk_str = sdisk;
214  bool mask = side_str.find("HalfCylinder_1")!=string::npos||
215  side_str.find("HalfCylinder_2")!=string::npos||
216  side_str.find("HalfCylinder_4")!=string::npos||
217  disk_str.find("Disk_2")!=string::npos;
218  // clutch to take all of FPIX, but no BPIX:
219  mask = false;
220  if(isPIB && mask) continue;
221 
222  thePixelStructure.insert(pair<uint32_t,SiPixelRawDataErrorModule*> (id,theModule));
223  }
224  }
225  }
226  LogDebug ("PixelDQM") << " ---> Adding Module for Additional Errors " << endl;
228  fedIds.first = 0;
229  fedIds.second = 39;
230  for (int fedId = fedIds.first; fedId <= fedIds.second; fedId++) {
231  //std::cout<<"Adding FED module: "<<fedId<<std::endl;
232  uint32_t id = static_cast<uint32_t> (fedId);
234  theFEDStructure.insert(pair<uint32_t,SiPixelRawDataErrorModule*> (id,theModule));
235  }
236 
237  LogInfo ("PixelDQM") << " *** Pixel Structure Size " << thePixelStructure.size() << endl;
238 }
239 //------------------------------------------------------------------
240 // Book MEs
241 //------------------------------------------------------------------
243  //cout<<"Entering SiPixelRawDataErrorSource::bookMEs now: "<<endl;
244  // Get DQM interface
246  theDMBE->setCurrentFolder("Pixel/AdditionalPixelErrors");
247  char title[80]; sprintf(title, "By-LumiSection Error counters");
248  byLumiErrors = theDMBE->book1D("byLumiErrors",title,2,0.,2.);
250  char title1[80]; sprintf(title1, "Errors per LumiSection;LumiSection;NErrors");
251  errorRate = theDMBE->book1D("errorRate",title1,5000,0.,5000.);
252 
253  std::map<uint32_t,SiPixelRawDataErrorModule*>::iterator struct_iter;
254  std::map<uint32_t,SiPixelRawDataErrorModule*>::iterator struct_iter2;
255 
256  SiPixelFolderOrganizer theSiPixelFolder;
257 
258  for(struct_iter = thePixelStructure.begin(); struct_iter != thePixelStructure.end(); struct_iter++){
260 
261  if(modOn){
262  if(theSiPixelFolder.setModuleFolder((*struct_iter).first)) {
263  (*struct_iter).second->book( conf_, 0 );
264  }
265  else {
266  //std::cout<<"PIB! not booking histograms for non-PIB modules!"<<std::endl;
267  if(!isPIB) throw cms::Exception("LogicError")
268  << "[SiPixelRawDataErrorSource::bookMEs] Creation of DQM folder failed";
269  }
270  }
271 
272  if(ladOn){
273  if(theSiPixelFolder.setModuleFolder((*struct_iter).first,1)) {
274  (*struct_iter).second->book( conf_, 1 );
275  }
276  else {
277  LogDebug ("PixelDQM") << "PROBLEM WITH LADDER-FOLDER\n";
278  }
279  }
280 
281  if(bladeOn){
282  if(theSiPixelFolder.setModuleFolder((*struct_iter).first,4)) {
283  (*struct_iter).second->book( conf_, 4 );
284  }
285  else {
286  LogDebug ("PixelDQM") << "PROBLEM WITH BLADE-FOLDER\n";
287  }
288  }
289 
290  }//for loop
291 
292  for(struct_iter2 = theFEDStructure.begin(); struct_iter2 != theFEDStructure.end(); struct_iter2++){
294  if(theSiPixelFolder.setFedFolder((*struct_iter2).first)) {
295  (*struct_iter2).second->bookFED( conf_ );
296  }
297  else {
298  throw cms::Exception("LogicError")
299  << "[SiPixelRawDataErrorSource::bookMEs] Creation of DQM folder failed";
300  }
301 
302  }
303  //cout<<"...leaving SiPixelRawDataErrorSource::bookMEs now! "<<endl;
304 
305 }
306 
#define LogDebug(id)
int plaquetteName() const
plaquetteId (in pannel)
T getParameter(std::string const &) const
void setBinContent(int binx, double content)
set content of bin (1-D)
std::map< uint32_t, SiPixelRawDataErrorModule * > theFEDStructure
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:873
edm::EDGetTokenT< edm::DetSetVector< SiPixelRawDataError > > src_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
virtual void analyze(const edm::Event &, const edm::EventSetup &)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual void beginRun(const edm::Run &, edm::EventSetup const &)
edm::LuminosityBlockNumber_t luminosityBlock() const
Definition: EventBase.h:59
static std::string const input
Definition: EdmProvDump.cc:44
void Fill(long long x)
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
int bladeName() const
blade id
bool setFedFolder(const uint32_t FedId)
Set folder name for a FED (used in the case of errors without detId)
int iEvent
Definition: GenABIO.cc:243
std::map< uint32_t, SiPixelRawDataErrorModule * > thePixelStructure
void save(const std::string &filename, const std::string &path="", const std::string &pattern="", const std::string &rewrite="", const uint32_t run=0, SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, const std::string &fileupdate="RECREATE")
Definition: DQMStore.cc:2297
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
bool isValid() const
Definition: HandleBase.h:76
SiPixelRawDataErrorSource(const edm::ParameterSet &conf)
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
Definition: DetId.h:18
const T & get() const
Definition: EventSetup.h:55
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
virtual void buildStructure(edm::EventSetup const &)
int pannelName() const
pannel id
int diskName() const
disk id
volatile std::atomic< bool > shutdown_flag false
void setLumiFlag(void)
this ME is meant to be stored for each luminosity section
HalfCylinder halfCylinder() const
Definition: vlib.h:208
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:585
Pixel error – collection of errors and error information.
Definition: Run.h:41
bool setModuleFolder(const uint32_t &rawdetid=0, int type=0)
Set folder name for a module or plaquette.