CMS 3D CMS Logo

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