CMS 3D CMS Logo

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