#include <L1Trigger/TextToDigi/src/RctTextToRctDigi.h>
Public Member Functions | |
RctTextToRctDigi (const edm::ParameterSet &) | |
~RctTextToRctDigi () | |
Private Member Functions | |
void | bxSynchro (int &, int) |
Synchronize bunch crossing. | |
virtual void | produce (edm::Event &, const edm::EventSetup &) |
void | putEmptyDigi (edm::Event &) |
Create empty digi collection. | |
Private Attributes | |
std::ifstream | m_file [18] |
file handle | |
int | m_fileEventOffset |
Number of events to be offset wrt input. | |
int | m_nevt |
Event counter. | |
std::string | m_textFileName |
Name out input file. |
Description: Makes RCT digis from the file format specified by Pam Klabbers
Definition at line 40 of file RctTextToRctDigi.h.
RctTextToRctDigi::RctTextToRctDigi | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 42 of file RctTextToRctDigi.cc.
References convertXMLtoSQLite_cfg::fileName, i, recoMuon::in, LogDebug, m_file, m_textFileName, and NUM_RCT_CRATES.
: m_textFileName(iConfig.getParameter<std::string>("TextFileName")), m_fileEventOffset(iConfig.getParameter<int>("FileEventOffset")), m_nevt(0) { // Produces collections produces<L1CaloEmCollection>(); produces<L1CaloRegionCollection>(); // Open the input files for (unsigned i=0; i<NUM_RCT_CRATES; i++){ std::stringstream fileStream; fileStream << m_textFileName << std::setw(2) << std::setfill('0') << i << ".txt"; std::string fileName(fileStream.str()); m_file[i].open(fileName.c_str(),std::ios::in); if(!m_file[i].good()) { //throw cms::Exception("RctTextToRctDigiTextFileOpenError") LogDebug("RctTextToRctDigi") << "RctTextToRctDigi::RctTextToRctDigi : " << " couldn't open the file " << fileName << "...skipping!" << std::endl; } } }
RctTextToRctDigi::~RctTextToRctDigi | ( | ) |
Definition at line 68 of file RctTextToRctDigi.cc.
References i, m_file, and NUM_RCT_CRATES.
void RctTextToRctDigi::bxSynchro | ( | int & | bx, |
int | crate | ||
) | [private] |
Synchronize bunch crossing.
Syncronize bunch crossing number/n.
Definition at line 95 of file RctTextToRctDigi.cc.
References Exception, j, m_file, m_fileEventOffset, m_nevt, and tmp.
Referenced by produce().
{ std::string tmp; // bypass bx input until correct bx is reached while(bx<m_nevt+m_fileEventOffset) { for (int j=0; j<6; j++){ getline(m_file[crate],tmp); } m_file[crate] >> tmp >> bx; if(tmp!="Crossing") throw cms::Exception("RctTextToRctDigiTextFileReadError") << "RctTextToRctDigi::bxSynchro : " << " something screwy happened Crossing!=" << tmp << std::endl; } }
void RctTextToRctDigi::produce | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [private, virtual] |
Synchronize bunch crossing
Implements edm::EDProducer.
Definition at line 111 of file RctTextToRctDigi.cc.
References bxSynchro(), ExpressReco_HICollisions_FallBack::et, Exception, i, j, LogDebug, m_file, m_fileEventOffset, m_nevt, m_textFileName, NUM_RCT_CRATES, edm::Event::put(), putEmptyDigi(), and tmp.
{ // Skip event if required if (m_nevt < m_fileEventOffset){ //string tmp; //for (int i=0; i<6; i++){ // getline(m_file[i],tmp); //} putEmptyDigi(iEvent); m_nevt++; return; } // New collections std::auto_ptr<L1CaloEmCollection> em (new L1CaloEmCollection); std::auto_ptr<L1CaloRegionCollection> rgn (new L1CaloRegionCollection); // Loop over RCT crates for (unsigned i=0; i<NUM_RCT_CRATES; i++){ if(!m_file[i].good()) { continue; } // Check we're not at the end of the file if(m_file[i].eof()) { //throw cms::Exception("RctTextToRctDigiTextFileReadError") LogDebug("RctTextToRctDigi") << "RctTextToRctDigi::produce : " << " unexpected end of file " << m_textFileName << i << " adding empty collection for event !" << std::endl; putEmptyDigi(iEvent); continue; } // Check we're at the start of an event std::string tmp; m_file[i]>> tmp; if(tmp!="Crossing") { throw cms::Exception("RctTextToRctDigiTextFileReadError") << "RctTextToRctDigi::produce : " << " something screwy happened Crossing!=" << tmp << std::endl; } // Read BX number dec(m_file[i]); int BXNum; m_file[i]>>BXNum; bxSynchro(BXNum,i); if(BXNum!=m_nevt+m_fileEventOffset) throw cms::Exception("RctTextToRctDigiTextSyncError") << "RctTextToRctDigi::produce : " << " something screwy happened " << "evt:" << m_nevt << " != bx:" << BXNum << " + " << m_fileEventOffset << std::endl; // Buffers unsigned long int uLongBuffer; bool mipBitBuffer[14],qBitBuffer[14]; // All in hex from now on hex(m_file[i]); // Isolated electrons for (unsigned j=0; j<4; j++){ m_file[i] >> uLongBuffer; em->push_back(L1CaloEmCand(uLongBuffer, i, true, j,BXNum,0)); } // Non-isolated electrons for (unsigned j=0; j<4; j++){ m_file[i] >> uLongBuffer; em->push_back(L1CaloEmCand(uLongBuffer, i, false, j,BXNum,0)); } // MIP bits for (unsigned j=0; j<14; j++){ m_file[i] >> mipBitBuffer[j]; } // Quiet bits for (unsigned j=0; j<14; j++){ m_file[i] >> qBitBuffer[j]; } // Barrel and endcap regions for (unsigned j=0; j<14; j++){ m_file[i] >> uLongBuffer; unsigned et = uLongBuffer & 0x3ff; // put the first 10 bits of rawData into the Et uLongBuffer >>= 10; // shift the remaining bits down to remove the 10 bits of Et bool overFlow = ((uLongBuffer & 0x1) != 0); //LSB is now overflow bit bool tauVeto = (((uLongBuffer & 0x2) >> 1) != 0); //2nd bit is tauveto rgn->push_back(L1CaloRegion(et,overFlow,tauVeto,mipBitBuffer[j],qBitBuffer[j],i,j/2,j%2)); } // HF for (unsigned j=0; j<8; j++){ m_file[i] >> uLongBuffer; unsigned et = uLongBuffer & 0xff; // put the first 8 bits into the Et rgn->push_back(L1CaloRegion(et,true,i,j)); } dec(m_file[i]); } iEvent.put(em); iEvent.put(rgn); m_nevt++; }
void RctTextToRctDigi::putEmptyDigi | ( | edm::Event & | iEvent | ) | [private] |
Create empty digi collection.
Append empty digi collection/n.
Definition at line 77 of file RctTextToRctDigi.cc.
References i, j, NUM_RCT_CRATES, and edm::Event::put().
Referenced by produce().
{ std::auto_ptr<L1CaloEmCollection> em (new L1CaloEmCollection); std::auto_ptr<L1CaloRegionCollection> rgn (new L1CaloRegionCollection); for (unsigned i=0; i<NUM_RCT_CRATES; i++){ for (unsigned j=0; j<4; j++) { em->push_back(L1CaloEmCand(0, i, true)); em->push_back(L1CaloEmCand(0, i, false)); } for (unsigned j=0; j<14; j++) rgn->push_back(L1CaloRegion(0,false,false,false,false,i,j/2,j%2)); for (unsigned j=0; j<8; j++) rgn->push_back(L1CaloRegion(0,true,i,j)); } iEvent.put(em); iEvent.put(rgn); }
std::ifstream RctTextToRctDigi::m_file[18] [private] |
file handle
Definition at line 64 of file RctTextToRctDigi.h.
Referenced by bxSynchro(), produce(), RctTextToRctDigi(), and ~RctTextToRctDigi().
int RctTextToRctDigi::m_fileEventOffset [private] |
Number of events to be offset wrt input.
Definition at line 58 of file RctTextToRctDigi.h.
Referenced by bxSynchro(), and produce().
int RctTextToRctDigi::m_nevt [private] |
Event counter.
Definition at line 61 of file RctTextToRctDigi.h.
Referenced by bxSynchro(), and produce().
std::string RctTextToRctDigi::m_textFileName [private] |
Name out input file.
Definition at line 55 of file RctTextToRctDigi.h.
Referenced by produce(), and RctTextToRctDigi().