Main Page
Namespaces
Classes
Package Documentation
PhysicsTools
JetMCAlgos
plugins
GenJetBCEnergyRatio.cc
Go to the documentation of this file.
1
//
2
// Plugin to store B and C ratio for a GenJet in the event
3
// Author: Attilio
4
// Date: 05.10.2007
5
//
6
7
//=======================================================================
8
9
// user include files
10
#include "
FWCore/Framework/interface/global/EDProducer.h
"
11
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
12
#include "
FWCore/ParameterSet/interface/ParameterSetfwd.h
"
13
#include "
FWCore/Utilities/interface/InputTag.h
"
14
15
#include "
FWCore/Framework/interface/Event.h
"
16
#include "
FWCore/Framework/interface/EventSetup.h
"
17
#include "
FWCore/Framework/interface/MakerMacros.h
"
18
#include "
FWCore/Framework/interface/ESHandle.h
"
19
#include "
FWCore/Framework/interface/makeRefToBaseProdFrom.h
"
20
#include "
FWCore/MessageLogger/interface/MessageLogger.h
"
21
22
#include "
DataFormats/JetReco/interface/Jet.h
"
23
#include "
DataFormats/JetReco/interface/GenJetCollection.h
"
24
#include "
SimDataFormats/JetMatching/interface/JetFlavour.h
"
25
26
#include "
DataFormats/Common/interface/Ref.h
"
27
#include "
DataFormats/Candidate/interface/Candidate.h
"
28
#include "
DataFormats/Candidate/interface/CandidateFwd.h
"
29
#include "
DataFormats/JetReco/interface/JetFloatAssociation.h
"
30
#include "
DataFormats/Math/interface/Point3D.h
"
31
#include "
DataFormats/Math/interface/LorentzVector.h
"
32
33
#include "
DataFormats/JetReco/interface/GenJet.h
"
34
#include "
PhysicsTools/JetMCUtils/interface/JetMCTag.h
"
35
#include "
PhysicsTools/JetMCUtils/interface/CandMCTag.h
"
36
37
#include <memory>
38
#include <string>
39
#include <iostream>
40
#include <vector>
41
#include <Math/VectorUtil.h>
42
#include <TMath.h>
43
44
using namespace
std
;
45
using namespace
reco
;
46
using namespace
edm
;
47
using namespace
ROOT::Math::VectorUtil
;
48
using namespace
JetMCTagUtils
;
49
using namespace
CandMCTagUtils
;
50
51
class
GenJetBCEnergyRatio
:
public
edm::global::EDProducer
<> {
52
public
:
53
GenJetBCEnergyRatio
(
const
edm::ParameterSet
&);
54
~
GenJetBCEnergyRatio
()
override
;
55
56
typedef
reco::JetFloatAssociation::Container
JetBCEnergyRatioCollection
;
57
58
private
:
59
void
produce(
StreamID
,
edm::Event
&,
const
edm::EventSetup
&)
const override
;
60
edm::EDGetTokenT<View<Jet>
>
m_genjetsSrcToken
;
61
};
62
63
//=========================================================================
64
65
GenJetBCEnergyRatio::GenJetBCEnergyRatio
(
const
edm::ParameterSet
& iConfig) {
66
produces<JetBCEnergyRatioCollection>(
"bRatioCollection"
);
67
produces<JetBCEnergyRatioCollection>(
"cRatioCollection"
);
68
m_genjetsSrcToken = consumes<View<Jet> >(iConfig.
getParameter
<
edm::InputTag
>(
"genJets"
));
69
}
70
71
//=========================================================================
72
73
GenJetBCEnergyRatio::~GenJetBCEnergyRatio
() {}
74
75
// ------------ method called to produce the data ------------
76
77
void
GenJetBCEnergyRatio::produce
(
StreamID
,
Event
&
iEvent
,
const
EventSetup
& iEs)
const
{
78
Handle<View<Jet>
> genjets;
79
iEvent.
getByToken
(m_genjetsSrcToken, genjets);
80
81
typedef
edm::RefToBase<reco::Jet>
JetRef
;
82
83
JetBCEnergyRatioCollection
* jtc1;
84
JetBCEnergyRatioCollection
* jtc2;
85
86
if
(!genjets.
product
()->empty()) {
87
const
JetRef
jj
= genjets->refAt(0);
88
jtc1 =
new
JetBCEnergyRatioCollection
(
edm::makeRefToBaseProdFrom
(jj, iEvent));
89
jtc2 =
new
JetBCEnergyRatioCollection
(
edm::makeRefToBaseProdFrom
(jj, iEvent));
90
}
else
{
91
jtc1 =
new
JetBCEnergyRatioCollection
();
92
jtc2 =
new
JetBCEnergyRatioCollection
();
93
}
94
95
std::unique_ptr<JetBCEnergyRatioCollection> bRatioColl(jtc1);
96
std::unique_ptr<JetBCEnergyRatioCollection> cRatioColl(jtc2);
97
98
for
(
size_t
j
= 0;
j
!= genjets->size(); ++
j
) {
99
float
bRatio =
EnergyRatioFromBHadrons
((*genjets)[
j
]);
100
float
cRatio =
EnergyRatioFromCHadrons
((*genjets)[j]);
101
102
const
JetRef& aJet = genjets->refAt(j);
103
104
JetFloatAssociation::setValue
(*bRatioColl, aJet, bRatio);
105
JetFloatAssociation::setValue
(*cRatioColl, aJet, cRatio);
106
}
107
108
iEvent.
put
(
std::move
(bRatioColl),
"bRatioCollection"
);
109
iEvent.
put
(
std::move
(cRatioColl),
"cRatioCollection"
);
110
}
111
112
//define this as a plug-in
113
DEFINE_FWK_MODULE
(
GenJetBCEnergyRatio
);
JetMCTag.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
makeRefToBaseProdFrom.h
EDProducer.h
MessageLogger.h
edm::Event::put
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition:
Event.h:131
JetFloatAssociation.h
dqmiolumiharvest.j
j
Definition:
dqmiolumiharvest.py:66
edm::Event::getByToken
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition:
Event.h:525
edm::StreamID
Definition:
StreamID.h:30
Event.h
GenJetBCEnergyRatio
Definition:
GenJetBCEnergyRatio.cc:51
MakerMacros.h
edm::Handle
Definition:
AssociativeIterator.h:49
EventSetup.h
l1t::JetRef
edm::Ref< JetBxCollection > JetRef
Definition:
Jet.h:12
CandMCTag.h
JetFlavour.h
JetMCTagUtils
Definition:
JetMCTag.h:7
std
Definition:
JetResolutionObject.h:76
edm::EDGetTokenT
Definition:
EDGetToken.h:33
reco::JetExtendedAssociation::setValue
bool setValue(Container &, const reco::JetBaseRef &, const JetExtendedData &)
associate jet with value. Returns false and associate nothing if jet is already associated ...
Definition:
JetExtendedAssociation.cc:44
Point3D.h
edm::RefToBase< reco::Jet >
ParameterSetfwd.h
ParameterSet.h
Candidate.h
GenJetBCEnergyRatio::~GenJetBCEnergyRatio
~GenJetBCEnergyRatio() override
Definition:
GenJetBCEnergyRatio.cc:73
iEvent
int iEvent
Definition:
GenABIO.cc:224
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition:
MakerMacros.h:16
JetMCTagUtils::EnergyRatioFromBHadrons
double EnergyRatioFromBHadrons(const reco::Candidate &c)
Definition:
JetMCTag.cc:12
GenJetCollection.h
ESHandle.h
VectorUtil
edm::EventSetup
Definition:
EventSetup.h:57
Jet.h
edm::global::EDProducer
Definition:
EDProducer.h:32
LorentzVector.h
GenJetBCEnergyRatio::GenJetBCEnergyRatio
GenJetBCEnergyRatio(const edm::ParameterSet &)
Definition:
GenJetBCEnergyRatio.cc:65
edm::AssociationVector
Definition:
AssociationVector.h:67
edm::Handle::product
T const * product() const
Definition:
Handle.h:69
Ref.h
findQualityFiles.jj
string jj
Definition:
findQualityFiles.py:188
edm::makeRefToBaseProdFrom
RefToBaseProd< T > makeRefToBaseProdFrom(RefToBase< T > const &iRef, Event const &iEvent)
Definition:
makeRefToBaseProdFrom.h:34
reco
fixed size matrix
Definition:
AlignmentAlgorithmBase.h:45
edm
HLT enums.
Definition:
AlignableModifier.h:19
edm::InputTag
Definition:
InputTag.h:15
InputTag.h
JetMCTagUtils::EnergyRatioFromCHadrons
double EnergyRatioFromCHadrons(const reco::Candidate &c)
Definition:
JetMCTag.cc:24
edm::ParameterSet
Definition:
ParameterSet.h:36
CandidateFwd.h
edm::Event
Definition:
Event.h:72
CandMCTagUtils
Definition:
CandMCTag.h:6
eostools.move
def move(src, dest)
Definition:
eostools.py:511
GenJetBCEnergyRatio::produce
void produce(StreamID, edm::Event &, const edm::EventSetup &) const override
Definition:
GenJetBCEnergyRatio.cc:77
GenJetBCEnergyRatio::m_genjetsSrcToken
edm::EDGetTokenT< View< Jet > > m_genjetsSrcToken
Definition:
GenJetBCEnergyRatio.cc:60
GenJet.h
GenJetBCEnergyRatio::JetBCEnergyRatioCollection
reco::JetFloatAssociation::Container JetBCEnergyRatioCollection
Definition:
GenJetBCEnergyRatio.cc:56
Generated for CMSSW Reference Manual by
1.8.11