CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
BJetEnergyRegressionVarProducer< T > Class Template Reference

#include <PhysicsTools/NanoAOD/plugins/BJetEnergyRegressionVarProducer.cc>

Inheritance diagram for BJetEnergyRegressionVarProducer< T >:
edm::global::EDProducer<> edm::global::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 BJetEnergyRegressionVarProducer (const edm::ParameterSet &iConfig)
 
 ~BJetEnergyRegressionVarProducer () override
 
- Public Member Functions inherited from edm::global::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
EDProduceroperator= (const EDProducer &)=delete
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
bool wantsInputProcessBlocks () const final
 
bool wantsProcessBlocks () const final
 
bool wantsStreamLuminosityBlocks () const final
 
bool wantsStreamRuns () const final
 
- Public Member Functions inherited from edm::global::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () const
 
 ~EDProducerBase () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
std::vector< bool > const & recordProvenanceList () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
TypeLabelList const & typeLabelList () const
 used by the fwk to register the list of products of this module More...
 
 ~ProducerBase () noexcept(false) override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
std::vector< ESProxyIndex > const & esGetTokenIndicesVector (edm::Transition iTrans) const
 
std::vector< ESRecordIndex > const & esGetTokenRecordIndicesVector (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::array< std::vector< ModuleDescription const *> *, NumBranchTypes > &modulesAll, std::vector< ModuleProcessName > &modulesInPreviousProcesses, ProductRegistry const &preg, std::map< std::string, ModuleDescription const *> const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void selectInputProcessBlocks (ProductRegistry const &productRegistry, ProcessBlockHelperBase const &processBlockHelperBase)
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 
- Static Public Member Functions inherited from edm::global::EDProducerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 

Private Member Functions

std::tuple< float, float, float > calculatePtRatioRelSimple (edm::Ptr< reco::Candidate > lep, edm::Ptr< pat::Jet > jet) const
 
void produce (edm::StreamID, edm::Event &, const edm::EventSetup &) const override
 

Private Attributes

edm::EDGetTokenT< edm::View< pat::Jet > > srcJet_
 
edm::EDGetTokenT< edm::View< reco::VertexCompositePtrCandidate > > srcSV_
 
edm::EDGetTokenT< std::vector< reco::Vertex > > srcVtx_
 

Additional Inherited Members

- Public Types inherited from edm::global::EDProducerBase
typedef EDProducerBase ModuleType
 
- Public Types inherited from edm::ProducerBase
template<typename T >
using BranchAliasSetterT = ProductRegistryHelper::BranchAliasSetterT< T >
 
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex > >
 
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Protected Member Functions inherited from edm::ProducerBase
template<Transition Tr = Transition::Event>
auto produces (std::string instanceName) noexcept
 declare what type of product will make and with which optional label More...
 
template<Transition B>
BranchAliasSetter produces (const TypeID &id, std::string instanceName=std::string(), bool recordProvenance=true)
 
template<BranchType B>
BranchAliasSetter produces (const TypeID &id, std::string instanceName=std::string(), bool recordProvenance=true)
 
BranchAliasSetter produces (const TypeID &id, std::string instanceName=std::string(), bool recordProvenance=true)
 
template<typename ProductType , Transition B>
BranchAliasSetterT< ProductType > produces (std::string instanceName)
 
template<class ProductType >
BranchAliasSetterT< ProductType > produces ()
 
template<typename ProductType , BranchType B>
BranchAliasSetterT< ProductType > produces (std::string instanceName)
 
template<typename ProductType , BranchType B>
BranchAliasSetterT< ProductType > produces ()
 
template<class ProductType >
BranchAliasSetterT< ProductType > produces (std::string instanceName)
 
template<typename ProductType , Transition B>
BranchAliasSetterT< ProductType > produces ()
 
template<Transition Tr = Transition::Event>
auto produces () noexcept
 
ProducesCollector producesCollector ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
template<BranchType B = InEvent>
EDConsumerBaseAdaptor< Bconsumes (edm::InputTag tag) noexcept
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<Transition Tr = Transition::Event>
constexpr auto esConsumes ()
 
template<Transition Tr = Transition::Event>
auto esConsumes (ESInputTag tag)
 
template<Transition Tr = Transition::Event>
ESGetTokenGeneric esConsumes (eventsetup::EventSetupRecordKey const &iRecord, eventsetup::DataKey const &iKey)
 Used with EventSetupRecord::doGet. More...
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
void resetItemsToGetFrom (BranchType iType)
 

Detailed Description

template<typename T>
class BJetEnergyRegressionVarProducer< T >

Description: [one line class summary]

Implementation: [Notes on implementation]

Definition at line 48 of file BJetEnergyRegressionVarProducer.cc.

Constructor & Destructor Documentation

◆ BJetEnergyRegressionVarProducer()

template<typename T >
BJetEnergyRegressionVarProducer< T >::BJetEnergyRegressionVarProducer ( const edm::ParameterSet iConfig)
inlineexplicit

Definition at line 50 of file BJetEnergyRegressionVarProducer.cc.

52  srcVtx_(consumes<std::vector<reco::Vertex>>(iConfig.getParameter<edm::InputTag>("pvsrc"))),
54  //un prodotto da copiare
55  produces<edm::ValueMap<float>>("leptonPtRelv0");
56  produces<edm::ValueMap<float>>("leptonPtRelInvv0"); //v0 ~ heppy?
57  produces<edm::ValueMap<int>>("leptonPdgId");
58  produces<edm::ValueMap<float>>("leptonDeltaR");
59  produces<edm::ValueMap<float>>("leadTrackPt");
60  produces<edm::ValueMap<float>>("vtxPt");
61  produces<edm::ValueMap<float>>("vtxMass");
62  produces<edm::ValueMap<float>>("vtx3dL");
63  produces<edm::ValueMap<float>>("vtx3deL");
64  produces<edm::ValueMap<int>>("vtxNtrk");
65  produces<edm::ValueMap<float>>("ptD");
66  }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDGetTokenT< edm::View< reco::VertexCompositePtrCandidate > > srcSV_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< edm::View< pat::Jet > > srcJet_
edm::EDGetTokenT< std::vector< reco::Vertex > > srcVtx_

