#include <HcalTBWriter.h>
Public Member Functions | |
virtual void | analyze (const edm::Event &e, const edm::EventSetup &es) |
virtual void | endJob () |
HcalTBWriter (const edm::ParameterSet &pset) | |
Private Member Functions | |
void | buildTree (const FEDRawDataCollection &raw) |
void | extractEventInfo (const FEDRawDataCollection &raw, const edm::EventID &id) |
Private Attributes | |
std::map< int, std::string > | blockToName_ |
CDFChunk * | chunkList_ [1024] |
std::map< int, int > | chunkMap_ |
CDFEventInfo * | eventInfo_ |
TFile * | file_ |
std::string | namePattern_ |
CDFRunInfo | ri_ |
TTree * | tree_ |
int | trigChunk_ |
Writes HCAL-style ROOT files from the RawData block
Definition at line 26 of file HcalTBWriter.h.
HcalTBWriter::HcalTBWriter | ( | const edm::ParameterSet & | pset | ) |
Definition at line 13 of file HcalTBWriter.cc.
References blockToName_, eventInfo_, file_, edm::ParameterSet::getUntrackedParameter(), j, AlCaRecoCosmics_cfg::name, h::names, and tree_.
: namePattern_(pset.getUntrackedParameter<std::string>("FilenamePattern","/tmp/HTB_%06d.root")) { std::vector<edm::ParameterSet> names=pset.getUntrackedParameter<std::vector<edm::ParameterSet> >("ChunkNames"); std::vector<edm::ParameterSet>::iterator j; for (j=names.begin(); j!=names.end(); j++) { std::string name=j->getUntrackedParameter<std::string>("Name"); int num=j->getUntrackedParameter<int>("Number"); blockToName_[num]=name; } file_=0; tree_=0; eventInfo_=0; }
void HcalTBWriter::analyze | ( | const edm::Event & | e, |
const edm::EventSetup & | es | ||
) | [virtual] |
Implements edm::EDAnalyzer.
Definition at line 48 of file HcalTBWriter.cc.
References CDFChunk::adoptBuffer(), buildTree(), trackerHits::c, chunkList_, chunkMap_, FEDRawData::data(), extractEventInfo(), file_, alignmentValidation::fname, edm::Event::getByType(), i, edm::EventBase::id(), namePattern_, CDFChunk::releaseBuffer(), ri_, edm::EventID::run(), CDFRunInfo::setInfo(), FEDRawData::size(), and tree_.
{ edm::Handle<FEDRawDataCollection> raw; e.getByType(raw); // assume just one! if (file_==0) { char fname[4096]; snprintf(fname,4096, namePattern_.c_str(),e.id().run()); edm::LogInfo("HCAL") << "Opening " << fname << " for writing HCAL-format file."; file_=new TFile(fname,"RECREATE"); ri_.setInfo("OriginalFile",fname); buildTree(*raw); } // adopt the buffers for writing for (std::map<int,int>::const_iterator i=chunkMap_.begin(); i!=chunkMap_.end(); i++) { CDFChunk* c=chunkList_[i->second]; const FEDRawData& frd=raw->FEDData(i->first); c->adoptBuffer((ULong64_t*)frd.data(),frd.size()/8); } // copy the event info bits extractEventInfo(*raw,e.id()); // fill the tree tree_->Fill(); // release all the buffers for (std::map<int,int>::const_iterator i=chunkMap_.begin(); i!=chunkMap_.end(); i++) { CDFChunk* c=chunkList_[i->second]; c->releaseBuffer(); } }
void HcalTBWriter::buildTree | ( | const FEDRawDataCollection & | raw | ) | [private] |
Definition at line 81 of file HcalTBWriter.cc.
References blockToName_, trackerHits::c, chunkList_, chunkMap_, eventInfo_, FEDRawDataCollection::FEDData(), i, j, AlCaRecoCosmics_cfg::name, FEDRawData::size(), tree_, and trigChunk_.
Referenced by analyze().
{ tree_=new TTree("CMSRAW","CMS Common Data Format Tree"); chunkMap_.clear(); trigChunk_=-1; int j=0; for (int i=0; i<2048; i++) { const FEDRawData& frd=raw.FEDData(i); if (frd.size()<16) continue; // it's empty... like std::string name; if (blockToName_.find(i)!=blockToName_.end()) name=blockToName_[i]; else { char sname[64]; snprintf(sname,64,"Chunk%03d",i); name=sname; } CDFChunk* c=new CDFChunk(name.c_str()); chunkList_[j]=c; tree_->Branch(name.c_str(),"CDFChunk",&(chunkList_[j])); chunkMap_[i]=j; if (name=="HCAL_Trigger" || name=="SliceTest_Trigger") trigChunk_=j; j++; } eventInfo_=new CDFEventInfo(); tree_->Branch("CDFEventInfo","CDFEventInfo",&eventInfo_,16000,2); }
void HcalTBWriter::endJob | ( | void | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 30 of file HcalTBWriter.cc.
References chunkMap_, eventInfo_, file_, ri_, CDFRunInfo::setInfo(), CDFRunInfo::store(), and tree_.
void HcalTBWriter::extractEventInfo | ( | const FEDRawDataCollection & | raw, |
const edm::EventID & | id | ||
) | [private] |
Definition at line 135 of file HcalTBWriter.cc.
References FEDHeader::bxID(), chunkList_, eventInfo_, FEDHeader::lvl1ID(), StandardTrgMsgBlkStruct::orbitNumber, newExtendedTrgMsgBlkStruct::runNumberSequenceId, CDFEventInfo::Set(), newExtendedTrgMsgBlkStruct::stdBlock, and trigChunk_.
Referenced by analyze().
{ int runno=id.run(); const char* seqid=""; int eventNo=id.event(); int l1aNo=eventNo; int orbitNo=0; int bunchNo=0; if (trigChunk_>=0) { const newExtendedTrgMsgBlk* tinfo=(const newExtendedTrgMsgBlk*)(chunkList_[trigChunk_]->getData()+2); // assume 2 64-bit words for the CDF header orbitNo=tinfo->stdBlock.orbitNumber; seqid=tinfo->runNumberSequenceId; FEDHeader head((const unsigned char*)chunkList_[trigChunk_]->getData()); bunchNo=head.bxID(); l1aNo=head.lvl1ID(); } eventInfo_->Set(runno,seqid,eventNo,l1aNo,orbitNo,bunchNo); }
std::map<int, std::string> HcalTBWriter::blockToName_ [private] |
Definition at line 34 of file HcalTBWriter.h.
Referenced by buildTree(), and HcalTBWriter().
CDFChunk* HcalTBWriter::chunkList_[1024] [private] |
Definition at line 42 of file HcalTBWriter.h.
Referenced by analyze(), buildTree(), and extractEventInfo().
std::map<int,int> HcalTBWriter::chunkMap_ [private] |
Definition at line 41 of file HcalTBWriter.h.
Referenced by analyze(), buildTree(), and endJob().
CDFEventInfo* HcalTBWriter::eventInfo_ [private] |
Definition at line 39 of file HcalTBWriter.h.
Referenced by buildTree(), endJob(), extractEventInfo(), and HcalTBWriter().
TFile* HcalTBWriter::file_ [private] |
Definition at line 37 of file HcalTBWriter.h.
Referenced by analyze(), endJob(), and HcalTBWriter().
std::string HcalTBWriter::namePattern_ [private] |
Definition at line 32 of file HcalTBWriter.h.
Referenced by analyze().
CDFRunInfo HcalTBWriter::ri_ [private] |
Definition at line 40 of file HcalTBWriter.h.
TTree* HcalTBWriter::tree_ [private] |
Definition at line 38 of file HcalTBWriter.h.
Referenced by analyze(), buildTree(), endJob(), and HcalTBWriter().
int HcalTBWriter::trigChunk_ [private] |
Definition at line 43 of file HcalTBWriter.h.
Referenced by buildTree(), and extractEventInfo().