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 > calculatePtRatioRel (edm::Ptr< reco::Candidate > lep, edm::Ptr< pat::Jet > jet) const
 
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< std::vector< reco::GenParticle > > srcGP_
 
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
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 73 of file BJetEnergyRegressionVarProducer.cc.

Constructor & Destructor Documentation

◆ BJetEnergyRegressionVarProducer()

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

Definition at line 75 of file BJetEnergyRegressionVarProducer.cc.

77  srcVtx_(consumes<std::vector<reco::Vertex>>(iConfig.getParameter<edm::InputTag>("pvsrc"))),
79  srcGP_(consumes<std::vector<reco::GenParticle>>(iConfig.getParameter<edm::InputTag>("gpsrc"))) {
80  //un prodotto da copiare
81  produces<edm::ValueMap<float>>("leptonPtRel");
82  produces<edm::ValueMap<float>>("leptonPtRatio");
83  produces<edm::ValueMap<float>>("leptonPtRelInv"); //wrong variable?
84  produces<edm::ValueMap<float>>("leptonPtRelv0");
85  produces<edm::ValueMap<float>>("leptonPtRatiov0");
86  produces<edm::ValueMap<float>>("leptonPtRelInvv0"); //v0 ~ heppy?
87  produces<edm::ValueMap<float>>("leptonPt");
88  produces<edm::ValueMap<int>>("leptonPdgId");
89  produces<edm::ValueMap<float>>("leptonDeltaR");
90  produces<edm::ValueMap<float>>("leadTrackPt");
91  produces<edm::ValueMap<float>>("vtxPt");
92  produces<edm::ValueMap<float>>("vtxMass");
93  produces<edm::ValueMap<float>>("vtx3dL");
94  produces<edm::ValueMap<float>>("vtx3deL");
95  produces<edm::ValueMap<int>>("vtxNtrk");
96  produces<edm::ValueMap<float>>("ptD");
97  produces<edm::ValueMap<float>>("genPtwNu");
98  }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::EDGetTokenT< edm::View< reco::VertexCompositePtrCandidate > > srcSV_
edm::EDGetTokenT< std::vector< reco::GenParticle > > srcGP_
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 99 of file BJetEnergyRegressionVarProducer.cc.

99 {};

Member Function Documentation

◆ calculatePtRatioRel()

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

Definition at line 353 of file BJetEnergyRegressionVarProducer.cc.

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

354  {
355  auto rawp4 = jet->correctedP4("Uncorrected");
356  auto lepp4 = lep->p4();
357 
358  if ((rawp4 - lepp4).R() < 1e-4)
359  return std::tuple<float, float, float>(1.0, 0.0, 0.0);
360 
361  auto jetp4 = (rawp4 - lepp4 * (1.0 / jet->jecFactor("L1FastJet"))) * (jet->pt() / rawp4.pt()) + lepp4;
362  auto ptratio = lepp4.pt() / jetp4.pt();
363  auto ptrel = lepp4.Vect().Cross((jetp4 - lepp4).Vect().Unit()).R();
364  auto ptrelinv = (jetp4 - lepp4).Vect().Cross((lepp4).Vect().Unit()).R();
365 
366  return std::tuple<float, float, float>(ptratio, ptrel, ptrelinv);
367 }
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector

◆ 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 370 of file BJetEnergyRegressionVarProducer.cc.

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

371  {
372  auto lepp4 = lep->p4();
373  auto rawp4 = jet->correctedP4("Uncorrected");
374 
375  if ((rawp4 - lepp4).R() < 1e-4)
376  return std::tuple<float, float, float>(1.0, 0.0, 0.0);
377 
378  auto jetp4 = jet->p4(); //(rawp4 - lepp4*(1.0/jet->jecFactor("L1FastJet")))*(jet->pt()/rawp4.pt())+lepp4;
379  auto ptratio = lepp4.pt() / jetp4.pt();
380  auto ptrel = lepp4.Vect().Cross((jetp4).Vect().Unit()).R();
381  auto ptrelinv = jetp4.Vect().Cross((lepp4).Vect().Unit()).R();
382 
383  return std::tuple<float, float, float>(ptratio, ptrel, ptrelinv);
384 }
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 388 of file BJetEnergyRegressionVarProducer.cc.

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