◆ ~BJetEnergyRegressionVarProducer()

template<typename T >
BJetEnergyRegressionVarProducer< T >::~BJetEnergyRegressionVarProducer ( )
inlineoverride

Definition at line 67 of file BJetEnergyRegressionVarProducer.cc.

67 {};

Member Function Documentation

◆ calculatePtRatioRelSimple()

template<typename T >
std::tuple< float, float, float > BJetEnergyRegressionVarProducer< T >::calculatePtRatioRelSimple ( edm::Ptr< reco::Candidate lep,
edm::Ptr< pat::Jet jet 
) const
private

Definition at line 248 of file BJetEnergyRegressionVarProducer.cc.

References MillePedeFileConverter_cfg::e, metsig::jet, reco::Candidate::p4(), and dttmaxenums::R.

249  {
250  auto lepp4 = lep->p4();
251  auto rawp4 = jet->correctedP4("Uncorrected");
252 
253  if ((rawp4 - lepp4).R() < 1e-4)
254  return std::tuple<float, float, float>(1.0, 0.0, 0.0);
255 
256  auto jetp4 = jet->p4(); //(rawp4 - lepp4*(1.0/jet->jecFactor("L1FastJet")))*(jet->pt()/rawp4.pt())+lepp4;
257  auto ptratio = lepp4.pt() / jetp4.pt();
258  auto ptrel = lepp4.Vect().Cross((jetp4).Vect().Unit()).R();
259  auto ptrelinv = jetp4.Vect().Cross((lepp4).Vect().Unit()).R();
260 
261  return std::tuple<float, float, float>(ptratio, ptrel, ptrelinv);
262 }
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector

◆ fillDescriptions()

template<typename T >
void BJetEnergyRegressionVarProducer< T >::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 266 of file BJetEnergyRegressionVarProducer.cc.

