CMS 3D CMS Logo

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