CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
PFMEtSignInterfaceBase.h
Go to the documentation of this file.
1 #ifndef RecoMET_METPUSubtraction_PFMEtSignInterfaceBase_h
2 #define RecoMET_METPUSubtraction_PFMEtSignInterfaceBase_h
3 
18 
20 
36 
37 #include <TFile.h>
38 #include <TH2.h>
39 
41 public:
44 
45  template <typename T>
46  metsig::SigInputObj compResolution(const T* particle) const {
47  double pt = particle->pt();
48  double eta = particle->eta();
49  double phi = particle->phi();
50 
51  if (dynamic_cast<const reco::GsfElectron*>(particle) != nullptr ||
52  dynamic_cast<const pat::Electron*>(particle) != nullptr) {
53  std::string particleType = "electron";
54  // WARNING: SignAlgoResolutions::PFtype2 needs to be kept in sync with reco::PFCandidate::e !!
55  double dpt = pfMEtResolution_->eval(metsig::PFtype2, metsig::ET, pt, phi, eta);
56  double dphi = pfMEtResolution_->eval(metsig::PFtype2, metsig::PHI, pt, phi, eta);
57  //std::cout << "electron: pt = " << pt << ", eta = " << eta << ", phi = " << phi
58  // << " --> dpt = " << dpt << ", dphi = " << dphi << std::endl;
59  return metsig::SigInputObj(particleType, pt, phi, dpt, dphi);
60  } else if (dynamic_cast<const reco::Photon*>(particle) != nullptr ||
61  dynamic_cast<const pat::Photon*>(particle) != nullptr) {
62  // CV: assume resolutions for photons to be identical to electron resolutions
63  std::string particleType = "electron";
64  // WARNING: SignAlgoResolutions::PFtype2 needs to be kept in sync with reco::PFCandidate::e !!
65  double dpt = pfMEtResolution_->eval(metsig::PFtype2, metsig::ET, pt, phi, eta);
66  double dphi = pfMEtResolution_->eval(metsig::PFtype2, metsig::PHI, pt, phi, eta);
67  //std::cout << "photon: pt = " << pt << ", eta = " << eta << ", phi = " << phi
68  // << " --> dpt = " << dpt << ", dphi = " << dphi << std::endl;
69  return metsig::SigInputObj(particleType, pt, phi, dpt, dphi);
70  } else if (dynamic_cast<const reco::Muon*>(particle) != nullptr ||
71  dynamic_cast<const pat::Muon*>(particle) != nullptr) {
72  std::string particleType = "muon";
73  double dpt, dphi;
74  const reco::Track* muonTrack = nullptr;
75  if (dynamic_cast<const pat::Muon*>(particle) != nullptr) {
76  const pat::Muon* muon = dynamic_cast<const pat::Muon*>(particle);
77  if (muon->track().isNonnull() && muon->track().isAvailable())
78  muonTrack = muon->track().get();
79  } else if (dynamic_cast<const reco::Muon*>(particle) != nullptr) {
80  const reco::Muon* muon = dynamic_cast<const reco::Muon*>(particle);
81  if (muon->track().isNonnull() && muon->track().isAvailable())
82  muonTrack = muon->track().get();
83  } else
84  assert(0);
85  if (muonTrack) {
86  dpt = muonTrack->ptError();
87  dphi = pt * muonTrack->phiError(); // CV: pt*dphi is indeed correct
88  } else {
89  // WARNING: SignAlgoResolutions::PFtype3 needs to be kept in sync with reco::PFCandidate::mu !!
90  dpt = pfMEtResolution_->eval(metsig::PFtype3, metsig::ET, pt, phi, eta);
91  dphi = pfMEtResolution_->eval(metsig::PFtype3, metsig::PHI, pt, phi, eta);
92  }
93  //std::cout << "muon: pt = " << pt << ", eta = " << eta << ", phi = " << phi
94  // << " --> dpt = " << dpt << ", dphi = " << dphi << std::endl;
95  return metsig::SigInputObj(particleType, pt, phi, dpt, dphi);
96  } else if (dynamic_cast<const reco::PFTau*>(particle) != nullptr ||
97  dynamic_cast<const pat::Tau*>(particle) != nullptr) {
98  // CV: use PFJet resolutions for PFTaus for now...
99  // (until PFTau specific resolutions are available)
100  if (dynamic_cast<const pat::Tau*>(particle) != nullptr) {
101  const pat::Tau* pfTau = dynamic_cast<const pat::Tau*>(particle);
102  return pfMEtResolution_->evalPFJet(pfTau->pfJetRef().get());
103  } else if (dynamic_cast<const reco::PFTau*>(particle) != nullptr) {
104  const reco::PFTau* pfTau = dynamic_cast<const reco::PFTau*>(particle);
105  return pfMEtResolution_->evalPFJet(pfTau->jetRef().get());
106  } else
107  assert(0);
108  } else if (dynamic_cast<const reco::PFJet*>(particle) != nullptr ||
109  dynamic_cast<const pat::Jet*>(particle) != nullptr) {
110  metsig::SigInputObj pfJetResolution;
111  if (dynamic_cast<const reco::PFJet*>(particle) != nullptr) {
112  const reco::PFJet* pfJet = dynamic_cast<const reco::PFJet*>(particle);
113  pfJetResolution = pfMEtResolution_->evalPFJet(pfJet);
114  } else if (dynamic_cast<const pat::Jet*>(particle) != nullptr) {
115  const pat::Jet* jet = dynamic_cast<const pat::Jet*>(particle);
116  if (jet->isPFJet()) {
117  reco::PFJet pfJet(jet->p4(), jet->vertex(), jet->pfSpecific(), jet->getJetConstituents());
118  pfJetResolution = pfMEtResolution_->evalPFJet(&pfJet);
119  } else
120  throw cms::Exception("addPFMEtSignObjects") << "PAT jet not of PF-type !!\n";
121  } else
122  assert(0);
123  //std::cout << "pfJet: pt = " << pt << ", eta = " << eta << ", phi = " << phi << std::endl;
124  // CV: apply additional jet energy resolution corrections
125  // not included in (PF)MEt significance algorithm yet
126  // (cf. CMS AN-11/400 vs. CMS AN-11/330)
127  if (lut_ && pt > 10.) {
128  double x = std::abs(eta);
129  double y = pt;
130  if (x > lut_->GetXaxis()->GetXmin() && x < lut_->GetXaxis()->GetXmax() && y > lut_->GetYaxis()->GetXmin() &&
131  y < lut_->GetYaxis()->GetXmax()) {
132  int binIndex = lut_->FindBin(x, y);
133  double addJERcorrFactor = lut_->GetBinContent(binIndex);
134  //std::cout << " addJERcorrFactor = " << addJERcorrFactor << std::endl;
135  pfJetResolution.set(pfJetResolution.get_type(),
136  pfJetResolution.get_energy(),
137  pfJetResolution.get_phi(),
138  addJERcorrFactor * pfJetResolution.get_sigma_e(),
139  pfJetResolution.get_sigma_tan());
140  }
141  }
142  return pfJetResolution;
143  } else if (dynamic_cast<const reco::PFCandidate*>(particle) != nullptr) {
144  const reco::PFCandidate* pfCandidate = dynamic_cast<const reco::PFCandidate*>(particle);
145  //std::cout << "pfCandidate: pt = " << pt << ", eta = " << eta << ", phi = " << phi << std::endl;
146  return pfMEtResolution_->evalPF(pfCandidate);
147  } else
148  throw cms::Exception("addPFMEtSignObjects")
149  << "Invalid type of particle:"
150  << " valid types = { reco::GsfElectron/pat::Electron, reco::Photon/pat::Photon, reco::Muon/pat::Muon, "
151  "reco::PFTau/pat::Tau,"
152  << " reco::PFJet/pat::Jet, reco::PFCandidate } !!\n";
153  }
154 
155  template <typename T>
156  std::vector<metsig::SigInputObj> compResolution(const std::list<T*>& particles) const {
157  LogDebug("compResolution") << " particles: entries = " << particles.size() << std::endl;
158 
159  std::vector<metsig::SigInputObj> pfMEtSignObjects;
160  addPFMEtSignObjects(pfMEtSignObjects, particles);
161 
162  return pfMEtSignObjects;
163  }
164 
165  reco::METCovMatrix operator()(const std::vector<metsig::SigInputObj>&) const;
166 
167  template <typename T>
168  reco::METCovMatrix operator()(const std::list<T*>& particles) const {
169  std::vector<metsig::SigInputObj> pfMEtSignObjects = compResolution(particles);
170  return this->operator()(pfMEtSignObjects);
171  }
172 
173 protected:
174  template <typename T>
175  void addPFMEtSignObjects(std::vector<metsig::SigInputObj>& metSignObjects, const std::list<T*>& particles) const {
176  for (typename std::list<T*>::const_iterator particle = particles.begin(); particle != particles.end(); ++particle) {
177  metSignObjects.push_back(this->compResolution(*particle));
178  }
179  }
180 
181 private:
183 
184  // CV: look-up table for additional jet energy resolution corrections
185  // not included in (PF)MEt significance algorithm yet
186  // (cf. CMS AN-11/400 vs. CMS AN-11/330)
187  TFile* inputFile_;
188  TH2* lut_;
189 
191 };
192 
193 #endif
bool isAvailable() const
Definition: Ref.h:537
value_type const * get() const
Definition: RefToBase.h:209
const reco::JetBaseRef & pfJetRef() const
Definition: Tau.h:127
double get_sigma_tan() const
Definition: SigInputObj.h:40
metsig::SignAlgoResolutions * pfMEtResolution_
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
PFMEtSignInterfaceBase(const edm::ParameterSet &)
double get_phi() const
Definition: SigInputObj.h:38
const PFSpecific & pfSpecific() const
retrieve the pf specific part of the jet
Definition: Jet.h:293
ROOT::Math::SMatrix< double, 2 > METCovMatrix
Definition: MET.h:39
metsig::SigInputObj evalPF(const reco::PFCandidate *candidate) const
void addPFMEtSignObjects(std::vector< metsig::SigInputObj > &metSignObjects, const std::list< T * > &particles) const
const Point & vertex() const override
vertex position (overwritten by PF...)
virtual Constituents getJetConstituents() const
list of constituents
TrackRef track() const override
reference to a Track
Definition: Muon.h:46
double eval(const resolutionType &type, const resolutionFunc &func, const double &et, const double &phi, const double &eta, const double &p) const
assert(be >=bs)
Jets made from PFObjects.
Definition: PFJet.h:20
const LorentzVector & p4() const final
four-momentum Lorentz vector
void set(const std::string &m_type, const double &m_energy, const double &m_phi, const double &m_sigma_e, const double &m_sigma_tan)
Definition: SigInputObj.h:42
metsig::SigInputObj evalPFJet(const reco::Jet *jet) const
std::vector< metsig::SigInputObj > compResolution(const std::list< T * > &particles) const
const JetBaseRef & jetRef() const
Definition: PFTau.cc:56
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:754
double phiError() const
error on phi
Definition: TrackBase.h:766
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:232
Analysis-level tau class.
Definition: Tau.h:53
bool isPFJet() const
check to see if the jet is a reco::PFJet
Definition: Jet.h:275
reco::TrackRef track() const override
reference to Track reconstructed in the tracker only (reimplemented from reco::Muon) ...
double get_energy() const
Definition: SigInputObj.h:37
Analysis-level calorimeter jet class.
Definition: Jet.h:77
reco::METCovMatrix operator()(const std::list< T * > &particles) const
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
reco::METCovMatrix operator()(const std::vector< metsig::SigInputObj > &) const
std::string get_type() const
Definition: SigInputObj.h:36
metsig::SigInputObj compResolution(const T *particle) const
long double T
Analysis-level muon class.
Definition: Muon.h:51
double get_sigma_e() const
Definition: SigInputObj.h:39
#define LogDebug(id)