HLTrigger
HLTanalyzers
plugins
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
9
#include "
EventHeader.h
"
10
11
EventHeader::EventHeader
()
12
:
fEvent
(0), fLumiBlock(-1), fRun(-1), fBx(-1), fOrbit(-1), fAvgInstDelLumi(-999.), _Debug(
false
) {}
13
14
EventHeader::~EventHeader
() =
default
;
15
16
/* Setup the analysis to put the branch-variables into the tree. */
17
void
EventHeader::setup
(
edm::ConsumesCollector
&& iC, TTree* HltTree) {
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
32
lumi_Token
= iC.consumes<
LumiSummary
,
edm::InLumi
>(
edm::InputTag
(
"lumiProducer"
));
33
}
34
35
/* **Analyze the event** */
36
void
EventHeader::analyze
(
edm::Event
const
&
iEvent
, TTree* HltTree) {
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();
45
edm::Handle<LumiSummary>
lumiSummary
;
46
try
{
47
iLumi.
getByToken
(
lumi_Token
,
lumiSummary
);
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
}
EventHeader::analyze
void analyze(edm::Event const &iEvent, TTree *tree)
Definition:
EventHeader.cc:36
EventHeader::setup
void setup(edm::ConsumesCollector &&iC, TTree *tree)
Definition:
EventHeader.cc:17
bphysicsOniaDQM_cfi.lumiSummary
lumiSummary
Definition:
bphysicsOniaDQM_cfi.py:8
LumiSummary
Definition:
LumiSummary.h:20
EventHeader::_Debug
bool _Debug
Definition:
EventHeader.h:38
funct::false
false
Definition:
Factorize.h:29
edm::LuminosityBlock
Definition:
LuminosityBlock.h:50
edm::Handle
Definition:
AssociativeIterator.h:50
HLT_2022v15_cff.InputTag
InputTag
Definition:
HLT_2022v15_cff.py:60424
EventHeader::~EventHeader
~EventHeader()
EventHeader::fBx
int fBx
Definition:
EventHeader.h:33
edm::InLumi
Definition:
BranchType.h:11
iEvent
int iEvent
Definition:
GenABIO.cc:224
EventHeader::fLumiBlock
int fLumiBlock
Definition:
EventHeader.h:31
hcaldqm::fEvent
Definition:
DQTask.h:32
EventHeader::lumi_Token
edm::EDGetTokenT< LumiSummary > lumi_Token
Definition:
EventHeader.h:39
EventHeader.h
EventHeader::fRun
int fRun
Definition:
EventHeader.h:32
EventHeader::EventHeader
EventHeader()
Definition:
EventHeader.cc:11
EventHeader::fOrbit
int fOrbit
Definition:
EventHeader.h:34
EventHeader::fAvgInstDelLumi
double fAvgInstDelLumi
Definition:
EventHeader.h:35
edm::LuminosityBlock::getByToken
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition:
LuminosityBlock.h:321
cms::Exception
Definition:
Exception.h:70
EventHeader::fEvent
unsigned long long fEvent
Definition:
EventHeader.h:30
gather_cfg.cout
cout
Definition:
gather_cfg.py:144
edm::Event
Definition:
Event.h:73
edm::ConsumesCollector
Definition:
ConsumesCollector.h:45
Generated for CMSSW Reference Manual by
1.8.14