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 //____________________________________________________________________________||
13 namespace cms
14 {
15 
16 //____________________________________________________________________________||
18  {
19 
20  std::vector<edm::InputTag> srcLeptonsTags = iConfig.getParameter< std::vector<edm::InputTag> >("srcLeptons");
21  for(std::vector<edm::InputTag>::const_iterator it=srcLeptonsTags.begin();it!=srcLeptonsTags.end();it++) {
22  lepTokens_.push_back( consumes<edm::View<reco::Candidate> >( *it ) );
23  }
24 
25  pfjetsToken_ = consumes<edm::View<reco::Jet> >(iConfig.getParameter<edm::InputTag>("srcPfJets"));
26 
27  metToken_ = consumes<edm::View<reco::MET> >(iConfig.getParameter<edm::InputTag>("srcMet"));
28  pfCandidatesToken_ = consumes<edm::View<reco::Candidate> >(iConfig.getParameter<edm::InputTag>("srcPFCandidates"));
29 
30  jetSFType_ = iConfig.getParameter<std::string>("srcJetSF");
31  jetResPtType_ = iConfig.getParameter<std::string>("srcJetResPt");
32  jetResPhiType_ = iConfig.getParameter<std::string>("srcJetResPhi");
33  rhoToken_ = consumes<double>(iConfig.getParameter<edm::InputTag>("srcRho"));
34 
35  metSigAlgo_ = new metsig::METSignificance(iConfig);
36 
37  produces<double>("METSignificance");
38  produces<math::Error<2>::type>("METCovariance");
39 
40  }
41 
43  {
44  delete metSigAlgo_;
45  }
46 
47 //____________________________________________________________________________||
49  {
50  //
51  // met
52  //
54  event.getByToken(metToken_, metHandle);
55  const reco::MET& met = (*metHandle)[0];
56 
57  //
58  // candidates
59  //
61  event.getByToken( pfCandidatesToken_, pfCandidates );
62 
63  //
64  // leptons
65  //
66  std::vector< edm::Handle<reco::CandidateView> > leptons;
67  for ( std::vector<edm::EDGetTokenT<edm::View<reco::Candidate> > >::const_iterator srcLeptons_i = lepTokens_.begin();
68  srcLeptons_i != lepTokens_.end(); ++srcLeptons_i ) {
69 
71  event.getByToken(*srcLeptons_i, leptons_i);
72  leptons.push_back( leptons_i );
73 
74  }
75 
76  //
77  // jets
78  //
80  event.getByToken( pfjetsToken_, jets );
81 
83  event.getByToken(rhoToken_, rho);
84 
88 
89  //
90  // compute the significance
91  //
92  const reco::METCovMatrix cov = metSigAlgo_->getCovariance( *jets, leptons, pfCandidates, *rho, resPtObj, resPhiObj, resSFObj, event.isRealData() );
93  double sig = metSigAlgo_->getSignificance(cov, met);
94 
95  auto significance = std::make_unique<double>();
96  (*significance) = sig;
97 
99  (*covPtr)(0,0) = cov(0,0);
100  (*covPtr)(1,0) = cov(1,0);
101  (*covPtr)(1,1) = cov(1,1);
102 
103  event.put(std::move(covPtr), "METCovariance" );
104  event.put(std::move(significance), "METSignificance" );
105 
106  }
107 
108 //____________________________________________________________________________||
110 }
111 
112 //____________________________________________________________________________||
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 &)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< edm::View< reco::Candidate > > pfCandidatesToken_
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:2
ROOT::Math::SMatrix< double, 2 > METCovMatrix
Definition: MET.h:40
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)
std::vector< edm::EDGetTokenT< edm::View< reco::Candidate > > > lepTokens_
bool isRealData() const
Definition: EventBase.h:64
edm::EDGetTokenT< double > rhoToken_
Definition: MET.h:42
vector< PseudoJet > jets
Namespace of DDCMS conversion namespace.
METSignificanceProducer(const edm::ParameterSet &)
METSignificance
____________________________________________________________________________||
met
===> hadronic RAZOR
static const JetResolutionScaleFactor get(const edm::EventSetup &, const std::string &)
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