CMS 3D CMS Logo

EventHeader.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <sstream>
3 #include <istream>
4 #include <fstream>
5 #include <iomanip>
6 #include <cstdlib>
7 #include <cstring>
8 
10 
12  fEvent( 0 ),
13  fLumiBlock( -1 ),
14  fRun( -1 ),
15  fBx( -1 ),
16  fOrbit( -1 ),
17  fAvgInstDelLumi( -999. ),
18  _Debug( false )
19 { }
20 
21 EventHeader::~EventHeader() = default;
22 
23 /* Setup the analysis to put the branch-variables into the tree. */
24 void EventHeader::setup(edm::ConsumesCollector && iC, TTree* HltTree) {
25 
26  fEvent = 0;
27  fLumiBlock=-1;
28  fRun = -1;
29  fBx = -1;
30  fOrbit = -1;
31  fAvgInstDelLumi = -999.;
32 
33  HltTree->Branch("Event", &fEvent, "Event/l");
34  HltTree->Branch("LumiBlock", &fLumiBlock, "LumiBlock/I");
35  HltTree->Branch("Run", &fRun, "Run/I");
36  HltTree->Branch("Bx", &fBx, "Bx/I");
37  HltTree->Branch("Orbit", &fOrbit, "Orbit/I");
38  HltTree->Branch("AvgInstDelLumi", &fAvgInstDelLumi, "AvgInstDelLumi/D");
39 
40  lumi_Token = iC.consumes<LumiSummary,edm::InLumi>(edm::InputTag("lumiProducer"));
41 }
42 
43 /* **Analyze the event** */
44 void EventHeader::analyze(edm::Event const& iEvent, TTree* HltTree) {
45 
46  fEvent = iEvent.id().event();
47  fLumiBlock = iEvent.luminosityBlock();
48  fRun = iEvent.id().run();
49  fBx = iEvent.bunchCrossing();
50  fOrbit = iEvent.orbitNumber();
51 
52 
53  bool lumiException = false;
54  const edm::LuminosityBlock& iLumi = iEvent.getLuminosityBlock();
56  try{
57  iLumi.getByToken(lumi_Token, lumiSummary);
58  lumiSummary->isValid();
59  }
60  catch(cms::Exception&){
61  lumiException = true;
62  }
63  if(!lumiException)
64  fAvgInstDelLumi = lumiSummary->avgInsDelLumi();
65  else
66  fAvgInstDelLumi = -999.;
67 
68 
69  if (_Debug) {
70  std::cout << "EventHeader -- event = " << fEvent << std::endl;
71  std::cout << "EventHeader -- lumisection = " << fLumiBlock << std::endl;
72  std::cout << "EventHeader -- run = " << fRun << std::endl;
73  std::cout << "EventHeader -- bunch crossing = " << fBx << std::endl;
74  std::cout << "EventHeader -- orbit number = " << fOrbit << std::endl;
75  }
76 
77 }
RunNumber_t run() const
Definition: EventID.h:39
EventNumber_t event() const
Definition: EventID.h:41
void analyze(edm::Event const &iEvent, TTree *tree)
Definition: EventHeader.cc:44
void setup(edm::ConsumesCollector &&iC, TTree *tree)
Definition: EventHeader.cc:24
bool getByToken(EDGetToken token, Handle< PROD > &result) const
int bunchCrossing() const
Definition: EventBase.h:66
edm::LuminosityBlockNumber_t luminosityBlock() const
Definition: EventBase.h:63
int iEvent
Definition: GenABIO.cc:230
int fLumiBlock
Definition: EventHeader.h:32
edm::EDGetTokenT< LumiSummary > lumi_Token
Definition: EventHeader.h:40
def lumiSummary(schema, nlumils)
LuminosityBlock const & getLuminosityBlock() const
Definition: Event.h:91
int orbitNumber() const
Definition: EventBase.h:67
double fAvgInstDelLumi
Definition: EventHeader.h:36
float avgInsDelLumi() const
Definition: LumiSummary.cc:9
unsigned long long fEvent
Definition: EventHeader.h:31
edm::EventID id() const
Definition: EventBase.h:60
bool isValid() const
Definition: LumiSummary.cc:78