CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PFMETProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: METProducers
4 // Class: PFMETProducer
5 //
6 //
7 
8 //____________________________________________________________________________||
12 
13 //____________________________________________________________________________||
14 namespace cms {
15 
16  //____________________________________________________________________________||
18  : src_(iConfig.getParameter<edm::InputTag>("src")),
19  inputToken_(consumes<edm::View<reco::Candidate>>(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.exists("alias") ? iConfig.getParameter<std::string>("alias") : "";
49 
50  produces<reco::PFMETCollection>().setBranchAlias(alias);
51  }
52 
53  //____________________________________________________________________________||
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 
68  PFSpecificAlgo pf;
70  specific = pf.run(*input.product(), weights_);
71 
72  reco::PFMET pfmet(specific, commonMETdata.sumet, p4, vtx);
74 
76  reco::METCovMatrix sigcov = getMETCovMatrix(event, setup, input);
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  }
85 
87  const edm::EventSetup& setup,
88  const edm::Handle<edm::View<reco::Candidate>>& candInput) const {
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;
118  reco::METCovMatrix cov = metSigAlgo_->getCovariance(*inputJets,
119  leptons,
120  candInput,
121  *rho,
122  resPtObj,
123  resPhiObj,
124  resSFObj,
125  event.isRealData(),
126  sumPtUnclustered,
127  weights_);
128 
129  return cov;
130  }
131 
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  }
156 
157  //____________________________________________________________________________||
159 } // namespace cms
160 
161 //____________________________________________________________________________||
void setIsWeighted(bool isWeighted)
Set boolean if weights were applied by algorithm (e.g. PUPPI weights)
Definition: MET.h:77
SpecificPFMETData run(const edm::View< reco::Candidate > &pfCands, edm::ValueMap< float > const *weights=nullptr)
dictionary specific
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)
ParameterDescriptionBase * addOptional(U const &iLabel, T const &value)
edm::ESGetToken< JME::JetResolutionObject, JetResolutionRcd > jetResPtToken_
Definition: PFMETProducer.h:81
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ROOT::Math::SMatrix< double, 2 > METCovMatrix
Definition: MET.h:39
void setSignificanceMatrix(const reco::METCovMatrix &matrix)
Definition: MET.cc:142
std::vector< edm::EDGetTokenT< edm::View< reco::Candidate > > > lepTokens_
Definition: PFMETProducer.h:78
bool exists(std::string const &parameterName) const
checks if a parameter exists
PFMETProducer(const edm::ParameterSet &)
edm::EDGetTokenT< edm::View< reco::Jet > > jetToken_
Definition: PFMETProducer.h:77
edm::InputTag src_
Definition: PFMETProducer.h:68
bool isRealData() const
Definition: EventBase.h:62
static std::string const input
Definition: EdmProvDump.cc:47
bool getData(T &iHolder) const
Definition: EventSetup.h:128
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
tuple METSignificance
____________________________________________________________________________||
edm::ValueMap< float > const * weights_
Definition: PFMETProducer.h:86
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void produce(edm::Event &, const edm::EventSetup &) override
def move
Definition: eostools.py:511
metsig::METSignificance * metSigAlgo_
Definition: PFMETProducer.h:72
Structure containing data common to all types of MET.
Definition: CommonMETData.h:12
ParameterDescriptionBase * add(U const &iLabel, T const &value)
const_reference_type get(ProductID id, size_t idx) const
Definition: ValueMap.h:144
MET made from Particle Flow Candidates.
edm::ESGetToken< JME::JetResolutionObject, JetResolutionScaleFactorRcd > jetSFToken_
Definition: PFMETProducer.h:80
reco::METCovMatrix getMETCovMatrix(const edm::Event &event, const edm::EventSetup &, const edm::Handle< edm::View< reco::Candidate >> &input) const
T const * product() const
Definition: Handle.h:70
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::string const & label() const
Definition: InputTag.h:36
edm::EDGetTokenT< edm::View< reco::Candidate > > inputToken_
Definition: PFMETProducer.h:69
CommonMETData run(const edm::View< reco::Candidate > &candidates, double globalThreshold=0.0, edm::ValueMap< float > const *weights=nullptr)
Definition: METAlgo.cc:16
edm::EDGetTokenT< double > rhoToken_
Definition: PFMETProducer.h:83
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
edm::ESGetToken< JME::JetResolutionObject, JetResolutionRcd > jetResPhiToken_
Definition: PFMETProducer.h:82
edm::EDGetTokenT< edm::ValueMap< float > > weightsToken_
Definition: PFMETProducer.h:85