CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiPixelHistoryDQMService.cc
Go to the documentation of this file.
6 
7 
9 : DQMHistoryServiceBase::DQMHistoryServiceBase(iConfig, aReg), iConfig_(iConfig)
10 {
11  //setSeperator("@@#@@"); // Change the seperator used in DB
12  edm::LogInfo("SiPixelHistoryDQMService") << "[SiPixelHistoryDQMService::SiPixelHistoryDQMService]";
13 }
14 
15 
17  edm::LogInfo("SiPixelHistoryDQMService") << "[SiPixelHistoryDQMService::~SiPixelHistoryDQMService]";
18 }
19 
20 
22  LogTrace("SiPixelHistoryDQMService") << "[SiPixelHistoryDQMService::returnDetComponent]";
23  std::string str=ME->getName();
24  size_t __key_length__=7;
25  size_t __detid_length__=9;
26 
27 
28  if(str.find("__det__")!= std::string::npos){
29  return atoi(str.substr(str.find("__det__")+__key_length__,__detid_length__).c_str());
30  }
31  else if(str.find("Barrel")!= std::string::npos)
32  {return sipixelsummary::Barrel;}
33  else if(str.find("Shell_mI")!= std::string::npos)
34  {return sipixelsummary::Shell_mI;}
35  else if(str.find("Shell_mO")!= std::string::npos)
36  {return sipixelsummary::Shell_mO;}
37  else if(str.find("Shell_pI")!= std::string::npos)
38  {return sipixelsummary::Shell_pI;}
39  else if(str.find("Shell_pO")!= std::string::npos)
40  {return sipixelsummary::Shell_pO;}
41  else if(str.find("Endcap")!= std::string::npos)
42  {return sipixelsummary::Endcap;}
43  else if(str.find("HalfCylinder_mI")!= std::string::npos)
45  else if(str.find("HalfCylinder_mO")!= std::string::npos)
47  else if(str.find("HalfCylinder_pI")!= std::string::npos)
49  else if(str.find("HalfCylinder_pO")!= std::string::npos)
51  else {
52  return sipixelsummary::TRACKER; //Full Tracker
53  }
54 
55 }
56 
57 
58 
59 
60 
61 //Example on how to define an user function for the statistic extraction
62 bool SiPixelHistoryDQMService::setDBLabelsForUser (std::string& keyName, std::vector<std::string>& userDBContent, std::string& quantity){
63  if (quantity == "user_ymean") {
64  userDBContent.push_back(keyName+fSep+std::string("yMean"));
65  userDBContent.push_back(keyName+fSep+std::string("yError"));
66  } else if (quantity == "user_A") {
67  userDBContent.push_back(keyName+fSep+std::string("NTracksPixOverAll"));
68  userDBContent.push_back(keyName+fSep+std::string("NTracksPixOverAllError"));
69  } else if (quantity == "user_B") {
70  userDBContent.push_back(keyName+fSep+std::string("NTracksFPixOverBPix"));
71  userDBContent.push_back(keyName+fSep+std::string("NTracksFPixOverBPixError"));
72  } else {
73  edm::LogError("SiPixelHistoryDQMService") << "ERROR: quantity does not exist in SiPixelHistoryDQMService::setDBValuesForUser(): " << quantity;
74  return false;
75  }
76  return true;
77 }
78 bool SiPixelHistoryDQMService::setDBValuesForUser(std::vector<MonitorElement*>::const_iterator iterMes, HDQMSummary::InputVector& values, std::string& quantity ){
79 
80  if (quantity == "user_ymean") {
81  TH1F* Hist = (TH1F*) (*iterMes)->getTH1F()->Clone();
82  // if( Hist == 0 || Hist->Integral() == 0 ) {
83  // std::cout << "Error: histogram not found or empty!!" << std::endl;
84  // values.push_back( 0. );
85  // values.push_back( 0. );
86  // }
87  // else
88  // if( ) {
89  Hist->Fit("pol0");
90  TF1* Fit = Hist->GetFunction("pol0");
91  float FitValue = Fit ? Fit->GetParameter(0) : 0;
92  float FitError = Fit ? Fit->GetParError(0) : 0;
93  std::cout << "FITERROR: " << FitError << std::endl;
94 
95  values.push_back( FitValue );
96  values.push_back( FitError );
97  // }
98  } else if (quantity == "user_A") {
99  TH1F* Hist = (TH1F*) (*iterMes)->getTH1F();
100  if( Hist->GetBinContent(1) != 0 && Hist->GetBinContent(2) != 0 ) {
101  values.push_back( Hist->GetBinContent(2) / Hist->GetBinContent(1) );
102  values.push_back( TMath::Abs(Hist->GetBinContent(2) / Hist->GetBinContent(1)) * TMath::Sqrt( ( TMath::Power( Hist->GetBinError(1)/Hist->GetBinContent(1), 2) + TMath::Power( Hist->GetBinError(2)/Hist->GetBinContent(2), 2) )) );
103  }
104  else {
105  values.push_back( 0. );
106  values.push_back( 0. );
107  }
108  } else if (quantity == "user_B") {
109  TH1F* Hist = (TH1F*) (*iterMes)->getTH1F();
110  if( Hist->GetBinContent(3) != 0 && Hist->GetBinContent(4) != 0 ) {
111  values.push_back( Hist->GetBinContent(4) / Hist->GetBinContent(3) );
112  values.push_back( TMath::Abs(Hist->GetBinContent(4) / Hist->GetBinContent(3)) * TMath::Sqrt( ( TMath::Power( Hist->GetBinError(3)/Hist->GetBinContent(3), 2) + TMath::Power( Hist->GetBinError(4)/Hist->GetBinContent(4), 2) )) );
113  }
114  else {
115  values.push_back( 0. );
116  values.push_back( 0. );
117  }
118  } else {
119  edm::LogError("SiPixelHistoryDQMService") << "ERROR: quantity does not exist in SiPixelHistoryDQMService::setDBValuesForUser(): " << quantity;
120  return false;
121  }
122 
123  return true;
124 }
125 
const std::string & getName(void) const
get name of ME
SiPixelHistoryDQMService(const edm::ParameterSet &, const edm::ActivityRegistry &)
Definition: ME.h:11
Definition: Fit.h:34
std::vector< float > InputVector
Definition: HDQMSummary.h:58
#define LogTrace(id)
uint32_t returnDetComponent(const MonitorElement *ME)
bool setDBLabelsForUser(std::string &keyName, std::vector< std::string > &userDBContent, std::string &quantity)
bool setDBValuesForUser(std::vector< MonitorElement * >::const_iterator iterMes, HDQMSummary::InputVector &values, std::string &quantity)
tuple cout
Definition: gather_cfg.py:41