CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiStripPedestalsDQM.cc
Go to the documentation of this file.
3 #include "TCanvas.h"
4 
5 // -----
7  edm::ParameterSet const& hPSet,
8  edm::ParameterSet const& fPSet):SiStripBaseCondObjDQM(eSetup, hPSet, fPSet){
9 
10  // Build the Histo_TkMap:
11  if(HistoMaps_On_ ) Tk_HM_ = new TkHistoMap("SiStrip/Histo_Map","MeanPed_TkMap",0.);
12 }
13 // -----
14 
15 
16 
17 // -----
19 // -----
20 
21 
22 // -----
24 
25  getConditionObject(eSetup);
26  pedestalHandle_->getDetIds(activeDetIds);
27 
28 }
29 // -----
30 
31 
32 // -----
33 void SiStripPedestalsDQM::fillModMEs(const std::vector<uint32_t> & selectedDetIds){
34 
35  ModMEs CondObj_ME;
36 
37 
38  for(std::vector<uint32_t>::const_iterator detIter_ = selectedDetIds.begin();
39  detIter_!= selectedDetIds.end();detIter_++){
40 
41  fillMEsForDet(CondObj_ME,*detIter_);
42 
43  }
44 }
45 // -----
46 
47 
48 
49 
50 // -----
51 void SiStripPedestalsDQM::fillMEsForDet(ModMEs selModME_, uint32_t selDetId_){
52 
53  getModMEs(selModME_,selDetId_);
54 
55  SiStripPedestals::Range pedRange = pedestalHandle_->getRange(selDetId_);
56  int nStrip = reader->getNumberOfApvsAndStripLength(selDetId_).first*128;
57 
58  for( int istrip=0;istrip<nStrip;++istrip){
59  if( CondObj_fillId_ =="onlyProfile" || CondObj_fillId_ =="ProfileAndCumul"){
60  selModME_.ProfileDistr->Fill(istrip+1,pedestalHandle_->getPed(istrip,pedRange));
61  }
62  }// istrip
63 
64 }
65 // -----
66 
67 
68 
69 // -----
70 void SiStripPedestalsDQM::fillSummaryMEs(const std::vector<uint32_t> & selectedDetIds){
71 
72  for(std::vector<uint32_t>::const_iterator detIter_ = selectedDetIds.begin();
73  detIter_!= selectedDetIds.end();detIter_++){
74  fillMEsForLayer(/*SummaryMEsMap_,*/ *detIter_);
75  }
76 
77  for (std::map<uint32_t, ModMEs>::iterator iter=SummaryMEsMap_.begin(); iter!=SummaryMEsMap_.end(); iter++){
78 
79  ModMEs selME;
80  selME = iter->second;
81 
82  if(hPSet_.getParameter<bool>("FillSummaryProfileAtLayerLevel") && fPSet_.getParameter<bool>("OutputSummaryProfileAtLayerLevelAsImage")){
83 
84  if( CondObj_fillId_ =="onlyProfile" || CondObj_fillId_ =="ProfileAndCumul"){
85 
86  TCanvas c1("c1");
87  selME.SummaryOfProfileDistr->getTProfile()->Draw();
88  std::string name (selME.SummaryOfProfileDistr->getTProfile()->GetTitle());
89  name+=".png";
90  c1.Print(name.c_str());
91  }
92  }
93  if(hPSet_.getParameter<bool>("FillSummaryAtLayerLevel") && fPSet_.getParameter<bool>("OutputSummaryAtLayerLevelAsImage")){
94 
95  TCanvas c1("c1");
96  selME.SummaryDistr->getTH1()->Draw();
97  std::string name (selME.SummaryDistr->getTH1()->GetTitle());
98  name+=".png";
99  c1.Print(name.c_str());
100  }
101 
102  }
103 
104 }
105 
106 // -----
107 
108 
109 
110 // -----
111 void SiStripPedestalsDQM::fillMEsForLayer( /*std::map<uint32_t, ModMEs> selMEsMap_,*/ uint32_t selDetId_){
112 
113  // ----
114  int subdetectorId_ = ((selDetId_>>25)&0x7);
115 
116  if( subdetectorId_<3 || subdetectorId_>6 ){
117  edm::LogError("SiStripPedestalsDQM")
118  << "[SiStripPedestalsDQM::fillMEsForLayer] WRONG INPUT : no such subdetector type : "
119  << subdetectorId_ << " no folder set!"
120  << std::endl;
121  return;
122  }
123  // ----
124 
125 // // Cumulative distribution with average Ped value on a layer (not needed):
126 
127  std::map<uint32_t, ModMEs>::iterator selMEsMapIter_ = SummaryMEsMap_.find(getLayerNameAndId(selDetId_).second);
128  ModMEs selME_;
129  if ( selMEsMapIter_ != SummaryMEsMap_.end())
130  selME_ =selMEsMapIter_->second;
131  getSummaryMEs(selME_,selDetId_);
132 
133  SiStripPedestals::Range pedRange = pedestalHandle_->getRange(selDetId_);
134 
135  int nStrip = reader->getNumberOfApvsAndStripLength(selDetId_).first*128;
136 
138 
139  if(hPSet_.getParameter<bool>("FillSummaryProfileAtLayerLevel")){
140 
141  // --> profile summary
142 
143  std::string hSummaryOfProfile_description;
144  hSummaryOfProfile_description = hPSet_.getParameter<std::string>("SummaryOfProfile_description");
145 
146  std::string hSummaryOfProfile_name;
147 
148 
149  hSummaryOfProfile_name = hidmanager.createHistoLayer(hSummaryOfProfile_description,
150  "layer",
151  getLayerNameAndId(selDetId_).first,
152  "") ;
153 
154  for( int istrip=0;istrip<nStrip;++istrip){
155 
156  if( CondObj_fillId_ =="onlyProfile" || CondObj_fillId_ =="ProfileAndCumul"){
157  selME_.SummaryOfProfileDistr->Fill(istrip+1,pedestalHandle_->getPed(istrip,pedRange));
158  }
159 
160  //fill the TkMap
161  if(fPSet_.getParameter<bool>("TkMap_On") || hPSet_.getParameter<bool>("TkMap_On")){
162  fillTkMap(selDetId_, pedestalHandle_->getPed(istrip,pedRange));
163  }
164 
165  }// istrip
166  }//if Fill ...
167 
168  if(hPSet_.getParameter<bool>("FillSummaryAtLayerLevel")){
169 
170  // --> summary
171 
172  std::string hSummary_description;
173  hSummary_description = hPSet_.getParameter<std::string>("Summary_description");
174 
175  std::string hSummary_name;
176  hSummary_name = hidmanager.createHistoLayer(hSummary_description,
177  "layer",
178  getLayerNameAndId(selDetId_).first,
179  "") ;
180  float meanPedestal=0;
181 
182  for( int istrip=0;istrip<nStrip;++istrip){
183 
184  meanPedestal = meanPedestal + pedestalHandle_->getPed(istrip,pedRange);
185 
186  }//istrip
187 
188  meanPedestal = meanPedestal/nStrip;
189 
190 
191  // -----
192  // get detIds belonging to same layer to fill X-axis with detId-number
193 
194  std::vector<uint32_t> sameLayerDetIds_;
195 
196  sameLayerDetIds_.clear();
197 
198  sameLayerDetIds_=GetSameLayerDetId(activeDetIds,selDetId_);
199 
200 
201  unsigned int iBin=0;
202  for(unsigned int i=0;i<sameLayerDetIds_.size();i++){
203  if(sameLayerDetIds_[i]==selDetId_){iBin=i+1;}
204  }
205 
206  selME_.SummaryDistr->Fill(iBin,meanPedestal);
207 
208  // Fill the Histo_TkMap with the mean Pedestal:
209  if(HistoMaps_On_ ) Tk_HM_->fill(selDetId_, meanPedestal);
210 
211 
212  }//if Fill ...
213 
214 
215 
216 }
217 // -----
218 
void getModMEs(ModMEs &CondObj_ME, const uint32_t &detId_)
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
edm::ESHandle< SiStripPedestals > pedestalHandle_
ADDITON OF TK_HISTO_MAP.
const std::pair< unsigned short, double > getNumberOfApvsAndStripLength(uint32_t detId) const
std::vector< uint32_t > GetSameLayerDetId(std::vector< uint32_t > activeDetIds, uint32_t selDetId)
void fillMEsForLayer(uint32_t selDetId_)
std::pair< ContainerIterator, ContainerIterator > Range
void fillTkMap(const uint32_t &detid, const float &value)
SiStripPedestalsDQM(const edm::EventSetup &eSetup, edm::ParameterSet const &hPSet, edm::ParameterSet const &fPSet)
void Fill(long long x)
U second(std::pair< T, U > const &p)
std::map< uint32_t, ModMEs > SummaryMEsMap_
SiStripDetInfoFileReader * reader
void fillSummaryMEs(const std::vector< uint32_t > &selectedDetIds)
void fill(uint32_t &detid, float value)
Definition: TkHistoMap.cc:130
void fillMEsForDet(ModMEs selModME_, uint32_t selDetId_)
void fillModMEs(const std::vector< uint32_t > &selectedDetIds)
TH1 * getTH1(void) const
bool first
Definition: L1TdeRCT.cc:94
void getConditionObject(const edm::EventSetup &eSetup)
std::pair< std::string, uint32_t > getLayerNameAndId(const uint32_t &detId_)
void getActiveDetIds(const edm::EventSetup &eSetup)
TProfile * getTProfile(void) const
std::vector< uint32_t > activeDetIds
std::string createHistoLayer(std::string description, std::string id_type, std::string path, std::string flag)
void getSummaryMEs(ModMEs &CondObj_ME, const uint32_t &detId_)