CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
SiStripPedestalsDQM.cc
Go to the documentation of this file.
2 #include "TCanvas.h"
3 
5  edm::RunNumber_t iRun,
6  edm::ParameterSet const &hPSet,
7  edm::ParameterSet const &fPSet,
8  const TrackerTopology *tTopo,
9  const TkDetMap *tkDetMap)
10  : SiStripBaseCondObjDQMGet<SiStripPedestals, SiStripPedestalsRcd>{token, iRun, hPSet, fPSet, tTopo} {
11  if (HistoMaps_On_) {
12  Tk_HM_ = std::make_unique<TkHistoMap>(tkDetMap, "SiStrip/Histo_Map", "MeanPed_TkMap", 0.);
13  }
14 }
15 
17 
19  getConditionObject(eSetup);
21 }
22 
23 void SiStripPedestalsDQM::fillModMEs(const std::vector<uint32_t> &selectedDetIds) {
24  ModMEs CondObj_ME;
25  for (const auto det : selectedDetIds) {
26  fillMEsForDet(CondObj_ME, det);
27  }
28 }
29 
30 void SiStripPedestalsDQM::fillMEsForDet(const ModMEs &_selModME_, uint32_t selDetId_) {
31  ModMEs selModME_ = _selModME_;
32  getModMEs(selModME_, selDetId_);
33 
34  const auto pedRange = condObj_->getRange(selDetId_);
35  int nStrip = detInfo_.getNumberOfApvsAndStripLength(selDetId_).first * 128;
36 
37  for (int istrip = 0; istrip < nStrip; ++istrip) {
38  if (CondObj_fillId_ == "onlyProfile" || CondObj_fillId_ == "ProfileAndCumul") {
39  selModME_.ProfileDistr->Fill(istrip + 1, condObj_->getPed(istrip, pedRange));
40  }
41  } // istrip
42 }
43 
44 void SiStripPedestalsDQM::fillSummaryMEs(const std::vector<uint32_t> &selectedDetIds) {
45  for (const auto det : selectedDetIds) {
46  fillMEsForLayer(/*SummaryMEsMap_,*/ det);
47  }
48 
49  for (const auto &itm : SummaryMEsMap_) {
50  ModMEs selME;
51  selME = itm.second;
52 
53  if (hPSet_.getParameter<bool>("FillSummaryProfileAtLayerLevel") &&
54  fPSet_.getParameter<bool>("OutputSummaryProfileAtLayerLevelAsImage")) {
55  if (CondObj_fillId_ == "onlyProfile" || CondObj_fillId_ == "ProfileAndCumul") {
56  TCanvas c1("c1");
57  selME.SummaryOfProfileDistr->getTProfile()->Draw();
58  std::string name(selME.SummaryOfProfileDistr->getTProfile()->GetTitle());
59  name += ".png";
60  c1.Print(name.c_str());
61  }
62  }
63  if (hPSet_.getParameter<bool>("FillSummaryAtLayerLevel") &&
64  fPSet_.getParameter<bool>("OutputSummaryAtLayerLevelAsImage")) {
65  TCanvas c1("c1");
66  selME.SummaryDistr->getTH1()->Draw();
67  std::string name(selME.SummaryDistr->getTitle());
68  name += ".png";
69  c1.Print(name.c_str());
70  }
71  }
72 }
73 
75  /*std::map<uint32_t, ModMEs> selMEsMap_,*/ uint32_t selDetId_) {
76  // ----
77  int subdetectorId_ = ((selDetId_ >> 25) & 0x7);
78 
79  if (subdetectorId_ < 3 || subdetectorId_ > 6) {
80  edm::LogError("SiStripPedestalsDQM") << "[SiStripPedestalsDQM::fillMEsForLayer] WRONG INPUT : no such "
81  "subdetector type : "
82  << subdetectorId_ << " no folder set!" << std::endl;
83  return;
84  }
85  // ----
86 
87  // // Cumulative distribution with average Ped value on a layer (not
88  // needed):
89 
90  const auto selMEsMapIter_ = SummaryMEsMap_.find(getLayerNameAndId(selDetId_).second);
91  ModMEs selME_;
92  if (selMEsMapIter_ != SummaryMEsMap_.end())
93  selME_ = selMEsMapIter_->second;
94  getSummaryMEs(selME_, selDetId_);
95 
96  const auto pedRange = condObj_->getRange(selDetId_);
97 
98  int nStrip = detInfo_.getNumberOfApvsAndStripLength(selDetId_).first * 128;
99 
101 
102  if (hPSet_.getParameter<bool>("FillSummaryProfileAtLayerLevel")) {
103  // --> profile summary
104 
105  std::string hSummaryOfProfile_description;
106  hSummaryOfProfile_description = hPSet_.getParameter<std::string>("SummaryOfProfile_description");
107 
108  std::string hSummaryOfProfile_name;
109 
110  hSummaryOfProfile_name =
111  hidmanager.createHistoLayer(hSummaryOfProfile_description, "layer", getLayerNameAndId(selDetId_).first, "");
112 
113  for (int istrip = 0; istrip < nStrip; ++istrip) {
114  if (CondObj_fillId_ == "onlyProfile" || CondObj_fillId_ == "ProfileAndCumul") {
115  selME_.SummaryOfProfileDistr->Fill(istrip + 1, condObj_->getPed(istrip, pedRange));
116  }
117 
118  // fill the TkMap
119  if (fPSet_.getParameter<bool>("TkMap_On") || hPSet_.getParameter<bool>("TkMap_On")) {
120  fillTkMap(selDetId_, condObj_->getPed(istrip, pedRange));
121  }
122 
123  } // istrip
124  } // if Fill ...
125 
126  if (hPSet_.getParameter<bool>("FillSummaryAtLayerLevel")) {
127  // --> summary
128 
129  std::string hSummary_description;
130  hSummary_description = hPSet_.getParameter<std::string>("Summary_description");
131 
132  std::string hSummary_name;
133  hSummary_name = hidmanager.createHistoLayer(hSummary_description, "layer", getLayerNameAndId(selDetId_).first, "");
134  float meanPedestal = 0;
135 
136  for (int istrip = 0; istrip < nStrip; ++istrip) {
137  meanPedestal = meanPedestal + condObj_->getPed(istrip, pedRange);
138 
139  } // istrip
140 
141  meanPedestal = meanPedestal / nStrip;
142 
143  // -----
144  // get detIds belonging to same layer to fill X-axis with detId-number
145 
146  std::vector<uint32_t> sameLayerDetIds_;
147 
148  sameLayerDetIds_.clear();
149 
150  sameLayerDetIds_ = GetSameLayerDetId(activeDetIds, selDetId_);
151 
152  unsigned int iBin = 0;
153  for (unsigned int i = 0; i < sameLayerDetIds_.size(); i++) {
154  if (sameLayerDetIds_[i] == selDetId_) {
155  iBin = i + 1;
156  }
157  }
158 
159  selME_.SummaryDistr->Fill(iBin, meanPedestal);
160 
161  // Fill the Histo_TkMap with the mean Pedestal:
162  if (HistoMaps_On_)
163  Tk_HM_->fill(selDetId_, meanPedestal);
164 
165  } // if Fill ...
166 }
void getModMEs(ModMEs &CondObj_ME, const uint32_t &detId_)
void fillMEsForLayer(uint32_t selDetId_) override
void fillSummaryMEs(const std::vector< uint32_t > &selectedDetIds) override
SiStripPedestalsDQM(edm::ESGetToken< SiStripPedestals, SiStripPedestalsRcd > token, edm::RunNumber_t iRun, edm::ParameterSet const &hPSet, edm::ParameterSet const &fPSet, const TrackerTopology *tTopo, const TkDetMap *tkDetMap)
Log< level::Error, false > LogError
void fillTkMap(const uint32_t &detid, const float &value)
float getPed(const uint16_t &strip, const Range &range) const
U second(std::pair< T, U > const &p)
~SiStripPedestalsDQM() override
std::map< uint32_t, ModMEs > SummaryMEsMap_
void fillModMEs(const std::vector< uint32_t > &selectedDetIds) override
const std::pair< unsigned short, double > getNumberOfApvsAndStripLength(uint32_t detId) const
void getDetIds(std::vector< uint32_t > &DetIds_) const
void getConditionObject(const edm::EventSetup &eSetup) override
void fillMEsForDet(const ModMEs &selModME_, uint32_t selDetId_) override
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::pair< std::string, uint32_t > getLayerNameAndId(const uint32_t &detId_)
std::vector< uint32_t > GetSameLayerDetId(const std::vector< uint32_t > &activeDetIds, uint32_t selDetId)
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)
void getSummaryMEs(ModMEs &CondObj_ME, const uint32_t &detId_)
std::unique_ptr< TkHistoMap > Tk_HM_