CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Attributes
BJetEnergyRegressionMVA Class Reference
Inheritance diagram for BJetEnergyRegressionMVA:
BaseMVAValueMapProducer< pat::Jet > edm::stream::EDProducer<>

Public Member Functions

 BJetEnergyRegressionMVA (const edm::ParameterSet &iConfig)
 
void fillAdditionalVariables (const pat::Jet &j) override
 
void readAdditionalCollections (edm::Event &iEvent, const edm::EventSetup &) override
 to be implemented in derived classes, filling values for additional variables More...
 
- Public Member Functions inherited from BaseMVAValueMapProducer< pat::Jet >
 BaseMVAValueMapProducer (const edm::ParameterSet &iConfig)
 
void setValue (const std::string var, float val)
 
 ~BaseMVAValueMapProducer () override
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
bool hasAbilityToProduceInLumis () const final
 
bool hasAbilityToProduceInRuns () const final
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 
- Static Public Member Functions inherited from BaseMVAValueMapProducer< pat::Jet >
static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 
static edm::ParameterSetDescription getDescription ()
 

Private Attributes

edm::Handle< std::vector< reco::Vertex > > pvs_
 
const edm::EDGetTokenT< std::vector< reco::Vertex > > pvsrc_
 
edm::Handle< double > rho_
 
edm::EDGetTokenT< double > rhosrc_
 
edm::Handle< edm::View< reco::VertexCompositePtrCandidate > > svs_
 
const edm::EDGetTokenT< edm::View< reco::VertexCompositePtrCandidate > > svsrc_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
typedef CacheContexts< T... > CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T... > HasAbility
 
typedef CacheTypes::LuminosityBlockCache LuminosityBlockCache
 
typedef LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCacheLuminosityBlockContext
 
typedef CacheTypes::LuminosityBlockSummaryCache LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache, GlobalCacheRunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 

Detailed Description

Definition at line 31 of file BJetEnergyRegressionMVA.cc.

Constructor & Destructor Documentation

BJetEnergyRegressionMVA::BJetEnergyRegressionMVA ( const edm::ParameterSet iConfig)
inlineexplicit

Definition at line 33 of file BJetEnergyRegressionMVA.cc.

33  :
35  pvsrc_(edm::stream::EDProducer<>::consumes<std::vector<reco::Vertex>>(iConfig.getParameter<edm::InputTag>("pvsrc"))),
38 
39  {
40 
41 
42  }
T getParameter(std::string const &) const
const edm::EDGetTokenT< edm::View< reco::VertexCompositePtrCandidate > > svsrc_
edm::EDGetTokenT< double > rhosrc_
const edm::EDGetTokenT< std::vector< reco::Vertex > > pvsrc_

Member Function Documentation

void BJetEnergyRegressionMVA::fillAdditionalVariables ( const pat::Jet j)
inlineoverridevirtual

Reimplemented from BaseMVAValueMapProducer< pat::Jet >.

Definition at line 50 of file BJetEnergyRegressionMVA.cc.

References funct::abs(), edmIntegrityCheck::d, pat::Jet::daughterPtrVector(), deltaR(), pat::Jet::jecFactor(), reco::LeafCandidate::p4(), BPhysicsValidation_cfi::pdgid, edm::Handle< T >::product(), pvs_, rho_, and BaseMVAValueMapProducer< pat::Jet >::setValue().

