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.
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, const edm::EventSetup& es){
37 
38  //Retrieve tracker topology from geometry
40  es.get<TrackerTopologyRcd>().get(tTopoHandle);
41  const TrackerTopology* const tTopo = tTopoHandle.product();
42 
43  ModMEs CondObj_ME;
44 
45  for(std::vector<uint32_t>::const_iterator detIter_ =selectedDetIds.begin();
46  detIter_!=selectedDetIds.end();++detIter_){
47  fillMEsForDet(CondObj_ME,*detIter_,tTopo);
48  }
49 }
50 
51 
52 // -----
53 void SiStripApvGainsDQM::fillMEsForDet(const ModMEs& _selModME_, uint32_t selDetId_, const TrackerTopology* tTopo){
54  ModMEs selModME_ = _selModME_;
55  std::vector<uint32_t> DetIds;
56  gainHandle_->getDetIds(DetIds);
57 
58  SiStripApvGain::Range gainRange = gainHandle_->getRange(selDetId_);
59 
60  int nApv = reader->getNumberOfApvsAndStripLength(selDetId_).first;
61 
62  getModMEs(selModME_,selDetId_, tTopo);
63 
64  for( int iapv=0;iapv<nApv;++iapv){
65  if( CondObj_fillId_ =="onlyProfile" || CondObj_fillId_ =="ProfileAndCumul"){
66  selModME_.ProfileDistr->Fill(iapv+1,gainHandle_->getApvGain(iapv,gainRange));
67  }
68  if( CondObj_fillId_ =="onlyCumul" || CondObj_fillId_ =="ProfileAndCumul"){
69  selModME_.CumulDistr ->Fill(gainHandle_->getApvGain(iapv,gainRange));
70  }
71 
72  // Fill the TkMap
73  if(fPSet_.getParameter<bool>("TkMap_On") || hPSet_.getParameter<bool>("TkMap_On")){
74  fillTkMap(selDetId_, gainHandle_->getApvGain(iapv,gainRange));
75  }
76 
77  }
78 }
79 
80 // -----
81 void SiStripApvGainsDQM::fillSummaryMEs(const std::vector<uint32_t> & selectedDetIds, const edm::EventSetup& es){
82 
83  //Retrieve tracker topology from geometry
85  es.get<TrackerTopologyRcd>().get(tTopoHandle);
86  const TrackerTopology* const tTopo = tTopoHandle.product();
87 
88  for(std::vector<uint32_t>::const_iterator detIter_ = selectedDetIds.begin();
89  detIter_!= selectedDetIds.end();detIter_++){
90  fillMEsForLayer(/*SummaryMEsMap_,*/ *detIter_, tTopo);
91  }
92 
93  for (std::map<uint32_t, ModMEs>::iterator iter=SummaryMEsMap_.begin(); iter!=SummaryMEsMap_.end(); iter++){
94 
95  ModMEs selME;
96  selME = iter->second;
97 
98  if(hPSet_.getParameter<bool>("FillSummaryProfileAtLayerLevel") && fPSet_.getParameter<bool>("OutputSummaryProfileAtLayerLevelAsImage")){
99 
100  TCanvas c1("c1");
101  selME.SummaryOfProfileDistr->getTProfile()->Draw();
102  std::string name (selME.SummaryOfProfileDistr->getTProfile()->GetTitle());
103  name+=".png";
104  c1.Print(name.c_str());
105  }
106  if(hPSet_.getParameter<bool>("FillSummaryAtLayerLevel") && fPSet_.getParameter<bool>("OutputSummaryAtLayerLevelAsImage")){
107 
108  TCanvas c1("c1");
109  selME.SummaryDistr->getTH1()->Draw();
110  std::string name (selME.SummaryDistr->getTH1()->GetTitle());
111  name+=".png";
112  c1.Print(name.c_str());
113  }
114  }
115 
116 }
117 
118 // -----
119 void SiStripApvGainsDQM::fillMEsForLayer( /*std::map<uint32_t, ModMEs> selMEsMap_, */uint32_t selDetId_, const TrackerTopology* tTopo){
120 
121  int subdetectorId_ = ((selDetId_>>25)&0x7);
122 
123  if( subdetectorId_<3 ||subdetectorId_>6 ){
124  edm::LogError("SiStripApvGainsDQM")
125  << "[SiStripApvGainsDQM::fillMEsForLayer] WRONG INPUT : no such subdetector type : "
126  << subdetectorId_ << " no folder set!"
127  << std::endl;
128  return;
129  }
130  // ----
131 
132  std::map<uint32_t, ModMEs>::iterator selMEsMapIter_ = SummaryMEsMap_.find(getLayerNameAndId(selDetId_,tTopo).second);
133  ModMEs selME_;
134  if ( selMEsMapIter_ != SummaryMEsMap_.end())
135  selME_ =selMEsMapIter_->second;
136  getSummaryMEs(selME_,selDetId_,tTopo);
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_,tTopo).first, "") ;
155 
156  for( int iapv=0;iapv<nApv;++iapv){
157 
158  meanApvGain = meanApvGain +gainHandle_ ->getApvGain(iapv,gainRange);
159  selME_.SummaryOfProfileDistr->Fill(iapv+1,gainHandle_->getApvGain(iapv,gainRange));
160 
161  // Fill the TkMap
162  if(fPSet_.getParameter<bool>("TkMap_On") || hPSet_.getParameter<bool>("TkMap_On")){
163  fillTkMap(selDetId_, gainHandle_->getApvGain(iapv,gainRange));
164  }
165 
166  }// iapv
167 
168  meanApvGain = meanApvGain/nApv;
169 
170  // Fill the TkHistoMap with meanApvgain:
171  if(HistoMaps_On_ ) Tk_HM_->setBinContent(selDetId_, meanApvGain);
172 
173  }//if Fill ...
174 
175 
176  if(hPSet_.getParameter<bool>("FillSummaryAtLayerLevel")){
177 
178  // --> summary
179  std::string hSummary_description;
180  hSummary_description = hPSet_.getParameter<std::string>("Summary_description");
181 
182  std::string hSummary_name;
183  hSummary_name = hidmanager.createHistoLayer(hSummary_description,
184  "layer",
185  getLayerNameAndId(selDetId_,tTopo).first,
186  "") ;
187 
188 
189  // get detIds belonging to same layer to fill X-axis with detId-number
190 
191  std::vector<uint32_t> sameLayerDetIds_;
192 
193  sameLayerDetIds_.clear();
194 
195  sameLayerDetIds_=GetSameLayerDetId(activeDetIds,selDetId_,tTopo);
196 
197  unsigned int iBin=0;
198  for(unsigned int i=0;i<sameLayerDetIds_.size();i++){
199  if(sameLayerDetIds_[i]==selDetId_){iBin=i+1;}
200  }
201 
202  for( int iapv=0;iapv<nApv;++iapv){
203  meanApvGain = meanApvGain +gainHandle_ ->getApvGain(iapv,gainRange);
204  selME_.SummaryDistr->Fill(iBin,gainHandle_->getApvGain(iapv,gainRange));
205  }//iapv
206  meanApvGain = meanApvGain/nApv;
207 
208  // Fill the TkHistoMap with meanApvgain:
209  // if(HistoMaps_On_ ) Tk_HM_->setBinContent(selDetId_, meanApvGain);
210 
211  }//if Fill ...
212 }
213 // -----
214 
215 
216 
217 
218 
219 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
SiStripApvGainsDQM(const edm::EventSetup &eSetup, edm::ParameterSet const &hPSet, edm::ParameterSet const &fPSet)
const std::pair< unsigned short, double > getNumberOfApvsAndStripLength(uint32_t detId) const
void fillMEsForLayer(uint32_t selDetId_, const TrackerTopology *tTopo)
void getActiveDetIds(const edm::EventSetup &eSetup)
void fillSummaryMEs(const std::vector< uint32_t > &selectedDetIds, const edm::EventSetup &es)
void getSummaryMEs(ModMEs &CondObj_ME, const uint32_t &detId_, const TrackerTopology *tTopo)
void fillModMEs(const std::vector< uint32_t > &selectedDetIds, const edm::EventSetup &es)
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_
SiStripDetInfoFileReader * reader
std::vector< uint32_t > GetSameLayerDetId(const std::vector< uint32_t > &activeDetIds, uint32_t selDetId, const TrackerTopology *tTopo)
std::pair< ContainerIterator, ContainerIterator > Range
TH1 * getTH1(void) const
void getConditionObject(const edm::EventSetup &eSetup)
edm::ESHandle< SiStripApvGain > gainHandle_
void setBinContent(uint32_t &detid, float value)
Definition: TkHistoMap.cc:196
void getModMEs(ModMEs &CondObj_ME, const uint32_t &detId_, const TrackerTopology *tTopo)
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
std::pair< std::string, uint32_t > getLayerNameAndId(const uint32_t &detId_, const TrackerTopology *tTopo)
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 fillMEsForDet(const ModMEs &selModME_, uint32_t selDetId_, const TrackerTopology *tTopo)