CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
cms::PFMETProducer Class Reference

#include <PFMETProducer.h>

Inheritance diagram for cms::PFMETProducer:
edm::stream::EDProducer<>

Public Member Functions

 PFMETProducer (const edm::ParameterSet &)
 
void produce (edm::Event &, const edm::EventSetup &) override
 
 ~PFMETProducer () override
 
- Public Member Functions inherited from edm::stream::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
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Member Functions

reco::METCovMatrix getMETCovMatrix (const edm::Event &event, const edm::EventSetup &, const edm::Handle< edm::View< reco::Candidate >> &input) const
 

Private Attributes

bool applyWeight_
 
bool calculateSignificance_
 
double globalThreshold_
 
edm::EDGetTokenT< edm::View< reco::Candidate > > inputToken_
 
edm::ESGetToken< JME::JetResolutionObject, JetResolutionRcdjetResPhiToken_
 
edm::ESGetToken< JME::JetResolutionObject, JetResolutionRcdjetResPtToken_
 
edm::ESGetToken< JME::JetResolutionObject, JetResolutionScaleFactorRcdjetSFToken_
 
double jetThreshold_
 
edm::EDGetTokenT< edm::View< reco::Jet > > jetToken_
 
std::vector< edm::EDGetTokenT< edm::View< reco::Candidate > > > lepTokens_
 
metsig::METSignificancemetSigAlgo_
 
edm::EDGetTokenT< double > rhoToken_
 
edm::InputTag src_
 
edm::ValueMap< float > const * weights_
 
edm::EDGetTokenT< edm::ValueMap< float > > weightsToken_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

Definition at line 57 of file PFMETProducer.h.

Constructor & Destructor Documentation

◆ PFMETProducer()

PFMETProducer::PFMETProducer ( const edm::ParameterSet iConfig)
explicit

Definition at line 17 of file PFMETProducer.cc.

References SiStripOfflineCRack_cfg::alias, applyWeight_, calculateSignificance_, deDxTools::esConsumes(), edm::ParameterSet::getParameter(), jetResPhiToken_, jetResPtToken_, jetSFToken_, jetToken_, edm::InputTag::label(), lepTokens_, metSigAlgo_, METSignificance_cfi::METSignificance, rhoToken_, src_, HLT_2022v12_cff::srcWeights, AlCaHLTBitMon_QueryRunRegistry::string, and weightsToken_.

18  : src_(iConfig.getParameter<edm::InputTag>("src")),
20  calculateSignificance_(iConfig.getParameter<bool>("calculateSignificance")),
21  globalThreshold_(iConfig.getParameter<double>("globalThreshold")),
22  applyWeight_(iConfig.getParameter<bool>("applyWeight")),
23  weights_(nullptr) {
24  if (applyWeight_) {
25  edm::InputTag srcWeights = iConfig.getParameter<edm::InputTag>("srcWeights");
26  if (srcWeights.label().empty())
27  throw cms::Exception("InvalidInput") << "applyWeight set to True, but no weights given in PFMETProducer\n";
28  if (srcWeights.label() == src_.label())
29  edm::LogWarning("PFMETProducer")
30  << "Particle and weights collection have the same label. You may be applying the same weights twice.\n";
31  weightsToken_ = consumes<edm::ValueMap<float>>(srcWeights);
32  }
34  metSigAlgo_ = new metsig::METSignificance(iConfig);
35 
36  jetToken_ = mayConsume<edm::View<reco::Jet>>(iConfig.getParameter<edm::InputTag>("srcJets"));
37  std::vector<edm::InputTag> srcLeptonsTags = iConfig.getParameter<std::vector<edm::InputTag>>("srcLeptons");
38  for (std::vector<edm::InputTag>::const_iterator it = srcLeptonsTags.begin(); it != srcLeptonsTags.end(); it++) {
39  lepTokens_.push_back(mayConsume<edm::View<reco::Candidate>>(*it));
40  }
41 
42  jetSFToken_ = esConsumes(edm::ESInputTag("", iConfig.getParameter<std::string>("srcJetSF")));
43  jetResPtToken_ = esConsumes(edm::ESInputTag("", iConfig.getParameter<std::string>("srcJetResPt")));
44  jetResPhiToken_ = esConsumes(edm::ESInputTag("", iConfig.getParameter<std::string>("srcJetResPhi")));
45  rhoToken_ = consumes<double>(iConfig.getParameter<edm::InputTag>("srcRho"));
46  }
47 
48  std::string alias = iConfig.getParameter<std::string>("alias");
49 
50  produces<reco::PFMETCollection>().setBranchAlias(alias);
51  }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::ESGetToken< JME::JetResolutionObject, JetResolutionRcd > jetResPtToken_
Definition: PFMETProducer.h:81
std::vector< edm::EDGetTokenT< edm::View< reco::Candidate > > > lepTokens_
Definition: PFMETProducer.h:78
std::string const & label() const
Definition: InputTag.h:36
edm::EDGetTokenT< edm::View< reco::Jet > > jetToken_
Definition: PFMETProducer.h:77
edm::InputTag src_
Definition: PFMETProducer.h:68
edm::ValueMap< float > const * weights_
Definition: PFMETProducer.h:86
metsig::METSignificance * metSigAlgo_
Definition: PFMETProducer.h:72
edm::ESGetToken< JME::JetResolutionObject, JetResolutionScaleFactorRcd > jetSFToken_
Definition: PFMETProducer.h:80
METSignificance
____________________________________________________________________________||
edm::EDGetTokenT< edm::View< reco::Candidate > > inputToken_
Definition: PFMETProducer.h:69
edm::EDGetTokenT< double > rhoToken_
Definition: PFMETProducer.h:83
edm::ESGetToken< JME::JetResolutionObject, JetResolutionRcd > jetResPhiToken_
Definition: PFMETProducer.h:82
edm::EDGetTokenT< edm::ValueMap< float > > weightsToken_
Definition: PFMETProducer.h:85