50  {
51 
52  this->setValue("nPVs",pvs_->size());
53  this->setValue("rho",*(rho_.product()));
54 
55 
56  float cone_boundaries[] = { 0.05, 0.1, 0.2, 0.3,0.4 };
57  size_t ncone_boundaries = sizeof(cone_boundaries)/sizeof(float);
58  std::vector<float> emFractionEnergyRings(ncone_boundaries+1);
59  std::vector<float> chFractionEnergyRings(ncone_boundaries+1);
60  std::vector<float> neFractionEnergyRings(ncone_boundaries+1);
61  std::vector<float> muFractionEnergyRings(ncone_boundaries+1);
62  float jetRawEnergy=j.p4().E()*j.jecFactor("Uncorrected");
63  int numDaughtersPt03=0;
64  for (unsigned int ijcone = 0; ijcone<ncone_boundaries; ijcone++){
65  emFractionEnergyRings[ijcone] = 0;
66  muFractionEnergyRings[ijcone] = 0;
67  chFractionEnergyRings[ijcone] = 0;
68  neFractionEnergyRings[ijcone] = 0;
69  }
70  for(const auto & d : j.daughterPtrVector()){
71  float candDr = Geom::deltaR(d->p4(),j.p4());
72  size_t icone = std::lower_bound(&cone_boundaries[0],&cone_boundaries[ncone_boundaries],candDr) - &cone_boundaries[0];
73  float candEnergy = d->energy()/jetRawEnergy;
74  int pdgid = abs(d->pdgId()) ;
75  if( pdgid == 22 || pdgid == 11 ) {
76  emFractionEnergyRings[icone] += candEnergy;
77  } else if ( pdgid == 13 ) {
78  muFractionEnergyRings[icone] += candEnergy;
79  } else if ( d->charge() != 0 ) {
80  chFractionEnergyRings[icone] += candEnergy;
81  } else {
82  neFractionEnergyRings[icone] += candEnergy;
83  }
84  if(d->pt()>0.3) numDaughtersPt03+=1;
85  } // end of jet daughters loop
86 
87  this->setValue("Jet_energyRing_dR0_em_Jet_rawEnergy", emFractionEnergyRings[0]);
88  this->setValue("Jet_energyRing_dR1_em_Jet_rawEnergy", emFractionEnergyRings[1]);
89  this->setValue("Jet_energyRing_dR2_em_Jet_rawEnergy", emFractionEnergyRings[2]);
90  this->setValue("Jet_energyRing_dR3_em_Jet_rawEnergy", emFractionEnergyRings[3]);
91  this->setValue("Jet_energyRing_dR4_em_Jet_rawEnergy", emFractionEnergyRings[4]);
92  // this->setValue("Jet_energyRing_dR5_em_Jet_rawEnergy", emFractionEnergyRings[5]);
93 
94  this->setValue("Jet_energyRing_dR0_ch_Jet_rawEnergy", chFractionEnergyRings[0]);
95  this->setValue("Jet_energyRing_dR1_ch_Jet_rawEnergy", chFractionEnergyRings[1]);
96  this->setValue("Jet_energyRing_dR2_ch_Jet_rawEnergy", chFractionEnergyRings[2]);
97  this->setValue("Jet_energyRing_dR3_ch_Jet_rawEnergy", chFractionEnergyRings[3]);
98  this->setValue("Jet_energyRing_dR4_ch_Jet_rawEnergy", chFractionEnergyRings[4]);
99  // this->setValue("Jet_energyRing_dR5_ch_Jet_rawEnergy", chFractionEnergyRings[5]);
100 
101  this->setValue("Jet_energyRing_dR0_mu_Jet_rawEnergy", muFractionEnergyRings[0]);
102  this->setValue("Jet_energyRing_dR1_mu_Jet_rawEnergy", muFractionEnergyRings[1]);
103  this->setValue("Jet_energyRing_dR2_mu_Jet_rawEnergy", muFractionEnergyRings[2]);
104  this->setValue("Jet_energyRing_dR3_mu_Jet_rawEnergy", muFractionEnergyRings[3]);
105  this->setValue("Jet_energyRing_dR4_mu_Jet_rawEnergy", muFractionEnergyRings[4]);
106  // this->setValue("Jet_energyRing_dR5_mu_Jet_rawEnergy", muFractionEnergyRings[5]);
107 
108  this->setValue("Jet_energyRing_dR0_neut_Jet_rawEnergy", neFractionEnergyRings[0]);
109  this->setValue("Jet_energyRing_dR1_neut_Jet_rawEnergy", neFractionEnergyRings[1]);
110  this->setValue("Jet_energyRing_dR2_neut_Jet_rawEnergy", neFractionEnergyRings[2]);
111  this->setValue("Jet_energyRing_dR3_neut_Jet_rawEnergy", neFractionEnergyRings[3]);
112  this->setValue("Jet_energyRing_dR4_neut_Jet_rawEnergy", neFractionEnergyRings[4]);
113 // this->setValue("Jet_energyRing_dR5_neut_Jet_rawEnergy", neFractionEnergyRings[5]);
114 
115  this->setValue("Jet_numDaughters_pt03",numDaughtersPt03);
116 
117  }
void setValue(const std::string var, float val)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
edm::Handle< std::vector< reco::Vertex > > pvs_
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
T const * product() const
Definition: Handle.h:81
float jecFactor(const std::string &level, const std::string &flavor="none", const std::string &set="") const
const reco::CompositePtrCandidate::daughters & daughterPtrVector() const override
references to daughtes
static void BJetEnergyRegressionMVA::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
inlinestatic

