CMS 3D CMS Logo

METSignificanceProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: METProducers
4 // Class: METSignificanceProducer
5 //
6 //
7 
8 //____________________________________________________________________________||
10 
11 //____________________________________________________________________________||
12 namespace cms {
13 
14  //____________________________________________________________________________||
16  std::vector<edm::InputTag> srcLeptonsTags = iConfig.getParameter<std::vector<edm::InputTag> >("srcLeptons");
17  for (std::vector<edm::InputTag>::const_iterator it = srcLeptonsTags.begin(); it != srcLeptonsTags.end(); it++) {
18  lepTokens_.push_back(consumes<edm::View<reco::Candidate> >(*it));
19  }
20 
21  pfjetsToken_ = consumes<edm::View<reco::Jet> >(iConfig.getParameter<edm::InputTag>("srcPfJets"));
22 
23  metToken_ = consumes<edm::View<reco::MET> >(iConfig.getParameter<edm::InputTag>("srcMet"));
24  pfCandidatesToken_ = consumes<edm::View<reco::Candidate> >(iConfig.getParameter<edm::InputTag>("srcPFCandidates"));
25 
26  jetSFType_ = iConfig.getParameter<std::string>("srcJetSF");
27  jetResPtType_ = iConfig.getParameter<std::string>("srcJetResPt");
28  jetResPhiType_ = iConfig.getParameter<std::string>("srcJetResPhi");
29  rhoToken_ = consumes<double>(iConfig.getParameter<edm::InputTag>("srcRho"));
30 
31  metSigAlgo_ = new metsig::METSignificance(iConfig);
32 
33  produces<double>("METSignificance");
34  produces<math::Error<2>::type>("METCovariance");
35  }
36 
38 
39  //____________________________________________________________________________||
41  //
42  // met
43  //
45  event.getByToken(metToken_, metHandle);
46  const reco::MET& met = (*metHandle)[0];
47 
48  //
49  // candidates
50  //
52  event.getByToken(pfCandidatesToken_, pfCandidates);
53 
54  //
55  // leptons
56  //
57  std::vector<edm::Handle<reco::CandidateView> > leptons;
58  for (std::vector<edm::EDGetTokenT<edm::View<reco::Candidate> > >::const_iterator srcLeptons_i = lepTokens_.begin();
59  srcLeptons_i != lepTokens_.end();
60  ++srcLeptons_i) {
62  event.getByToken(*srcLeptons_i, leptons_i);
63  leptons.push_back(leptons_i);
64  }
65 
66  //
67  // jets
68  //
70  event.getByToken(pfjetsToken_, jets);
71 
73  event.getByToken(rhoToken_, rho);
74 
78 
79  //
80  // compute the significance
81  //
82  double sumPtUnclustered = 0;
84  *jets, leptons, pfCandidates, *rho, resPtObj, resPhiObj, resSFObj, event.isRealData(), sumPtUnclustered);
85  double sig = metSigAlgo_->getSignificance(cov, met);
86 
87  auto significance = std::make_unique<double>();
88  (*significance) = sig;
89 
91  (*covPtr)(0, 0) = cov(0, 0);
92  (*covPtr)(1, 0) = cov(1, 0);
93  (*covPtr)(1, 1) = cov(1, 1);
94 
95  event.put(std::move(covPtr), "METCovariance");
96  event.put(std::move(significance), "METSignificance");
97  }
98 
99  //____________________________________________________________________________||
101 } // namespace cms
102 
103 //____________________________________________________________________________||
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
void produce(edm::Event &, const edm::EventSetup &) override
static double getSignificance(const reco::METCovMatrix &cov, const reco::MET &met)
static const JetResolution get(const edm::EventSetup &, const std::string &)
ROOT::Math::SMatrix< double, 2 > METCovMatrix
Definition: MET.h:39
edm::EDGetTokenT< edm::View< reco::Candidate > > pfCandidatesToken_
std::vector< edm::EDGetTokenT< edm::View< reco::Candidate > > > lepTokens_
bool isRealData() const
Definition: EventBase.h:62
edm::EDGetTokenT< double > rhoToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Definition: MET.h:41
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)
Namespace of DDCMS conversion namespace.
METSignificanceProducer(const edm::ParameterSet &)
METSignificance
____________________________________________________________________________||
static const JetResolutionScaleFactor get(const edm::EventSetup &, const std::string &)
significance
Definition: met_cff.py:19
metsig::METSignificance * metSigAlgo_
edm::EDGetTokenT< edm::View< reco::Jet > > pfjetsToken_
edm::EDGetTokenT< edm::View< reco::MET > > metToken_
def move(src, dest)
Definition: eostools.py:511
Definition: event.py:1