◆ ~PFMETProducer()

cms::PFMETProducer::~PFMETProducer ( )
inlineoverride

Definition at line 60 of file PFMETProducer.h.

60 {}

Member Function Documentation

◆ fillDescriptions()

void PFMETProducer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 132 of file PFMETProducer.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), submitPVResolutionJobs::desc, HLT_2022v12_cff::InputTag, submitPVValidationJobs::params, and AlCaHLTBitMon_QueryRunRegistry::string.

132  {
134  desc.add<edm::InputTag>("src", edm::InputTag("particleFlow"));
135  desc.add<double>("globalThreshold", 0.);
136  desc.add<std::string>("alias", "@module_label");
137  desc.add<bool>("calculateSignificance", false);
138  desc.addOptional<edm::InputTag>("srcJets");
139  desc.addOptional<std::vector<edm::InputTag>>("srcLeptons");
140  desc.addOptional<std::string>("srcJetSF");
141  desc.addOptional<std::string>("srcJetResPt");
142  desc.addOptional<std::string>("srcJetResPhi");
143  desc.addOptional<edm::InputTag>("srcRho");
145  params.setAllowAnything(); // FIXME: This still needs to be defined in METSignficance
146  desc.addOptional<edm::ParameterSetDescription>("parameters", params);
149  desc1.add<bool>("applyWeight", false);
150  desc1.add<edm::InputTag>("srcWeights", edm::InputTag(""));
151  descriptions.add("pfMet", desc1);
152  desc2.add<bool>("applyWeight", true);
153  desc2.add<edm::InputTag>("srcWeights", edm::InputTag("puppiNoLep"));
154  descriptions.add("pfMetPuppi", desc2);
155  }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)

◆ getMETCovMatrix()

reco::METCovMatrix PFMETProducer::getMETCovMatrix ( const edm::Event event,
const edm::EventSetup setup,
const edm::Handle< edm::View< reco::Candidate >> &  input 
) const
private

Definition at line 86 of file PFMETProducer.cc.

References metsig::METSignificance::getCovariance(), HLT_2022v12_cff::inputJets, jetResPhiToken_, jetResPtToken_, jetSFToken_, jetToken_, lepTokens_, HLT_2022v12_cff::leptons, metSigAlgo_, rhoToken_, singleTopDQM_cfi::setup, met_cff::sumPtUnclustered, trackerHitRTTI::vector, and weights_.

Referenced by produce().