Definition at line 119 of file BJetEnergyRegressionMVA.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), and BaseMVAValueMapProducer< T >::getDescription().

119  {
121  desc.add<edm::InputTag>("pvsrc")->setComment("primary vertices input collection");
122  desc.add<edm::InputTag>("svsrc")->setComment("secondary vertices input collection");
123  desc.add<edm::InputTag>("rhosrc")->setComment("rho input collection");
124  descriptions.add("BJetEnergyRegressionMVA",desc);
125  }
static edm::ParameterSetDescription getDescription()
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void BJetEnergyRegressionMVA::readAdditionalCollections ( edm::Event ,
const edm::EventSetup  
)
inlineoverridevirtual

to be implemented in derived classes, filling values for additional variables

Reimplemented from BaseMVAValueMapProducer< pat::Jet >.

Definition at line 43 of file BJetEnergyRegressionMVA.cc.

References edm::Event::getByToken(), pvs_, pvsrc_, rho_, rhosrc_, svs_, and svsrc_.

43  {
44  iEvent.getByToken(pvsrc_, pvs_);
45  iEvent.getByToken(svsrc_, svs_);
46  iEvent.getByToken(rhosrc_,rho_);
47 
48  }
int iEvent
Definition: GenABIO.cc:230
const edm::EDGetTokenT< edm::View< reco::VertexCompositePtrCandidate > > svsrc_
edm::Handle< std::vector< reco::Vertex > > pvs_
edm::Handle< edm::View< reco::VertexCompositePtrCandidate > > svs_
edm::EDGetTokenT< double > rhosrc_
const edm::EDGetTokenT< std::vector< reco::Vertex > > pvsrc_

Member Data Documentation

edm::Handle<std::vector<reco::Vertex> > BJetEnergyRegressionMVA::pvs_
private
const edm::EDGetTokenT<std::vector<reco::Vertex> > BJetEnergyRegressionMVA::pvsrc_
private

Definition at line 128 of file BJetEnergyRegressionMVA.cc.

Referenced by readAdditionalCollections().

edm::Handle<double> BJetEnergyRegressionMVA::rho_
private
edm::EDGetTokenT<double> BJetEnergyRegressionMVA::rhosrc_
private

Definition at line 132 of file BJetEnergyRegressionMVA.cc.

Referenced by readAdditionalCollections().

edm::Handle<edm::View<reco::VertexCompositePtrCandidate> > BJetEnergyRegressionMVA::svs_
private

Definition at line 131 of file BJetEnergyRegressionMVA.cc.

Referenced by readAdditionalCollections().

const edm::EDGetTokenT<edm::View<reco::VertexCompositePtrCandidate> > BJetEnergyRegressionMVA::svsrc_
private

Definition at line 130 of file BJetEnergyRegressionMVA.cc.

Referenced by readAdditionalCollections().