CMS 3D CMS Logo

SiStripBackPlaneCorrectionDQM.cc
Go to the documentation of this file.
5 #include "TCanvas.h"
6 
7 // -----
9  edm::RunNumber_t iRun,
10  edm::ParameterSet const &hPSet,
11  edm::ParameterSet const &fPSet)
12  : SiStripBaseCondObjDQM(eSetup, iRun, hPSet, fPSet) {
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", "BP_TkMap", 0.);
18  }
19 }
20 // -----
21 
22 // -----
24 // -----
25 
26 // -----
28  getConditionObject(eSetup);
29 
30  std::map<uint32_t, float>::const_iterator BPMapIter_;
31  std::map<uint32_t, float> BPMap_ = bpcorrectionHandle_->getBackPlaneCorrections();
32 
33  for (BPMapIter_ = BPMap_.begin(); BPMapIter_ != BPMap_.end(); BPMapIter_++) {
34  activeDetIds.push_back((*BPMapIter_).first);
35  }
36 }
37 // -----
38 
39 // -----
40 void SiStripBackPlaneCorrectionDQM::fillSummaryMEs(const std::vector<uint32_t> &selectedDetIds,
41  const edm::EventSetup &es) {
42  // Retrieve tracker topology from geometry
44  es.get<TrackerTopologyRcd>().get(tTopoHandle);
45  const TrackerTopology *const tTopo = tTopoHandle.product();
46 
47  // -----
48  // BP on layer-level : fill at once all detIds belonging to same layer when
49  // encountering first detID in the layer
50 
51  bool fillNext = true;
52  for (unsigned int i = 0; i < selectedDetIds.size(); i++) {
53  int subDetId_ = ((selectedDetIds[i] >> 25) & 0x7);
54  if (subDetId_ < 3 || subDetId_ > 6) {
55  edm::LogError("SiStripBackPlaneCorrection")
56  << "[SiStripBackPlaneCorrection::fillSummaryMEs] WRONG INPUT : no "
57  "such subdetector type : "
58  << subDetId_ << " and detId " << selectedDetIds[i] << " therefore no filling!" << std::endl;
59  } else if (SummaryOnLayerLevel_On_) {
60  if (fillNext) {
61  fillMEsForLayer(/*SummaryMEsMap_,*/ selectedDetIds[i], tTopo);
62  }
63  if (getLayerNameAndId(selectedDetIds[i + 1], tTopo) == getLayerNameAndId(selectedDetIds[i], tTopo)) {
64  fillNext = false;
65  } else {
66  fillNext = true;
67  }
68  } else if (SummaryOnStringLevel_On_) {
69  if (fillNext) {
70  fillMEsForLayer(/*SummaryMEsMap_,*/ selectedDetIds[i], tTopo);
71  }
72  if (getStringNameAndId(selectedDetIds[i + 1], tTopo) == getStringNameAndId(selectedDetIds[i], tTopo)) {
73  fillNext = false;
74  } else {
75  fillNext = true;
76  }
77  }
78  }
79 
80  for (std::map<uint32_t, ModMEs>::iterator iter = SummaryMEsMap_.begin(); iter != SummaryMEsMap_.end(); iter++) {
81  ModMEs selME;
82  selME = iter->second;
83 
85  if (fPSet_.getParameter<bool>("OutputSummaryProfileAtLayerLevelAsImage")) {
86  TCanvas c1("c1");
87  selME.SummaryOfProfileDistr->getTProfile()->Draw();
88  std::string name(selME.SummaryOfProfileDistr->getTProfile()->GetTitle());
89  name += ".png";
90  c1.Print(name.c_str());
91  }
92 
93  if (fPSet_.getParameter<bool>("OutputCumulativeSummaryAtLayerLevelAsImage")) {
94  TCanvas c2("c2");
95  selME.SummaryOfCumulDistr->getTH1()->Draw();
96  std::string name2(selME.SummaryOfCumulDistr->getTH1()->GetTitle());
97  name2 += ".png";
98  c2.Print(name2.c_str());
99  }
100 
101  } else {
102  if (hPSet_.getParameter<bool>("FillSummaryProfileAtLayerLevel") &&
103  fPSet_.getParameter<bool>("OutputSummaryProfileAtLayerLevelAsImage")) {
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 
111  if (hPSet_.getParameter<bool>("FillCumulativeSummaryAtLayerLevel") &&
112  fPSet_.getParameter<bool>("OutputCumulativeSummaryAtLayerLevelAsImage")) {
113  TCanvas c1("c1");
114  selME.SummaryOfCumulDistr->getTH1()->Draw();
115  std::string name(selME.SummaryOfCumulDistr->getTH1()->GetTitle());
116  name += ".png";
117  c1.Print(name.c_str());
118  }
119  }
120  }
121 }
122 // -----
123 
124 // -----
126  /*std::map<uint32_t, ModMEs> selMEsMap_,*/ uint32_t selDetId_, const TrackerTopology *tTopo) {
128 
129  std::string hSummaryOfProfile_description;
130  hSummaryOfProfile_description = hPSet_.getParameter<std::string>("SummaryOfProfile_description");
131 
132  std::string hSummary_name;
133 
134  int subDetId_ = ((selDetId_ >> 25) & 0x7);
135 
136  if (subDetId_ < 3 || subDetId_ > 6) {
137  edm::LogError("SiStripBackPlaneCorrectionDQM")
138  << "[SiStripBackPlaneCorrectionDQM::fillMEsForLayer] WRONG INPUT : no "
139  "such subdetector type : "
140  << subDetId_ << " no folder set!" << std::endl;
141  return;
142  }
143 
144  uint32_t selSubDetId_ = ((selDetId_ >> 25) & 0x7);
145 
146  std::vector<uint32_t> sameLayerDetIds_;
147  sameLayerDetIds_.clear();
148 
149  if (SummaryOnStringLevel_On_) { // FILLING FOR STRING LEVEL
150 
151  hSummary_name = hidmanager.createHistoLayer(
152  hSummaryOfProfile_description, "layer", getStringNameAndId(selDetId_, tTopo).first, "");
153  std::map<uint32_t, ModMEs>::iterator selMEsMapIter_ =
154  SummaryMEsMap_.find(getStringNameAndId(selDetId_, tTopo).second);
155 
156  ModMEs selME_;
157  if (selMEsMapIter_ != SummaryMEsMap_.end())
158  selME_ = selMEsMapIter_->second;
159 
160  getSummaryMEs(selME_, selDetId_, tTopo);
161 
162  // -----
163  sameLayerDetIds_.clear();
164 
165  if (selSubDetId_ == 3) { // TIB
166  if (tTopo->tibIsInternalString(selDetId_)) {
168  activeDetIds, sameLayerDetIds_, tTopo, tTopo->tibLayer(selDetId_), 0, 1, tTopo->tibString(selDetId_));
169  }
170  if (tTopo->tibIsExternalString(selDetId_)) {
172  activeDetIds, sameLayerDetIds_, tTopo, tTopo->tibLayer(selDetId_), 0, 2, tTopo->tibString(selDetId_));
173  }
174  } else if (selSubDetId_ == 4) { // TID
175  SiStripSubStructure::getTIDDetectors(activeDetIds, sameLayerDetIds_, tTopo, 0, 0, 0, 0);
176  } else if (selSubDetId_ == 5) { // TOB
178  activeDetIds, sameLayerDetIds_, tTopo, tTopo->tobLayer(selDetId_), 0, tTopo->tobRod(selDetId_));
179  } else if (selSubDetId_ == 6) { // TEC
180  SiStripSubStructure::getTECDetectors(activeDetIds, sameLayerDetIds_, tTopo, 0, 0, 0, 0, 0, 0);
181  }
182 
183  // -----
184 
185  for (unsigned int i = 0; i < sameLayerDetIds_.size(); i++) {
186  selME_.SummaryOfProfileDistr->Fill(i + 1, bpcorrectionHandle_->getBackPlaneCorrection(sameLayerDetIds_[i]));
187 
188  // Fill the Histo_TkMap+TkMap with the BP:
189  if (HistoMaps_On_)
190  Tk_HM_->fill(sameLayerDetIds_[i], bpcorrectionHandle_->getBackPlaneCorrection(sameLayerDetIds_[i]));
191 
192  std::cout << sameLayerDetIds_[i] << "\t" << bpcorrectionHandle_->getBackPlaneCorrection(sameLayerDetIds_[i])
193  << std::endl;
194 
195  if (fPSet_.getParameter<bool>("TkMap_On") || hPSet_.getParameter<bool>("TkMap_On")) {
196  fillTkMap(sameLayerDetIds_[i], bpcorrectionHandle_->getBackPlaneCorrection(sameLayerDetIds_[i]));
197  }
198  }
199 
200  std::string hSummaryOfCumul_description;
201  hSummaryOfCumul_description = hPSet_.getParameter<std::string>("SummaryOfCumul_description");
202 
203  std::string hSummaryOfCumul_name;
204 
205  if (subDetId_ < 3 || subDetId_ > 6) {
206  edm::LogError("SiStripBackPlaneCorrectionDQM")
207  << "[SiStripBackPlaneCorrectionDQM::fillMEsForLayer] WRONG INPUT : "
208  "no such subdetector type : "
209  << subDetId_ << " no folder set!" << std::endl;
210  return;
211  }
212 
213  hSummaryOfCumul_name = hidmanager.createHistoLayer(
214  hSummaryOfCumul_description, "layer", getStringNameAndId(selDetId_, tTopo).first, "");
215 
216  for (unsigned int i = 0; i < sameLayerDetIds_.size(); i++) {
218  }
219  } // FILLING FOR STRING LEVEL
220 
221  else { // FILLING FOR LAYER LEVEL
222 
223  std::map<uint32_t, ModMEs>::iterator selMEsMapIter_ =
224  SummaryMEsMap_.find(getLayerNameAndId(selDetId_, tTopo).second);
225 
226  ModMEs selME_;
227  if (selMEsMapIter_ != SummaryMEsMap_.end())
228  selME_ = selMEsMapIter_->second;
229 
230  getSummaryMEs(selME_, selDetId_, tTopo);
231 
232  if (hPSet_.getParameter<bool>("FillSummaryProfileAtLayerLevel")) {
233  hSummary_name = hidmanager.createHistoLayer(
234  hSummaryOfProfile_description, "layer", getLayerNameAndId(selDetId_, tTopo).first, "");
235 
236  // -----
237  sameLayerDetIds_.clear();
238 
239  sameLayerDetIds_ = GetSameLayerDetId(activeDetIds, selDetId_, tTopo);
240 
241  for (unsigned int i = 0; i < sameLayerDetIds_.size(); i++) {
242  selME_.SummaryOfProfileDistr->Fill(i + 1, bpcorrectionHandle_->getBackPlaneCorrection(sameLayerDetIds_[i]));
243 
244  // Fill the Histo_TkMap with BP:
245  if (HistoMaps_On_)
246  Tk_HM_->fill(sameLayerDetIds_[i], bpcorrectionHandle_->getBackPlaneCorrection(sameLayerDetIds_[i]));
247 
248  if (fPSet_.getParameter<bool>("TkMap_On") || hPSet_.getParameter<bool>("TkMap_On")) {
249  fillTkMap(sameLayerDetIds_[i], bpcorrectionHandle_->getBackPlaneCorrection(sameLayerDetIds_[i]));
250  }
251  }
252  } // if Fill ...
253 
254  if (hPSet_.getParameter<bool>("FillCumulativeSummaryAtLayerLevel")) {
255  std::string hSummaryOfCumul_description;
256  hSummaryOfCumul_description = hPSet_.getParameter<std::string>("SummaryOfCumul_description");
257 
258  std::string hSummaryOfCumul_name;
259 
260  if (subDetId_ < 3 || subDetId_ > 6) {
261  edm::LogError("SiStripBackPlaneCorrectionDQM")
262  << "[SiStripBackPlaneCorrectionDQM::fillMEsForLayer] WRONG INPUT : "
263  "no such subdetector type : "
264  << subDetId_ << " no folder set!" << std::endl;
265  return;
266  }
267 
268  hSummaryOfCumul_name = hidmanager.createHistoLayer(
269  hSummaryOfCumul_description, "layer", getLayerNameAndId(selDetId_, tTopo).first, "");
270 
271  for (unsigned int i = 0; i < sameLayerDetIds_.size(); i++) {
273  }
274  } // if Fill ...
275  } // FILLING FOR LAYER LEVEL
276 }
277 // -----
SiStripBackPlaneCorrectionDQM(const edm::EventSetup &eSetup, edm::RunNumber_t iRun, edm::ParameterSet const &hPSet, edm::ParameterSet const &fPSet)
TProfile * getTProfile() const
T getParameter(std::string const &) const
edm::ESHandle< SiStripBackPlaneCorrection > bpcorrectionHandle_
void getTIBDetectors(const std::vector< uint32_t > &inputDetRawIds, std::vector< uint32_t > &tibDetRawIds, const TrackerTopology *trackerTopology, uint32_t layer=0, uint32_t bkw_frw=0, uint32_t int_ext=0, uint32_t string=0)
void fillSummaryMEs(const std::vector< uint32_t > &selectedDetIds, const edm::EventSetup &es) override
unsigned int tibLayer(const DetId &id) const
unsigned int tibString(const DetId &id) const
TH1 * getTH1() const
void getTIDDetectors(const std::vector< uint32_t > &inputDetRawIds, std::vector< uint32_t > &tidDetRawIds, const TrackerTopology *trackerTopology, uint32_t side=0, uint32_t wheel=0, uint32_t ring=0, uint32_t ster=0)
void getActiveDetIds(const edm::EventSetup &eSetup) override
void getTECDetectors(const std::vector< uint32_t > &inputDetRawIds, std::vector< uint32_t > &tecDetRawIds, const TrackerTopology *trackerTopology, uint32_t side=0, uint32_t wheel=0, uint32_t petal_bkw_frw=0, uint32_t petal=0, uint32_t ring=0, uint32_t ster=0)
void getSummaryMEs(ModMEs &CondObj_ME, const uint32_t &detId_, const TrackerTopology *tTopo)
void fillMEsForLayer(uint32_t selDetId_, const TrackerTopology *tTopo) override
void fillTkMap(const uint32_t &detid, const float &value)
void Fill(long long x)
U second(std::pair< T, U > const &p)
const std::map< unsigned int, float > & getBackPlaneCorrections() const
std::map< uint32_t, ModMEs > SummaryMEsMap_
void getConditionObject(const edm::EventSetup &eSetup) override
bool tibIsExternalString(const DetId &id) const
std::vector< uint32_t > GetSameLayerDetId(const std::vector< uint32_t > &activeDetIds, uint32_t selDetId, const TrackerTopology *tTopo)
float getBackPlaneCorrection(const uint32_t &) const
std::pair< std::string, uint32_t > getStringNameAndId(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:71
std::vector< uint32_t > activeDetIds
unsigned int RunNumber_t
void getTOBDetectors(const std::vector< uint32_t > &inputDetRawIds, std::vector< uint32_t > &tobDetRawIds, const TrackerTopology *trackerTopology, uint32_t layer=0, uint32_t bkw_frw=0, uint32_t rod=0)
std::string createHistoLayer(std::string description, std::string id_type, std::string path, std::string flag)
unsigned int tobRod(const DetId &id) const
bool tibIsInternalString(const DetId &id) const
T const * product() const
Definition: ESHandle.h:86
std::unique_ptr< TkHistoMap > Tk_HM_
unsigned int tobLayer(const DetId &id) const