Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:06:36

0001 // -*- C++ -*-
0002 //
0003 // Package:    SiStripTools
0004 // Class:      SiPixelQualityHistory
0005 //
0006 /**\class SiPixelQualityHistory SiPixelQualityHistory.cc DPGAnalysis/SiStripTools/plugins/SiPixelQualityHistory.cc
0007 
0008  Description: <one line class summary>
0009 
0010  Implementation:
0011      <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Andrea Venturi
0015 //         Created:  Tue Sep 18 17:52:00 CEST 2009
0016 //
0017 //
0018 
0019 // system include files
0020 #include <memory>
0021 
0022 // user include files
0023 
0024 #include <vector>
0025 #include <map>
0026 
0027 //#include "TGraph.h"
0028 #include "TH1F.h"
0029 
0030 #include "FWCore/Framework/interface/Frameworkfwd.h"
0031 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0032 
0033 #include "FWCore/Framework/interface/Event.h"
0034 #include "FWCore/Framework/interface/MakerMacros.h"
0035 #include "FWCore/Framework/interface/Run.h"
0036 
0037 #include "FWCore/Framework/interface/ESHandle.h"
0038 #include "FWCore/Framework/interface/EventSetup.h"
0039 
0040 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0041 
0042 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0043 
0044 #include "FWCore/ServiceRegistry/interface/Service.h"
0045 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0046 
0047 #include "FWCore/Utilities/interface/InputTag.h"
0048 
0049 #include "CondFormats/SiPixelObjects/interface/SiPixelQuality.h"
0050 #include "CondFormats/DataRecord/interface/SiPixelQualityRcd.h"
0051 
0052 #include "DPGAnalysis/SiStripTools/interface/RunHistogramManager.h"
0053 //
0054 // class decleration
0055 //
0056 
0057 class SiPixelQualityHistory : public edm::one::EDAnalyzer<edm::one::WatchRuns, edm::one::SharedResources> {
0058 public:
0059   explicit SiPixelQualityHistory(const edm::ParameterSet&);
0060   ~SiPixelQualityHistory() override;
0061 
0062   enum { Summary, Module, ROC };
0063 
0064 private:
0065   void beginJob() override;
0066   void beginRun(const edm::Run&, const edm::EventSetup&) override;
0067   void endRun(const edm::Run&, const edm::EventSetup&) override {}
0068   void analyze(const edm::Event&, const edm::EventSetup&) override;
0069   void endJob() override;
0070 
0071   // ----------member data ---------------------------
0072 
0073   RunHistogramManager m_rhm;
0074   const std::vector<edm::ParameterSet> m_monitoredspq;
0075   std::vector<edm::ESGetToken<SiPixelQuality, SiPixelQualityRcd>> m_spqTokens;
0076   const unsigned int m_mode;
0077   const bool m_run;
0078   const unsigned int m_maxLS;
0079   const unsigned int m_LSfrac;
0080   //  std::map<std::string,TGraph*> m_history;
0081   std::map<std::string, TH1F*> m_history;
0082   std::map<std::string, TProfile**> m_badmodrun;
0083 };
0084 
0085 //
0086 // constants, enums and typedefs
0087 //
0088 
0089 //
0090 // static data member definitions
0091 //
0092 
0093 //
0094 // constructors and destructor
0095 //
0096 SiPixelQualityHistory::SiPixelQualityHistory(const edm::ParameterSet& iConfig)
0097     : m_rhm(consumesCollector()),
0098       m_monitoredspq(iConfig.getParameter<std::vector<edm::ParameterSet>>("monitoredSiPixelQuality")),
0099       m_mode(iConfig.getUntrackedParameter<unsigned int>("granularityMode", Module)),
0100       m_run(iConfig.getParameter<bool>("runProcess")),
0101       m_maxLS(iConfig.getUntrackedParameter<unsigned int>("maxLSBeforeRebin", 100)),
0102       m_LSfrac(iConfig.getUntrackedParameter<unsigned int>("startingLSFraction", 4)),
0103       m_history(),
0104       m_badmodrun() {
0105   //now do what ever initialization is needed
0106   usesResource(TFileService::kSharedResource);
0107 
0108   edm::Service<TFileService> tfserv;
0109 
0110   for (const auto& ps : m_monitoredspq) {
0111     m_spqTokens.emplace_back(
0112         esConsumes<edm::Transition::BeginRun>(edm::ESInputTag{"", ps.getParameter<std::string>("spqLabel")}));
0113 
0114     std::string name = ps.getParameter<std::string>("name");
0115 
0116     if (m_run)
0117       m_history[name] = tfserv->make<TH1F>(name.c_str(), name.c_str(), 10, 0, 10);
0118 
0119     char hrunname[400];
0120     sprintf(hrunname, "badmodrun_%s", name.c_str());
0121     char hruntitle[400];
0122     sprintf(hruntitle, "Number of bad modules %s", name.c_str());
0123     m_badmodrun[name] = m_rhm.makeTProfile(hrunname, hruntitle, m_LSfrac * m_maxLS, 0, m_maxLS * 262144);
0124   }
0125 }
0126 
0127 SiPixelQualityHistory::~SiPixelQualityHistory() {
0128   // do anything here that needs to be done at desctruction time
0129   // (e.g. close files, deallocate resources etc.)
0130 }
0131 
0132 //
0133 // member functions
0134 //
0135 
0136 // ------------ method called to for each event  ------------
0137 void SiPixelQualityHistory::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0138   //  edm::LogInfo("EventProcessing") << "event being processed";
0139 
0140   for (std::size_t iMon = 0; iMon != m_monitoredspq.size(); ++iMon) {
0141     std::string name = m_monitoredspq[iMon].getParameter<std::string>("name");
0142     const auto& spq = iSetup.getData(m_spqTokens[iMon]);
0143 
0144     int nbad = 0;
0145 
0146     if (m_mode == Summary) {
0147       //      nbad = spq.BadModuleNumber();
0148 
0149     } else {
0150       std::vector<SiPixelQuality::disabledModuleType> bads = spq.getBadComponentList();
0151 
0152       LogDebug("BadComponents") << bads.size() << " bad components found";
0153 
0154       for (const auto& bc : bads) {
0155         if (m_mode == Module) {
0156           if (spq.IsModuleBad(bc.DetID))
0157             ++nbad;
0158           //      if(bc.errorType==0) ++nbad;
0159         } else if (m_mode == ROC) {
0160           for (int roc = 1; roc < 2 * 2 * 2 * 2 * 2 * 2 * 2 + 1; roc *= 2) {
0161             if ((bc.BadRocs & roc) > 0)
0162               ++nbad;
0163           }
0164         }
0165       }
0166     }
0167     if (m_badmodrun.find(name) != m_badmodrun.end() && m_badmodrun[name] && *m_badmodrun[name]) {
0168       (*m_badmodrun[name])->Fill(iEvent.orbitNumber(), nbad);
0169     }
0170   }
0171 }
0172 
0173 void SiPixelQualityHistory::beginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {
0174   m_rhm.beginRun(iRun);
0175 
0176   // loop on all the SiPixelQuality objects to be monitored
0177   for (std::size_t iMon = 0; iMon != m_monitoredspq.size(); ++iMon) {
0178     const auto& ps = m_monitoredspq[iMon];
0179     std::string name = ps.getParameter<std::string>("name");
0180 
0181     if (m_badmodrun.find(name) != m_badmodrun.end()) {
0182       if (m_badmodrun[name] && *m_badmodrun[name]) {
0183         (*m_badmodrun[name])->SetCanExtend(TH1::kXaxis);
0184         (*m_badmodrun[name])->GetXaxis()->SetTitle("time [Orb#]");
0185         (*m_badmodrun[name])->GetYaxis()->SetTitle("bad components");
0186       }
0187     }
0188 
0189     if (m_run) {
0190       const auto& spq = iSetup.getData(m_spqTokens[iMon]);
0191 
0192       int nbad = 0;
0193 
0194       if (m_mode == Summary) {
0195         //  nbad = spq.BadModuleNumber();
0196 
0197       } else {
0198         std::vector<SiPixelQuality::disabledModuleType> bads = spq.getBadComponentList();
0199 
0200         LogDebug("BadComponents") << bads.size() << " bad components found";
0201 
0202         for (const auto& bc : bads) {
0203           if (m_mode == Module) {
0204             if (spq.IsModuleBad(bc.DetID))
0205               ++nbad;
0206             //    if(bc.errorType==0) ++nbad;
0207           } else if (m_mode == ROC) {
0208             for (int roc = 1; roc < 2 * 2 * 2 * 2 * 2 * 2 * 2 + 1; roc *= 2) {
0209               if ((bc.BadRocs & roc) > 0)
0210                 ++nbad;
0211             }
0212           }
0213         }
0214       }
0215       char runname[100];
0216       sprintf(runname, "%d", iRun.run());
0217       LogDebug("AnalyzedRun") << name << " " << runname << " " << nbad;
0218       m_history[name]->Fill(runname, nbad);
0219     }
0220   }
0221 }
0222 
0223 // ------------ method called once each job just before starting event loop  ------------
0224 void SiPixelQualityHistory::beginJob() {}
0225 
0226 // ------------ method called once each job just after ending the event loop  ------------
0227 void SiPixelQualityHistory::endJob() {
0228   /*
0229   for(std::vector<edm::ParameterSet>::const_iterator ps=m_monitoredspq.begin();ps!=m_monitoredspq.end();++ps) {
0230 
0231     std::string name = ps->getParameter<std::string>("name");
0232     m_history[name]->Write();
0233 
0234   }
0235   */
0236 }
0237 
0238 //define this as a plug-in
0239 DEFINE_FWK_MODULE(SiPixelQualityHistory);