CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/src/DQM/SiPixelMonitorRawData/src/SiPixelRawDataErrorSource.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    SiPixelMonitorRawData
00004 // Class:      SiPixelRawDataErrorSource
00005 // 
00021 //
00022 // Original Author:  Andrew York
00023 //
00024 #include "DQM/SiPixelMonitorRawData/interface/SiPixelRawDataErrorSource.h"
00025 // Framework
00026 #include "FWCore/ServiceRegistry/interface/Service.h"
00027 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00028 // DQM Framework
00029 #include "DQM/SiPixelCommon/interface/SiPixelFolderOrganizer.h"
00030 #include "DQMServices/Core/interface/DQMStore.h"
00031 // Geometry
00032 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00033 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00034 #include "Geometry/TrackerGeometryBuilder/interface/PixelGeomDetUnit.h"
00035 #include "Geometry/CommonTopologies/interface/PixelTopology.h"
00036 // DataFormats
00037 #include "DataFormats/DetId/interface/DetId.h"
00038 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
00039 #include "DataFormats/FEDRawData/interface/FEDNumbering.h"
00040 #include "DataFormats/SiPixelDetId/interface/PixelBarrelName.h"
00041 #include "DataFormats/SiPixelDetId/interface/PixelBarrelNameUpgrade.h"
00042 #include "DataFormats/SiPixelDetId/interface/PixelEndcapName.h"
00043 #include "DataFormats/SiPixelDetId/interface/PixelEndcapNameUpgrade.h"
00044 //
00045 #include <string>
00046 #include <stdlib.h>
00047 
00048 using namespace std;
00049 using namespace edm;
00050 
00051 SiPixelRawDataErrorSource::SiPixelRawDataErrorSource(const edm::ParameterSet& iConfig) :
00052   conf_(iConfig),
00053   src_( conf_.getParameter<edm::InputTag>( "src" ) ),
00054   saveFile( conf_.getUntrackedParameter<bool>("saveFile",false) ),
00055   isPIB( conf_.getUntrackedParameter<bool>("isPIB",false) ),
00056   slowDown( conf_.getUntrackedParameter<bool>("slowDown",false) ),
00057   reducedSet( conf_.getUntrackedParameter<bool>("reducedSet",false) ),
00058   modOn( conf_.getUntrackedParameter<bool>("modOn",true) ),
00059   ladOn( conf_.getUntrackedParameter<bool>("ladOn",false) ), 
00060   bladeOn( conf_.getUntrackedParameter<bool>("bladeOn",false) ),
00061   isUpgrade( conf_.getUntrackedParameter<bool>("isUpgrade",false) )
00062 {
00063    theDMBE = edm::Service<DQMStore>().operator->();
00064    LogInfo ("PixelDQM") << "SiPixelRawDataErrorSource::SiPixelRawDataErrorSource: Got DQM BackEnd interface"<<endl;
00065 }
00066 
00067 
00068 SiPixelRawDataErrorSource::~SiPixelRawDataErrorSource()
00069 {
00070    // do anything here that needs to be done at desctruction time
00071    // (e.g. close files, deallocate resources etc.)
00072   LogInfo ("PixelDQM") << "SiPixelRawDataErrorSource::~SiPixelRawDataErrorSource: Destructor"<<endl;
00073 }
00074 
00075 
00076 void SiPixelRawDataErrorSource::beginJob(){
00077   firstRun = true;
00078 }
00079 
00080 void SiPixelRawDataErrorSource::beginRun(const edm::Run& r, const edm::EventSetup& iSetup){
00081 
00082   LogInfo ("PixelDQM") << " SiPixelRawDataErrorSource::beginRun - Initialisation ... " << std::endl;
00083   LogInfo ("PixelDQM") << "Mod/Lad/Blade " << modOn << "/" << ladOn << "/" << bladeOn << std::endl;
00084 
00085   if(firstRun){
00086     eventNo = 0;
00087     // Build map
00088     buildStructure(iSetup);
00089     // Book Monitoring Elements
00090     bookMEs();
00091     
00092     firstRun = false;
00093   }
00094 }
00095 
00096 
00097 void SiPixelRawDataErrorSource::endJob(void){
00098 
00099   if(saveFile) {
00100     LogInfo ("PixelDQM") << " SiPixelRawDataErrorSource::endJob - Saving Root File " << std::endl;
00101     std::string outputFile = conf_.getParameter<std::string>("outputFile");
00102     theDMBE->save( outputFile.c_str() );
00103   }
00104 
00105 }
00106 
00107 //------------------------------------------------------------------
00108 // Method called for every event
00109 //------------------------------------------------------------------
00110 void SiPixelRawDataErrorSource::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00111 {
00112   eventNo++;
00113   //std::cout<<"Event number: "<<eventNo<<std::endl;
00114   // get input data
00115   edm::Handle< edm::DetSetVector<SiPixelRawDataError> >  input;
00116   iEvent.getByLabel( src_, input );
00117   if (!input.isValid()) return; 
00118 
00119   // Get DQM interface
00120   DQMStore* theDMBE = edm::Service<DQMStore>().operator->();
00121    
00122   //float iOrbitSec = iEvent.orbitNumber()/11223.;
00123   //int bx = iEvent.bunchCrossing();
00124   //long long tbx = (long long)iEvent.orbitNumber() * 3564 + bx;
00125   int lumiSection = (int)iEvent.luminosityBlock();
00126   
00127   int nEventBPIXModuleErrors = 0; int nEventFPIXModuleErrors = 0; int nEventBPIXFEDErrors = 0; int nEventFPIXFEDErrors = 0;
00128   int nErrors = 0;
00129 
00130   std::map<uint32_t,SiPixelRawDataErrorModule*>::iterator struct_iter;
00131   std::map<uint32_t,SiPixelRawDataErrorModule*>::iterator struct_iter2;
00132   for (struct_iter = thePixelStructure.begin() ; struct_iter != thePixelStructure.end() ; struct_iter++) {
00133     
00134     int numberOfModuleErrors = (*struct_iter).second->fill(*input, modOn, ladOn, bladeOn);
00135     if(DetId( (*struct_iter).first ).subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel)) nEventBPIXModuleErrors = nEventBPIXModuleErrors + numberOfModuleErrors;
00136     if(DetId( (*struct_iter).first ).subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap)) nEventFPIXModuleErrors = nEventFPIXModuleErrors + numberOfModuleErrors;
00137     //cout<<"NErrors: "<<nEventBPIXModuleErrors<<" "<<nEventFPIXModuleErrors<<endl;
00138     nErrors = nErrors + numberOfModuleErrors;
00139     //if(nErrors>0) cout<<"MODULES: nErrors: "<<nErrors<<endl;
00140   }
00141   for (struct_iter2 = theFEDStructure.begin() ; struct_iter2 != theFEDStructure.end() ; struct_iter2++) {
00142     
00143     int numberOfFEDErrors = (*struct_iter2).second->fillFED(*input);
00144     if((*struct_iter2).first <= 31) nEventBPIXFEDErrors = nEventBPIXFEDErrors + numberOfFEDErrors; // (*struct_iter2).first >= 0, since (*struct_iter2).first is unsigned
00145     if((*struct_iter2).first >= 32 && (*struct_iter2).first <= 39) nEventFPIXFEDErrors = nEventFPIXFEDErrors + numberOfFEDErrors;    
00146     //cout<<"NFEDErrors: "<<nEventBPIXFEDErrors<<" "<<nEventFPIXFEDErrors<<endl;
00147     nErrors = nErrors + numberOfFEDErrors;
00148     //if(nErrors>0) cout<<"FEDS: nErrors: "<<nErrors<<endl;
00149   }
00150   MonitorElement* me = theDMBE->get("Pixel/AdditionalPixelErrors/byLumiErrors");
00151   if(me){
00152     me->setBinContent(0,eventNo);
00153     //cout<<"NErrors: "<<nEventBPIXModuleErrors<<" "<<nEventFPIXModuleErrors<<" "<<nEventBPIXFEDErrors<<" "<<nEventFPIXFEDErrors<<endl;
00154     //if(nEventBPIXModuleErrors+nEventBPIXFEDErrors>0) me->setBinContent(1,me->getBinContent(1)+nEventBPIXModuleErrors+nEventBPIXFEDErrors);
00155     //if(nEventFPIXModuleErrors+nEventFPIXFEDErrors>0) me->setBinContent(2,me->getBinContent(2)+nEventFPIXModuleErrors+nEventFPIXFEDErrors);
00156     if(nEventBPIXModuleErrors+nEventBPIXFEDErrors>0) me->Fill(0,1.);
00157     if(nEventFPIXModuleErrors+nEventFPIXFEDErrors>0) me->Fill(1,1.);
00158     //cout<<"histo: "<<me->getBinContent(0)<<" "<<me->getBinContent(1)<<" "<<me->getBinContent(2)<<endl;
00159   }  
00160   
00161   // Rate of errors per lumi section:
00162   MonitorElement* me1 = theDMBE->get("Pixel/AdditionalPixelErrors/errorRate");
00163   if(me1) me1->Fill(lumiSection, nErrors);
00164 
00165   // slow down...
00166   if(slowDown) usleep(100000);
00167   
00168 }
00169 
00170 //------------------------------------------------------------------
00171 // Build data structure
00172 //------------------------------------------------------------------
00173 void SiPixelRawDataErrorSource::buildStructure(const edm::EventSetup& iSetup){
00174 
00175   LogInfo ("PixelDQM") <<" SiPixelRawDataErrorSource::buildStructure" ;
00176   edm::ESHandle<TrackerGeometry> pDD;
00177   iSetup.get<TrackerDigiGeometryRecord>().get( pDD );
00178 
00179   LogVerbatim ("PixelDQM") << " *** Geometry node for TrackerGeom is  "<<&(*pDD)<<std::endl;
00180   LogVerbatim ("PixelDQM") << " *** I have " << pDD->dets().size() <<" detectors"<<std::endl;
00181   LogVerbatim ("PixelDQM") << " *** I have " << pDD->detTypes().size() <<" types"<<std::endl;
00182 
00183   for(TrackerGeometry::DetContainer::const_iterator it = pDD->dets().begin(); it != pDD->dets().end(); it++){
00184 
00185     if( ((*it)->subDetector()==GeomDetEnumerators::PixelBarrel) || ((*it)->subDetector()==GeomDetEnumerators::PixelEndcap) ){
00186       DetId detId = (*it)->geographicalId();
00187       const GeomDetUnit      * geoUnit = pDD->idToDetUnit( detId );
00188       const PixelGeomDetUnit * pixDet  = dynamic_cast<const PixelGeomDetUnit*>(geoUnit);
00189       int nrows = (pixDet->specificTopology()).nrows();
00190       int ncols = (pixDet->specificTopology()).ncolumns();
00191 
00192       if(detId.subdetId() == static_cast<int>(PixelSubdetector::PixelBarrel)) {
00193         if(isPIB) continue;
00194         LogDebug ("PixelDQM") << " ---> Adding Barrel Module " <<  detId.rawId() << endl;
00195         uint32_t id = detId();
00196         SiPixelRawDataErrorModule* theModule = new SiPixelRawDataErrorModule(id, ncols, nrows);
00197         thePixelStructure.insert(pair<uint32_t,SiPixelRawDataErrorModule*> (id,theModule));
00198 
00199       } else if( (detId.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap)) && (!isUpgrade)) {
00200         LogDebug ("PixelDQM") << " ---> Adding Endcap Module " <<  detId.rawId() << endl;
00201         uint32_t id = detId();
00202         SiPixelRawDataErrorModule* theModule = new SiPixelRawDataErrorModule(id, ncols, nrows);
00203         
00204         PixelEndcapName::HalfCylinder side = PixelEndcapName(DetId(id)).halfCylinder();
00205         int disk   = PixelEndcapName(DetId(id)).diskName();
00206         int blade  = PixelEndcapName(DetId(id)).bladeName();
00207         int panel  = PixelEndcapName(DetId(id)).pannelName();
00208         int module = PixelEndcapName(DetId(id)).plaquetteName();
00209 
00210         char sside[80];  sprintf(sside,  "HalfCylinder_%i",side);
00211         char sdisk[80];  sprintf(sdisk,  "Disk_%i",disk);
00212         char sblade[80]; sprintf(sblade, "Blade_%02i",blade);
00213         char spanel[80]; sprintf(spanel, "Panel_%i",panel);
00214         char smodule[80];sprintf(smodule,"Module_%i",module);
00215         std::string side_str = sside;
00216         std::string disk_str = sdisk;
00217         bool mask = side_str.find("HalfCylinder_1")!=string::npos||
00218                     side_str.find("HalfCylinder_2")!=string::npos||
00219                     side_str.find("HalfCylinder_4")!=string::npos||
00220                     disk_str.find("Disk_2")!=string::npos;
00221         // clutch to take all of FPIX, but no BPIX:
00222         mask = false;
00223         if(isPIB && mask) continue;
00224                 
00225         thePixelStructure.insert(pair<uint32_t,SiPixelRawDataErrorModule*> (id,theModule));
00226       } else if( (detId.subdetId() == static_cast<int>(PixelSubdetector::PixelEndcap)) && (isUpgrade)) {
00227         LogDebug ("PixelDQM") << " ---> Adding Endcap Module " <<  detId.rawId() << endl;
00228         uint32_t id = detId();
00229         SiPixelRawDataErrorModule* theModule = new SiPixelRawDataErrorModule(id, ncols, nrows);
00230         
00231         PixelEndcapNameUpgrade::HalfCylinder side = PixelEndcapNameUpgrade(DetId(id)).halfCylinder();
00232         int disk   = PixelEndcapNameUpgrade(DetId(id)).diskName();
00233         int blade  = PixelEndcapNameUpgrade(DetId(id)).bladeName();
00234         int panel  = PixelEndcapNameUpgrade(DetId(id)).pannelName();
00235         int module = PixelEndcapNameUpgrade(DetId(id)).plaquetteName();
00236 
00237         char sside[80];  sprintf(sside,  "HalfCylinder_%i",side);
00238         char sdisk[80];  sprintf(sdisk,  "Disk_%i",disk);
00239         char sblade[80]; sprintf(sblade, "Blade_%02i",blade);
00240         char spanel[80]; sprintf(spanel, "Panel_%i",panel);
00241         char smodule[80];sprintf(smodule,"Module_%i",module);
00242         std::string side_str = sside;
00243         std::string disk_str = sdisk;
00244         bool mask = side_str.find("HalfCylinder_1")!=string::npos||
00245                     side_str.find("HalfCylinder_2")!=string::npos||
00246                     side_str.find("HalfCylinder_4")!=string::npos||
00247                     disk_str.find("Disk_2")!=string::npos;
00248         // clutch to take all of FPIX, but no BPIX:
00249         mask = false;
00250         if(isPIB && mask) continue;
00251                 
00252         thePixelStructure.insert(pair<uint32_t,SiPixelRawDataErrorModule*> (id,theModule));
00253       }//endif(isUpgrade)
00254     }
00255   }
00256   LogDebug ("PixelDQM") << " ---> Adding Module for Additional Errors " << endl;
00257   pair<int,int> fedIds (FEDNumbering::MINSiPixelFEDID, FEDNumbering::MAXSiPixelFEDID);
00258   fedIds.first = 0;
00259   fedIds.second = 39;
00260   for (int fedId = fedIds.first; fedId <= fedIds.second; fedId++) {
00261     //std::cout<<"Adding FED module: "<<fedId<<std::endl;
00262     uint32_t id = static_cast<uint32_t> (fedId);
00263     SiPixelRawDataErrorModule* theModule = new SiPixelRawDataErrorModule(id);
00264     theFEDStructure.insert(pair<uint32_t,SiPixelRawDataErrorModule*> (id,theModule));
00265   }
00266   
00267   LogInfo ("PixelDQM") << " *** Pixel Structure Size " << thePixelStructure.size() << endl;
00268 }
00269 //------------------------------------------------------------------
00270 // Book MEs
00271 //------------------------------------------------------------------
00272 void SiPixelRawDataErrorSource::bookMEs(){
00273   //cout<<"Entering SiPixelRawDataErrorSource::bookMEs now: "<<endl;
00274   // Get DQM interface
00275   DQMStore* theDMBE = edm::Service<DQMStore>().operator->();
00276   theDMBE->setCurrentFolder("Pixel/AdditionalPixelErrors");
00277   char title[80]; sprintf(title, "By-LumiSection Error counters");
00278   byLumiErrors = theDMBE->book1D("byLumiErrors",title,2,0.,2.);
00279   byLumiErrors->setLumiFlag();
00280   char title1[80]; sprintf(title1, "Errors per LumiSection;LumiSection;NErrors");
00281   errorRate = theDMBE->book1D("errorRate",title1,5000,0.,5000.);
00282   
00283   std::map<uint32_t,SiPixelRawDataErrorModule*>::iterator struct_iter;
00284   std::map<uint32_t,SiPixelRawDataErrorModule*>::iterator struct_iter2;
00285   
00286   SiPixelFolderOrganizer theSiPixelFolder;
00287   
00288   for(struct_iter = thePixelStructure.begin(); struct_iter != thePixelStructure.end(); struct_iter++){
00290 
00291     if(modOn){
00292       if(theSiPixelFolder.setModuleFolder((*struct_iter).first,0,isUpgrade)) {
00293         (*struct_iter).second->book( conf_, 0 ,isUpgrade );
00294       }
00295       else {
00296         //std::cout<<"PIB! not booking histograms for non-PIB modules!"<<std::endl;
00297         if(!isPIB) throw cms::Exception("LogicError")
00298                        << "[SiPixelRawDataErrorSource::bookMEs] Creation of DQM folder failed";
00299       }
00300     }
00301     
00302     if(ladOn){
00303       if(theSiPixelFolder.setModuleFolder((*struct_iter).first,1,isUpgrade)) {
00304         (*struct_iter).second->book( conf_, 1 , isUpgrade );
00305       }
00306       else {
00307         LogDebug ("PixelDQM") << "PROBLEM WITH LADDER-FOLDER\n";
00308       }
00309     }
00310     
00311     if(bladeOn){
00312       if(theSiPixelFolder.setModuleFolder((*struct_iter).first,4,isUpgrade)) {
00313         (*struct_iter).second->book( conf_, 4 ,isUpgrade );
00314       }
00315       else {
00316         LogDebug ("PixelDQM") << "PROBLEM WITH BLADE-FOLDER\n";
00317       }
00318     }
00319     
00320   }//for loop
00321 
00322   for(struct_iter2 = theFEDStructure.begin(); struct_iter2 != theFEDStructure.end(); struct_iter2++){
00324     if(theSiPixelFolder.setFedFolder((*struct_iter2).first)) {
00325       (*struct_iter2).second->bookFED( conf_ );
00326     }
00327     else {
00328       throw cms::Exception("LogicError")
00329         << "[SiPixelRawDataErrorSource::bookMEs] Creation of DQM folder failed";
00330     }
00331 
00332   }
00333   //cout<<"...leaving SiPixelRawDataErrorSource::bookMEs now! "<<endl;
00334 
00335 }
00336 
00337 //define this as a plug-in
00338 DEFINE_FWK_MODULE(SiPixelRawDataErrorSource);