Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

reading event id from EVM FED (71X backport) #4886

Merged
merged 2 commits into from Aug 7, 2014
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
56 changes: 50 additions & 6 deletions EventFilter/Utilities/plugins/DaqFakeReader.cc
Expand Up @@ -8,6 +8,10 @@
#include "DataFormats/FEDRawData/interface/FEDTrailer.h"
#include "DataFormats/FEDRawData/interface/FEDNumbering.h"

#include "EventFilter/FEDInterface/interface/GlobalEventNumber.h"
#include "EventFilter/FEDInterface/interface/fed_header.h"
#include "EventFilter/FEDInterface/interface/fed_trailer.h"

#include "DataFormats/FEDRawData/interface/FEDRawData.h"

#include "FWCore/ParameterSet/interface/ParameterSet.h"
Expand Down Expand Up @@ -89,7 +93,11 @@ int DaqFakeReader::fillRawData(Event& e,
fillFEDs(FEDNumbering::MINHCALFEDID,
FEDNumbering::MAXHCALFEDID,
eID, *data, meansize, width);
fillFED1023(eID,*data);

timeval now;
gettimeofday(&now,0);
fillGTPFED(eID, *data,&now);
fillFED1023(eID,*data,&now);
}
return 1;
}
Expand Down Expand Up @@ -144,21 +152,57 @@ void DaqFakeReader::fillFEDs(const int fedmin, const int fedmax,
}
}

void DaqFakeReader::fillGTPFED(EventID& eID,
FEDRawDataCollection& data, timeval *now)
{
uint32_t fedId = FEDNumbering::MINTriggerGTPFEDID;
FEDRawData& feddata = data.FEDData(fedId);
uint32_t size = evf::evtn::SLINK_WORD_SIZE*37-16;//BST52_3BX
feddata.resize(size+16);

FEDHeader::set(feddata.data(),
1, // Trigger type
eID.event(), // LV1_id (24 bits)
0, // BX_id
fedId); // source_id

int crc = 0; // FIXME : get CRC
FEDTrailer::set(feddata.data()+8+size,
size/8+2, // in 64 bit words!!!
crc,
0, // Evt_stat
0); // TTS bits

unsigned char * pOffset = feddata.data() + sizeof(fedh_t);
//fill in event ID
*( (uint32_t*)(pOffset + evf::evtn::EVM_BOARDID_OFFSET * evf::evtn::SLINK_WORD_SIZE / 2)) = evf::evtn::EVM_BOARDID_VALUE << evf::evtn::EVM_BOARDID_SHIFT;
*( (uint32_t*)(pOffset + sizeof(fedh_t) + (9*2 + evf::evtn::EVM_TCS_TRIGNR_OFFSET) * evf::evtn::SLINK_WORD_SIZE / 2)) = eID.event();
//fill in timestamp
*( (uint32_t*) (pOffset + evf::evtn::EVM_GTFE_BSTGPS_OFFSET * evf::evtn::SLINK_WORD_SIZE / 2)) = now->tv_sec;
*( (uint32_t*) (pOffset + sizeof(fedh_t) + evf::evtn::EVM_GTFE_BSTGPS_OFFSET * evf::evtn::SLINK_WORD_SIZE / 2 + evf::evtn::SLINK_HALFWORD_SIZE)) = now->tv_usec;

//*( (uint16_t*) (pOffset + (evtn::EVM_GTFE_BLOCK*2 + evtn::EVM_TCS_LSBLNR_OFFSET)*evtn::SLINK_HALFWORD_SIZE)) = (unsigned short)fakeLs_-1;

//we could also generate lumiblock, bcr, orbit,... but they are not currently used by the FRD input source

}


void DaqFakeReader::fillFED1023(EventID& eID,
FEDRawDataCollection& data)
FEDRawDataCollection& data,timeval *now)
{
FEDRawData& feddata = data.FEDData(frb.fedId());
// Allocate space for header+trailer+payload
feddata.resize(frb.size());
frb.putHeader(eID.event(),0);
if(eID.event()%modulo_error_events==0) frb.setDAQDiaWord1(1ll); else frb.setDAQDiaWord1(0ll);
timeval now;
gettimeofday(&now, 0);
frb.setRBTimeStamp(((uint64_t) (now.tv_sec) << 32) + (uint64_t) (now.tv_usec));
frb.setRBTimeStamp(((uint64_t) (now->tv_sec) << 32) + (uint64_t) (now->tv_usec));
frb.putTrailer();
memcpy(feddata.data(),frb.getPayload(),frb.size());
}

