CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_4/src/DPGAnalysis/SiStripTools/plugins/APVCyclePhaseMonitor.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    SiStripTools
00004 // Class:      APVCyclePhaseMonitor
00005 // 
00013 //
00014 // Original Author:  Andrea Venturi
00015 //         Created:  Tue Jul 19 11:56:00 CEST 2009
00016 //
00017 //
00018 
00019 
00020 // system include files
00021 #include <memory>
00022 
00023 // user include files
00024 #include "TH1F.h"
00025 #include "TProfile.h"
00026 
00027 #include "FWCore/Framework/interface/Frameworkfwd.h"
00028 #include "FWCore/Framework/interface/EDAnalyzer.h"
00029 
00030 #include "FWCore/Framework/interface/Event.h"
00031 #include "FWCore/Framework/interface/MakerMacros.h"
00032 #include "FWCore/Framework/interface/Run.h"
00033 
00034 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00035 
00036 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00037 
00038 #include "FWCore/ServiceRegistry/interface/Service.h"
00039 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00040 
00041 #include "FWCore/Utilities/interface/InputTag.h"
00042 
00043 #include "DPGAnalysis/SiStripTools/interface/APVCyclePhaseCollection.h"
00044 
00045 //
00046 // class decleration
00047 //
00048 
00049 class APVCyclePhaseMonitor : public edm::EDAnalyzer {
00050  public:
00051     explicit APVCyclePhaseMonitor(const edm::ParameterSet&);
00052     ~APVCyclePhaseMonitor();
00053 
00054 
00055    private:
00056       virtual void beginJob() ;
00057       virtual void beginRun(const edm::Run&, const edm::EventSetup&) ;
00058       virtual void endRun(const edm::Run&, const edm::EventSetup&) ;
00059       virtual void analyze(const edm::Event&, const edm::EventSetup&);
00060       virtual void endJob() ;
00061 
00062       // ----------member data ---------------------------
00063 
00064   const edm::InputTag _apvphasecollection;
00065   std::vector<std::string> _selectedparts;
00066   std::vector<std::string> _selectedvectorparts;
00067   std::map<std::string,TH1F*> _hphases;
00068   std::map<std::string,TH1F*> _hselectedphases;
00069   std::map<std::string,TH1F*> _hselectedphasesvector;
00070   std::map<std::string,TH1F*> _hselectedphasessize;
00071   std::map<std::string,TProfile*> _hphasevsorbit;
00072   std::map<std::string,TProfile*> _hselectedphasevsorbit;
00073   std::map<std::string,TProfile*> _hselectedphasevectorvsorbit;
00074   unsigned int _nevents;
00075 };
00076 
00077 //
00078 // constants, enums and typedefs
00079 //
00080 
00081 //
00082 // static data member definitions
00083 //
00084 
00085 //
00086 // constructors and destructor
00087 //
00088 APVCyclePhaseMonitor::APVCyclePhaseMonitor(const edm::ParameterSet& iConfig):
00089   _apvphasecollection(iConfig.getParameter<edm::InputTag>("apvCyclePhaseCollection")),
00090   _selectedparts(iConfig.getUntrackedParameter<std::vector<std::string> >("selectedPartitions",std::vector<std::string>())),
00091   _selectedvectorparts(iConfig.getUntrackedParameter<std::vector<std::string> >("selectedVectorPartitions",std::vector<std::string>())),
00092   _hphases(),_hselectedphases(),_hselectedphasesvector(),_hselectedphasessize(),
00093   _hphasevsorbit(),_hselectedphasevsorbit(),_hselectedphasevectorvsorbit(),
00094   _nevents(0)
00095 {
00096    //now do what ever initialization is needed
00097 
00098   edm::LogInfo("UsedAPVCyclePhaseCollection") << " APVCyclePhaseCollection " << _apvphasecollection << " used";
00099 
00100 }
00101 
00102 
00103 APVCyclePhaseMonitor::~APVCyclePhaseMonitor()
00104 {
00105  
00106    // do anything here that needs to be done at desctruction time
00107    // (e.g. close files, deallocate resources etc.)
00108 
00109 }
00110 
00111 
00112 //
00113 // member functions
00114 //
00115 
00116 // ------------ method called to for each event  ------------
00117 void
00118 APVCyclePhaseMonitor::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00119 {
00120    using namespace edm;
00121 
00122    _nevents++;
00123 
00124    edm::Handle<APVCyclePhaseCollection> apvphases;
00125    iEvent.getByLabel(_apvphasecollection,apvphases);
00126 
00127    // improve the matchin between default and actual partitions
00128    
00129    edm::Service<TFileService> tfserv;
00130 
00131    for(std::map<std::string,int>::const_iterator phase = apvphases->get().begin(); phase != apvphases->get().end(); ++phase) {
00132 
00133      if(_hphases.find(phase->first)==_hphases.end()) {
00134        char dirname[300];
00135        sprintf(dirname,"run_%d",iEvent.run());
00136        TFileDirectory subrun = tfserv->mkdir(dirname);
00137 
00138        char hname[300];
00139 
00140        sprintf(hname,"phase_%s",phase->first.c_str());
00141        edm::LogInfo("TH1FBeingBooked") << "TH1F " << hname << " being booked" ;
00142        _hphases[phase->first] = subrun.make<TH1F>(hname,hname,70,-0.5,69.5);
00143        _hphases[phase->first]->GetXaxis()->SetTitle("BX mod 70"); _hphases[phase->first]->GetYaxis()->SetTitle("Events");
00144 
00145        sprintf(hname,"phasevsorbit_%s",phase->first.c_str());
00146        edm::LogInfo("TProfileBeingBooked") << "TProfile " << hname << " being booked" ;
00147        _hphasevsorbit[phase->first] = subrun.make<TProfile>(hname,hname,1800,0.,1800*11223);
00148        _hphasevsorbit[phase->first]->SetBit(TH1::kCanRebin);
00149        _hphasevsorbit[phase->first]->GetXaxis()->SetTitle("time [orbit#]"); _hphasevsorbit[phase->first]->GetYaxis()->SetTitle("Phase");
00150        
00151      }
00152      _hphases[phase->first]->Fill(phase->second);
00153      _hphasevsorbit[phase->first]->Fill(iEvent.orbitNumber(),phase->second);
00154    }
00155 
00156    // selected partitions
00157    
00158    for(std::vector<std::string>::const_iterator part=_selectedparts.begin();part!=_selectedparts.end();++part) {
00159      if(_hselectedphases.find(*part)!=_hselectedphases.end()) {
00160        _hselectedphases[*part]->Fill(apvphases->getPhase(*part));
00161      }
00162      if(_hselectedphasevsorbit.find(*part)!=_hselectedphasevsorbit.end()) {
00163        _hselectedphasevsorbit[*part]->Fill(iEvent.orbitNumber(),apvphases->getPhase(*part));
00164      }
00165    }
00166 
00167    for(std::vector<std::string>::const_iterator part=_selectedvectorparts.begin();
00168        part!=_selectedvectorparts.end();++part) {
00169 
00170      const std::vector<int> phases = apvphases->getPhases(*part);
00171 
00172       if(_hselectedphasessize.find(*part)!=_hselectedphasessize.end()) {
00173         _hselectedphasessize[*part]->Fill(phases.size());
00174       }
00175 
00176      for(std::vector<int>::const_iterator phase=phases.begin();phase!=phases.end();++phase) {
00177        if(_hselectedphasesvector.find(*part)!=_hselectedphasesvector.end()) {
00178          _hselectedphasesvector[*part]->Fill(*phase);
00179        }
00180        if(_hselectedphasevectorvsorbit.find(*part)!=_hselectedphasevectorvsorbit.end()) {
00181          _hselectedphasevectorvsorbit[*part]->Fill(iEvent.orbitNumber(),*phase);
00182        }
00183      }
00184 
00185    }
00186 }
00187 
00188 void 
00189 APVCyclePhaseMonitor::beginRun(const edm::Run& iRun, const edm::EventSetup&)
00190 {
00191   
00192   _hphases.clear();
00193   _hselectedphases.clear();
00194   _hselectedphasesvector.clear();
00195   _hselectedphasessize.clear();
00196   _hphasevsorbit.clear();
00197   _hselectedphasevsorbit.clear();
00198   _hselectedphasevectorvsorbit.clear();
00199 
00200   edm::Service<TFileService> tfserv;
00201 
00202   char dirname[300];
00203   sprintf(dirname,"run_%d",iRun.run());
00204   TFileDirectory subrun = tfserv->mkdir(dirname);
00205   
00206   for(std::vector<std::string>::const_iterator part=_selectedparts.begin();part!=_selectedparts.end();++part) {
00207 
00208     char hname[300];
00209     
00210     sprintf(hname,"selected_phase_%s",part->c_str());
00211     edm::LogInfo("SelectedTH1FBeingBooked") << "TH1F " << hname << " being booked" ;
00212     _hselectedphases[*part] = subrun.make<TH1F>(hname,hname,70,-0.5,69.5);
00213     _hselectedphases[*part]->GetXaxis()->SetTitle("BX mod 70"); _hselectedphases[*part]->GetYaxis()->SetTitle("Events");
00214     
00215     sprintf(hname,"selected_phasevsorbit_%s",part->c_str());
00216     edm::LogInfo("SelectedTProfileBeingBooked") << "TProfile " << hname << " being booked" ;
00217     _hselectedphasevsorbit[*part] = subrun.make<TProfile>(hname,hname,1800,0.,1800*11223);
00218     _hselectedphasevsorbit[*part]->SetBit(TH1::kCanRebin);
00219     _hselectedphasevsorbit[*part]->GetXaxis()->SetTitle("time [orbit#]"); 
00220     _hselectedphasevsorbit[*part]->GetYaxis()->SetTitle("Phase");
00221   }
00222 
00223   for(std::vector<std::string>::const_iterator part=_selectedvectorparts.begin();
00224       part!=_selectedvectorparts.end();++part) {
00225 
00226     char hname[300];
00227     
00228     sprintf(hname,"selected_phase_vector_%s",part->c_str());
00229     edm::LogInfo("SelectedVectTH1FBeingBooked") << "TH1F " << hname << " being booked" ;
00230     _hselectedphasesvector[*part] = subrun.make<TH1F>(hname,hname,70,-0.5,69.5);
00231     _hselectedphasesvector[*part]->GetXaxis()->SetTitle("BX mod 70"); 
00232     _hselectedphasesvector[*part]->GetYaxis()->SetTitle("Events");
00233     
00234     sprintf(hname,"selected_phase_size_%s",part->c_str());
00235     edm::LogInfo("SelectedVectSizeTH1FBeingBooked") << "TH1F " << hname << " being booked" ;
00236     _hselectedphasessize[*part] = subrun.make<TH1F>(hname,hname,10,-0.5,9.5);
00237     _hselectedphasessize[*part]->GetXaxis()->SetTitle("Number of Phases"); 
00238     _hselectedphasessize[*part]->GetYaxis()->SetTitle("Events");
00239     
00240     sprintf(hname,"selected_phasevectorvsorbit_%s",part->c_str());
00241     edm::LogInfo("SelectedVectTProfileBeingBooked") << "TProfile " << hname << " being booked" ;
00242     _hselectedphasevectorvsorbit[*part] = subrun.make<TProfile>(hname,hname,1800,0.,1800*11223);
00243     _hselectedphasevectorvsorbit[*part]->SetBit(TH1::kCanRebin);
00244     _hselectedphasevectorvsorbit[*part]->GetXaxis()->SetTitle("time [orbit#]"); 
00245     _hselectedphasevectorvsorbit[*part]->GetYaxis()->SetTitle("Phase");
00246   }
00247 }
00248 
00249 void 
00250 APVCyclePhaseMonitor::endRun(const edm::Run& iRun, const edm::EventSetup&)
00251 {
00252 }
00253 
00254 
00255 // ------------ method called once each job just before starting event loop  ------------
00256 void 
00257 APVCyclePhaseMonitor::beginJob()
00258 {
00259 
00260 }
00261 
00262 // ------------ method called once each job just after ending the event loop  ------------
00263 void 
00264 APVCyclePhaseMonitor::endJob() {
00265 
00266   edm::LogInfo("EndOfJob") << _nevents << " analyzed events";
00267 
00268 }
00269 
00270 
00271 //define this as a plug-in
00272 DEFINE_FWK_MODULE(APVCyclePhaseMonitor);