Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:08:53

0001 #include "DQM/SiStripMonitorSummary/plugins/SiStripPlotGain.h"
0002 
0003 SiStripPlotGain::SiStripPlotGain(const edm::ParameterSet &iConfig)
0004     : gainToken_{esConsumes<edm::Transition::BeginRun>()}, tTopoToken_{esConsumes<edm::Transition::BeginRun>()} {
0005   // now do what ever initialization is needed
0006   file = new TFile("correlTest.root", "RECREATE");
0007   tkmap = new TrackerMap();
0008 }
0009 
0010 void SiStripPlotGain::beginRun(const edm::Run &run, const edm::EventSetup &es) {
0011   if (gainWatcher_.check(es)) {
0012     edm::LogInfo("") << "[SiStripPlotGain::beginRun] cacheID " << es.get<SiStripApvGainRcd>().cacheIdentifier()
0013                      << std::endl;
0014     DoAnalysis(es.getData(tTopoToken_), es.getData(gainToken_));
0015   }
0016 }
0017 
0018 void SiStripPlotGain::DoAnalysis(const TrackerTopology &tTopo, const SiStripApvGain &gain) {
0019   edm::LogInfo("") << "[Doanalysis]";
0020 
0021   std::vector<TH1F *> histos;
0022 
0023   SiStripApvGain::RegistryPointers p = gain.getRegistryPointers();
0024   SiStripApvGain::RegistryConstIterator iter, iterE;
0025   iter = p.detid_begin;
0026   iterE = p.detid_end;
0027 
0028   float value;
0029 
0030   // Divide result by d
0031   for (; iter != iterE; ++iter) {
0032     getHistos(*iter, tTopo, histos);
0033     SiStripApvGain::Range range = SiStripApvGain::Range(p.getFirstElement(iter), p.getLastElement(iter));
0034 
0035     edm::LogInfo("") << "[Doanalysis] detid " << *iter << " range " << range.second - range.first;
0036     size_t apv = 0, apvE = (range.second - range.first);
0037     for (; apv < apvE; apv += 2) {
0038       value = gain.getApvGain(apv, range);
0039       tkmap->fill(*iter, value);
0040       for (size_t i = 0; i < histos.size(); ++i)
0041         histos[i]->Fill(value);
0042     }
0043   }
0044 }
0045 
0046 void SiStripPlotGain::getHistos(DetId detid, const TrackerTopology &tTopo, std::vector<TH1F *> &histos) {
0047   histos.clear();
0048 
0049   int subdet = -999;
0050   int component = -999;
0051   if (detid.subdetId() == 3) {
0052     subdet = 0;
0053     component = tTopo.tibLayer(detid);
0054   } else if (detid.subdetId() == 4) {
0055     subdet = 1;
0056     component = tTopo.tidSide(detid) == 2 ? tTopo.tidWheel(detid) : tTopo.tidWheel(detid) + 3;
0057   } else if (detid.subdetId() == 5) {
0058     subdet = 2;
0059     component = tTopo.tobLayer(detid);
0060   } else if (detid.subdetId() == 6) {
0061     subdet = 3;
0062     component = tTopo.tecSide(detid) == 2 ? tTopo.tecWheel(detid) : tTopo.tecWheel(detid) + 9;
0063   }
0064 
0065   int index = 100 + subdet * 100 + component;
0066 
0067   histos.push_back(getHisto(subdet + 1));
0068   histos.push_back(getHisto(index));
0069 }
0070 
0071 TH1F *SiStripPlotGain::getHisto(const long unsigned int &index) {
0072   if (vTH1.size() < index + 1)
0073     vTH1.resize(index + 1, nullptr);
0074 
0075   if (vTH1[index] == nullptr) {
0076     char name[128];
0077     sprintf(name, "%lu", index);
0078     edm::LogInfo("") << "[getHisto] creating index " << index << std::endl;
0079     vTH1[index] = new TH1F(name, name, 150, 0., 5.);
0080   }
0081 
0082   return vTH1[index];
0083 }
0084 
0085 void SiStripPlotGain::endJob() {
0086   for (size_t i = 0; i < vTH1.size(); i++)
0087     if (vTH1[i] != nullptr)
0088       vTH1[i]->Write();
0089 
0090   file->Write();
0091   file->Close();
0092 
0093   tkmap->save(false, 0, 0, "testTkMap.png");
0094 }