#include <DPGAnalysis/SiStripTools/plugins/APVCyclePhaseProducerFromL1TS.cc>
Public Member Functions | |
APVCyclePhaseProducerFromL1TS (const edm::ParameterSet &) | |
~APVCyclePhaseProducerFromL1TS () | |
Private Member Functions | |
virtual void | beginJob () |
virtual void | beginRun (edm::Run &, const edm::EventSetup &) |
virtual void | endJob () |
virtual void | endRun (edm::Run &, const edm::EventSetup &) |
virtual void | produce (edm::Event &, const edm::EventSetup &) |
Private Attributes | |
const std::vector< std::string > | _defpartnames |
const std::vector< int > | _defphases |
const unsigned int | _firstgoodrun |
TH1F * | _hdlec0lresync |
TH1F * | _hdlresynclHR |
TH1F * | _hlEC0 |
TH1F * | _hlHR |
TH1F * | _hlOC0 |
TH1F * | _hlresync |
TH1F * | _hlstart |
TH1F * | _hlTE |
TH1F * | _hsize |
const edm::InputTag | _l1tscollection |
long long | _lastEventCounter0 |
long long | _lastHardReset |
long long | _lastOrbitCounter0 |
long long | _lastResync |
long long | _lastStart |
long long | _lastTestEnable |
const int | _magicOffset |
const bool | _useEC0 |
const bool | _wantHistos |
Description: EDproducer for APVCyclePhaseCollection which uses the configuration file to assign a phase to the run
Implementation: <Notes on="" implementation>="">
Definition at line 55 of file APVCyclePhaseProducerFromL1TS.cc.
APVCyclePhaseProducerFromL1TS::APVCyclePhaseProducerFromL1TS | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 111 of file APVCyclePhaseProducerFromL1TS.cc.
: _l1tscollection(iConfig.getParameter<edm::InputTag>("l1TSCollection")), _defpartnames(iConfig.getParameter<std::vector<std::string> >("defaultPartitionNames")), _defphases(iConfig.getParameter<std::vector<int> >("defaultPhases")), _wantHistos(iConfig.getUntrackedParameter<bool>("wantHistos",false)), _useEC0(iConfig.getUntrackedParameter<bool>("useEC0",false)), _magicOffset(iConfig.getUntrackedParameter<int>("magicOffset",8)), _firstgoodrun(131768), _lastResync(-1),_lastHardReset(-1),_lastStart(-1), _lastEventCounter0(-1),_lastOrbitCounter0(-1),_lastTestEnable(-1) { produces<APVCyclePhaseCollection,edm::InEvent>(); //now do what ever other initialization is needed }
APVCyclePhaseProducerFromL1TS::~APVCyclePhaseProducerFromL1TS | ( | ) |
Definition at line 130 of file APVCyclePhaseProducerFromL1TS.cc.
{ // do anything here that needs to be done at desctruction time // (e.g. close files, deallocate resources etc.) }
void APVCyclePhaseProducerFromL1TS::beginJob | ( | void | ) | [private, virtual] |
Reimplemented from edm::EDProducer.
Definition at line 301 of file APVCyclePhaseProducerFromL1TS.cc.
{ }
void APVCyclePhaseProducerFromL1TS::beginRun | ( | edm::Run & | iRun, |
const edm::EventSetup & | iSetup | ||
) | [private, virtual] |
Reimplemented from edm::EDProducer.
Definition at line 145 of file APVCyclePhaseProducerFromL1TS.cc.
References _firstgoodrun, _hdlec0lresync, _hdlresynclHR, _hlEC0, _hlHR, _hlOC0, _hlresync, _hlstart, _hlTE, _hsize, _wantHistos, LogDebug, and edm::RunBase::run().
{ // reset offset vector if(_wantHistos) { edm::Service<TFileService> tfserv; char dirname[300]; sprintf(dirname,"run_%d",iRun.run()); TFileDirectory subrun = tfserv->mkdir(dirname); _hsize = subrun.make<TH1F>("size","Level1TriggerScalers Collection size",20,-0.5,19.5); _hlresync = subrun.make<TH1F>("lresync","Orbit of last resync",3600,0.,3600*11223); _hlresync->GetXaxis()->SetTitle("Orbit"); _hlresync->GetYaxis()->SetTitle("Events"); _hlresync->SetBit(TH1::kCanRebin); _hlOC0 = subrun.make<TH1F>("lOC0","Orbit of last OC0",3600,0.,3600*11223); _hlOC0->GetXaxis()->SetTitle("Orbit"); _hlOC0->GetYaxis()->SetTitle("Events"); _hlOC0->SetBit(TH1::kCanRebin); _hlTE = subrun.make<TH1F>("lTE","Orbit of last TestEnable",3600,0.,3600*11223); _hlTE->GetXaxis()->SetTitle("Orbit"); _hlTE->GetYaxis()->SetTitle("Events"); _hlTE->SetBit(TH1::kCanRebin); _hlstart = subrun.make<TH1F>("lstart","Orbit of last Start",3600,0.,3600*11223); _hlstart->GetXaxis()->SetTitle("Orbit"); _hlstart->GetYaxis()->SetTitle("Events"); _hlstart->SetBit(TH1::kCanRebin); _hlEC0 = subrun.make<TH1F>("lEC0","Orbit of last EC0",3600,0.,3600*11223); _hlEC0->GetXaxis()->SetTitle("Orbit"); _hlEC0->GetYaxis()->SetTitle("Events"); _hlEC0->SetBit(TH1::kCanRebin); _hlHR = subrun.make<TH1F>("lHR","Orbit of last HardReset",3600,0.,3600*11223); _hlHR->GetXaxis()->SetTitle("Orbit"); _hlHR->GetYaxis()->SetTitle("Events"); _hlHR->SetBit(TH1::kCanRebin); _hdlec0lresync = subrun.make<TH1F>("dlec0lresync","Orbit difference EC0-Resync",4000,-1999.5,2000.5); _hdlec0lresync->GetXaxis()->SetTitle("lastEC0-lastResync"); _hdlresynclHR = subrun.make<TH1F>("dlresynclHR","Orbit difference Resync-HR",4000,-1999.5,2000.5); _hdlresynclHR->GetXaxis()->SetTitle("lastEC0-lastResync"); } if(iRun.run() < _firstgoodrun) { LogDebug("UnreliableMissingL1TriggerScalers") << "In this run L1TriggerScalers is missing or unreliable for phase determination: default phases will be used"; } }
void APVCyclePhaseProducerFromL1TS::endJob | ( | void | ) | [private, virtual] |
Reimplemented from edm::EDProducer.
Definition at line 307 of file APVCyclePhaseProducerFromL1TS.cc.
{ }
void APVCyclePhaseProducerFromL1TS::endRun | ( | edm::Run & | , |
const edm::EventSetup & | |||
) | [private, virtual] |
Reimplemented from edm::EDProducer.
Definition at line 201 of file APVCyclePhaseProducerFromL1TS.cc.
{
// summary of absolute bx offset vector
}
void APVCyclePhaseProducerFromL1TS::produce | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [private, virtual] |
Implements edm::EDProducer.
Definition at line 209 of file APVCyclePhaseProducerFromL1TS.cc.
References _defpartnames, _defphases, _firstgoodrun, _hdlec0lresync, _hdlresynclHR, _hlEC0, _hlHR, _hlOC0, _hlresync, _hlstart, _hlTE, _hsize, _l1tscollection, _lastEventCounter0, _lastHardReset, _lastOrbitCounter0, _lastResync, _lastStart, _lastTestEnable, _magicOffset, _useEC0, _wantHistos, Exception, edm::Event::getByLabel(), LogDebug, edm::Event::put(), and edm::Event::run().
{ using namespace edm; std::auto_ptr<APVCyclePhaseCollection> apvphases(new APVCyclePhaseCollection() ); const std::vector<int>& phases = _defphases; const std::vector<std::string>& partnames = _defpartnames; int phasechange = 0; if(iEvent.run() >= _firstgoodrun ) { Handle<Level1TriggerScalersCollection> l1ts; iEvent.getByLabel(_l1tscollection,l1ts); if(_wantHistos) _hsize->Fill(l1ts->size()); // offset computation long long orbitoffset = 0; if(l1ts->size()>0) { if((*l1ts)[0].lastResync()!=0) { orbitoffset = _useEC0 ? (*l1ts)[0].lastEventCounter0() + _magicOffset : (*l1ts)[0].lastResync() + _magicOffset; } if(_wantHistos) { _hlresync->Fill((*l1ts)[0].lastResync()); _hlOC0->Fill((*l1ts)[0].lastOrbitCounter0()); _hlTE->Fill((*l1ts)[0].lastTestEnable()); _hlstart->Fill((*l1ts)[0].lastStart()); _hlEC0->Fill((*l1ts)[0].lastEventCounter0()); _hlHR->Fill((*l1ts)[0].lastHardReset()); } if(_lastResync != (*l1ts)[0].lastResync()) { _lastResync = (*l1ts)[0].lastResync(); if(_wantHistos) _hdlec0lresync->Fill((*l1ts)[0].lastEventCounter0()-(*l1ts)[0].lastResync()); LogDebug("TTCSignalReceived") << "New Resync at orbit " << _lastResync ; } if(_lastHardReset != (*l1ts)[0].lastHardReset()) { _lastHardReset = (*l1ts)[0].lastHardReset(); if(_wantHistos) _hdlresynclHR->Fill((*l1ts)[0].lastResync()-(*l1ts)[0].lastHardReset()); LogDebug("TTCSignalReceived") << "New HardReset at orbit " << _lastHardReset ; } if(_lastTestEnable != (*l1ts)[0].lastTestEnable()) { _lastTestEnable = (*l1ts)[0].lastTestEnable(); // LogDebug("TTCSignalReceived") << "New TestEnable at orbit " << _lastTestEnable ; } if(_lastOrbitCounter0 != (*l1ts)[0].lastOrbitCounter0()) { _lastOrbitCounter0 = (*l1ts)[0].lastOrbitCounter0(); LogDebug("TTCSignalReceived") << "New OrbitCounter0 at orbit " << _lastOrbitCounter0 ; } if(_lastEventCounter0 != (*l1ts)[0].lastEventCounter0()) { _lastEventCounter0 = (*l1ts)[0].lastEventCounter0(); LogDebug("TTCSignalReceived") << "New EventCounter0 at orbit " << _lastEventCounter0 ; } if(_lastStart != (*l1ts)[0].lastStart()) { _lastStart = (*l1ts)[0].lastStart(); LogDebug("TTCSignalReceived") << "New Start at orbit " << _lastStart ; } phasechange = ((long long)(orbitoffset*3564))%70; } } if(phases.size() < partnames.size() ) { // throw exception throw cms::Exception("InvalidAPVCyclePhases") << " Inconsistent phases/partitions vector sizes: " << phases.size() << " " << partnames.size(); } for(unsigned int ipart=0;ipart<partnames.size();++ipart) { if(phases[ipart]>=0) { apvphases->get()[partnames[ipart]] = (phases[ipart]+phasechange)%70; } } iEvent.put(apvphases); }
const std::vector<std::string> APVCyclePhaseProducerFromL1TS::_defpartnames [private] |
Definition at line 70 of file APVCyclePhaseProducerFromL1TS.cc.
Referenced by produce().
const std::vector<int> APVCyclePhaseProducerFromL1TS::_defphases [private] |
Definition at line 71 of file APVCyclePhaseProducerFromL1TS.cc.
Referenced by produce().
const unsigned int APVCyclePhaseProducerFromL1TS::_firstgoodrun [private] |
Definition at line 87 of file APVCyclePhaseProducerFromL1TS.cc.
Referenced by beginRun(), and produce().
TH1F* APVCyclePhaseProducerFromL1TS::_hdlec0lresync [private] |
Definition at line 84 of file APVCyclePhaseProducerFromL1TS.cc.
Referenced by beginRun(), and produce().
TH1F* APVCyclePhaseProducerFromL1TS::_hdlresynclHR [private] |
Definition at line 85 of file APVCyclePhaseProducerFromL1TS.cc.
Referenced by beginRun(), and produce().
TH1F* APVCyclePhaseProducerFromL1TS::_hlEC0 [private] |
Definition at line 80 of file APVCyclePhaseProducerFromL1TS.cc.
Referenced by beginRun(), and produce().
TH1F* APVCyclePhaseProducerFromL1TS::_hlHR [private] |
Definition at line 82 of file APVCyclePhaseProducerFromL1TS.cc.
Referenced by beginRun(), and produce().
TH1F* APVCyclePhaseProducerFromL1TS::_hlOC0 [private] |
Definition at line 78 of file APVCyclePhaseProducerFromL1TS.cc.
Referenced by beginRun(), and produce().
TH1F* APVCyclePhaseProducerFromL1TS::_hlresync [private] |
Definition at line 77 of file APVCyclePhaseProducerFromL1TS.cc.
Referenced by beginRun(), and produce().
TH1F* APVCyclePhaseProducerFromL1TS::_hlstart [private] |
Definition at line 81 of file APVCyclePhaseProducerFromL1TS.cc.
Referenced by beginRun(), and produce().
TH1F* APVCyclePhaseProducerFromL1TS::_hlTE [private] |
Definition at line 79 of file APVCyclePhaseProducerFromL1TS.cc.
Referenced by beginRun(), and produce().
TH1F* APVCyclePhaseProducerFromL1TS::_hsize [private] |
Definition at line 76 of file APVCyclePhaseProducerFromL1TS.cc.
Referenced by beginRun(), and produce().
const edm::InputTag APVCyclePhaseProducerFromL1TS::_l1tscollection [private] |
Definition at line 69 of file APVCyclePhaseProducerFromL1TS.cc.
Referenced by produce().
long long APVCyclePhaseProducerFromL1TS::_lastEventCounter0 [private] |
Definition at line 92 of file APVCyclePhaseProducerFromL1TS.cc.
Referenced by produce().
long long APVCyclePhaseProducerFromL1TS::_lastHardReset [private] |
Definition at line 90 of file APVCyclePhaseProducerFromL1TS.cc.
Referenced by produce().
long long APVCyclePhaseProducerFromL1TS::_lastOrbitCounter0 [private] |
Definition at line 93 of file APVCyclePhaseProducerFromL1TS.cc.
Referenced by produce().
long long APVCyclePhaseProducerFromL1TS::_lastResync [private] |
Definition at line 89 of file APVCyclePhaseProducerFromL1TS.cc.
Referenced by produce().
long long APVCyclePhaseProducerFromL1TS::_lastStart [private] |
Definition at line 91 of file APVCyclePhaseProducerFromL1TS.cc.
Referenced by produce().
long long APVCyclePhaseProducerFromL1TS::_lastTestEnable [private] |
Definition at line 94 of file APVCyclePhaseProducerFromL1TS.cc.
Referenced by produce().
const int APVCyclePhaseProducerFromL1TS::_magicOffset [private] |
Definition at line 74 of file APVCyclePhaseProducerFromL1TS.cc.
Referenced by produce().
const bool APVCyclePhaseProducerFromL1TS::_useEC0 [private] |
Definition at line 73 of file APVCyclePhaseProducerFromL1TS.cc.
Referenced by produce().
const bool APVCyclePhaseProducerFromL1TS::_wantHistos [private] |
Definition at line 72 of file APVCyclePhaseProducerFromL1TS.cc.
Referenced by beginRun(), and produce().