88  {
89  // leptons
90  std::vector<edm::Handle<reco::CandidateView>> leptons;
91  for (std::vector<edm::EDGetTokenT<edm::View<reco::Candidate>>>::const_iterator srcLeptons_i = lepTokens_.begin();
92  srcLeptons_i != lepTokens_.end();
93  ++srcLeptons_i) {
95  event.getByToken(*srcLeptons_i, leptons_i);
96  leptons.push_back(leptons_i);
97  /*
98  for ( reco::CandidateView::const_iterator lepton = leptons_i->begin();
99  lepton != leptons_i->end(); ++lepton ) {
100  leptons.push_back(*lepton);
101  }
102  */
103  }
104 
105  // jets
107  event.getByToken(jetToken_, inputJets);
108 
109  JME::JetResolution resPtObj(setup.getData(jetResPtToken_));
110  JME::JetResolution resPhiObj(setup.getData(jetResPhiToken_));
112 
114  event.getByToken(rhoToken_, rho);
115 
116  //Compute the covariance matrix and fill it
117  double sumPtUnclustered = 0;
119  leptons,
120  candInput,
121  *rho,
122  resPtObj,
123  resPhiObj,
124  resSFObj,
125  event.isRealData(),
127  weights_);
128 
129  return cov;
130  }
reco::METCovMatrix getCovariance(const edm::View< reco::Jet > &jets, const std::vector< edm::Handle< reco::CandidateView > > &leptons, const edm::Handle< edm::View< reco::Candidate > > &pfCandidates, double rho, JME::JetResolution &resPtObj, JME::JetResolution &resPhiObj, JME::JetResolutionScaleFactor &resSFObj, bool isRealData, double &sumPtUnclustered, edm::ValueMap< float > const *weights=nullptr)
edm::ESGetToken< JME::JetResolutionObject, JetResolutionRcd > jetResPtToken_
Definition: PFMETProducer.h:81
ROOT::Math::SMatrix< double, 2 > METCovMatrix
Definition: MET.h:39
std::vector< edm::EDGetTokenT< edm::View< reco::Candidate > > > lepTokens_
Definition: PFMETProducer.h:78
edm::EDGetTokenT< edm::View< reco::Jet > > jetToken_
Definition: PFMETProducer.h:77
edm::ValueMap< float > const * weights_
Definition: PFMETProducer.h:86
metsig::METSignificance * metSigAlgo_
Definition: PFMETProducer.h:72
edm::ESGetToken< JME::JetResolutionObject, JetResolutionScaleFactorRcd > jetSFToken_
Definition: PFMETProducer.h:80
sumPtUnclustered
Definition: met_cff.py:19
edm::EDGetTokenT< double > rhoToken_
Definition: PFMETProducer.h:83
edm::ESGetToken< JME::JetResolutionObject, JetResolutionRcd > jetResPhiToken_
Definition: PFMETProducer.h:82
Definition: event.py:1

◆ produce()

void PFMETProducer::produce ( edm::Event event,
const edm::EventSetup setup 
)
override

Definition at line 54 of file PFMETProducer.cc.

References applyWeight_, calculateSignificance_, edm::ValueMap< T >::get(), getMETCovMatrix(), globalThreshold_, input, inputToken_, CommonMETData::met, CommonMETData::mex, CommonMETData::mey, eostools::move(), packedPFCandidateRefMixer_cfi::pf, reco::MET::setIsWeighted(), reco::MET::setSignificanceMatrix(), singleTopDQM_cfi::setup, timingPdfMaker::specific, CommonMETData::sumet, extraflags_cff::vtx, weights_, and weightsToken_.

