CMS 3D CMS Logo

Public Member Functions | Public Attributes | Private Attributes

EventHeader Class Reference

#include <EventHeader.h>

List of all members.

Public Member Functions

void analyze (edm::Event const &iEvent, TTree *tree)
 EventHeader ()
void setup (TTree *tree)
 ~EventHeader ()

Public Attributes

char_uint32 droppedEventsCount_
char_uint32 event_
Header header_
char_uint32 lumi_
char_uint32 origDataSize_
char_uint32 outModId_
uint8 protocolVersion_
char_uint32 run_

Private Attributes

bool _Debug
double fAvgInstDelLumi
int fBx
int fEvent
int fLumiBlock
int fOrbit
int fRun

Detailed Description

Event Message Represented here

Protocol Versions 1-4: code 1 | size 4 | run 4 | event 4 | lumi 4 | reserved 4 | l1_count 4| l1bits l1_count/8 | hlt_count 4| hltbits hlt_count/4 | eventdatalength 4 | eventdata blob {variable}

Protocol Version 5: code 1 | size 4 | protocol version 1 | run 4 | event 4 | lumi 4 | origDataSize 4 | outModId 4 | l1_count 4| l1bits l1_count/8 | hlt_count 4| hltbits hlt_count/4 | eventdatalength 4 | eventdata blob {variable}

Protocol Version 6: // no change here, only INIT msg change code 1 | size 4 | protocol version 1 | run 4 | event 4 | lumi 4 | origDataSize 4 | outModId 4 | l1_count 4| l1bits l1_count/8 | hlt_count 4| hltbits hlt_count/4 | eventdatalength 4 | eventdata blob {variable}

Protocol Version 7: // no change here, only data blob changes code 1 | size 4 | protocol version 1 | run 4 | event 4 | lumi 4 | origDataSize 4 | outModId 4 | l1_count 4| l1bits l1_count/8 | hlt_count 4| hltbits hlt_count/4 | eventdatalength 4 | eventdata blob {variable}

Protocol Version 8: // add in checksum of data blob changes code 1 | size 4 | protocol version 1 | run 4 | event 4 | lumi 4 | origDataSize 4 | outModId 4 | l1_count 4| l1bits l1_count/8 | hlt_count 4| hltbits hlt_count/4 | adler32_chksum 4 | host name length 1 | host name {Fixed size} eventdatalength 4 | eventdata blob {variable}

Protocol Version 9: // add dropped events counter code 1 | size 4 | protocol version 1 | run 4 | event 4 | lumi 4 | origDataSize 4 | outModId 4 | droppedEventsCount 4 | l1_count 4 | l1bits l1_count/8 | hlt_count 4 | hltbits hlt_count/4 | adler32_chksum 4 | host name length 1 | host name {Fixed size} eventdatalength 4 | eventdata blob {variable}

$Date: November 2006 $Revision:

Author:
V. Rekovic - UMinn

Definition at line 17 of file EventHeader.h.


Constructor & Destructor Documentation

EventHeader::EventHeader ( )

Definition at line 11 of file EventHeader.cc.

                         :
  fRun( -1 ),
  fEvent( -1 ),
  fLumiBlock( -1 ),
  fBx( -1 ),
  fOrbit( -1 ),
  fAvgInstDelLumi( -999. ),
  _Debug( false )
{ }
EventHeader::~EventHeader ( )

Definition at line 21 of file EventHeader.cc.

                          {

}

Member Function Documentation

void EventHeader::analyze ( edm::Event const &  iEvent,
TTree *  tree 
)

Analyze the Data

Definition at line 44 of file EventHeader.cc.

References _Debug, edm::EventBase::bunchCrossing(), gather_cfg::cout, edm::EventID::event(), fAvgInstDelLumi, fBx, fEvent, fLumiBlock, fOrbit, fRun, edm::LuminosityBlock::getByLabel(), edm::Event::getLuminosityBlock(), edm::EventBase::id(), edm::HandleBase::isValid(), edm::EventBase::luminosityBlock(), runregparse::lumiSummary, edm::EventBase::orbitNumber(), and edm::EventID::run().

