CMS 3D CMS Logo

SiStripPlotGain.cc
Go to the documentation of this file.
2 
6 
7 SiStripPlotGain::SiStripPlotGain(const edm::ParameterSet &iConfig) : cacheID(0xFFFFFFFF) {
8  // now do what ever initialization is needed
9  file = new TFile("correlTest.root", "RECREATE");
10  tkmap = new TrackerMap();
11 }
12 
14 
15 //
16 
18  if (getCache(es) == cacheID)
19  return;
20  cacheID = getCache(es);
21 
22  edm::LogInfo("") << "[SiStripPlotGain::beginRun] cacheID " << cacheID << std::endl;
23 
24  es.get<SiStripApvGainRcd>().get(Handle_);
25  DoAnalysis(es, *Handle_.product());
26 }
27 
29  edm::LogInfo("") << "[Doanalysis]";
30 
31  // Retrieve tracker topology from geometry
33  es.get<TrackerTopologyRcd>().get(tTopoHandle);
34  const TrackerTopology *const tTopo = tTopoHandle.product();
35 
36  std::vector<TH1F *> histos;
37 
40  iter = p.detid_begin;
41  iterE = p.detid_end;
42 
43  float value;
44 
45  // Divide result by d
46  for (; iter != iterE; ++iter) {
47  getHistos(*iter, tTopo, histos);
49 
50  edm::LogInfo("") << "[Doanalysis] detid " << *iter << " range " << range.second - range.first;
51  size_t apv = 0, apvE = (range.second - range.first);
52  for (; apv < apvE; apv += 2) {
53  value = gain.getApvGain(apv, range);
54  tkmap->fill(*iter, value);
55  for (size_t i = 0; i < histos.size(); ++i)
56  histos[i]->Fill(value);
57  }
58  }
59 }
60 
61 void SiStripPlotGain::getHistos(const uint32_t &detid, const TrackerTopology *tTopo, std::vector<TH1F *> &histos) {
62  histos.clear();
63 
64  int subdet = -999;
65  int component = -999;
66  SiStripDetId a(detid);
67  if (a.subdetId() == 3) {
68  subdet = 0;
69  component = tTopo->tibLayer(detid);
70  } else if (a.subdetId() == 4) {
71  subdet = 1;
72  component = tTopo->tidSide(detid) == 2 ? tTopo->tidWheel(detid) : tTopo->tidWheel(detid) + 3;
73  } else if (a.subdetId() == 5) {
74  subdet = 2;
75  component = tTopo->tobLayer(detid);
76  } else if (a.subdetId() == 6) {
77  subdet = 3;
78  component = tTopo->tecSide(detid) == 2 ? tTopo->tecWheel(detid) : tTopo->tecWheel(detid) + 9;
79  }
80 
81  int index = 100 + subdet * 100 + component;
82 
83  histos.push_back(getHisto(subdet + 1));
84  histos.push_back(getHisto(index));
85 }
86 
87 TH1F *SiStripPlotGain::getHisto(const long unsigned int &index) {
88  if (vTH1.size() < index + 1)
89  vTH1.resize(index + 1, nullptr);
90 
91  if (vTH1[index] == nullptr) {
92  char name[128];
93  sprintf(name, "%lu", index);
94  edm::LogInfo("") << "[getHisto] creating index " << index << std::endl;
95  vTH1[index] = new TH1F(name, name, 150, 0., 5.);
96  }
97 
98  return vTH1[index];
99 }
100 
102  for (size_t i = 0; i < vTH1.size(); i++)
103  if (vTH1[i] != nullptr)
104  vTH1[i]->Write();
105 
106  file->Write();
107  file->Close();
108 
109  tkmap->save(false, 0, 0, "testTkMap.png");
110 }
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:811
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
std::pair< ContainerIterator, ContainerIterator > Range
Detector identifier class for the strip tracker.
Definition: SiStripDetId.h:18
ContainerIterator getLastElement(RegistryConstIterator &idet)
TH1F * getHisto(const long unsigned int &index)
TrackerMap * tkmap
histos
Definition: combine.py:4
double a
Definition: hdecay.h:119
T get() const
Definition: EventSetup.h:73
SiStripPlotGain(const edm::ParameterSet &)
std::vector< TH1F * > vTH1
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:3289
unsigned int tobLayer(const DetId &id) const
Definition: Run.h:45
~SiStripPlotGain() override
unsigned int tecSide(const DetId &id) const