void DaqFakeReader::beginLuminosityBlock(LuminosityBlock const& iL, EventSetup const& iE)
{
std::cout << "DaqFakeReader begin Lumi " << iL.luminosityBlock() << std::endl;}
std::cout << "DaqFakeReader begin Lumi " << iL.luminosityBlock() << std::endl;
fakeLs_=iL.luminosityBlock();
}
5 changes: 4 additions & 1 deletion EventFilter/Utilities/plugins/DaqFakeReader.h
Expand Up @@ -48,8 +48,10 @@ class DaqFakeReader : public edm::EDProducer
FEDRawDataCollection& data,
float meansize,
float width);
void fillGTPFED(edm::EventID& eID,
FEDRawDataCollection& data,timeval * now);
void fillFED1023(edm::EventID& eID,
FEDRawDataCollection& data);
FEDRawDataCollection& data,timeval * now);
virtual void beginLuminosityBlock(edm::LuminosityBlock const& iL, edm::EventSetup const& iE);
private:
//
Expand All @@ -62,6 +64,7 @@ class DaqFakeReader : public edm::EDProducer
unsigned int width;
unsigned int injected_errors_per_million_events;
unsigned int modulo_error_events;
unsigned int fakeLs_=0;
evf::EvffedFillerRB frb;
};

