CMS 3D CMS Logo

BJetEnergyRegressionMVA.cc
Go to the documentation of this file.
1 //
2 // Original Author: Andre Rizzi
3 // Created: Mon, 07 Sep 2017 09:18:03 GMT
4 //
5 //
6 
9 
12 
15 
19 
22 
26 
28 #include <vector>
29 
31 public:
32  explicit BJetEnergyRegressionMVA(const edm::ParameterSet& iConfig)
33  : BaseMVAValueMapProducer<pat::Jet>(iConfig),
34  pvsrc_(edm::stream::EDProducer<>::consumes<std::vector<reco::Vertex>>(
35  iConfig.getParameter<edm::InputTag>("pvsrc"))),
36  svsrc_(edm::stream::EDProducer<>::consumes<edm::View<reco::VertexCompositePtrCandidate>>(
37  iConfig.getParameter<edm::InputTag>("svsrc"))),
38  rhosrc_(edm::stream::EDProducer<>::consumes<double>(iConfig.getParameter<edm::InputTag>("rhosrc")))
39 
40  {}
42  iEvent.getByToken(pvsrc_, pvs_);
43  iEvent.getByToken(svsrc_, svs_);
44  iEvent.getByToken(rhosrc_, rho_);
45  }
46 
47  void fillAdditionalVariables(const pat::Jet& j) override {
48  this->setValue("nPVs", pvs_->size());
49  this->setValue("rho", *(rho_.product()));
50 
51  float cone_boundaries[] = {0.05, 0.1, 0.2, 0.3, 0.4};
52  size_t ncone_boundaries = sizeof(cone_boundaries) / sizeof(float);
53  std::vector<float> emFractionEnergyRings(ncone_boundaries + 1);
54  std::vector<float> chFractionEnergyRings(ncone_boundaries + 1);
55  std::vector<float> neFractionEnergyRings(ncone_boundaries + 1);
56  std::vector<float> muFractionEnergyRings(ncone_boundaries + 1);
57  float jetRawEnergy = j.p4().E() * j.jecFactor("Uncorrected");
58  int numDaughtersPt03 = 0;
59  for (unsigned int ijcone = 0; ijcone < ncone_boundaries; ijcone++) {
60  emFractionEnergyRings[ijcone] = 0;
61  muFractionEnergyRings[ijcone] = 0;
62  chFractionEnergyRings[ijcone] = 0;
63  neFractionEnergyRings[ijcone] = 0;
64  }
65  for (const auto& d : j.daughterPtrVector()) {
66  float candDr = Geom::deltaR(d->p4(), j.p4());
67  size_t icone =
68  std::lower_bound(&cone_boundaries[0], &cone_boundaries[ncone_boundaries], candDr) - &cone_boundaries[0];
69  float candEnergy = d->energy() / jetRawEnergy;
70  int pdgid = abs(d->pdgId());
71  if (pdgid == 22 || pdgid == 11) {
72  emFractionEnergyRings[icone] += candEnergy;
73  } else if (pdgid == 13) {
74  muFractionEnergyRings[icone] += candEnergy;
75  } else if (d->charge() != 0) {
76  chFractionEnergyRings[icone] += candEnergy;
77  } else {
78  neFractionEnergyRings[icone] += candEnergy;
79  }
80  if (d->pt() > 0.3)
81  numDaughtersPt03 += 1;
82  } // end of jet daughters loop
83 
84  this->setValue("Jet_energyRing_dR0_em_Jet_rawEnergy", emFractionEnergyRings[0]);
85  this->setValue("Jet_energyRing_dR1_em_Jet_rawEnergy", emFractionEnergyRings[1]);
86  this->setValue("Jet_energyRing_dR2_em_Jet_rawEnergy", emFractionEnergyRings[2]);
87  this->setValue("Jet_energyRing_dR3_em_Jet_rawEnergy", emFractionEnergyRings[3]);
88  this->setValue("Jet_energyRing_dR4_em_Jet_rawEnergy", emFractionEnergyRings[4]);
89  // this->setValue("Jet_energyRing_dR5_em_Jet_rawEnergy", emFractionEnergyRings[5]);
90 
91  this->setValue("Jet_energyRing_dR0_ch_Jet_rawEnergy", chFractionEnergyRings[0]);
92  this->setValue("Jet_energyRing_dR1_ch_Jet_rawEnergy", chFractionEnergyRings[1]);
93  this->setValue("Jet_energyRing_dR2_ch_Jet_rawEnergy", chFractionEnergyRings[2]);
94  this->setValue("Jet_energyRing_dR3_ch_Jet_rawEnergy", chFractionEnergyRings[3]);
95  this->setValue("Jet_energyRing_dR4_ch_Jet_rawEnergy", chFractionEnergyRings[4]);
96  // this->setValue("Jet_energyRing_dR5_ch_Jet_rawEnergy", chFractionEnergyRings[5]);
97 
98  this->setValue("Jet_energyRing_dR0_mu_Jet_rawEnergy", muFractionEnergyRings[0]);
99  this->setValue("Jet_energyRing_dR1_mu_Jet_rawEnergy", muFractionEnergyRings[1]);
100  this->setValue("Jet_energyRing_dR2_mu_Jet_rawEnergy", muFractionEnergyRings[2]);
101  this->setValue("Jet_energyRing_dR3_mu_Jet_rawEnergy", muFractionEnergyRings[3]);
102  this->setValue("Jet_energyRing_dR4_mu_Jet_rawEnergy", muFractionEnergyRings[4]);
103  // this->setValue("Jet_energyRing_dR5_mu_Jet_rawEnergy", muFractionEnergyRings[5]);
104 
105  this->setValue("Jet_energyRing_dR0_neut_Jet_rawEnergy", neFractionEnergyRings[0]);
106  this->setValue("Jet_energyRing_dR1_neut_Jet_rawEnergy", neFractionEnergyRings[1]);
107  this->setValue("Jet_energyRing_dR2_neut_Jet_rawEnergy", neFractionEnergyRings[2]);
108  this->setValue("Jet_energyRing_dR3_neut_Jet_rawEnergy", neFractionEnergyRings[3]);
109  this->setValue("Jet_energyRing_dR4_neut_Jet_rawEnergy", neFractionEnergyRings[4]);
110  // this->setValue("Jet_energyRing_dR5_neut_Jet_rawEnergy", neFractionEnergyRings[5]);
111 
112  this->setValue("Jet_numDaughters_pt03", numDaughtersPt03);
113  }
114 
117  desc.add<edm::InputTag>("pvsrc")->setComment("primary vertices input collection");
118  desc.add<edm::InputTag>("svsrc")->setComment("secondary vertices input collection");
119  desc.add<edm::InputTag>("rhosrc")->setComment("rho input collection");
120  descriptions.add("BJetEnergyRegressionMVA", desc);
121  }
122 
123 private:
130 };
131 
132 //define this as a plug-in
BJetEnergyRegressionMVA
Definition: BJetEnergyRegressionMVA.cc:30
edm::Handle::product
T const * product() const
Definition: Handle.h:70
sistrip::View
View
Definition: ConstantsForView.h:26
edm::EDGetTokenT
Definition: EDGetToken.h:33
edm
HLT enums.
Definition: AlignableModifier.h:19
Muon.h
BaseMVAValueMapProducer< pat::Jet >::setValue
void setValue(const std::string var, float val)
Definition: BaseMVAValueMapProducer.h:112
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89301
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
cms::cuda::stream
uint32_t const T *__restrict__ const uint32_t *__restrict__ int32_t int Histo::index_type cudaStream_t stream
Definition: HistoContainer.h:51
EDProducer.h
BJetEnergyRegressionMVA::svsrc_
const edm::EDGetTokenT< edm::View< reco::VertexCompositePtrCandidate > > svsrc_
Definition: BJetEnergyRegressionMVA.cc:126
ConvertToFromReco.h
VertexDistance3D.h
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:46
edm::Handle
Definition: AssociativeIterator.h:50
BJetEnergyRegressionMVA::pvs_
edm::Handle< std::vector< reco::Vertex > > pvs_
Definition: BJetEnergyRegressionMVA.cc:125
MakerMacros.h
pat::Jet
Analysis-level calorimeter jet class.
Definition: Jet.h:77
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Jet
Definition: Jet.py:1
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
BJetEnergyRegressionMVA::svs_
edm::Handle< edm::View< reco::VertexCompositePtrCandidate > > svs_
Definition: BJetEnergyRegressionMVA.cc:127
BJetEnergyRegressionMVA::rho_
edm::Handle< double > rho_
Definition: BJetEnergyRegressionMVA.cc:129
VertexCompositePtrCandidate.h
BJetEnergyRegressionMVA::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: BJetEnergyRegressionMVA.cc:115
VertexState.h
BaseMVAValueMapProducer
Definition: BaseMVAValueMapProducer.h:56
PbPb_ZMuSkimMuonDPG_cff.deltaR
deltaR
Definition: PbPb_ZMuSkimMuonDPG_cff.py:63
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
BaseMVAValueMapProducer.h
BJetEnergyRegressionMVA::BJetEnergyRegressionMVA
BJetEnergyRegressionMVA(const edm::ParameterSet &iConfig)
Definition: BJetEnergyRegressionMVA.cc:32
edm::stream::EDProducer<>::EDProducer
EDProducer()=default
BJetEnergyRegressionMVA::pvsrc_
const edm::EDGetTokenT< std::vector< reco::Vertex > > pvsrc_
Definition: BJetEnergyRegressionMVA.cc:124
pfDeepBoostedJetPreprocessParams_cfi.lower_bound
lower_bound
Definition: pfDeepBoostedJetPreprocessParams_cfi.py:15
edm::ParameterSet
Definition: ParameterSet.h:47
Event.h
iEvent
int iEvent
Definition: GenABIO.cc:224
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
edm::EventSetup
Definition: EventSetup.h:58
pat
Definition: HeavyIon.h:7
Jet.h
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
std
Definition: JetResolutionObject.h:76
HltBtagValidation_cff.Vertex
Vertex
Definition: HltBtagValidation_cff.py:32
Vertex.h
Frameworkfwd.h
BJetEnergyRegressionMVA::rhosrc_
edm::EDGetTokenT< double > rhosrc_
Definition: BJetEnergyRegressionMVA.cc:128
BJetEnergyRegressionMVA::fillAdditionalVariables
void fillAdditionalVariables(const pat::Jet &j) override
Definition: BJetEnergyRegressionMVA.cc:47
Electron.h
BaseMVAValueMapProducer::getDescription
static edm::ParameterSetDescription getDescription()
Definition: BaseMVAValueMapProducer.h:201
ztail.d
d
Definition: ztail.py:151
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterSet.h
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
edm::Event
Definition: Event.h:73
EgammaValidation_cff.pdgid
pdgid
Definition: EgammaValidation_cff.py:29
BJetEnergyRegressionMVA::readAdditionalCollections
void readAdditionalCollections(edm::Event &iEvent, const edm::EventSetup &) override
to be implemented in derived classes, filling values for additional variables
Definition: BJetEnergyRegressionMVA.cc:41
StreamID.h
edm::InputTag
Definition: InputTag.h:15