388  {
389  //The following says we do not know what parameters are allowed so do no validation
390  // Please change this to state exactly what you do use, even if it is no parameters
392  desc.add<edm::InputTag>("src")->setComment("jet input collection");
393  // desc.add<edm::InputTag>("musrc")->setComment("muons input collection");
394  // desc.add<edm::InputTag>("elesrc")->setComment("electrons input collection");
395  desc.add<edm::InputTag>("pvsrc")->setComment("primary vertex input collection");
396  desc.add<edm::InputTag>("svsrc")->setComment("secondary vertex input collection");
397  desc.add<edm::InputTag>("gpsrc")->setComment("genparticles for nu recovery");
399  if (typeid(T) == typeid(pat::Jet))
400  modname += "Jet";
401  modname += "RegressionVarProducer";
402  descriptions.add(modname, desc);
403 }
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 132 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(), jetsAK4_CHS_cff::genPtwNu, runTauDisplay::gp, iEvent, edm::helper::Filler< Map >::insert(), metsig::jet, jetsAK4_CHS_cff::leadTrackPt, reco::btau::leptonDeltaR, jetsAK4_CHS_cff::leptonPdgId, jetsAK4_CHS_cff::leptonPt, jetsAK4_CHS_cff::leptonPtRatio, reco::btau::leptonPtRel, jetsAK4_CHS_cff::leptonPtRelInv, eostools::move(), jetsAK4_CHS_cff::ptD, AlignmentTrackSelector_cfi::ptMax, MetAnalyzer::pv(), Measurement1D::significance(), mathSSE::sqrt(), electrons_cff::srcJet, jetAnalyzer_cff::srcVtx, 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.

