CMS 3D CMS Logo

SiStripPlotGain.cc
Go to the documentation of this file.
2 
5 
6 SiStripPlotGain::SiStripPlotGain(const edm::ParameterSet &iConfig) : cacheID(0xFFFFFFFF) {
7  // now do what ever initialization is needed
8  if (!edm::Service<SiStripDetInfoFileReader>().isAvailable()) {
9  edm::LogError("TkLayerMap") << "\n------------------------------------------"
10  "\nUnAvailable Service SiStripDetInfoFileReader: please insert in "
11  "the configuration file an instance like"
12  "\n\tprocess.SiStripDetInfoFileReader = "
13  "cms.Service(\"SiStripDetInfoFileReader\")"
14  "\n------------------------------------------";
15  }
16 
18  file = new TFile("correlTest.root", "RECREATE");
19  tkmap = new TrackerMap();
20 }
21 
23 
24 //
25 
27  if (getCache(es) == cacheID)
28  return;
29  cacheID = getCache(es);
30 
31  edm::LogInfo("") << "[SiStripPlotGain::beginRun] cacheID " << cacheID << std::endl;
32 
33  es.get<SiStripApvGainRcd>().get(Handle_);
34  DoAnalysis(es, *Handle_.product());
35 }
36 
38  edm::LogInfo("") << "[Doanalysis]";
39 
40  // Retrieve tracker topology from geometry
42  es.get<TrackerTopologyRcd>().get(tTopoHandle);
43  const TrackerTopology *const tTopo = tTopoHandle.product();
44 
45  std::vector<TH1F *> histos;
46 
49  iter = p.detid_begin;
50  iterE = p.detid_end;
51 
52  float value;
53 
54  // Divide result by d
55  for (; iter != iterE; ++iter) {
56  getHistos(*iter, tTopo, histos);
58 
59  edm::LogInfo("") << "[Doanalysis] detid " << *iter << " range " << range.second - range.first;
60  size_t apv = 0, apvE = (range.second - range.first);
61  for (; apv < apvE; apv += 2) {
62  value = gain.getApvGain(apv, range);
63  tkmap->fill(*iter, value);
64  for (size_t i = 0; i < histos.size(); ++i)
65  histos[i]->Fill(value);
66  }
67  }
68 }
69 
70 void SiStripPlotGain::getHistos(const uint32_t &detid, const TrackerTopology *tTopo, std::vector<TH1F *> &histos) {
71  histos.clear();
72 
73  int subdet = -999;
74  int component = -999;
75  SiStripDetId a(detid);
76  if (a.subdetId() == 3) {
77  subdet = 0;
78  component = tTopo->tibLayer(detid);
79  } else if (a.subdetId() == 4) {
80  subdet = 1;
81  component = tTopo->tidSide(detid) == 2 ? tTopo->tidWheel(detid) : tTopo->tidWheel(detid) + 3;
82  } else if (a.subdetId() == 5) {
83  subdet = 2;
84  component = tTopo->tobLayer(detid);
85  } else if (a.subdetId() == 6) {
86  subdet = 3;
87  component = tTopo->tecSide(detid) == 2 ? tTopo->tecWheel(detid) : tTopo->tecWheel(detid) + 9;
88  }
89 
90  int index = 100 + subdet * 100 + component;
91 
92  histos.push_back(getHisto(subdet + 1));
93  histos.push_back(getHisto(index));
94 }
95 
96 TH1F *SiStripPlotGain::getHisto(const long unsigned int &index) {
97  if (vTH1.size() < index + 1)
98  vTH1.resize(index + 1, nullptr);
99 
100  if (vTH1[index] == nullptr) {
101  char name[128];
102  sprintf(name, "%lu", index);
103  edm::LogInfo("") << "[getHisto] creating index " << index << std::endl;
104  vTH1[index] = new TH1F(name, name, 150, 0., 5.);
105  }
106 
107  return vTH1[index];
108 }
109 
111  for (size_t i = 0; i < vTH1.size(); i++)
112  if (vTH1[i] != nullptr)
113  vTH1[i]->Write();
114 
115  file->Write();
116  file->Close();
117 
118  tkmap->save(false, 0, 0, "testTkMap.png");
119 }
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
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:71
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:2786
unsigned int tobLayer(const DetId &id) const
Definition: Run.h:45
~SiStripPlotGain() override
unsigned int tecSide(const DetId &id) const