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< 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"))),
37 
38  {
39 
40 
41  }
T getParameter(std::string const &) const
const edm::EDGetTokenT< edm::View< reco::VertexCompositePtrCandidate > > svsrc_
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 47 of file BJetEnergyRegressionMVA.cc.

References RecoVertex::convertError(), RecoVertex::convertPos(), edmIntegrityCheck::d, pat::Jet::daughterPtrVector(), reco::deltaR2(), VertexDistance3D::distance(), edm::PtrVectorBase::empty(), Measurement1D::error(), pat::PATObject< ObjectType >::overlaps(), AlignmentTrackSelector_cfi::ptMax, MetAnalyzer::pv(), pvs_, reco::LeafCandidate::px(), reco::LeafCandidate::py(), reco::LeafCandidate::pz(), BaseMVAValueMapProducer< pat::Jet >::setValue(), BaseMVAValueMapProducer< T >::setValue(), Measurement1D::significance(), svs_, and Measurement1D::value().

47  {
48  this->setValue("nPVs",pvs_->size());
50 
51  if(!j.overlaps("muons").empty()) {
52  const auto *lep=dynamic_cast<const pat::Muon *>(&*j.overlaps("muons")[0]);
53  if(lep!=nullptr) {BaseMVAValueMapProducer<pat::Jet>::setValue("Jet_leptonPtRel",lep->userFloat("ptRel"));}
54  }
55  else if(!j.overlaps("electrons").empty()) {
56  const auto *lep=dynamic_cast<const pat::Electron *>(&*j.overlaps("electrons")[0]);
57  if(lep!=nullptr) {BaseMVAValueMapProducer<pat::Jet>::setValue("Jet_leptonPtRel",lep->userFloat("ptRel"));}
58  }
59 
60  float ptMax=0;
61  for(const auto & d : j.daughterPtrVector()){if(d->pt()>ptMax) ptMax=d->pt();}
62  BaseMVAValueMapProducer<pat::Jet>::setValue("Jet_leadTrackPt",ptMax);
63 
64  //Fill vertex properties
65  VertexDistance3D vdist;
66  float maxFoundSignificance=0;
67  const auto & pv = (*pvs_)[0];
68  this->setValue("Jet_vtxPt",0);
69  this->setValue("Jet_vtxMass",0);
70  this->setValue("Jet_vtx3dL",0);
71  this->setValue("Jet_vtx3deL",0);
72  this->setValue("Jet_vtxNtrk",0);
73 
74  for(const auto &sv: *svs_){
75  GlobalVector flightDir(sv.vertex().x() - pv.x(), sv.vertex().y() - pv.y(),sv.vertex().z() - pv.z());
76  GlobalVector jetDir(j.px(),j.py(),j.pz());
77  if( Geom::deltaR2( flightDir, jetDir ) < 0.09 ){
79  if(dl.significance() > maxFoundSignificance){
80  maxFoundSignificance=dl.significance();
81  this->setValue("Jet_vtxPt",sv.pt());
82  this->setValue("Jet_vtxMass",sv.p4().M());
83  this->setValue("Jet_vtx3dL",dl.value());
84  this->setValue("Jet_vtx3deL",dl.error());
85  this->setValue("Jet_vtxNtrk",sv.numberOfSourceCandidatePtrs());
86  }
87  }
88  }
89 
90  }
reco::Vertex::Point convertPos(const GlobalPoint &p)
double px() const final
x coordinate of momentum vector
bool empty() const
Is the RefVector empty.
Definition: PtrVectorBase.h:71
reco::Vertex::Error convertError(const GlobalError &ge)
Definition: ConvertError.h:8
double error() const
Definition: Measurement1D.h:30
void setValue(const std::string var, float val)
const reco::CandidatePtrVector & overlaps(const std::string &label) const
Definition: PATObject.h:742
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
double pz() const final
z coordinate of momentum vector
def pv(vc)
Definition: MetAnalyzer.py:6
edm::Handle< std::vector< reco::Vertex > > pvs_
edm::Handle< edm::View< reco::VertexCompositePtrCandidate > > svs_
double significance() const
Definition: Measurement1D.h:32
Analysis-level electron class.
Definition: Electron.h:52
double py() const final
y coordinate of momentum vector
double value() const
Definition: Measurement1D.h:28
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
const reco::CompositePtrCandidate::daughters & daughterPtrVector() const override
references to daughtes
Analysis-level muon class.
Definition: Muon.h:50
static void BJetEnergyRegressionMVA::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
inlinestatic

Definition at line 92 of file BJetEnergyRegressionMVA.cc.

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

92  {
94  desc.add<edm::InputTag>("pvsrc")->setComment("primary vertices input collection");
95  desc.add<edm::InputTag>("svsrc")->setComment("secondary vertices input collection");
96  descriptions.add("BJetEnergyRegressionMVA",desc);
97  }
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 42 of file BJetEnergyRegressionMVA.cc.

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

42  {
43  iEvent.getByToken(pvsrc_, pvs_);
44  iEvent.getByToken(svsrc_, svs_);
45  }
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_
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 100 of file BJetEnergyRegressionMVA.cc.

Referenced by readAdditionalCollections().

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

Definition at line 102 of file BJetEnergyRegressionMVA.cc.

Referenced by readAdditionalCollections().