134  {
136  iEvent.getByToken(srcJet_, srcJet);
138  iEvent.getByToken(srcVtx_, srcVtx);
140  iEvent.getByToken(srcSV_, srcSV);
142  iEvent.getByToken(srcGP_, srcGP);
143 
144  unsigned int nJet = srcJet->size();
145  // unsigned int nLep = srcLep->size();
146 
147  std::vector<float> leptonPtRel(nJet, 0);
148  std::vector<float> leptonPtRatio(nJet, 0);
149  std::vector<float> leptonPtRelInv(nJet, 0);
150  std::vector<float> leptonPtRel_v0(nJet, 0);
151  std::vector<float> leptonPtRatio_v0(nJet, 0);
152  std::vector<float> leptonPtRelInv_v0(nJet, 0);
153  std::vector<int> leptonPdgId(nJet, 0);
154  std::vector<float> leptonPt(nJet, 0);
155  std::vector<float> leptonDeltaR(nJet, 0);
156  std::vector<float> leadTrackPt(nJet, 0);
157  std::vector<float> vtxPt(nJet, 0);
158  std::vector<float> vtxMass(nJet, 0);
159  std::vector<float> vtx3dL(nJet, 0);
160  std::vector<float> vtx3deL(nJet, 0);
161  std::vector<int> vtxNtrk(nJet, 0);
162  std::vector<float> ptD(nJet, 0);
163  std::vector<float> genPtwNu(nJet, 0);
164 
165  const auto& pv = (*srcVtx)[0];
166  for (unsigned int ij = 0; ij < nJet; ij++) {
167  auto jet = srcJet->ptrAt(ij);
168 
169  if (jet->genJet() != nullptr) {
170  auto genp4 = jet->genJet()->p4();
171  auto gep4wNu = genp4;
172  for (const auto& gp : *srcGP) {
173  if ((abs(gp.pdgId()) == 12 || abs(gp.pdgId()) == 14 || abs(gp.pdgId()) == 16) && gp.status() == 1) {
174  if (reco::deltaR(genp4, gp.p4()) < 0.4) {
175  // std::cout<<" from "<<gep4wNu.pt()<<std::endl;
176  gep4wNu = gep4wNu + gp.p4();
177  // std::cout<<" to "<<gep4wNu.pt()<<std::endl;
178  }
179  }
180  }
181 
182  genPtwNu[ij] = gep4wNu.pt();
183  }
184 
185  float ptMax = 0;
186  float sumWeight = 0;
187  float sumPt = 0;
188 
189  for (const auto& d : jet->daughterPtrVector()) {
190  sumWeight += (d->pt()) * (d->pt());
191  sumPt += d->pt();
192  if (d->pt() > ptMax)
193  ptMax = d->pt();
194  }
195  leadTrackPt[ij] = ptMax;
196  ptD[ij] = (sumWeight > 0 ? sqrt(sumWeight) / sumPt : 0);
197 
198  //lepton properties
199  float maxLepPt = 0;
200  leptonPtRel[ij] = 0;
201 
202  for (const auto& d : jet->daughterPtrVector()) {
203  if (abs(d->pdgId()) == 11 || abs(d->pdgId()) == 13) {
204  if (d->pt() < maxLepPt)
205  continue;
206  auto res = calculatePtRatioRel(d, jet);
207  leptonPtRatio[ij] = std::get<0>(res);
208  leptonPtRel[ij] = std::get<1>(res);
209  leptonPtRelInv[ij] = std::get<2>(res);
210  auto res2 = calculatePtRatioRelSimple(d, jet);
211  leptonPtRatio_v0[ij] = std::get<0>(res2);
212  leptonPtRel_v0[ij] = std::get<1>(res2);
213  leptonPtRelInv_v0[ij] = std::get<2>(res2);
214  leptonPdgId[ij] = d->pdgId();
215  leptonDeltaR[ij] = reco::deltaR(jet->p4(), d->p4());
216  leptonPt[ij] = d->pt();
217  maxLepPt = d->pt();
218  }
219  }
220 
221  //Fill vertex properties
222  VertexDistance3D vdist;
223  float maxFoundSignificance = 0;
224 
225  vtxPt[ij] = 0;
226  vtxMass[ij] = 0;
227  vtx3dL[ij] = 0;
228  vtx3deL[ij] = 0;
229  vtxNtrk[ij] = 0;
230 
231  for (const auto& sv : *srcSV) {
232  GlobalVector flightDir(sv.vertex().x() - pv.x(), sv.vertex().y() - pv.y(), sv.vertex().z() - pv.z());
233  GlobalVector jetDir(jet->px(), jet->py(), jet->pz());
234  if (reco::deltaR2(flightDir, jetDir) < 0.09) {
235  Measurement1D dl = vdist.distance(
237  if (dl.significance() > maxFoundSignificance) {
238  maxFoundSignificance = dl.significance();
239  vtxPt[ij] = sv.pt();
240  vtxMass[ij] = sv.p4().M();
241  vtx3dL[ij] = dl.value();
242  vtx3deL[ij] = dl.error();
243  vtxNtrk[ij] = sv.numberOfSourceCandidatePtrs();
244  }
245  }
246  }
247  }
248 
249  std::unique_ptr<edm::ValueMap<float>> leptonPtRelV(new edm::ValueMap<float>());
250  edm::ValueMap<float>::Filler fillerRel(*leptonPtRelV);
251  fillerRel.insert(srcJet, leptonPtRel.begin(), leptonPtRel.end());
252  fillerRel.fill();
253  iEvent.put(std::move(leptonPtRelV), "leptonPtRel");
254 
255  std::unique_ptr<edm::ValueMap<float>> leptonPtRatioV(new edm::ValueMap<float>());
256  edm::ValueMap<float>::Filler fillerRatio(*leptonPtRatioV);
257  fillerRatio.insert(srcJet, leptonPtRatio.begin(), leptonPtRatio.end());
258  fillerRatio.fill();
259  iEvent.put(std::move(leptonPtRatioV), "leptonPtRatio");
260 
261  std::unique_ptr<edm::ValueMap<float>> leptonPtRelInvV(new edm::ValueMap<float>());
262  edm::ValueMap<float>::Filler fillerRelInv(*leptonPtRelInvV);
263  fillerRelInv.insert(srcJet, leptonPtRelInv.begin(), leptonPtRelInv.end());
264  fillerRelInv.fill();
265  iEvent.put(std::move(leptonPtRelInvV), "leptonPtRelInv");
266 
267  std::unique_ptr<edm::ValueMap<float>> leptonPtRelV_v0(new edm::ValueMap<float>());
268  edm::ValueMap<float>::Filler fillerRel_v0(*leptonPtRelV_v0);
269  fillerRel_v0.insert(srcJet, leptonPtRel_v0.begin(), leptonPtRel_v0.end());
270  fillerRel_v0.fill();
271  iEvent.put(std::move(leptonPtRelV_v0), "leptonPtRelv0");
272 
273  std::unique_ptr<edm::ValueMap<float>> leptonPtRatioV_v0(new edm::ValueMap<float>());
274  edm::ValueMap<float>::Filler fillerRatio_v0(*leptonPtRatioV_v0);
275  fillerRatio_v0.insert(srcJet, leptonPtRatio_v0.begin(), leptonPtRatio_v0.end());
276  fillerRatio_v0.fill();
277  iEvent.put(std::move(leptonPtRatioV_v0), "leptonPtRatiov0");
278 
279  std::unique_ptr<edm::ValueMap<float>> leptonPtRelInvV_v0(new edm::ValueMap<float>());
280  edm::ValueMap<float>::Filler fillerRelInv_v0(*leptonPtRelInvV_v0);
281  fillerRelInv_v0.insert(srcJet, leptonPtRelInv_v0.begin(), leptonPtRelInv_v0.end());
282  fillerRelInv_v0.fill();
283  iEvent.put(std::move(leptonPtRelInvV_v0), "leptonPtRelInvv0");
284 
285  std::unique_ptr<edm::ValueMap<float>> leptonPtV(new edm::ValueMap<float>());
286  edm::ValueMap<float>::Filler fillerLpt(*leptonPtV);
287  fillerLpt.insert(srcJet, leptonPt.begin(), leptonPt.end());
288  fillerLpt.fill();
289  iEvent.put(std::move(leptonPtV), "leptonPt");
290 
291  std::unique_ptr<edm::ValueMap<float>> leptonDeltaRV(new edm::ValueMap<float>());
292  edm::ValueMap<float>::Filler fillerLdR(*leptonDeltaRV);
293  fillerLdR.insert(srcJet, leptonDeltaR.begin(), leptonDeltaR.end());
294  fillerLdR.fill();
295  iEvent.put(std::move(leptonDeltaRV), "leptonDeltaR");
296 
297  std::unique_ptr<edm::ValueMap<int>> leptonPtPdgIdV(new edm::ValueMap<int>());
298  edm::ValueMap<int>::Filler fillerId(*leptonPtPdgIdV);
299  fillerId.insert(srcJet, leptonPdgId.begin(), leptonPdgId.end());
300  fillerId.fill();
301  iEvent.put(std::move(leptonPtPdgIdV), "leptonPdgId");
302 
303  std::unique_ptr<edm::ValueMap<float>> leadTrackPtV(new edm::ValueMap<float>());
304  edm::ValueMap<float>::Filler fillerLT(*leadTrackPtV);
305  fillerLT.insert(srcJet, leadTrackPt.begin(), leadTrackPt.end());
306  fillerLT.fill();
307  iEvent.put(std::move(leadTrackPtV), "leadTrackPt");
308 
309  std::unique_ptr<edm::ValueMap<float>> vtxPtV(new edm::ValueMap<float>());
310  edm::ValueMap<float>::Filler fillerVtxPt(*vtxPtV);
311  fillerVtxPt.insert(srcJet, vtxPt.begin(), vtxPt.end());
312  fillerVtxPt.fill();
313  iEvent.put(std::move(vtxPtV), "vtxPt");
314 
315  std::unique_ptr<edm::ValueMap<float>> vtxMassV(new edm::ValueMap<float>());
316  edm::ValueMap<float>::Filler fillerVtxMass(*vtxMassV);
317  fillerVtxMass.insert(srcJet, vtxMass.begin(), vtxMass.end());
318  fillerVtxMass.fill();
319  iEvent.put(std::move(vtxMassV), "vtxMass");
320 
321  std::unique_ptr<edm::ValueMap<float>> vtx3dLV(new edm::ValueMap<float>());
322  edm::ValueMap<float>::Filler fillerVtx3dL(*vtx3dLV);
323  fillerVtx3dL.insert(srcJet, vtx3dL.begin(), vtx3dL.end());
324  fillerVtx3dL.fill();
325  iEvent.put(std::move(vtx3dLV), "vtx3dL");
326 
327  std::unique_ptr<edm::ValueMap<float>> vtx3deLV(new edm::ValueMap<float>());
328  edm::ValueMap<float>::Filler fillerVtx3deL(*vtx3deLV);
329  fillerVtx3deL.insert(srcJet, vtx3deL.begin(), vtx3deL.end());
330  fillerVtx3deL.fill();
331  iEvent.put(std::move(vtx3deLV), "vtx3deL");
332 
333  std::unique_ptr<edm::ValueMap<int>> vtxNtrkV(new edm::ValueMap<int>());
334  edm::ValueMap<int>::Filler fillerVtxNT(*vtxNtrkV);
335  fillerVtxNT.insert(srcJet, vtxNtrk.begin(), vtxNtrk.end());
336  fillerVtxNT.fill();
337  iEvent.put(std::move(vtxNtrkV), "vtxNtrk");
338 
339  std::unique_ptr<edm::ValueMap<float>> ptDV(new edm::ValueMap<float>());
340  edm::ValueMap<float>::Filler fillerPtD(*ptDV);
341  fillerPtD.insert(srcJet, ptD.begin(), ptD.end());
342  fillerPtD.fill();
343  iEvent.put(std::move(ptDV), "ptD");
344 
345  std::unique_ptr<edm::ValueMap<float>> genptV(new edm::ValueMap<float>());
346  edm::ValueMap<float>::Filler fillergenpt(*genptV);
347  fillergenpt.insert(srcJet, genPtwNu.begin(), genPtwNu.end());
348  fillergenpt.fill();
349  iEvent.put(std::move(genptV), "genPtwNu");
350 }
reco::Vertex::Point convertPos(const GlobalPoint &p)
edm::EDGetTokenT< edm::View< reco::VertexCompositePtrCandidate > > srcSV_
edm::EDGetTokenT< std::vector< reco::GenParticle > > srcGP_
reco::Vertex::Error convertError(const GlobalError &ge)
Definition: ConvertError.h:8
Definition: Electron.h:6
edm::EDGetTokenT< edm::View< pat::Jet > > srcJet_
int iEvent
Definition: GenABIO.cc:224
T sqrt(T t)
Definition: SSEVec.h:19
def pv(vc)
Definition: MetAnalyzer.py:7
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_
std::tuple< float, float, float > calculatePtRatioRel(edm::Ptr< reco::Candidate > lep, edm::Ptr< pat::Jet > jet) const
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

◆ srcGP_

template<typename T >
edm::EDGetTokenT<std::vector<reco::GenParticle> > BJetEnergyRegressionVarProducer< T >::srcGP_
private

Definition at line 115 of file BJetEnergyRegressionVarProducer.cc.

◆ srcJet_

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

Definition at line 112 of file BJetEnergyRegressionVarProducer.cc.

◆ srcSV_

Definition at line 114 of file BJetEnergyRegressionVarProducer.cc.

◆ srcVtx_

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

Definition at line 113 of file BJetEnergyRegressionVarProducer.cc.