CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/DQM/SiStripHistoricInfoClient/plugins/SiStripHistoryDQMService.cc

Go to the documentation of this file.
00001 #include "DQM/SiStripHistoricInfoClient/plugins/SiStripHistoryDQMService.h"
00002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00003 #include "DQMServices/Core/interface/MonitorElement.h"
00004 #include "DQMServices/Diagnostic/interface/HDQMfitUtilities.h"
00005 #include "DataFormats/SiStripDetId/interface/TIBDetId.h"
00006 #include "DataFormats/SiStripDetId/interface/TIDDetId.h"
00007 #include "DataFormats/SiStripDetId/interface/TOBDetId.h"
00008 #include "DataFormats/SiStripDetId/interface/TECDetId.h"
00009 
00010 
00011 SiStripHistoryDQMService::SiStripHistoryDQMService(const edm::ParameterSet& iConfig,const edm::ActivityRegistry& aReg)
00012 : DQMHistoryServiceBase::DQMHistoryServiceBase(iConfig, aReg), iConfig_(iConfig)
00013 {
00014   edm::LogInfo("SiStripHistoryDQMService") <<  "[SiStripHistoryDQMService::SiStripHistoryDQMService]";
00015 }
00016 
00017 
00018 SiStripHistoryDQMService::~SiStripHistoryDQMService() { 
00019   edm::LogInfo("SiStripHistoryDQMService") <<  "[SiStripHistoryDQMService::~SiStripHistoryDQMService]";
00020 }
00021 
00022 
00023 uint32_t SiStripHistoryDQMService::returnDetComponent(const MonitorElement* ME){
00024   LogTrace("SiStripHistoryDQMService") <<  "[SiStripHistoryDQMService::returnDetComponent]";
00025   std::string str=ME->getName();
00026   size_t __key_length__=7;
00027   size_t __detid_length__=9;
00028 
00029   uint32_t layer=0,side=0;
00030 
00031   if(str.find("__det__")!= std::string::npos){
00032     return atoi(str.substr(str.find("__det__")+__key_length__,__detid_length__).c_str());
00033   }
00034   //TIB
00035   else if(str.find("TIB")!= std::string::npos){
00036     if (str.find("layer")!= std::string::npos) 
00037       layer=atoi(str.substr(str.find("layer__")+__key_length__,1).c_str());
00038     return TIBDetId(layer,0,0,0,0,0).rawId();
00039   }
00040   //TOB
00041   else if(str.find("TOB")!= std::string::npos){
00042     if (str.find("layer")!= std::string::npos) 
00043       layer=atoi(str.substr(str.find("layer__")+__key_length__,1).c_str());
00044     return TOBDetId(layer,0,0,0,0).rawId();
00045   }
00046   //TID
00047   else if(str.find("TID")!= std::string::npos){  
00048     if (str.find("side")!= std::string::npos){
00049       side=atoi(str.substr(str.find("_side__")+__key_length__,1).c_str());
00050       if (str.find("wheel")!= std::string::npos){
00051         layer=atoi(str.substr(str.find("wheel__")+__key_length__,1).c_str());
00052       }
00053     }
00054     return TIDDetId(side,layer,0,0,0,0).rawId();
00055   } 
00056   //TEC
00057   else if(str.find("TEC")!= std::string::npos){  
00058     if (str.find("side")!= std::string::npos){
00059       side=atoi(str.substr(str.find("_side__")+__key_length__,1).c_str());
00060       if (str.find("wheel")!= std::string::npos){
00061         layer=atoi(str.substr(str.find("wheel__")+__key_length__,1).c_str());
00062       }
00063     }
00064     return TECDetId(side,layer,0,0,0,0,0).rawId();
00065   } 
00066   else 
00067     return SiStripDetId(DetId::Tracker,0).rawId(); //Full Tracker
00068 }
00069 
00070 //Example on how to define an user function for the statistic extraction
00071 bool SiStripHistoryDQMService::setDBLabelsForUser  (std::string& keyName, std::vector<std::string>& userDBContent, std::string& quantity){
00072   if (quantity == "user_2DYmean") {
00073     userDBContent.push_back(keyName+fSep+std::string("yMean"));
00074     userDBContent.push_back(keyName+fSep+std::string("yError"));
00075   } else {
00076     edm::LogError("SiStripHistoryDQMService") << "ERROR: quantity does not exist in SiStripHistoryDQMService::setDBValuesForUser(): " << quantity;
00077     return false;
00078   }
00079   return true;
00080 }
00081 bool SiStripHistoryDQMService::setDBValuesForUser(std::vector<MonitorElement*>::const_iterator iterMes, HDQMSummary::InputVector& values, std::string& quantity  ){
00082   if (quantity == "user_2DYmean") {
00083     TH2F* Hist = (TH2F*) (*iterMes)->getTH2F();
00084     values.push_back( Hist->GetMean(2) );
00085     values.push_back( Hist->GetRMS(2) );
00086   } else {
00087     edm::LogError("SiStripHistoryDQMService") << "ERROR: quantity does not exist in SiStripHistoryDQMService::setDBValuesForUser(): " << quantity;
00088     return false;
00089   }
00090   return true;
00091 }
00092