References edm::ConfigurationDescriptions::add(), submitPVResolutionJobs::desc, timingPdfMaker::modname, and AlCaHLTBitMon_QueryRunRegistry::string.

266  {
267  //The following says we do not know what parameters are allowed so do no validation
268  // Please change this to state exactly what you do use, even if it is no parameters
270  desc.add<edm::InputTag>("src")->setComment("jet input collection");
271  desc.add<edm::InputTag>("pvsrc")->setComment("primary vertex input collection");
272  desc.add<edm::InputTag>("svsrc")->setComment("secondary vertex input collection");
274  if (typeid(T) == typeid(pat::Jet))
275  modname += "Jet";
276  modname += "RegressionVarProducer";
277  descriptions.add(modname, desc);
278 }
Analysis-level calorimeter jet class.
Definition: Jet.h:77
void add(std::string const &label, ParameterSetDescription const &psetDescription)
long double T

◆ produce()

template<typename T >
void BJetEnergyRegressionVarProducer< T >::produce ( edm::StreamID  streamID,
edm::Event iEvent,
const edm::EventSetup iSetup 
) const
overrideprivatevirtual

Implements edm::global::EDProducerBase.

Definition at line 98 of file BJetEnergyRegressionVarProducer.cc.

References funct::abs(), RecoVertex::convertError(), RecoVertex::convertPos(), ztail::d, reco::deltaR(), reco::deltaR2(), VertexDistance3D::distance(), Measurement1D::error(), edm::helper::Filler< Map >::fill(), iEvent, edm::helper::Filler< Map >::insert(), metsig::jet, jetsAK4_CHS_cff::leadTrackPt, reco::btau::leptonDeltaR, jetsAK4_CHS_cff::leptonPdgId, eostools::move(), jetsAK4_CHS_cff::ptD, AlignmentTrackSelector_cfi::ptMax, Measurement1D::significance(), mathSSE::sqrt(), electrons_cff::srcJet, TtFullHadEvtBuilder_cfi::sumPt, pfDeepBoostedJetPreprocessParams_cfi::sv, Measurement1D::value(), jetsAK4_CHS_cff::vtx3deL, jetsAK4_CHS_cff::vtx3dL, jetsAK4_CHS_cff::vtxMass, jetsAK4_CHS_cff::vtxNtrk, and jetsAK4_CHS_cff::vtxPt.