54  {
56  event.getByToken(inputToken_, input);
57 
58  if (applyWeight_)
59  weights_ = &event.get(weightsToken_);
60 
61  METAlgo algo;
62  CommonMETData commonMETdata;
63  commonMETdata = algo.run(*input.product(), globalThreshold_, weights_);
64 
65  const math::XYZTLorentzVector p4(commonMETdata.mex, commonMETdata.mey, 0.0, commonMETdata.met);
66  const math::XYZPoint vtx(0.0, 0.0, 0.0);
67 
70  specific = pf.run(*input.product(), weights_);
71 
72  reco::PFMET pfmet(specific, commonMETdata.sumet, p4, vtx);
73  pfmet.setIsWeighted(applyWeight_);
74 
77  pfmet.setSignificanceMatrix(sigcov);
78  }
79 
80  auto pfmetcoll = std::make_unique<reco::PFMETCollection>();
81 
82  pfmetcoll->push_back(pfmet);
83  event.put(std::move(pfmetcoll));
84  }
const_reference_type get(ProductID id, size_t idx) const
Definition: ValueMap.h:144
ROOT::Math::SMatrix< double, 2 > METCovMatrix
Definition: MET.h:39
static std::string const input
Definition: EdmProvDump.cc:47
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
reco::METCovMatrix getMETCovMatrix(const edm::Event &event, const edm::EventSetup &, const edm::Handle< edm::View< reco::Candidate >> &input) const
edm::ValueMap< float > const * weights_
Definition: PFMETProducer.h:86
Structure containing data common to all types of MET.
Definition: CommonMETData.h:12
MET made from Particle Flow Candidates.
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
edm::EDGetTokenT< edm::View< reco::Candidate > > inputToken_
Definition: PFMETProducer.h:69
def move(src, dest)
Definition: eostools.py:511
Definition: event.py:1
edm::EDGetTokenT< edm::ValueMap< float > > weightsToken_
Definition: PFMETProducer.h:85

Member Data Documentation

◆ applyWeight_

bool cms::PFMETProducer::applyWeight_
private

Definition at line 84 of file PFMETProducer.h.

Referenced by PFMETProducer(), and produce().

◆ calculateSignificance_

bool cms::PFMETProducer::calculateSignificance_
private

Definition at line 71 of file PFMETProducer.h.

Referenced by PFMETProducer(), and produce().

◆ globalThreshold_

double cms::PFMETProducer::globalThreshold_
private

Definition at line 74 of file PFMETProducer.h.

Referenced by produce().

◆ inputToken_

edm::EDGetTokenT<edm::View<reco::Candidate> > cms::PFMETProducer::inputToken_
private

Definition at line 69 of file PFMETProducer.h.

Referenced by produce().

◆ jetResPhiToken_

edm::ESGetToken<JME::JetResolutionObject, JetResolutionRcd> cms::PFMETProducer::jetResPhiToken_
private

Definition at line 82 of file PFMETProducer.h.

Referenced by getMETCovMatrix(), and PFMETProducer().

◆ jetResPtToken_

edm::ESGetToken<JME::JetResolutionObject, JetResolutionRcd> cms::PFMETProducer::jetResPtToken_
private

Definition at line 81 of file PFMETProducer.h.

Referenced by getMETCovMatrix(), and PFMETProducer().

◆ jetSFToken_

edm::ESGetToken<JME::JetResolutionObject, JetResolutionScaleFactorRcd> cms::PFMETProducer::jetSFToken_
private

Definition at line 80 of file PFMETProducer.h.

Referenced by getMETCovMatrix(), and PFMETProducer().

◆ jetThreshold_

double cms::PFMETProducer::jetThreshold_
private

Definition at line 75 of file PFMETProducer.h.

◆ jetToken_

edm::EDGetTokenT<edm::View<reco::Jet> > cms::PFMETProducer::jetToken_
private

Definition at line 77 of file PFMETProducer.h.

Referenced by getMETCovMatrix(), and PFMETProducer().

◆ lepTokens_

std::vector<edm::EDGetTokenT<edm::View<reco::Candidate> > > cms::PFMETProducer::lepTokens_
private

Definition at line 78 of file PFMETProducer.h.

Referenced by getMETCovMatrix(), and PFMETProducer().

◆ metSigAlgo_

metsig::METSignificance* cms::PFMETProducer::metSigAlgo_
private

Definition at line 72 of file PFMETProducer.h.

Referenced by getMETCovMatrix(), and PFMETProducer().

◆ rhoToken_

edm::EDGetTokenT<double> cms::PFMETProducer::rhoToken_
private

Definition at line 83 of file PFMETProducer.h.

Referenced by getMETCovMatrix(), and PFMETProducer().

◆ src_

edm::InputTag cms::PFMETProducer::src_
private

Definition at line 68 of file PFMETProducer.h.

Referenced by PFMETProducer().

◆ weights_

edm::ValueMap<float> const* cms::PFMETProducer::weights_
private

Definition at line 86 of file PFMETProducer.h.

Referenced by getMETCovMatrix(), and produce().

◆ weightsToken_

edm::EDGetTokenT<edm::ValueMap<float> > cms::PFMETProducer::weightsToken_
private

Definition at line 85 of file PFMETProducer.h.

Referenced by PFMETProducer(), and produce().