Main Page
Namespaces
Classes
Package Documentation
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Pages
src
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/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/MessageLogger/interface/MessageLogger.h
"
20
21
#include "
DataFormats/JetReco/interface/Jet.h
"
22
#include "
DataFormats/JetReco/interface/GenJetCollection.h
"
23
#include "
SimDataFormats/JetMatching/interface/JetFlavour.h
"
24
25
#include "
DataFormats/Common/interface/Ref.h
"
26
#include "
DataFormats/Candidate/interface/Candidate.h
"
27
#include "
DataFormats/Candidate/interface/CandidateFwd.h
"
28
#include "
DataFormats/JetReco/interface/JetFloatAssociation.h
"
29
#include "
DataFormats/Math/interface/Point3D.h
"
30
#include "
DataFormats/Math/interface/LorentzVector.h
"
31
32
#include "
DataFormats/JetReco/interface/GenJet.h
"
33
#include "
PhysicsTools/JetMCUtils/interface/JetMCTag.h
"
34
#include "
PhysicsTools/JetMCUtils/interface/CandMCTag.h
"
35
36
#include <memory>
37
#include <string>
38
#include <iostream>
39
#include <vector>
40
#include <Math/VectorUtil.h>
41
#include <TMath.h>
42
43
using namespace
std;
44
using namespace
reco
;
45
using namespace
edm;
46
using namespace
ROOT::Math::VectorUtil;
47
using namespace
JetMCTagUtils;
48
using namespace
CandMCTagUtils;
49
50
class
GenJetBCEnergyRatio
:
public
edm::EDProducer
51
{
52
public
:
53
GenJetBCEnergyRatio
(
const
edm::ParameterSet
& );
54
~
GenJetBCEnergyRatio
();
55
56
typedef
reco::JetFloatAssociation::Container
JetBCEnergyRatioCollection
;
57
58
private
:
59
virtual
void
produce(
edm::Event
&,
const
edm::EventSetup
& )
override
;
60
Handle< View <Jet>
>
genjets
;
61
edm::EDGetTokenT< View <Jet>
>
m_genjetsSrcToken
;
62
63
};
64
65
//=========================================================================
66
67
GenJetBCEnergyRatio::GenJetBCEnergyRatio
(
const
edm::ParameterSet
& iConfig )
68
{
69
produces<JetBCEnergyRatioCollection>(
"bRatioCollection"
);
70
produces<JetBCEnergyRatioCollection>(
"cRatioCollection"
);
71
m_genjetsSrcToken = consumes< View <Jet> >(iConfig.
getParameter
<
edm::InputTag
>(
"genJets"
));
72
}
73
74
//=========================================================================
75
76
GenJetBCEnergyRatio::~GenJetBCEnergyRatio
()
77
{
78
}
79
80
// ------------ method called to produce the data ------------
81
82
void
GenJetBCEnergyRatio::produce
(
Event
&
iEvent
,
const
EventSetup
& iEs )
83
{
84
iEvent.
getByToken
(m_genjetsSrcToken, genjets);
85
86
typedef
edm::RefToBase<reco::Jet>
JetRef
;
87
88
JetBCEnergyRatioCollection
* jtc1;
89
JetBCEnergyRatioCollection
* jtc2;
90
91
if
(genjets.product()->size() > 0) {
92
const
JetRef
jj
= genjets->refAt(0);
93
jtc1 =
new
JetBCEnergyRatioCollection
(
RefToBaseProd<Jet>
(jj));
94
jtc2 =
new
JetBCEnergyRatioCollection
(
RefToBaseProd<Jet>
(jj));
95
}
else
{
96
jtc1 =
new
JetBCEnergyRatioCollection
();
97
jtc2 =
new
JetBCEnergyRatioCollection
();
98
}
99
100
std::auto_ptr<JetBCEnergyRatioCollection> bRatioColl(jtc1);
101
std::auto_ptr<JetBCEnergyRatioCollection> cRatioColl(jtc2);
102
103
for
(
size_t
j
= 0;
j
!= genjets->size(); ++
j
) {
104
105
float
bRatio =
EnergyRatioFromBHadrons
( (*genjets)[
j
] );
106
float
cRatio =
EnergyRatioFromCHadrons
( (*genjets)[j] );
107
108
const
JetRef & aJet = genjets->refAt(j) ;
109
110
JetFloatAssociation::setValue
(*bRatioColl, aJet, bRatio);
111
JetFloatAssociation::setValue
(*cRatioColl, aJet, cRatio);
112
113
}
114
115
116
iEvent.
put
(bRatioColl,
"bRatioCollection"
);
117
iEvent.
put
(cRatioColl,
"cRatioCollection"
);
118
119
}
120
121
//define this as a plug-in
122
DEFINE_FWK_MODULE
(
GenJetBCEnergyRatio
);
123
JetMCTag.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
MessageLogger.h
JetFloatAssociation.h
GenJetBCEnergyRatio::produce
virtual void produce(edm::Event &, const edm::EventSetup &) override
Definition:
GenJetBCEnergyRatio.cc:82
edm::Event::getByToken
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition:
Event.h:434
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition:
MakerMacros.h:17
Event.h
GenJetBCEnergyRatio
Definition:
GenJetBCEnergyRatio.cc:50
MakerMacros.h
EventSetup.h
CandMCTag.h
JetFlavour.h
pat::JetRef
edm::Ref< JetCollection > JetRef
Definition:
Jet.h:51
edm::RefToBase< reco::Jet >
edm::Handle
Definition:
AssociativeIterator.h:47
dt_dqm_sourceclient_common_cff.reco
tuple reco
Definition:
dt_dqm_sourceclient_common_cff.py:105
edm::EDGetTokenT
Definition:
EDGetToken.h:32
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
ParameterSetfwd.h
edm::EDProducer
Definition:
EDProducer.h:30
ParameterSet.h
Candidate.h
iEvent
int iEvent
Definition:
GenABIO.cc:243
edm::RefToBaseProd
Definition:
RefToBase.h:62
edm::Event::put
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition:
Event.h:116
JetMCTagUtils::EnergyRatioFromBHadrons
double EnergyRatioFromBHadrons(const reco::Candidate &c)
Definition:
JetMCTag.cc:12
GenJetCollection.h
ESHandle.h
j
int j
Definition:
DBlmapReader.cc:9
edm::EventSetup
Definition:
EventSetup.h:44
Jet.h
LorentzVector.h
GenJetBCEnergyRatio::GenJetBCEnergyRatio
GenJetBCEnergyRatio(const edm::ParameterSet &)
Definition:
GenJetBCEnergyRatio.cc:67
EDProducer.h
edm::AssociationVector< reco::JetRefBaseProd, Values >
GenJetBCEnergyRatio::m_genjetsSrcToken
edm::EDGetTokenT< View< Jet > > m_genjetsSrcToken
Definition:
GenJetBCEnergyRatio.cc:61
Ref.h
GenJetBCEnergyRatio::~GenJetBCEnergyRatio
~GenJetBCEnergyRatio()
Definition:
GenJetBCEnergyRatio.cc:76
edm::InputTag
Definition:
InputTag.h:17
InputTag.h
JetMCTagUtils::EnergyRatioFromCHadrons
double EnergyRatioFromCHadrons(const reco::Candidate &c)
Definition:
JetMCTag.cc:26
GenJetBCEnergyRatio::genjets
Handle< View< Jet > > genjets
Definition:
GenJetBCEnergyRatio.cc:60
edm::ParameterSet
Definition:
ParameterSet.h:35
CandidateFwd.h
edm::Event
Definition:
Event.h:62
findQualityFiles.jj
string jj
Definition:
findQualityFiles.py:186
GenJet.h
GenJetBCEnergyRatio::JetBCEnergyRatioCollection
reco::JetFloatAssociation::Container JetBCEnergyRatioCollection
Definition:
GenJetBCEnergyRatio.cc:56
Generated for CMSSW Reference Manual by
1.8.5