CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/DPGAnalysis/SiStripTools/plugins/EventTimeDistribution.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    SiStripTools
00004 // Class:      EventTimeDistribution
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 <string>
00025 
00026 #include "TH1F.h"
00027 #include "TH2F.h"
00028 
00029 #include "FWCore/Framework/interface/Frameworkfwd.h"
00030 #include "FWCore/Framework/interface/EDAnalyzer.h"
00031 
00032 #include "FWCore/Framework/interface/Event.h"
00033 #include "FWCore/Framework/interface/MakerMacros.h"
00034 #include "FWCore/Framework/interface/Run.h"
00035 
00036 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00037 
00038 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00039 
00040 #include "FWCore/ServiceRegistry/interface/Service.h"
00041 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00042 
00043 #include "FWCore/Utilities/interface/InputTag.h"
00044 
00045 #include "DPGAnalysis/SiStripTools/interface/EventWithHistory.h"
00046 #include "DPGAnalysis/SiStripTools/interface/APVCyclePhaseCollection.h"
00047 
00048 #include "DPGAnalysis/SiStripTools/interface/RunHistogramManager.h"
00049 //
00050 // class decleration
00051 //
00052 
00053 class EventTimeDistribution : public edm::EDAnalyzer {
00054  public:
00055     explicit EventTimeDistribution(const edm::ParameterSet&);
00056     ~EventTimeDistribution();
00057 
00058 
00059    private:
00060       virtual void beginJob() ;
00061       virtual void beginRun(const edm::Run&, const edm::EventSetup&) ;
00062       virtual void endRun(const edm::Run&, const edm::EventSetup&) ;
00063       virtual void analyze(const edm::Event&, const edm::EventSetup&);
00064       virtual void endJob() ;
00065 
00066       // ----------member data ---------------------------
00067 
00068   const edm::InputTag _historyProduct;
00069   const edm::InputTag _apvphasecoll;
00070   const std::string _phasepart;
00071   const bool _wantdbxvsbxincycle;
00072   const bool _wantdbxvsbx;
00073   const bool _wantbxincyclevsbx;
00074   const bool _wantorbitvsbxincycle;
00075   unsigned int _nevents;
00076   const double _binsize;
00077 
00078   RunHistogramManager _rhm;
00079 
00080   TH1F** _dbx;
00081   TH1F** _bx;
00082   TH1F** _bxincycle;
00083   TH1F** _orbit;
00084   TH2F** _dbxvsbxincycle;
00085   TH2F** _dbxvsbx;
00086   TH2F** _bxincyclevsbx;
00087   TH2F** _orbitvsbxincycle;
00088 
00089 };
00090 
00091 //
00092 // constants, enums and typedefs
00093 //
00094 
00095 //
00096 // static data member definitions
00097 //
00098 
00099 //
00100 // constructors and destructor
00101 //
00102 EventTimeDistribution::EventTimeDistribution(const edm::ParameterSet& iConfig):
00103   _historyProduct(iConfig.getParameter<edm::InputTag>("historyProduct")),
00104   _apvphasecoll(iConfig.getParameter<edm::InputTag>("apvPhaseCollection")),
00105   _phasepart(iConfig.getUntrackedParameter<std::string>("phasePartition","None")),
00106   _wantdbxvsbxincycle(iConfig.getUntrackedParameter<bool>("wantDBXvsBXincycle",false)),
00107   _wantdbxvsbx(iConfig.getUntrackedParameter<bool>("wantDBXvsBX",false)),
00108   _wantbxincyclevsbx(iConfig.getUntrackedParameter<bool>("wantBXincyclevsBX",false)),
00109   _wantorbitvsbxincycle(iConfig.getUntrackedParameter<bool>("wantOrbitvsBXincycle",false)),
00110   _nevents(0),
00111   _binsize(iConfig.getUntrackedParameter<double>("minBinSizeInSec",1.)),
00112   _rhm(),
00113   _dbxvsbxincycle(0),   _dbxvsbx(0),   _bxincyclevsbx(0),   _orbitvsbxincycle(0)
00114 {
00115    //now do what ever initialization is needed
00116 
00117   _dbx = _rhm.makeTH1F("dbx","dbx",1000,-0.5,999.5);
00118   _bx = _rhm.makeTH1F("bx","BX number",3564,-0.5,3563.5);
00119   _bxincycle = _rhm.makeTH1F("bxcycle","bxcycle",70,-0.5,69.5);
00120   _orbit = _rhm.makeTH1F("orbit","orbit",3600,0,11223*_binsize*3600);
00121   if(_wantdbxvsbxincycle) _dbxvsbxincycle = _rhm.makeTH2F("dbxvsbxincycle","dbxvsbxincycle",70,-0.5,69.5,1000,-0.5,999.5);
00122   if(_wantdbxvsbx) _dbxvsbx = _rhm.makeTH2F("dbxvsbx","dbxvsbx",3564,-0.5,3563.5,1000,-0.5,999.5);
00123   if(_wantbxincyclevsbx) _bxincyclevsbx = _rhm.makeTH2F("bxincyclevsbx","bxincyclevsbx",3564,-0.5,3563.5,70,-0.5,69.5);
00124   if(_wantorbitvsbxincycle) _orbitvsbxincycle = _rhm.makeTH2F("orbitvsbxincycle","orbitvsbxincycle",70,-0.5,69.5,3600,0,11223*_binsize*3600);
00125 
00126   edm::LogInfo("UsedAPVCyclePhaseCollection") << " APVCyclePhaseCollection " << _apvphasecoll << " used";
00127 
00128 }
00129 
00130 
00131 EventTimeDistribution::~EventTimeDistribution()
00132 {
00133  
00134    // do anything here that needs to be done at desctruction time
00135    // (e.g. close files, deallocate resources etc.)
00136 
00137 }
00138 
00139 
00140 //
00141 // member functions
00142 //
00143 
00144 // ------------ method called to for each event  ------------
00145 void
00146 EventTimeDistribution::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00147 {
00148    using namespace edm;
00149 
00150    _nevents++;
00151 
00152    edm::Handle<EventWithHistory> he;
00153    iEvent.getByLabel(_historyProduct,he);
00154 
00155    edm::Handle<APVCyclePhaseCollection> apvphase;
00156    iEvent.getByLabel(_apvphasecoll,apvphase);
00157 
00158    long long tbx = he->absoluteBX();
00159    if(apvphase.isValid() && !apvphase.failedToGet()) {
00160      const int thephase = apvphase->getPhase(_phasepart); 
00161      if(thephase!=APVCyclePhaseCollection::invalid &&
00162         thephase!=APVCyclePhaseCollection::multiphase &&
00163         thephase!=APVCyclePhaseCollection::nopartition)
00164        tbx -= thephase;
00165    }
00166 
00167 
00168 
00169    // improve the matchin between default and actual partitions
00170    
00171    (*_dbx)->Fill(he->deltaBX());
00172    (*_bx)->Fill(iEvent.bunchCrossing());
00173    (*_bxincycle)->Fill(tbx%70);
00174    (*_orbit)->Fill(iEvent.orbitNumber());
00175    if(_dbxvsbxincycle && *_dbxvsbxincycle) (*_dbxvsbxincycle)->Fill(tbx%70,he->deltaBX());
00176    if(_dbxvsbx && *_dbxvsbx) (*_dbxvsbx)->Fill(iEvent.bunchCrossing(),he->deltaBX());
00177    if(_bxincyclevsbx && *_bxincyclevsbx) (*_bxincyclevsbx)->Fill(iEvent.bunchCrossing(),tbx%70);
00178    if(_orbitvsbxincycle && *_orbitvsbxincycle) (*_orbitvsbxincycle)->Fill(tbx%70,iEvent.orbitNumber());
00179 
00180 
00181 }
00182 
00183 void 
00184 EventTimeDistribution::beginRun(const edm::Run& iRun, const edm::EventSetup&)
00185 {
00186 
00187   _rhm.beginRun(iRun);
00188   if(*_dbx) {    (*_dbx)->GetXaxis()->SetTitle("#DeltaBX"); }
00189 
00190   if(*_bx) { (*_bx)->GetXaxis()->SetTitle("BX");  }
00191 
00192   if(*_bxincycle) {  (*_bxincycle)->GetXaxis()->SetTitle("Event BX mod(70)"); }
00193 
00194   if(*_orbit) {
00195     (*_orbit)->SetBit(TH1::kCanRebin);
00196     (*_orbit)->GetXaxis()->SetTitle("time [Orb#]"); 
00197   }
00198 
00199   if(_dbxvsbxincycle && *_dbxvsbxincycle) {
00200     (*_dbxvsbxincycle)->GetXaxis()->SetTitle("Event BX mod(70)"); (*_dbxvsbxincycle)->GetYaxis()->SetTitle("#DeltaBX"); 
00201   }
00202 
00203   if(_dbxvsbx && *_dbxvsbx) { (*_dbxvsbx)->GetXaxis()->SetTitle("BX"); (*_dbxvsbx)->GetYaxis()->SetTitle("#DeltaBX"); }
00204 
00205   if(_bxincyclevsbx && *_bxincyclevsbx) {
00206     (*_bxincyclevsbx)->GetXaxis()->SetTitle("BX"); (*_bxincyclevsbx)->GetYaxis()->SetTitle("Event BX mod(70)");
00207   }
00208 
00209   if(_orbitvsbxincycle && *_orbitvsbxincycle) {
00210     (*_orbitvsbxincycle)->SetBit(TH1::kCanRebin);
00211     (*_orbitvsbxincycle)->GetXaxis()->SetTitle("Event BX mod(70)"); (*_orbitvsbxincycle)->GetYaxis()->SetTitle("time [Orb#]"); 
00212   }
00213 }
00214 
00215 void 
00216 EventTimeDistribution::endRun(const edm::Run& iRun, const edm::EventSetup&)
00217 {
00218 }
00219 
00220 
00221 // ------------ method called once each job just before starting event loop  ------------
00222 void 
00223 EventTimeDistribution::beginJob()
00224 {
00225 
00226 }
00227 
00228 // ------------ method called once each job just after ending the event loop  ------------
00229 void 
00230 EventTimeDistribution::endJob() {
00231 
00232   edm::LogInfo("EndOfJob") << _nevents << " analyzed events";
00233 
00234 }
00235 
00236 
00237 //define this as a plug-in
00238 DEFINE_FWK_MODULE(EventTimeDistribution);