Referenced by HLTAnalyzer::analyze(), and HLTBitAnalyzer::analyze().

                                                                {
                                        
  fRun          = iEvent.id().run();
  fEvent        = iEvent.id().event();
  fLumiBlock    = iEvent.luminosityBlock();
  fBx           = iEvent.bunchCrossing();
  fOrbit        = iEvent.orbitNumber();
  
 
  bool lumiException = false;
  const edm::LuminosityBlock& iLumi = iEvent.getLuminosityBlock(); 
  edm::Handle<LumiSummary> lumiSummary; 
  try{
    iLumi.getByLabel("lumiProducer", lumiSummary); 
    lumiSummary->isValid();
  }
  catch(cms::Exception&){
    lumiException = true;
  }
  if(!lumiException)
    fAvgInstDelLumi = lumiSummary->avgInsDelLumi(); 
  else 
    fAvgInstDelLumi = -999.; 
  
  
  if (_Debug) { 
    std::cout << "EventHeader -- run   = "          << fRun       << std::endl;
    std::cout << "EventHeader -- event = "          << fEvent     << std::endl;
    std::cout << "EventHeader -- lumisection = "    << fLumiBlock << std::endl; 
    std::cout << "EventHeader -- bunch crossing = " << fBx        << std::endl; 
    std::cout << "EventHeader -- orbit number = "   << fOrbit     << std::endl; 
  }

}
void EventHeader::setup ( TTree *  tree)

Definition at line 26 of file EventHeader.cc.

References fAvgInstDelLumi, fBx, fEvent, fLumiBlock, fOrbit, and fRun.

Referenced by HLTAnalyzer::HLTAnalyzer(), and HLTBitAnalyzer::HLTBitAnalyzer().

                                      {

  fRun = -1;
  fEvent = -1;
  fLumiBlock=-1;
  fBx = -1;
  fOrbit = -1;
  fAvgInstDelLumi = -999.; 

  HltTree->Branch("Run",       &fRun,         "Run/I");
  HltTree->Branch("Event",     &fEvent,       "Event/I");
  HltTree->Branch("LumiBlock", &fLumiBlock,   "LumiBlock/I"); 
  HltTree->Branch("Bx",        &fBx,          "Bx/I"); 
  HltTree->Branch("Orbit",     &fOrbit,       "Orbit/I"); 
  HltTree->Branch("AvgInstDelLumi", &fAvgInstDelLumi, "AvgInstDelLumi/D");
}

Member Data Documentation

bool EventHeader::_Debug [private]

Definition at line 38 of file EventHeader.h.

Referenced by analyze().

double EventHeader::fAvgInstDelLumi [private]

Definition at line 35 of file EventHeader.h.

Referenced by analyze(), and setup().

int EventHeader::fBx [private]

Definition at line 33 of file EventHeader.h.

Referenced by analyze(), and setup().

int EventHeader::fEvent [private]

Definition at line 31 of file EventHeader.h.

Referenced by analyze(), and setup().

int EventHeader::fLumiBlock [private]

Definition at line 32 of file EventHeader.h.

Referenced by analyze(), and setup().

int EventHeader::fOrbit [private]

Definition at line 34 of file EventHeader.h.

Referenced by analyze(), and setup().

int EventHeader::fRun [private]

Definition at line 30 of file EventHeader.h.

Referenced by analyze(), and setup().

Definition at line 59 of file EventMessage.h.

Referenced by EventMsgBuilder::setEventLength().

Definition at line 64 of file EventMessage.h.

Referenced by EventMsgView::origDataSize(), and EventMsgBuilder::setOrigDataSize().

Definition at line 65 of file EventMessage.h.

Referenced by EventMsgBuilder::EventMsgBuilder(), and EventMsgView::outModId().