CMS 3D CMS Logo

SiStripLorentzAngleDQM.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", "LA_TkMap", 0.);
18  }
19 }
20 // -----
21 
22 // -----
24 // -----
25 
26 // -----
28  getConditionObject(eSetup);
29 
30  std::map<uint32_t, float>::const_iterator LAMapIter_;
31  std::map<uint32_t, float> LAMap_ = lorentzangleHandle_->getLorentzAngles();
32 
33  for (LAMapIter_ = LAMap_.begin(); LAMapIter_ != LAMap_.end(); LAMapIter_++) {
34  activeDetIds.push_back((*LAMapIter_).first);
35  }
36 }
37 // -----
38 
39 // -----
40 void SiStripLorentzAngleDQM::fillSummaryMEs(const std::vector<uint32_t> &selectedDetIds, const edm::EventSetup &es) {
41  // Retrieve tracker topology from geometry
43  es.get<TrackerTopologyRcd>().get(tTopoHandle);
44  const TrackerTopology *const tTopo = tTopoHandle.product();
45 
46  // -----
47  // LA on layer-level : fill at once all detIds belonging to same layer when
48  // encountering first detID in the layer
49 
50  bool fillNext = true;
51  for (unsigned int i = 0; i < selectedDetIds.size(); i++) {
52  int subDetId_ = DetId(selectedDetIds[i]).subdetId();
53  if (subDetId_ < 3 || subDetId_ > 6) {
54  edm::LogError("SiStripLorentzAngle")
55  << "[SiStripLorentzAngle::fillSummaryMEs] WRONG INPUT : no such "
56  "subdetector type : "
57  << subDetId_ << " and detId " << selectedDetIds[i] << " therefore no filling!" << std::endl;
58  } else if (SummaryOnLayerLevel_On_) {
59  if (fillNext) {
60  fillMEsForLayer(/*SummaryMEsMap_,*/ selectedDetIds[i], tTopo);
61  }
62  if (getLayerNameAndId(selectedDetIds[i + 1], tTopo) == getLayerNameAndId(selectedDetIds[i], tTopo)) {
63  fillNext = false;
64  } else {
65  fillNext = true;
66  }
67  } else if (SummaryOnStringLevel_On_) {
68  if (fillNext) {
69  fillMEsForLayer(/*SummaryMEsMap_,*/ selectedDetIds[i], tTopo);
70  }
71  if (getStringNameAndId(selectedDetIds[i + 1], tTopo) == getStringNameAndId(selectedDetIds[i], tTopo)) {
72  fillNext = false;
73  } else {
74  fillNext = true;
75  }
76  }
77  }
78 
79  for (std::map<uint32_t, ModMEs>::iterator iter = SummaryMEsMap_.begin(); iter != SummaryMEsMap_.end(); iter++) {
80  ModMEs selME;
81  selME = iter->second;
82 
84  if (fPSet_.getParameter<bool>("OutputSummaryProfileAtLayerLevelAsImage")) {
85  TCanvas c1("c1");
86  selME.SummaryOfProfileDistr->getTProfile()->Draw();
87  std::string name(selME.SummaryOfProfileDistr->getTProfile()->GetTitle());
88  name += ".png";
89  c1.Print(name.c_str());
90  }
91 
92  if (fPSet_.getParameter<bool>("OutputCumulativeSummaryAtLayerLevelAsImage")) {
93  TCanvas c2("c2");
94  selME.SummaryOfCumulDistr->getTH1()->Draw();
95  std::string name2(selME.SummaryOfCumulDistr->getTH1()->GetTitle());
96  name2 += ".png";
97  c2.Print(name2.c_str());
98  }
99 
100  } else {
101  if (hPSet_.getParameter<bool>("FillSummaryProfileAtLayerLevel") &&
102  fPSet_.getParameter<bool>("OutputSummaryProfileAtLayerLevelAsImage")) {
103  TCanvas c1("c1");
104  selME.SummaryOfProfileDistr->getTProfile()->Draw();
105  std::string name(selME.SummaryOfProfileDistr->getTProfile()->GetTitle());
106  name += ".png";
107  c1.Print(name.c_str());
108  }
109 
110  if (hPSet_.getParameter<bool>("FillCumulativeSummaryAtLayerLevel") &&
111  fPSet_.getParameter<bool>("OutputCumulativeSummaryAtLayerLevelAsImage")) {
112  TCanvas c1("c1");
113  selME.SummaryOfCumulDistr->getTH1()->Draw();
114  std::string name(selME.SummaryOfCumulDistr->getTH1()->GetTitle());
115  name += ".png";
116  c1.Print(name.c_str());
117  }
118  }
119  }
120 }
121 // -----
122 
123 // -----
125  /*std::map<uint32_t, ModMEs> selMEsMap_,*/ uint32_t selDetId_, const TrackerTopology *tTopo) {
127 
128  std::string hSummaryOfProfile_description;
129  hSummaryOfProfile_description = hPSet_.getParameter<std::string>("SummaryOfProfile_description");
130 
131  std::string hSummary_name;
132 
133  int subDetId_ = DetId(selDetId_).subdetId();
134 
135  if (subDetId_ < 3 || subDetId_ > 6) {
136  edm::LogError("SiStripLorentzAngleDQM") << "[SiStripLorentzAngleDQM::fillMEsForLayer] WRONG INPUT : no such "
137  "subdetector type : "
138  << subDetId_ << " no folder set!" << std::endl;
139  return;
140  }
141 
142  std::vector<uint32_t> sameLayerDetIds_;
143  sameLayerDetIds_.clear();
144 
145  if (SummaryOnStringLevel_On_) { // FILLING FOR STRING LEVEL
146 
147  hSummary_name = hidmanager.createHistoLayer(
148  hSummaryOfProfile_description, "layer", getStringNameAndId(selDetId_, tTopo).first, "");
149  std::map<uint32_t, ModMEs>::iterator selMEsMapIter_ =
150  SummaryMEsMap_.find(getStringNameAndId(selDetId_, tTopo).second);
151 
152  ModMEs selME_;
153  if (selMEsMapIter_ != SummaryMEsMap_.end())
154  selME_ = selMEsMapIter_->second;
155 
156  getSummaryMEs(selME_, selDetId_, tTopo);
157 
158  // -----
159  sameLayerDetIds_.clear();
160 
161  switch (DetId(selDetId_).subdetId()) {
163  if (tTopo->tibIsInternalString(selDetId_)) {
165  activeDetIds, sameLayerDetIds_, tTopo, tTopo->tibLayer(selDetId_), 0, 1, tTopo->tibString(selDetId_));
166  }
167  if (tTopo->tibIsExternalString(selDetId_)) {
169  activeDetIds, sameLayerDetIds_, tTopo, tTopo->tibLayer(selDetId_), 0, 2, tTopo->tibString(selDetId_));
170  }
171  break;
173  SiStripSubStructure::getTIDDetectors(activeDetIds, sameLayerDetIds_, tTopo, 0, 0, 0, 0);
174  break;
177  activeDetIds, sameLayerDetIds_, tTopo, tTopo->tobLayer(selDetId_), 0, tTopo->tobRod(selDetId_));
178  break;
180  SiStripSubStructure::getTECDetectors(activeDetIds, sameLayerDetIds_, tTopo, 0, 0, 0, 0, 0, 0);
181  break;
182  }
183 
184  // -----
185 
186  for (unsigned int i = 0; i < sameLayerDetIds_.size(); i++) {
187  selME_.SummaryOfProfileDistr->Fill(i + 1, lorentzangleHandle_->getLorentzAngle(sameLayerDetIds_[i]));
188 
189  // Fill the Histo_TkMap+TkMap with the LA:
190  if (HistoMaps_On_)
191  Tk_HM_->fill(sameLayerDetIds_[i], lorentzangleHandle_->getLorentzAngle(sameLayerDetIds_[i]));
192 
193  std::cout << sameLayerDetIds_[i] << "\t" << lorentzangleHandle_->getLorentzAngle(sameLayerDetIds_[i])
194  << std::endl;
195 
196  if (fPSet_.getParameter<bool>("TkMap_On") || hPSet_.getParameter<bool>("TkMap_On")) {
197  fillTkMap(sameLayerDetIds_[i], lorentzangleHandle_->getLorentzAngle(sameLayerDetIds_[i]));
198  }
199  }
200 
201  std::string hSummaryOfCumul_description;
202  hSummaryOfCumul_description = hPSet_.getParameter<std::string>("SummaryOfCumul_description");
203 
204  std::string hSummaryOfCumul_name;
205 
206  if (subDetId_ < 3 || subDetId_ > 6) {
207  edm::LogError("SiStripLorentzAngleDQM") << "[SiStripLorentzAngleDQM::fillMEsForLayer] WRONG INPUT : no such "
208  "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++) {
217  selME_.SummaryOfCumulDistr->Fill(lorentzangleHandle_->getLorentzAngle(sameLayerDetIds_[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, lorentzangleHandle_->getLorentzAngle(sameLayerDetIds_[i]));
243 
244  // Fill the Histo_TkMap with LA:
245  if (HistoMaps_On_)
246  Tk_HM_->fill(sameLayerDetIds_[i], lorentzangleHandle_->getLorentzAngle(sameLayerDetIds_[i]));
247 
248  if (fPSet_.getParameter<bool>("TkMap_On") || hPSet_.getParameter<bool>("TkMap_On")) {
249  fillTkMap(sameLayerDetIds_[i], lorentzangleHandle_->getLorentzAngle(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("SiStripLorentzAngleDQM") << "[SiStripLorentzAngleDQM::fillMEsForLayer] WRONG INPUT : no "
262  "such subdetector type : "
263  << subDetId_ << " no folder set!" << std::endl;
264  return;
265  }
266 
267  hSummaryOfCumul_name = hidmanager.createHistoLayer(
268  hSummaryOfCumul_description, "layer", getLayerNameAndId(selDetId_, tTopo).first, "");
269 
270  for (unsigned int i = 0; i < sameLayerDetIds_.size(); i++) {
271  selME_.SummaryOfCumulDistr->Fill(lorentzangleHandle_->getLorentzAngle(sameLayerDetIds_[i]));
272  }
273  } // if Fill ...
274  } // FILLING FOR LAYER LEVEL
275 }
276 // -----
TProfile * getTProfile() const
T getParameter(std::string const &) const
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)
unsigned int tibLayer(const DetId &id) const
unsigned int tibString(const DetId &id) const
void fillSummaryMEs(const std::vector< uint32_t > &selectedDetIds, const edm::EventSetup &es) override
edm::ESHandle< SiStripLorentzAngle > lorentzangleHandle_
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 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)
SiStripLorentzAngleDQM(const edm::EventSetup &eSetup, edm::RunNumber_t iRun, edm::ParameterSet const &hPSet, edm::ParameterSet const &fPSet)
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 > & getLorentzAngles() const
std::map< uint32_t, ModMEs > SummaryMEsMap_
bool tibIsExternalString(const DetId &id) const
float getLorentzAngle(const uint32_t &) const
void getActiveDetIds(const edm::EventSetup &eSetup) override
std::vector< uint32_t > GetSameLayerDetId(const std::vector< uint32_t > &activeDetIds, uint32_t selDetId, const TrackerTopology *tTopo)
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:41
std::pair< std::string, uint32_t > getStringNameAndId(const uint32_t &detId_, const TrackerTopology *tTopo)
Definition: DetId.h:18
void fillMEsForLayer(uint32_t selDetId_, const TrackerTopology *tTopo) override
void getConditionObject(const edm::EventSetup &eSetup) override
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