100  {
101  const auto& srcJet = iEvent.getHandle(srcJet_);
102  const auto& vtxProd = iEvent.get(srcVtx_);
103  const auto& svProd = iEvent.get(srcSV_);
104 
105  auto nJet = srcJet->size();
106 
107  std::vector<float> leptonPtRel_v0(nJet, 0);
108  std::vector<float> leptonPtRelInv_v0(nJet, 0);
109  std::vector<int> leptonPdgId(nJet, 0);
110  std::vector<float> leptonDeltaR(nJet, 0);
111  std::vector<float> leadTrackPt(nJet, 0);
112  std::vector<float> vtxPt(nJet, 0);
113  std::vector<float> vtxMass(nJet, 0);
114  std::vector<float> vtx3dL(nJet, 0);
115  std::vector<float> vtx3deL(nJet, 0);
116  std::vector<int> vtxNtrk(nJet, 0);
117  std::vector<float> ptD(nJet, 0);
118 
119  const auto& pv = vtxProd.at(0);
120  for (unsigned int ij = 0; ij < nJet; ij++) {
121  auto jet = srcJet->ptrAt(ij);
122 
123  float ptMax = 0;
124  float sumWeight = 0;
125  float sumPt = 0;
126 
127  for (const auto& d : jet->daughterPtrVector()) {
128  sumWeight += (d->pt()) * (d->pt());
129  sumPt += d->pt();
130  if (d->pt() > ptMax)
131  ptMax = d->pt();
132  }
133  leadTrackPt[ij] = ptMax;
134  ptD[ij] = (sumWeight > 0 ? sqrt(sumWeight) / sumPt : 0);
135 
136  //lepton properties
137  float maxLepPt = 0;
138 
139  for (const auto& d : jet->daughterPtrVector()) {
140  if (abs(d->pdgId()) == 11 || abs(d->pdgId()) == 13) {
141  if (d->pt() < maxLepPt)
142  continue;
143  auto res2 = calculatePtRatioRelSimple(d, jet);
144  leptonPtRel_v0[ij] = std::get<1>(res2);
145  leptonPtRelInv_v0[ij] = std::get<2>(res2);
146  leptonPdgId[ij] = d->pdgId();
147  leptonDeltaR[ij] = reco::deltaR(jet->p4(), d->p4());
148  maxLepPt = d->pt();
149  }
150  }
151 
152  //Fill vertex properties
153  VertexDistance3D vdist;
154  float maxFoundSignificance = 0;
155 
156  vtxPt[ij] = 0;
157  vtxMass[ij] = 0;
158  vtx3dL[ij] = 0;
159  vtx3deL[ij] = 0;
160  vtxNtrk[ij] = 0;
161 
162  for (const auto& sv : svProd) {
163  GlobalVector flightDir(sv.vertex().x() - pv.x(), sv.vertex().y() - pv.y(), sv.vertex().z() - pv.z());
164  GlobalVector jetDir(jet->px(), jet->py(), jet->pz());
165  if (reco::deltaR2(flightDir, jetDir) < 0.09) {
166  Measurement1D dl = vdist.distance(
168  if (dl.significance() > maxFoundSignificance) {
169  maxFoundSignificance = dl.significance();
170  vtxPt[ij] = sv.pt();
171  vtxMass[ij] = sv.p4().M();
172  vtx3dL[ij] = dl.value();
173  vtx3deL[ij] = dl.error();
174  vtxNtrk[ij] = sv.numberOfSourceCandidatePtrs();
175  }
176  }
177  }
178  }
179 
180  std::unique_ptr<edm::ValueMap<float>> leptonPtRelV_v0(new edm::ValueMap<float>());
181  edm::ValueMap<float>::Filler fillerRel_v0(*leptonPtRelV_v0);
182  fillerRel_v0.insert(srcJet, leptonPtRel_v0.begin(), leptonPtRel_v0.end());
183  fillerRel_v0.fill();
184  iEvent.put(std::move(leptonPtRelV_v0), "leptonPtRelv0");
185 
186  std::unique_ptr<edm::ValueMap<float>> leptonPtRelInvV_v0(new edm::ValueMap<float>());
187  edm::ValueMap<float>::Filler fillerRelInv_v0(*leptonPtRelInvV_v0);
188  fillerRelInv_v0.insert(srcJet, leptonPtRelInv_v0.begin(), leptonPtRelInv_v0.end());
189  fillerRelInv_v0.fill();
190  iEvent.put(std::move(leptonPtRelInvV_v0), "leptonPtRelInvv0");
191 
192  std::unique_ptr<edm::ValueMap<float>> leptonDeltaRV(new edm::ValueMap<float>());
193  edm::ValueMap<float>::Filler fillerLdR(*leptonDeltaRV);
194  fillerLdR.insert(srcJet, leptonDeltaR.begin(), leptonDeltaR.end());
195  fillerLdR.fill();
196  iEvent.put(std::move(leptonDeltaRV), "leptonDeltaR");
197 
198  std::unique_ptr<edm::ValueMap<int>> leptonPtPdgIdV(new edm::ValueMap<int>());
199  edm::ValueMap<int>::Filler fillerId(*leptonPtPdgIdV);
200  fillerId.insert(srcJet, leptonPdgId.begin(), leptonPdgId.end());
201  fillerId.fill();
202  iEvent.put(std::move(leptonPtPdgIdV), "leptonPdgId");
203 
204  std::unique_ptr<edm::ValueMap<float>> leadTrackPtV(new edm::ValueMap<float>());
205  edm::ValueMap<float>::Filler fillerLT(*leadTrackPtV);
206  fillerLT.insert(srcJet, leadTrackPt.begin(), leadTrackPt.end());
207  fillerLT.fill();
208  iEvent.put(std::move(leadTrackPtV), "leadTrackPt");
209 
210  std::unique_ptr<edm::ValueMap<float>> vtxPtV(new edm::ValueMap<float>());
211  edm::ValueMap<float>::Filler fillerVtxPt(*vtxPtV);
212  fillerVtxPt.insert(srcJet, vtxPt.begin(), vtxPt.end());
213  fillerVtxPt.fill();
214  iEvent.put(std::move(vtxPtV), "vtxPt");
215 
216  std::unique_ptr<edm::ValueMap<float>> vtxMassV(new edm::ValueMap<float>());
217  edm::ValueMap<float>::Filler fillerVtxMass(*vtxMassV);
218  fillerVtxMass.insert(srcJet, vtxMass.begin(), vtxMass.end());
219  fillerVtxMass.fill();
220  iEvent.put(std::move(vtxMassV), "vtxMass");
221 
222  std::unique_ptr<edm::ValueMap<float>> vtx3dLV(new edm::ValueMap<float>());
223  edm::ValueMap<float>::Filler fillerVtx3dL(*vtx3dLV);
224  fillerVtx3dL.insert(srcJet, vtx3dL.begin(), vtx3dL.end());
225  fillerVtx3dL.fill();
226  iEvent.put(std::move(vtx3dLV), "vtx3dL");
227 
228  std::unique_ptr<edm::ValueMap<float>> vtx3deLV(new edm::ValueMap<float>());
229  edm::ValueMap<float>::Filler fillerVtx3deL(*vtx3deLV);
230  fillerVtx3deL.insert(srcJet, vtx3deL.begin(), vtx3deL.end());
231  fillerVtx3deL.fill();
232  iEvent.put(std::move(vtx3deLV), "vtx3deL");
233 
234  std::unique_ptr<edm::ValueMap<int>> vtxNtrkV(new edm::ValueMap<int>());
235  edm::ValueMap<int>::Filler fillerVtxNT(*vtxNtrkV);
236  fillerVtxNT.insert(srcJet, vtxNtrk.begin(), vtxNtrk.end());
237  fillerVtxNT.fill();
238  iEvent.put(std::move(vtxNtrkV), "vtxNtrk");
239 
240  std::unique_ptr<edm::ValueMap<float>> ptDV(new edm::ValueMap<float>());
241  edm::ValueMap<float>::Filler fillerPtD(*ptDV);
242  fillerPtD.insert(srcJet, ptD.begin(), ptD.end());
243  fillerPtD.fill();
244  iEvent.put(std::move(ptDV), "ptD");
245 }
reco::Vertex::Point convertPos(const GlobalPoint &p)
edm::EDGetTokenT< edm::View< reco::VertexCompositePtrCandidate > > srcSV_
reco::Vertex::Error convertError(const GlobalError &ge)
Definition: ConvertError.h:8
edm::EDGetTokenT< edm::View< pat::Jet > > srcJet_
int iEvent
Definition: GenABIO.cc:224
T sqrt(T t)
Definition: SSEVec.h:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
d
Definition: ztail.py:151
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
edm::EDGetTokenT< std::vector< reco::Vertex > > srcVtx_
double value() const
Definition: Measurement1D.h:25
double significance() const
Definition: Measurement1D.h:29
std::tuple< float, float, float > calculatePtRatioRelSimple(edm::Ptr< reco::Candidate > lep, edm::Ptr< pat::Jet > jet) const
double error() const
Definition: Measurement1D.h:27
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

◆ srcJet_

template<typename T >
edm::EDGetTokenT<edm::View<pat::Jet> > BJetEnergyRegressionVarProducer< T >::srcJet_
private

Definition at line 79 of file BJetEnergyRegressionVarProducer.cc.

◆ srcSV_

Definition at line 81 of file BJetEnergyRegressionVarProducer.cc.

◆ srcVtx_

template<typename T >
edm::EDGetTokenT<std::vector<reco::Vertex> > BJetEnergyRegressionVarProducer< T >::srcVtx_
private

Definition at line 80 of file BJetEnergyRegressionVarProducer.cc.