CMS 3D CMS Logo

SiStripPlotGain.cc
Go to the documentation of this file.
2 
4  : gainToken_{esConsumes<edm::Transition::BeginRun>()}, tTopoToken_{esConsumes<edm::Transition::BeginRun>()} {
5  // now do what ever initialization is needed
6  file = new TFile("correlTest.root", "RECREATE");
7  tkmap = new TrackerMap();
8 }
9 
11  if (gainWatcher_.check(es)) {
12  edm::LogInfo("") << "[SiStripPlotGain::beginRun] cacheID " << es.get<SiStripApvGainRcd>().cacheIdentifier()
13  << std::endl;
15  }
16 }
17 
19  edm::LogInfo("") << "[Doanalysis]";
20 
21  std::vector<TH1F *> histos;
22 
23  SiStripApvGain::RegistryPointers p = gain.getRegistryPointers();
25  iter = p.detid_begin;
26  iterE = p.detid_end;
27 
28  float value;
29 
30  // Divide result by d
31  for (; iter != iterE; ++iter) {
32  getHistos(*iter, tTopo, histos);
33  SiStripApvGain::Range range = SiStripApvGain::Range(p.getFirstElement(iter), p.getLastElement(iter));
34 
35  edm::LogInfo("") << "[Doanalysis] detid " << *iter << " range " << range.second - range.first;
36  size_t apv = 0, apvE = (range.second - range.first);
37  for (; apv < apvE; apv += 2) {
38  value = gain.getApvGain(apv, range);
39  tkmap->fill(*iter, value);
40  for (size_t i = 0; i < histos.size(); ++i)
41  histos[i]->Fill(value);
42  }
43  }
44 }
45 
46 void SiStripPlotGain::getHistos(DetId detid, const TrackerTopology &tTopo, std::vector<TH1F *> &histos) {
47  histos.clear();
48 
49  int subdet = -999;
50  int component = -999;
51  if (detid.subdetId() == 3) {
52  subdet = 0;
53  component = tTopo.tibLayer(detid);
54  } else if (detid.subdetId() == 4) {
55  subdet = 1;
56  component = tTopo.tidSide(detid) == 2 ? tTopo.tidWheel(detid) : tTopo.tidWheel(detid) + 3;
57  } else if (detid.subdetId() == 5) {
58  subdet = 2;
59  component = tTopo.tobLayer(detid);
60  } else if (detid.subdetId() == 6) {
61  subdet = 3;
62  component = tTopo.tecSide(detid) == 2 ? tTopo.tecWheel(detid) : tTopo.tecWheel(detid) + 9;
63  }
64 
65  int index = 100 + subdet * 100 + component;
66 
67  histos.push_back(getHisto(subdet + 1));
68  histos.push_back(getHisto(index));
69 }
70 
71 TH1F *SiStripPlotGain::getHisto(const long unsigned int &index) {
72  if (vTH1.size() < index + 1)
73  vTH1.resize(index + 1, nullptr);
74 
75  if (vTH1[index] == nullptr) {
76  char name[128];
77  sprintf(name, "%lu", index);
78  edm::LogInfo("") << "[getHisto] creating index " << index << std::endl;
79  vTH1[index] = new TH1F(name, name, 150, 0., 5.);
80  }
81 
82  return vTH1[index];
83 }
84 
86  for (size_t i = 0; i < vTH1.size(); i++)
87  if (vTH1[i] != nullptr)
88  vTH1[i]->Write();
89 
90  file->Write();
91  file->Close();
92 
93  tkmap->save(false, 0, 0, "testTkMap.png");
94 }
unsigned int tobLayer(const DetId &id) const
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > tTopoToken_
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
Registry::const_iterator RegistryConstIterator
unsigned int tidSide(const DetId &id) const
unsigned int tidWheel(const DetId &id) const
unsigned int tecWheel(const DetId &id) const
void DoAnalysis(const TrackerTopology &tTopo, const SiStripApvGain &)
unsigned int tecSide(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:844
std::pair< ContainerIterator, ContainerIterator > Range
T get() const
Definition: EventSetup.h:79
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
Definition: value.py:1
Log< level::Info, false > LogInfo
Definition: DetId.h:17
TH1F * getHisto(const long unsigned int &index)
edm::ESGetToken< SiStripApvGain, SiStripApvGainRcd > gainToken_
TrackerMap * tkmap
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
histos
Definition: combine.py:4
SiStripPlotGain(const edm::ParameterSet &)
std::vector< TH1F * > vTH1
unsigned int tibLayer(const DetId &id) const
edm::ESWatcher< SiStripApvGainRcd > gainWatcher_
void endJob() override
void beginRun(const edm::Run &run, const edm::EventSetup &es) override
void fill(int layer, int ring, int nmod, float x)
Definition: TrackerMap.cc:3443
Definition: Run.h:45
void getHistos(DetId detid, const TrackerTopology &tTopo, std::vector< TH1F *> &histos)