Expand Down
23 changes: 20 additions & 3 deletions EventFilter/Utilities/plugins/FRDStreamSource.cc
Expand Up @@ -15,7 +15,8 @@
FRDStreamSource::FRDStreamSource(edm::ParameterSet const& pset,
edm::InputSourceDescription const& desc)
: ProducerSourceFromFiles(pset,desc,true),
verifyAdler32_(pset.getUntrackedParameter<bool> ("verifyAdler32", true))
verifyAdler32_(pset.getUntrackedParameter<bool> ("verifyAdler32", true)),
useL1EventID_(pset.getUntrackedParameter<bool> ("useL1EventID", false))
{
itFileName_=fileNames().begin();
openFile(*itFileName_);
Expand Down Expand Up @@ -54,7 +55,8 @@ bool FRDStreamSource::setRunAndEventInfo(edm::EventID& id, edm::TimeValue_t& the
}

std::unique_ptr<FRDEventMsgView> frdEventMsg(new FRDEventMsgView(&buffer_[0]));
id = edm::EventID(frdEventMsg->run(), frdEventMsg->lumi(), frdEventMsg->event());
if (useL1EventID_)
id = edm::EventID(frdEventMsg->run(), frdEventMsg->lumi(), frdEventMsg->event());

const uint32_t totalSize = frdEventMsg->size();
if ( totalSize > buffer_.size() ) {
Expand Down Expand Up @@ -85,6 +87,7 @@ bool FRDStreamSource::setRunAndEventInfo(edm::EventID& id, edm::TimeValue_t& the

uint32_t eventSize = frdEventMsg->eventSize();
char* event = (char*)frdEventMsg->payload();
bool foundGTPFED=false;

while (eventSize > 0) {
eventSize -= sizeof(fedt_t);
Expand All @@ -94,11 +97,25 @@ bool FRDStreamSource::setRunAndEventInfo(edm::EventID& id, edm::TimeValue_t& the
const fedh_t* fedHeader = (fedh_t *) (event + eventSize);
const uint16_t fedId = FED_SOID_EXTRACT(fedHeader->sourceid);
if (fedId == FEDNumbering::MINTriggerGTPFEDID) {
evf::evtn::evm_board_setformat(fedSize);
foundGTPFED=true;
const bool GTPEvmBoardSense=evf::evtn::evm_board_sense((unsigned char*) fedHeader,fedSize);
if (!useL1EventID_) {
if (GTPEvmBoardSense)
id = edm::EventID(frdEventMsg->run(), frdEventMsg->lumi(), evf::evtn::get((unsigned char*) fedHeader,true));
else
id = edm::EventID(frdEventMsg->run(), frdEventMsg->lumi(), evf::evtn::get((unsigned char*) fedHeader,false));
}
//evf::evtn::evm_board_setformat(fedSize);
const uint64_t gpsl = evf::evtn::getgpslow((unsigned char*) fedHeader);
const uint64_t gpsh = evf::evtn::getgpshigh((unsigned char*) fedHeader);
theTime = static_cast<edm::TimeValue_t>((gpsh << 32) + gpsl);
}
//take event ID from GTPE FED
if (fedId == FEDNumbering::MINTriggerEGTPFEDID && !foundGTPFED && !useL1EventID_) {
if (evf::evtn::gtpe_board_sense((unsigned char*)fedHeader)) {
id = edm::EventID(frdEventMsg->run(), frdEventMsg->lumi(), evf::evtn::gtpe_get((unsigned char*) fedHeader));
}
}
FEDRawData& fedData = rawData_->FEDData(fedId);
fedData.resize(fedSize);
memcpy(fedData.data(), event + eventSize, fedSize);
Expand Down
1 change: 1 addition & 0 deletions EventFilter/Utilities/plugins/FRDStreamSource.h
Expand Up @@ -46,6 +46,7 @@ class FRDStreamSource : public edm::ProducerSourceFromFiles {
std::auto_ptr<FEDRawDataCollection> rawData_;
std::vector<char> buffer_;
const bool verifyAdler32_;
const bool useL1EventID_;
unsigned int detectedFRDversion_=0;
};

Expand Down
33 changes: 27 additions & 6 deletions EventFilter/Utilities/plugins/FedRawDataInputSource.cc
Expand Up @@ -55,6 +55,7 @@ FedRawDataInputSource::FedRawDataInputSource(edm::ParameterSet const& pset,
numBuffers_(pset.getUntrackedParameter<unsigned int> ("numBuffers",1)),
getLSFromFilename_(pset.getUntrackedParameter<bool> ("getLSFromFilename", true)),
verifyAdler32_(pset.getUntrackedParameter<bool> ("verifyAdler32", true)),
useL1EventID_(pset.getUntrackedParameter<bool> ("useL1EventID", false)),
testModeNoBuilderUnit_(edm::Service<evf::EvFDaqDirector>()->getTestModeNoBuilderUnit()),
runNumber_(edm::Service<evf::EvFDaqDirector>()->getRunNumber()),
fuOutputDir_(edm::Service<evf::EvFDaqDirector>()->baseRunDir()),
Expand Down Expand Up @@ -228,7 +229,8 @@ bool FedRawDataInputSource::checkNextEvent()
maybeOpenNewLumiSection( event_->lumi() );
}
}
eventID_ = edm::EventID(event_->run(), currentLumiSection_, event_->event());
eventRunNumber_=event_->run();
L1EventID_ = event_->event();

setEventCached();

Expand Down Expand Up @@ -406,8 +408,9 @@ inline evf::EvFDaqDirector::FileStatus FedRawDataInputSource::getNextEvent()

unsigned char *dataPosition = currentFile_->chunks_[0]->buf_+ currentFile_->chunkPosition_;


if (!bufferInputRead_ || bufferInputRead_ < headerSize[detectedFRDversion_])
//conditions when read amount is not sufficient for the header to fit
if (!bufferInputRead_ || bufferInputRead_ < headerSize[detectedFRDversion_]
|| eventChunkSize_ - currentFile_->chunkPosition_ < headerSize[detectedFRDversion_])
{
readNextChunkIntoBuffer(currentFile_);

Expand Down Expand Up @@ -558,11 +561,19 @@ void FedRawDataInputSource::deleteFile(std::string const& fileName)
}
}


void FedRawDataInputSource::read(edm::EventPrincipal& eventPrincipal)
{
std::auto_ptr<FEDRawDataCollection> rawData(new FEDRawDataCollection);
edm::Timestamp tstamp = fillFEDRawDataCollection(rawData);

if (useL1EventID_)
eventID_ = edm::EventID(eventRunNumber_, currentLumiSection_, L1EventID_);
else {
assert(GTPEventID_);
eventID_ = edm::EventID(eventRunNumber_, currentLumiSection_, GTPEventID_);
}

edm::EventAuxiliary aux(eventID_, processGUID(), tstamp, true,
edm::EventAuxiliary::PhysicsTrigger);
aux.setProcessHistoryID(processHistoryID_);
Expand Down Expand Up @@ -607,11 +618,12 @@ void FedRawDataInputSource::read(edm::EventPrincipal& eventPrincipal)
return;
}

edm::Timestamp FedRawDataInputSource::fillFEDRawDataCollection(std::auto_ptr<FEDRawDataCollection>& rawData) const
edm::Timestamp FedRawDataInputSource::fillFEDRawDataCollection(std::auto_ptr<FEDRawDataCollection>& rawData)
{
edm::Timestamp tstamp;
uint32_t eventSize = event_->eventSize();
char* event = (char*)event_->payload();
GTPEventID_=0;

while (eventSize > 0) {
eventSize -= sizeof(fedt_t);
Expand All @@ -621,11 +633,21 @@ edm::Timestamp FedRawDataInputSource::fillFEDRawDataCollection(std::auto_ptr<FED
const fedh_t* fedHeader = (fedh_t *) (event + eventSize);
const uint16_t fedId = FED_SOID_EXTRACT(fedHeader->sourceid);
if (fedId == FEDNumbering::MINTriggerGTPFEDID) {
evf::evtn::evm_board_setformat(fedSize);
if (evf::evtn::evm_board_sense((unsigned char*) fedHeader,fedSize))
GTPEventID_ = evf::evtn::get((unsigned char*) fedHeader,true);
else
GTPEventID_ = evf::evtn::get((unsigned char*) fedHeader,false);
//evf::evtn::evm_board_setformat(fedSize);
const uint64_t gpsl = evf::evtn::getgpslow((unsigned char*) fedHeader);
const uint64_t gpsh = evf::evtn::getgpshigh((unsigned char*) fedHeader);
tstamp = edm::Timestamp(static_cast<edm::TimeValue_t> ((gpsh << 32) + gpsl));
}
//take event ID from GTPE FED
if (fedId == FEDNumbering::MINTriggerEGTPFEDID && GTPEventID_==0) {
if (evf::evtn::gtpe_board_sense((unsigned char*)fedHeader)) {
GTPEventID_ = evf::evtn::gtpe_get((unsigned char*) fedHeader);
}
}
FEDRawData& fedData = rawData->FEDData(fedId);
fedData.resize(fedSize);
memcpy(fedData.data(), event + eventSize, fedSize);
Expand Down Expand Up @@ -839,7 +861,6 @@ void FedRawDataInputSource::readSupervisor()
int dbgcount=0;
if (status == evf::EvFDaqDirector::noFile) {
dbgcount++;
//if (!(dbgcount%20)) edm::LogInfo("FedRawDataInputSource") << "No file for me... sleep and try again...";
if (!(dbgcount%20)) LogDebug("FedRawDataInputSource") << "No file for me... sleep and try again...";
usleep(100000);
}
Expand Down
6 changes: 5 additions & 1 deletion EventFilter/Utilities/plugins/FedRawDataInputSource.h
Expand Up @@ -57,7 +57,7 @@ friend class InputChunk;
void maybeOpenNewLumiSection(const uint32_t lumiSection);
evf::EvFDaqDirector::FileStatus nextEvent();
evf::EvFDaqDirector::FileStatus getNextEvent();
edm::Timestamp fillFEDRawDataCollection(std::auto_ptr<FEDRawDataCollection>&) const;
edm::Timestamp fillFEDRawDataCollection(std::auto_ptr<FEDRawDataCollection>&);
void deleteFile(std::string const&);
int grabNextJsonFile(boost::filesystem::path const&);
void renameToNextFree(std::string const& fileName) const;
Expand Down Expand Up @@ -85,6 +85,7 @@ friend class InputChunk;
// get LS from filename instead of event header
const bool getLSFromFilename_;
const bool verifyAdler32_;
const bool useL1EventID_;
const bool testModeNoBuilderUnit_;

const edm::RunNumber_t runNumber_;
Expand All @@ -99,6 +100,9 @@ friend class InputChunk;
edm::ProcessHistoryID processHistoryID_;

unsigned int currentLumiSection_;
uint32_t eventRunNumber_=0;
uint32_t GTPEventID_ = 0;
uint32_t L1EventID_ = 0;
unsigned int eventsThisLumi_;
unsigned long eventsThisRun_ = 0;

Expand Down
9 changes: 5 additions & 4 deletions EventFilter/Utilities/python/FedRawDataInputSource_cfi.py
Expand Up @@ -2,9 +2,10 @@

source = cms.Source("FedRawDataInputSource",
getLSFromFilename = cms.untracked.bool(True),
eventChunkSize = cms.untracked.uint32(128),
eventChunkBlock = cms.untracked.uint32(128),
numBuffers = cms.untracked.uint32(1),
verifyAdler32 = cms.untracked.bool(True)
eventChunkSize = cms.untracked.uint32(32),
eventChunkBlock = cms.untracked.uint32(32),
numBuffers = cms.untracked.uint32(2),
verifyAdler32 = cms.untracked.bool(True),
useL1EventID = cms.untracked.bool(False)
)

1 change: 1 addition & 0 deletions EventFilter/Utilities/test/startFU.py
Expand Up @@ -88,6 +88,7 @@
getLSFromFilename = cms.untracked.bool(True),
testModeNoBuilderUnit = cms.untracked.bool(False),
verifyAdler32 = cms.untracked.bool(True),
useL1EventID = cms.untracked.bool(False),
eventChunkSize = cms.untracked.uint32(16),
numBuffers = cms.untracked.uint32(2),
eventChunkBlock = cms.untracked.uint32(1)
Expand Down