CMS 3D CMS Logo

List of all members | Public Member Functions | Public Attributes | Private Attributes
EventHeader Class Reference

#include <EventHeader.h>

Public Member Functions

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

Public Attributes

char_uint32 droppedEventsCount_
 
char_uint64 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
 
unsigned long long fEvent
 
int fLumiBlock
 
int fOrbit
 
int fRun
 
edm::EDGetTokenT< LumiSummarylumi_Token
 

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}

Protocol Version 10: identical to version 9, but incremented to keep in sync with init msg version

Protocol Version 11: identical to version 10, except event changed from 4 bytes to 8 bytes

$Date: November 2006 $Revision:

Author
V. Rekovic - UMinn

Definition at line 18 of file EventHeader.h.

Constructor & Destructor Documentation

◆ EventHeader()

EventHeader::EventHeader ( )

Definition at line 11 of file EventHeader.cc.

12  : fEvent(0), fLumiBlock(-1), fRun(-1), fBx(-1), fOrbit(-1), fAvgInstDelLumi(-999.), _Debug(false) {}
int fLumiBlock
Definition: EventHeader.h:31
double fAvgInstDelLumi
Definition: EventHeader.h:35
unsigned long long fEvent
Definition: EventHeader.h:30

◆ ~EventHeader()

EventHeader::~EventHeader ( )
default

Member Function Documentation

◆ analyze()

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

Analyze the Data

Definition at line 36 of file EventHeader.cc.

References _Debug, gather_cfg::cout, fAvgInstDelLumi, fBx, fEvent, fLumiBlock, fOrbit, fRun, edm::LuminosityBlock::getByToken(), iEvent, lumi_Token, and bphysicsOniaDQM_cfi::lumiSummary.

Referenced by HLTBitAnalyzer::analyze().

36  {
37  fEvent = iEvent.id().event();
38  fLumiBlock = iEvent.luminosityBlock();
39  fRun = iEvent.id().run();
40  fBx = iEvent.bunchCrossing();
41  fOrbit = iEvent.orbitNumber();
42 
43  bool lumiException = false;
44  const edm::LuminosityBlock& iLumi = iEvent.getLuminosityBlock();
46  try {
48  lumiSummary->isValid();
49  } catch (cms::Exception&) {
50  lumiException = true;
51  }
52  if (!lumiException)
53  fAvgInstDelLumi = lumiSummary->avgInsDelLumi();
54  else
55  fAvgInstDelLumi = -999.;
56 
57  if (_Debug) {
58  std::cout << "EventHeader -- event = " << fEvent << std::endl;
59  std::cout << "EventHeader -- lumisection = " << fLumiBlock << std::endl;
60  std::cout << "EventHeader -- run = " << fRun << std::endl;
61  std::cout << "EventHeader -- bunch crossing = " << fBx << std::endl;
62  std::cout << "EventHeader -- orbit number = " << fOrbit << std::endl;
63  }
64 }
int iEvent
Definition: GenABIO.cc:224
int fLumiBlock
Definition: EventHeader.h:31
edm::EDGetTokenT< LumiSummary > lumi_Token
Definition: EventHeader.h:39
double fAvgInstDelLumi
Definition: EventHeader.h:35
bool getByToken(EDGetToken token, Handle< PROD > &result) const
unsigned long long fEvent
Definition: EventHeader.h:30

◆ setup()

void EventHeader::setup ( edm::ConsumesCollector &&  iC,
TTree *  tree 
)

Definition at line 17 of file EventHeader.cc.

References fAvgInstDelLumi, fBx, fEvent, fLumiBlock, fOrbit, fRun, edm::InLumi, HLT_2022v12_cff::InputTag, and lumi_Token.

Referenced by HLTBitAnalyzer::HLTBitAnalyzer().

17  {
18  fEvent = 0;
19  fLumiBlock = -1;
20  fRun = -1;
21  fBx = -1;
22  fOrbit = -1;
23  fAvgInstDelLumi = -999.;
24 
25  HltTree->Branch("Event", &fEvent, "Event/l");
26  HltTree->Branch("LumiBlock", &fLumiBlock, "LumiBlock/I");
27  HltTree->Branch("Run", &fRun, "Run/I");
28  HltTree->Branch("Bx", &fBx, "Bx/I");
29  HltTree->Branch("Orbit", &fOrbit, "Orbit/I");
30  HltTree->Branch("AvgInstDelLumi", &fAvgInstDelLumi, "AvgInstDelLumi/D");
31 
33 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
int fLumiBlock
Definition: EventHeader.h:31
edm::EDGetTokenT< LumiSummary > lumi_Token
Definition: EventHeader.h:39
double fAvgInstDelLumi
Definition: EventHeader.h:35
unsigned long long fEvent
Definition: EventHeader.h:30

Member Data Documentation

◆ _Debug

bool EventHeader::_Debug
private

Definition at line 38 of file EventHeader.h.

Referenced by analyze().

◆ droppedEventsCount_

char_uint32 EventHeader::droppedEventsCount_

Definition at line 69 of file EventMessage.h.

◆ event_

char_uint64 EventHeader::event_

Definition at line 65 of file EventMessage.h.

Referenced by edm::StreamerInputFile::readEventMessage().

◆ fAvgInstDelLumi

double EventHeader::fAvgInstDelLumi
private

Definition at line 35 of file EventHeader.h.

Referenced by analyze(), and setup().

◆ fBx

int EventHeader::fBx
private

Definition at line 33 of file EventHeader.h.

Referenced by analyze(), and setup().

◆ fEvent

unsigned long long EventHeader::fEvent
private

Definition at line 30 of file EventHeader.h.

Referenced by analyze(), and setup().

◆ fLumiBlock

int EventHeader::fLumiBlock
private

Definition at line 31 of file EventHeader.h.

Referenced by analyze(), and setup().

◆ fOrbit

int EventHeader::fOrbit
private

Definition at line 34 of file EventHeader.h.

Referenced by analyze(), and setup().

◆ fRun

int EventHeader::fRun
private

Definition at line 32 of file EventHeader.h.

Referenced by analyze(), and setup().

◆ header_

Header EventHeader::header_

Definition at line 62 of file EventMessage.h.

◆ lumi_

char_uint32 EventHeader::lumi_

Definition at line 66 of file EventMessage.h.

Referenced by edm::StreamerInputFile::readEventMessage().

◆ lumi_Token

edm::EDGetTokenT<LumiSummary> EventHeader::lumi_Token
private

Definition at line 39 of file EventHeader.h.

Referenced by analyze(), and setup().

◆ origDataSize_

char_uint32 EventHeader::origDataSize_

Definition at line 67 of file EventMessage.h.

◆ outModId_

char_uint32 EventHeader::outModId_

Definition at line 68 of file EventMessage.h.

◆ protocolVersion_

uint8 EventHeader::protocolVersion_

Definition at line 63 of file EventMessage.h.

◆ run_

char_uint32 EventHeader::run_

Definition at line 64 of file EventMessage.h.

Referenced by edm::StreamerInputFile::readEventMessage().