CMS 3D CMS Logo

PFMETProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: METProducers
4 // Class: PFMETProducer
5 //
6 //
7 
8 //____________________________________________________________________________||
10 
11 //____________________________________________________________________________||
12 namespace cms
13 {
14 
15 //____________________________________________________________________________||
17  : inputToken_(consumes<edm::View<reco::Candidate> >(iConfig.getParameter<edm::InputTag>("src")))
18  , calculateSignificance_(iConfig.getParameter<bool>("calculateSignificance"))
19  , globalThreshold_(iConfig.getParameter<double>("globalThreshold"))
20  {
22  {
23  metSigAlgo_ = new metsig::METSignificance(iConfig);
24 
25  jetToken_ = mayConsume<edm::View<reco::Jet> >(iConfig.getParameter<edm::InputTag>("srcJets"));
26  std::vector<edm::InputTag> srcLeptonsTags = iConfig.getParameter< std::vector<edm::InputTag> >("srcLeptons");
27  for(std::vector<edm::InputTag>::const_iterator it=srcLeptonsTags.begin();it!=srcLeptonsTags.end();it++) {
28  lepTokens_.push_back( mayConsume<edm::View<reco::Candidate> >( *it ) );
29  }
30 
31  jetSFType_ = iConfig.getParameter<std::string>("srcJetSF");
32  jetResPtType_ = iConfig.getParameter<std::string>("srcJetResPt");
33  jetResPhiType_ = iConfig.getParameter<std::string>("srcJetResPhi");
34  rhoToken_ = consumes<double>(iConfig.getParameter<edm::InputTag>("srcRho"));
35  }
36 
37  std::string alias = iConfig.exists("alias") ? iConfig.getParameter<std::string>("alias") : "";
38 
39  produces<reco::PFMETCollection>().setBranchAlias(alias);
40  }
41 
42 //____________________________________________________________________________||
44  {
46  event.getByToken(inputToken_, input);
47 
48  METAlgo algo;
49  CommonMETData commonMETdata = algo.run(*input.product(), globalThreshold_);
50 
51  const math::XYZTLorentzVector p4(commonMETdata.mex, commonMETdata.mey, 0.0, commonMETdata.met);
52  const math::XYZPoint vtx(0.0, 0.0, 0.0);
53 
55  SpecificPFMETData specific = pf.run(*input.product());
56 
57  reco::PFMET pfmet(specific, commonMETdata.sumet, p4, vtx);
58 
60  {
61  reco::METCovMatrix sigcov = getMETCovMatrix(event, setup, input);
62  pfmet.setSignificanceMatrix(sigcov);
63  }
64 
65  auto pfmetcoll = std::make_unique<reco::PFMETCollection>();
66 
67  pfmetcoll->push_back(pfmet);
68  event.put(std::move(pfmetcoll));
69  }
70 
71 
72 
74 
75  // leptons
76  std::vector< edm::Handle<reco::CandidateView> > leptons;
77  for ( std::vector<edm::EDGetTokenT<edm::View<reco::Candidate> > >::const_iterator srcLeptons_i = lepTokens_.begin();
78  srcLeptons_i != lepTokens_.end(); ++srcLeptons_i ) {
80  event.getByToken(*srcLeptons_i, leptons_i);
81  leptons.push_back( leptons_i );
82  /*
83  for ( reco::CandidateView::const_iterator lepton = leptons_i->begin();
84  lepton != leptons_i->end(); ++lepton ) {
85  leptons.push_back(*lepton);
86  }
87  */
88  }
89 
90  // jets
92  event.getByToken( jetToken_, inputJets );
93 
97 
99  event.getByToken(rhoToken_, rho);
100 
101  //Compute the covariance matrix and fill it
102  reco::METCovMatrix cov = metSigAlgo_->getCovariance( *inputJets, leptons, candInput, *rho, resPtObj, resPhiObj, resSFObj, event.isRealData() );
103 
104  return cov;
105  }
106 
107 
108 
109 
110 
111 //____________________________________________________________________________||
113 }
114 
115 //____________________________________________________________________________||
reco::METCovMatrix getMETCovMatrix(const edm::Event &event, const edm::EventSetup &, const edm::Handle< edm::View< reco::Candidate > > &input) const
T getParameter(std::string const &) const
static const JetResolution get(const edm::EventSetup &, const std::string &)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void setSignificanceMatrix(const reco::METCovMatrix &matrix)
Definition: MET.cc:157
bool exists(std::string const &parameterName) const
checks if a parameter exists
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
PFMETProducer(const edm::ParameterSet &)
ROOT::Math::SMatrix< double, 2 > METCovMatrix
Definition: MET.h:40
std::string jetResPtType_
Definition: PFMETProducer.h:84
edm::EDGetTokenT< edm::View< reco::Jet > > jetToken_
Definition: PFMETProducer.h:81
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)
bool isRealData() const
Definition: EventBase.h:64
static std::string const input
Definition: EdmProvDump.cc:44
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
CommonMETData run(const edm::View< reco::Candidate > &candidates, double globalThreshold=0.0)
Definition: METAlgo.cc:16
double p4[4]
Definition: TauolaWrapper.h:92
void produce(edm::Event &, const edm::EventSetup &) override
metsig::METSignificance * metSigAlgo_
Definition: PFMETProducer.h:76
Structure containing data common to all types of MET.
Definition: CommonMETData.h:12
std::string jetSFType_
Definition: PFMETProducer.h:83
std::vector< edm::EDGetTokenT< edm::View< reco::Candidate > > > lepTokens_
Definition: PFMETProducer.h:82
MET made from Particle Flow Candidates.
SpecificPFMETData run(const edm::View< reco::Candidate > &pfCands)
Namespace of DDCMS conversion namespace.
std::string jetResPhiType_
Definition: PFMETProducer.h:85
METSignificance
____________________________________________________________________________||
T const * product() const
Definition: Handle.h:81
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
static const JetResolutionScaleFactor get(const edm::EventSetup &, const std::string &)
fixed size matrix
HLT enums.
edm::EDGetTokenT< edm::View< reco::Candidate > > inputToken_
Definition: PFMETProducer.h:73
edm::EDGetTokenT< double > rhoToken_
Definition: PFMETProducer.h:86
def move(src, dest)
Definition: eostools.py:510
Definition: event.py:1