CMS 3D CMS Logo

SiStripPlotGain.cc
Go to the documentation of this file.
2 
5 
6 
7 
9  cacheID(0xFFFFFFFF)
10 {
11  //now do what ever initialization is needed
12  if(!edm::Service<SiStripDetInfoFileReader>().isAvailable()){
13  edm::LogError("TkLayerMap") <<
14  "\n------------------------------------------"
15  "\nUnAvailable Service SiStripDetInfoFileReader: please insert in the configuration file an instance like"
16  "\n\tprocess.SiStripDetInfoFileReader = cms.Service(\"SiStripDetInfoFileReader\")"
17  "\n------------------------------------------";
18  }
19 
21  file = new TFile("correlTest.root","RECREATE");
22  tkmap = new TrackerMap();
23 }
24 
25 
27 {}
28 
29 //
30 
31 void
33 
34  if(getCache(es)==cacheID )
35  return;
36  cacheID=getCache(es);
37 
38  edm::LogInfo("") << "[SiStripPlotGain::beginRun] cacheID " << cacheID << std::endl;
39 
40  es.get<SiStripApvGainRcd>().get(Handle_);
41  DoAnalysis(es, *Handle_.product());
42 
43 
44 }
45 
46 void
48 
49  edm::LogInfo("") << "[Doanalysis]";
50 
51  //Retrieve tracker topology from geometry
53  es.get<TrackerTopologyRcd>().get(tTopoHandle);
54  const TrackerTopology* const tTopo = tTopoHandle.product();
55 
56  std::vector<TH1F *>histos;
57 
60  iter=p.detid_begin;
61  iterE=p.detid_end;
62 
63  float value;
64 
65  //Divide result by d
66  for(;iter!=iterE;++iter){
67  getHistos(*iter,tTopo,histos);
69 
70  edm::LogInfo("") << "[Doanalysis] detid " << *iter << " range " << range.second-range.first;
71  size_t apv=0, apvE= (range.second-range.first);
72  for (;apv<apvE;apv+=2){
73  value=gain.getApvGain(apv,range);
74  tkmap->fill(*iter,value);
75  for(size_t i=0;i<histos.size();++i)
76  histos[i]->Fill(value);
77  }
78 
79  }
80 }
81 
82 
83 void
84 SiStripPlotGain::getHistos(const uint32_t& detid, const TrackerTopology* tTopo, std::vector<TH1F*>& histos){
85 
86  histos.clear();
87 
88  int subdet=-999; int component=-999;
89  SiStripDetId a(detid);
90  if ( a.subdetId() == 3 ){
91  subdet=0;
92  component=tTopo->tibLayer(detid);
93  } else if ( a.subdetId() == 4 ) {
94  subdet=1;
95  component=tTopo->tidSide(detid)==2?tTopo->tidWheel(detid):tTopo->tidWheel(detid)+3;
96  } else if ( a.subdetId() == 5 ) {
97  subdet=2;
98  component=tTopo->tobLayer(detid);
99  } else if ( a.subdetId() == 6 ) {
100  subdet=3;
101  component=tTopo->tecSide(detid)==2?tTopo->tecWheel(detid):tTopo->tecWheel(detid)+9;
102  }
103 
104  int index=100+subdet*100+component;
105 
106 
107  histos.push_back(getHisto(subdet+1));
108  histos.push_back(getHisto(index));
109 
110 }
111 
112 TH1F*
113 SiStripPlotGain::getHisto(const long unsigned int& index){
114  if(vTH1.size()<index+1)
115  vTH1.resize(index+1,nullptr);
116 
117  if(vTH1[index]==nullptr){
118  char name[128];
119  sprintf(name,"%lu",index);
120  edm::LogInfo("")<<"[getHisto] creating index " << index << std::endl;
121  vTH1[index]=new TH1F(name,name,150,0.,5.);
122  }
123 
124  return vTH1[index];
125 }
126 
127 void
129  for(size_t i=0;i<vTH1.size();i++)
130  if(vTH1[i]!=nullptr)
131  vTH1[i]->Write();
132 
133  file->Write();
134  file->Close();
135 
136  tkmap->save(false,0,0,"testTkMap.png");
137 
138 }
139 
ContainerIterator getFirstElement(RegistryConstIterator &idet)
unsigned long long getCache(const edm::EventSetup &eSetup)
unsigned int tibLayer(const DetId &id) const
Registry::const_iterator RegistryConstIterator
void getHistos(const uint32_t &detid, const TrackerTopology *tTopo, std::vector< TH1F * > &histos)
static float getApvGain(uint16_t apv, const Range &range)
RegistryConstIterator detid_end
RegistryConstIterator detid_begin
unsigned int tidWheel(const DetId &id) const
void DoAnalysis(const edm::EventSetup &es, const SiStripApvGain &)
unsigned long long cacheID
unsigned int tidSide(const DetId &id) const
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
void save(bool print_total=true, float minval=0., float maxval=0., std::string s="svgmap.svg", int width=1500, int height=800)
Definition: TrackerMap.cc:699
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< ContainerIterator, ContainerIterator > Range
std::vector< TH1F * > vTH1
Detector identifier class for the strip tracker.
Definition: SiStripDetId.h:17
SiStripDetInfoFileReader * fr
ContainerIterator getLastElement(RegistryConstIterator &idet)
TH1F * getHisto(const long unsigned int &index)
TrackerMap * tkmap
double a
Definition: hdecay.h:121
T get() const
Definition: EventSetup.h:62
SiStripPlotGain(const edm::ParameterSet &)
edm::ESHandle< SiStripApvGain > Handle_
void endJob() override
void beginRun(const edm::Run &run, const edm::EventSetup &es) override
unsigned int tecWheel(const DetId &id) const
RegistryPointers getRegistryPointers() const
T const * product() const
Definition: ESHandle.h:86
void fill(int layer, int ring, int nmod, float x)
Definition: TrackerMap.cc:2786
unsigned int tobLayer(const DetId &id) const
Definition: Run.h:44
~SiStripPlotGain() override
unsigned int tecSide(const DetId &id) const