CMS 3D CMS Logo

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