Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 10:53:23

0001 // -*- C++ -*-
0002 //
0003 // Package:    SiStripTools
0004 // Class:      APVCyclePhaseMonitor
0005 //
0006 /**\class APVCyclePhaseMonitor APVCyclePhaseMonitor.cc DPGAnalysis/SiStripTools/plugins/APVCyclePhaseMonitor.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 Jul 19 11:56:00 CEST 2009
0016 //
0017 //
0018 
0019 // system include files
0020 #include <memory>
0021 
0022 // user include files
0023 #include "TH1F.h"
0024 #include "TProfile.h"
0025 
0026 #include "FWCore/Framework/interface/Frameworkfwd.h"
0027 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0028 
0029 #include "FWCore/Framework/interface/Event.h"
0030 #include "FWCore/Framework/interface/MakerMacros.h"
0031 #include "FWCore/Framework/interface/Run.h"
0032 
0033 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0034 
0035 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0036 
0037 #include "FWCore/ServiceRegistry/interface/Service.h"
0038 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0039 
0040 #include "FWCore/Utilities/interface/InputTag.h"
0041 
0042 #include "DPGAnalysis/SiStripTools/interface/APVCyclePhaseCollection.h"
0043 
0044 #include "DPGAnalysis/SiStripTools/interface/RunHistogramManager.h"
0045 
0046 //
0047 // class decleration
0048 //
0049 
0050 class APVCyclePhaseMonitor : public edm::one::EDAnalyzer<edm::one::WatchRuns, edm::one::SharedResources> {
0051 public:
0052   explicit APVCyclePhaseMonitor(const edm::ParameterSet&);
0053   ~APVCyclePhaseMonitor() override;
0054 
0055 private:
0056   void beginJob() override;
0057   void beginRun(const edm::Run&, const edm::EventSetup&) override;
0058   void endRun(const edm::Run&, const edm::EventSetup&) override;
0059   void analyze(const edm::Event&, const edm::EventSetup&) override;
0060   void endJob() override;
0061 
0062   // ----------member data ---------------------------
0063 
0064   edm::EDGetTokenT<APVCyclePhaseCollection> _apvphasecollectionToken;
0065   std::vector<std::string> _selectedparts;
0066   std::vector<std::string> _selectedvectorparts;
0067   const unsigned int m_maxLS;
0068   const unsigned int m_LSfrac;
0069   RunHistogramManager m_rhm;
0070   std::map<std::string, TH1F*> _hphases;
0071   std::map<std::string, TH1F**> _hselectedphases;
0072   std::map<std::string, TH1F**> _hselectedphasesvector;
0073   std::map<std::string, TH1F**> _hselectedphasessize;
0074   std::map<std::string, TProfile*> _hphasevsorbit;
0075   std::map<std::string, TProfile**> _hselectedphasevsorbit;
0076   std::map<std::string, TProfile**> _hselectedphasevectorvsorbit;
0077   unsigned int _nevents;
0078 };
0079 
0080 //
0081 // constants, enums and typedefs
0082 //
0083 
0084 //
0085 // static data member definitions
0086 //
0087 
0088 //
0089 // constructors and destructor
0090 //
0091 APVCyclePhaseMonitor::APVCyclePhaseMonitor(const edm::ParameterSet& iConfig)
0092     : _apvphasecollectionToken(
0093           consumes<APVCyclePhaseCollection>(iConfig.getParameter<edm::InputTag>("apvCyclePhaseCollection"))),
0094       _selectedparts(
0095           iConfig.getUntrackedParameter<std::vector<std::string> >("selectedPartitions", std::vector<std::string>())),
0096       _selectedvectorparts(iConfig.getUntrackedParameter<std::vector<std::string> >("selectedVectorPartitions",
0097                                                                                     std::vector<std::string>())),
0098       m_maxLS(iConfig.getUntrackedParameter<unsigned int>("maxLSBeforeRebin", 125)),
0099       m_LSfrac(iConfig.getUntrackedParameter<unsigned int>("startingLSFraction", 16)),
0100       m_rhm(consumesCollector()),
0101       _hphases(),
0102       _hselectedphases(),
0103       _hselectedphasesvector(),
0104       _hselectedphasessize(),
0105       _hphasevsorbit(),
0106       _hselectedphasevsorbit(),
0107       _hselectedphasevectorvsorbit(),
0108       _nevents(0) {
0109   //now do what ever initialization is needed
0110 
0111   usesResource(TFileService::kSharedResource);
0112 
0113   edm::LogInfo("UsedAPVCyclePhaseCollection")
0114       << " APVCyclePhaseCollection " << iConfig.getParameter<edm::InputTag>("apvCyclePhaseCollection") << " used";
0115 
0116   for (std::vector<std::string>::const_iterator part = _selectedparts.begin(); part != _selectedparts.end(); ++part) {
0117     char hname[300];
0118 
0119     sprintf(hname, "selected_phase_%s", part->c_str());
0120     edm::LogInfo("SelectedTH1FBeingBooked") << "TH1F " << hname << " being booked";
0121     _hselectedphases[*part] = m_rhm.makeTH1F(hname, hname, 70, -0.5, 69.5);
0122 
0123     sprintf(hname, "selected_phasevsorbit_%s", part->c_str());
0124     edm::LogInfo("SelectedTProfileBeingBooked") << "TProfile " << hname << " being booked";
0125     _hselectedphasevsorbit[*part] = m_rhm.makeTProfile(hname, hname, m_LSfrac * m_maxLS, 0, m_maxLS * 262144);
0126   }
0127 
0128   for (std::vector<std::string>::const_iterator part = _selectedvectorparts.begin(); part != _selectedvectorparts.end();
0129        ++part) {
0130     char hname[300];
0131 
0132     sprintf(hname, "selected_phase_vector_%s", part->c_str());
0133     edm::LogInfo("SelectedVectTH1FBeingBooked") << "TH1F " << hname << " being booked";
0134     _hselectedphasesvector[*part] = m_rhm.makeTH1F(hname, hname, 70, -0.5, 69.5);
0135 
0136     sprintf(hname, "selected_phase_size_%s", part->c_str());
0137     edm::LogInfo("SelectedVectSizeTH1FBeingBooked") << "TH1F " << hname << " being booked";
0138     _hselectedphasessize[*part] = m_rhm.makeTH1F(hname, hname, 10, -0.5, 9.5);
0139 
0140     sprintf(hname, "selected_phasevectorvsorbit_%s", part->c_str());
0141     edm::LogInfo("SelectedVectTProfileBeingBooked") << "TProfile " << hname << " being booked";
0142     _hselectedphasevectorvsorbit[*part] = m_rhm.makeTProfile(hname, hname, m_LSfrac * m_maxLS, 0, m_maxLS * 262144);
0143   }
0144 }
0145 
0146 APVCyclePhaseMonitor::~APVCyclePhaseMonitor() {
0147   // do anything here that needs to be done at desctruction time
0148   // (e.g. close files, deallocate resources etc.)
0149 }
0150 
0151 //
0152 // member functions
0153 //
0154 
0155 // ------------ method called to for each event  ------------
0156 void APVCyclePhaseMonitor::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0157   using namespace edm;
0158 
0159   _nevents++;
0160 
0161   edm::Handle<APVCyclePhaseCollection> apvphases;
0162   iEvent.getByToken(_apvphasecollectionToken, apvphases);
0163 
0164   // improve the matchin between default and actual partitions
0165 
0166   edm::Service<TFileService> tfserv;
0167 
0168   for (std::map<std::string, int>::const_iterator phase = apvphases->get().begin(); phase != apvphases->get().end();
0169        ++phase) {
0170     if (_hphases.find(phase->first) == _hphases.end()) {
0171       char dirname[300];
0172       sprintf(dirname, "run_%d", iEvent.run());
0173       TFileDirectory subrun = tfserv->mkdir(dirname);
0174 
0175       char hname[300];
0176 
0177       sprintf(hname, "phase_%s", phase->first.c_str());
0178       edm::LogInfo("TH1FBeingBooked") << "TH1F " << hname << " being booked";
0179       _hphases[phase->first] = subrun.make<TH1F>(hname, hname, 70, -0.5, 69.5);
0180       _hphases[phase->first]->GetXaxis()->SetTitle("BX mod 70");
0181       _hphases[phase->first]->GetYaxis()->SetTitle("Events");
0182 
0183       sprintf(hname, "phasevsorbit_%s", phase->first.c_str());
0184       edm::LogInfo("TProfileBeingBooked") << "TProfile " << hname << " being booked";
0185       _hphasevsorbit[phase->first] = subrun.make<TProfile>(hname, hname, m_LSfrac * m_maxLS, 0, m_maxLS * 262144);
0186       _hphasevsorbit[phase->first]->SetCanExtend(TH1::kXaxis);
0187       _hphasevsorbit[phase->first]->GetXaxis()->SetTitle("time [orbit#]");
0188       _hphasevsorbit[phase->first]->GetYaxis()->SetTitle("Phase");
0189     }
0190     _hphases[phase->first]->Fill(phase->second);
0191     _hphasevsorbit[phase->first]->Fill(iEvent.orbitNumber(), phase->second);
0192   }
0193 
0194   // selected partitions
0195 
0196   for (std::vector<std::string>::const_iterator part = _selectedparts.begin(); part != _selectedparts.end(); ++part) {
0197     if (_hselectedphases.find(*part) != _hselectedphases.end() && _hselectedphases[*part] && *_hselectedphases[*part]) {
0198       (*_hselectedphases[*part])->Fill(apvphases->getPhase(*part));
0199     }
0200     if (_hselectedphasevsorbit.find(*part) != _hselectedphasevsorbit.end() && _hselectedphasevsorbit[*part] &&
0201         *_hselectedphasevsorbit[*part]) {
0202       (*_hselectedphasevsorbit[*part])->Fill(iEvent.orbitNumber(), apvphases->getPhase(*part));
0203     }
0204   }
0205 
0206   for (std::vector<std::string>::const_iterator part = _selectedvectorparts.begin(); part != _selectedvectorparts.end();
0207        ++part) {
0208     const std::vector<int> phases = apvphases->getPhases(*part);
0209 
0210     if (_hselectedphasessize.find(*part) != _hselectedphasessize.end() && _hselectedphasessize[*part] &&
0211         *_hselectedphasessize[*part]) {
0212       (*_hselectedphasessize[*part])->Fill(phases.size());
0213     }
0214 
0215     for (std::vector<int>::const_iterator phase = phases.begin(); phase != phases.end(); ++phase) {
0216       if (_hselectedphasesvector.find(*part) != _hselectedphasesvector.end() && _hselectedphasesvector[*part] &&
0217           *_hselectedphasesvector[*part]) {
0218         (*_hselectedphasesvector[*part])->Fill(*phase);
0219       }
0220       if (_hselectedphasevectorvsorbit.find(*part) != _hselectedphasevectorvsorbit.end() &&
0221           _hselectedphasevectorvsorbit[*part] && *_hselectedphasevectorvsorbit[*part]) {
0222         (*_hselectedphasevectorvsorbit[*part])->Fill(iEvent.orbitNumber(), *phase);
0223       }
0224     }
0225   }
0226 }
0227 
0228 void APVCyclePhaseMonitor::beginRun(const edm::Run& iRun, const edm::EventSetup&) {
0229   _hphases.clear();
0230   _hphasevsorbit.clear();
0231 
0232   m_rhm.beginRun(iRun);
0233 
0234   for (std::map<std::string, TH1F**>::const_iterator hist = _hselectedphases.begin(); hist != _hselectedphases.end();
0235        ++hist) {
0236     if (*(hist->second)) {
0237       (*(hist->second))->GetXaxis()->SetTitle("BX mod 70");
0238       (*(hist->second))->GetYaxis()->SetTitle("Events");
0239     }
0240   }
0241   for (std::map<std::string, TProfile**>::const_iterator prof = _hselectedphasevsorbit.begin();
0242        prof != _hselectedphasevsorbit.end();
0243        ++prof) {
0244     if (*(prof->second)) {
0245       (*(prof->second))->SetCanExtend(TH1::kXaxis);
0246       (*(prof->second))->GetXaxis()->SetTitle("time [orbit#]");
0247       (*(prof->second))->GetYaxis()->SetTitle("Phase");
0248     }
0249   }
0250   for (std::map<std::string, TH1F**>::const_iterator hist = _hselectedphasesvector.begin();
0251        hist != _hselectedphasesvector.end();
0252        ++hist) {
0253     if (*(hist->second)) {
0254       (*(hist->second))->GetXaxis()->SetTitle("BX mod 70");
0255       (*(hist->second))->GetYaxis()->SetTitle("Events");
0256     }
0257   }
0258   for (std::map<std::string, TH1F**>::const_iterator hist = _hselectedphasessize.begin();
0259        hist != _hselectedphasessize.end();
0260        ++hist) {
0261     if (*(hist->second)) {
0262       (*(hist->second))->GetXaxis()->SetTitle("Number of Phases");
0263       (*(hist->second))->GetYaxis()->SetTitle("Events");
0264     }
0265   }
0266   for (std::map<std::string, TProfile**>::const_iterator prof = _hselectedphasevectorvsorbit.begin();
0267        prof != _hselectedphasevectorvsorbit.end();
0268        ++prof) {
0269     if (*(prof->second)) {
0270       (*(prof->second))->SetCanExtend(TH1::kXaxis);
0271       (*(prof->second))->GetXaxis()->SetTitle("time [orbit#]");
0272       (*(prof->second))->GetYaxis()->SetTitle("Phase");
0273     }
0274   }
0275 }
0276 
0277 void APVCyclePhaseMonitor::endRun(const edm::Run& iRun, const edm::EventSetup&) {}
0278 
0279 // ------------ method called once each job just before starting event loop  ------------
0280 void APVCyclePhaseMonitor::beginJob() {}
0281 
0282 // ------------ method called once each job just after ending the event loop  ------------
0283 void APVCyclePhaseMonitor::endJob() { edm::LogInfo("EndOfJob") << _nevents << " analyzed events"; }
0284 
0285 //define this as a plug-in
0286 DEFINE_FWK_MODULE(APVCyclePhaseMonitor);