#include <EventFilter/ScalersRawToDigi/src/ScalersRawToDigi.cc>
Public Member Functions | |
virtual void | produce (edm::Event &, const edm::EventSetup &) |
ScalersRawToDigi (const edm::ParameterSet &) | |
~ScalersRawToDigi () | |
Private Attributes | |
edm::InputTag | inputTag_ |
Description: Unpack FED data to Trigger and Lumi Scalers "bank" These Scalers are in FED id ScalersRaw::SCALERS_FED_ID
Definition at line 38 of file ScalersRawToDigi.cc.
ScalersRawToDigi::ScalersRawToDigi | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 51 of file ScalersRawToDigi.cc.
References edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), and inputTag_.
: inputTag_((char const *)"source") { produces<L1AcceptBunchCrossingCollection>(); produces<L1TriggerScalersCollection>(); produces<Level1TriggerScalersCollection>(); produces<LumiScalersCollection>(); produces<BeamSpotOnlineCollection>(); produces<DcsStatusCollection>(); if ( iConfig.exists("scalersInputTag") ) { inputTag_ = iConfig.getParameter<edm::InputTag>("scalersInputTag"); } }
ScalersRawToDigi::~ScalersRawToDigi | ( | ) |
Definition at line 67 of file ScalersRawToDigi.cc.
{}
void ScalersRawToDigi::produce | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
Take a reference to this FED's data
Implements edm::EDProducer.
Definition at line 70 of file ScalersRawToDigi.cc.
References runTheMatrix::data, edm::Event::getByLabel(), i, getHLTprescales::index, inputTag_, ScalersRaw::N_BX_v2, edm::Event::put(), runTheMatrix_dev::raw, ScalersRaw::SCALERS_FED_ID, FEDRawData::size(), and ScalersEventRecordRaw_v5::version.
{ using namespace edm; // Get a handle to the FED data collection edm::Handle<FEDRawDataCollection> rawdata; iEvent.getByLabel(inputTag_, rawdata); std::auto_ptr<LumiScalersCollection> pLumi(new LumiScalersCollection()); std::auto_ptr<L1TriggerScalersCollection> pOldTrigger(new L1TriggerScalersCollection()); std::auto_ptr<Level1TriggerScalersCollection> pTrigger(new Level1TriggerScalersCollection()); std::auto_ptr<L1AcceptBunchCrossingCollection> pBunch(new L1AcceptBunchCrossingCollection()); std::auto_ptr<BeamSpotOnlineCollection> pBeamSpotOnline(new BeamSpotOnlineCollection()); std::auto_ptr<DcsStatusCollection> pDcsStatus(new DcsStatusCollection()); const FEDRawData & fedData = rawdata->FEDData(ScalersRaw::SCALERS_FED_ID); unsigned short int length = fedData.size(); if ( length > 0 ) { int nWords = length / 8; int nBytesExtra = 0; const ScalersEventRecordRaw_v5 * raw = (struct ScalersEventRecordRaw_v5 *)fedData.data(); if ( ( raw->version == 1 ) || ( raw->version == 2 ) ) { L1TriggerScalers oldTriggerScalers(fedData.data()); pOldTrigger->push_back(oldTriggerScalers); nBytesExtra = length - sizeof(struct ScalersEventRecordRaw_v1); } else if ( raw->version >= 3 ) { Level1TriggerScalers triggerScalers(fedData.data()); pTrigger->push_back(triggerScalers); nBytesExtra = ScalersRaw::N_BX_v2 * sizeof(unsigned long long); } LumiScalers lumiScalers(fedData.data()); pLumi->push_back(lumiScalers); if (( nBytesExtra >= 8 ) && (( nBytesExtra % 8 ) == 0 )) { unsigned long long * data = (unsigned long long *)fedData.data(); int nWordsExtra = nBytesExtra / 8; for ( int i=0; i<nWordsExtra; i++) { int index = nWords - 5 + i; L1AcceptBunchCrossing bc(i,data[index]); pBunch->push_back(bc); } } if ( raw->version >= 4 ) { BeamSpotOnline beamSpotOnline(fedData.data()); pBeamSpotOnline->push_back(beamSpotOnline); DcsStatus dcsStatus(fedData.data()); pDcsStatus->push_back(dcsStatus); } } iEvent.put(pOldTrigger); iEvent.put(pTrigger); iEvent.put(pLumi); iEvent.put(pBunch); iEvent.put(pBeamSpotOnline); iEvent.put(pDcsStatus); }
edm::InputTag ScalersRawToDigi::inputTag_ [private] |
Definition at line 47 of file ScalersRawToDigi.cc.
Referenced by produce(), and ScalersRawToDigi().