#include <DPGAnalysis/SiStripTools/plugins/EventTimeDistribution.cc>
Public Member Functions | |
EventTimeDistribution (const edm::ParameterSet &) | |
~EventTimeDistribution () | |
Private Member Functions | |
virtual void | analyze (const edm::Event &, const edm::EventSetup &) |
virtual void | beginJob () |
virtual void | beginRun (const edm::Run &, const edm::EventSetup &) |
virtual void | endJob () |
virtual void | endRun (const edm::Run &, const edm::EventSetup &) |
Private Attributes | |
const edm::InputTag | _apvphasecoll |
const double | _binsize |
TH1F ** | _bx |
TH1F ** | _bxincycle |
TH2F ** | _bxincyclevsbx |
TH1F ** | _dbx |
TH2F ** | _dbxvsbx |
TH2F ** | _dbxvsbxincycle |
const edm::InputTag | _historyProduct |
unsigned int | _nevents |
TH1F ** | _orbit |
TH2F ** | _orbitvsbxincycle |
const std::string | _phasepart |
RunHistogramManager | _rhm |
const bool | _wantbxincyclevsbx |
const bool | _wantdbxvsbx |
const bool | _wantdbxvsbxincycle |
const bool | _wantorbitvsbxincycle |
Description: <one line="" class="" summary>="">
Implementation: <Notes on="" implementation>="">
Definition at line 53 of file EventTimeDistribution.cc.
EventTimeDistribution::EventTimeDistribution | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 102 of file EventTimeDistribution.cc.
References _apvphasecoll, _binsize, _bx, _bxincycle, _bxincyclevsbx, _dbx, _dbxvsbx, _dbxvsbxincycle, _orbit, _orbitvsbxincycle, _rhm, _wantbxincyclevsbx, _wantdbxvsbx, _wantdbxvsbxincycle, _wantorbitvsbxincycle, RunHistogramManager::makeTH1F(), and RunHistogramManager::makeTH2F().
: _historyProduct(iConfig.getParameter<edm::InputTag>("historyProduct")), _apvphasecoll(iConfig.getParameter<edm::InputTag>("apvPhaseCollection")), _phasepart(iConfig.getUntrackedParameter<std::string>("phasePartition","None")), _wantdbxvsbxincycle(iConfig.getUntrackedParameter<bool>("wantDBXvsBXincycle",false)), _wantdbxvsbx(iConfig.getUntrackedParameter<bool>("wantDBXvsBX",false)), _wantbxincyclevsbx(iConfig.getUntrackedParameter<bool>("wantBXincyclevsBX",false)), _wantorbitvsbxincycle(iConfig.getUntrackedParameter<bool>("wantOrbitvsBXincycle",false)), _nevents(0), _binsize(iConfig.getUntrackedParameter<double>("minBinSizeInSec",1.)), _rhm(), _dbxvsbxincycle(0), _dbxvsbx(0), _bxincyclevsbx(0), _orbitvsbxincycle(0) { //now do what ever initialization is needed _dbx = _rhm.makeTH1F("dbx","dbx",1000,-0.5,999.5); _bx = _rhm.makeTH1F("bx","BX number",3564,-0.5,3563.5); _bxincycle = _rhm.makeTH1F("bxcycle","bxcycle",70,-0.5,69.5); _orbit = _rhm.makeTH1F("orbit","orbit",3600,0,11223*_binsize*3600); if(_wantdbxvsbxincycle) _dbxvsbxincycle = _rhm.makeTH2F("dbxvsbxincycle","dbxvsbxincycle",70,-0.5,69.5,1000,-0.5,999.5); if(_wantdbxvsbx) _dbxvsbx = _rhm.makeTH2F("dbxvsbx","dbxvsbx",3564,-0.5,3563.5,1000,-0.5,999.5); if(_wantbxincyclevsbx) _bxincyclevsbx = _rhm.makeTH2F("bxincyclevsbx","bxincyclevsbx",3564,-0.5,3563.5,70,-0.5,69.5); if(_wantorbitvsbxincycle) _orbitvsbxincycle = _rhm.makeTH2F("orbitvsbxincycle","orbitvsbxincycle",70,-0.5,69.5,3600,0,11223*_binsize*3600); edm::LogInfo("UsedAPVCyclePhaseCollection") << " APVCyclePhaseCollection " << _apvphasecoll << " used"; }
EventTimeDistribution::~EventTimeDistribution | ( | ) |
Definition at line 131 of file EventTimeDistribution.cc.
{ // do anything here that needs to be done at desctruction time // (e.g. close files, deallocate resources etc.) }
void EventTimeDistribution::analyze | ( | const edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [private, virtual] |
Implements edm::EDAnalyzer.
Definition at line 146 of file EventTimeDistribution.cc.
References _apvphasecoll, _bxincyclevsbx, _dbxvsbx, _dbxvsbxincycle, _historyProduct, _nevents, _orbitvsbxincycle, _phasepart, edm::EventBase::bunchCrossing(), edm::Event::getByLabel(), APVCyclePhaseCollection::invalid, APVCyclePhaseCollection::multiphase, APVCyclePhaseCollection::nopartition, and edm::EventBase::orbitNumber().
{ using namespace edm; _nevents++; edm::Handle<EventWithHistory> he; iEvent.getByLabel(_historyProduct,he); edm::Handle<APVCyclePhaseCollection> apvphase; iEvent.getByLabel(_apvphasecoll,apvphase); long long tbx = he->absoluteBX(); if(apvphase.isValid() && !apvphase.failedToGet()) { const int thephase = apvphase->getPhase(_phasepart); if(thephase!=APVCyclePhaseCollection::invalid && thephase!=APVCyclePhaseCollection::multiphase && thephase!=APVCyclePhaseCollection::nopartition) tbx -= thephase; } // improve the matchin between default and actual partitions (*_dbx)->Fill(he->deltaBX()); (*_bx)->Fill(iEvent.bunchCrossing()); (*_bxincycle)->Fill(tbx%70); (*_orbit)->Fill(iEvent.orbitNumber()); if(_dbxvsbxincycle && *_dbxvsbxincycle) (*_dbxvsbxincycle)->Fill(tbx%70,he->deltaBX()); if(_dbxvsbx && *_dbxvsbx) (*_dbxvsbx)->Fill(iEvent.bunchCrossing(),he->deltaBX()); if(_bxincyclevsbx && *_bxincyclevsbx) (*_bxincyclevsbx)->Fill(iEvent.bunchCrossing(),tbx%70); if(_orbitvsbxincycle && *_orbitvsbxincycle) (*_orbitvsbxincycle)->Fill(tbx%70,iEvent.orbitNumber()); }
void EventTimeDistribution::beginJob | ( | void | ) | [private, virtual] |
void EventTimeDistribution::beginRun | ( | const edm::Run & | iRun, |
const edm::EventSetup & | |||
) | [private, virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 184 of file EventTimeDistribution.cc.
References _bx, _bxincycle, _bxincyclevsbx, _dbx, _dbxvsbx, _dbxvsbxincycle, _orbit, _orbitvsbxincycle, _rhm, and RunHistogramManager::beginRun().
{ _rhm.beginRun(iRun); if(*_dbx) { (*_dbx)->GetXaxis()->SetTitle("#DeltaBX"); } if(*_bx) { (*_bx)->GetXaxis()->SetTitle("BX"); } if(*_bxincycle) { (*_bxincycle)->GetXaxis()->SetTitle("Event BX mod(70)"); } if(*_orbit) { (*_orbit)->SetBit(TH1::kCanRebin); (*_orbit)->GetXaxis()->SetTitle("time [Orb#]"); } if(_dbxvsbxincycle && *_dbxvsbxincycle) { (*_dbxvsbxincycle)->GetXaxis()->SetTitle("Event BX mod(70)"); (*_dbxvsbxincycle)->GetYaxis()->SetTitle("#DeltaBX"); } if(_dbxvsbx && *_dbxvsbx) { (*_dbxvsbx)->GetXaxis()->SetTitle("BX"); (*_dbxvsbx)->GetYaxis()->SetTitle("#DeltaBX"); } if(_bxincyclevsbx && *_bxincyclevsbx) { (*_bxincyclevsbx)->GetXaxis()->SetTitle("BX"); (*_bxincyclevsbx)->GetYaxis()->SetTitle("Event BX mod(70)"); } if(_orbitvsbxincycle && *_orbitvsbxincycle) { (*_orbitvsbxincycle)->SetBit(TH1::kCanRebin); (*_orbitvsbxincycle)->GetXaxis()->SetTitle("Event BX mod(70)"); (*_orbitvsbxincycle)->GetYaxis()->SetTitle("time [Orb#]"); } }
void EventTimeDistribution::endJob | ( | void | ) | [private, virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 230 of file EventTimeDistribution.cc.
References _nevents.
{ edm::LogInfo("EndOfJob") << _nevents << " analyzed events"; }
void EventTimeDistribution::endRun | ( | const edm::Run & | iRun, |
const edm::EventSetup & | |||
) | [private, virtual] |
const edm::InputTag EventTimeDistribution::_apvphasecoll [private] |
Definition at line 69 of file EventTimeDistribution.cc.
Referenced by analyze(), and EventTimeDistribution().
const double EventTimeDistribution::_binsize [private] |
Definition at line 76 of file EventTimeDistribution.cc.
Referenced by EventTimeDistribution().
TH1F** EventTimeDistribution::_bx [private] |
Definition at line 81 of file EventTimeDistribution.cc.
Referenced by beginRun(), and EventTimeDistribution().
TH1F** EventTimeDistribution::_bxincycle [private] |
Definition at line 82 of file EventTimeDistribution.cc.
Referenced by beginRun(), and EventTimeDistribution().
TH2F** EventTimeDistribution::_bxincyclevsbx [private] |
Definition at line 86 of file EventTimeDistribution.cc.
Referenced by analyze(), beginRun(), and EventTimeDistribution().
TH1F** EventTimeDistribution::_dbx [private] |
Definition at line 80 of file EventTimeDistribution.cc.
Referenced by beginRun(), and EventTimeDistribution().
TH2F** EventTimeDistribution::_dbxvsbx [private] |
Definition at line 85 of file EventTimeDistribution.cc.
Referenced by analyze(), beginRun(), and EventTimeDistribution().
TH2F** EventTimeDistribution::_dbxvsbxincycle [private] |
Definition at line 84 of file EventTimeDistribution.cc.
Referenced by analyze(), beginRun(), and EventTimeDistribution().
const edm::InputTag EventTimeDistribution::_historyProduct [private] |
Definition at line 68 of file EventTimeDistribution.cc.
Referenced by analyze().
unsigned int EventTimeDistribution::_nevents [private] |
Definition at line 75 of file EventTimeDistribution.cc.
TH1F** EventTimeDistribution::_orbit [private] |
Definition at line 83 of file EventTimeDistribution.cc.
Referenced by beginRun(), and EventTimeDistribution().
TH2F** EventTimeDistribution::_orbitvsbxincycle [private] |
Definition at line 87 of file EventTimeDistribution.cc.
Referenced by analyze(), beginRun(), and EventTimeDistribution().
const std::string EventTimeDistribution::_phasepart [private] |
Definition at line 70 of file EventTimeDistribution.cc.
Referenced by analyze().
Definition at line 78 of file EventTimeDistribution.cc.
Referenced by beginRun(), and EventTimeDistribution().
const bool EventTimeDistribution::_wantbxincyclevsbx [private] |
Definition at line 73 of file EventTimeDistribution.cc.
Referenced by EventTimeDistribution().
const bool EventTimeDistribution::_wantdbxvsbx [private] |
Definition at line 72 of file EventTimeDistribution.cc.
Referenced by EventTimeDistribution().
const bool EventTimeDistribution::_wantdbxvsbxincycle [private] |
Definition at line 71 of file EventTimeDistribution.cc.
Referenced by EventTimeDistribution().
const bool EventTimeDistribution::_wantorbitvsbxincycle [private] |
Definition at line 74 of file EventTimeDistribution.cc.
Referenced by EventTimeDistribution().