CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Member Functions | Protected Member Functions | Protected Attributes | Friends
mvaMEtUtilities Class Reference

#include <mvaMEtUtilities.h>

Classes

struct  JetInfo
 
struct  leptonInfo
 
struct  pfCandInfo
 

Public Member Functions

std::vector< JetInfocleanJets (const std::vector< JetInfo > &, const std::vector< leptonInfo > &, double, double)
 
std::vector< pfCandInfocleanPFCands (const std::vector< pfCandInfo > &, const std::vector< leptonInfo > &, double, bool)
 
CommonMETData computeJetSum_neutral (const std::vector< JetInfo > &, bool)
 
CommonMETData computeNegNoPURecoil (const CommonMETData &, const std::vector< pfCandInfo > &, const std::vector< JetInfo > &, double)
 
CommonMETData computeNegPFRecoil (const CommonMETData &, const std::vector< pfCandInfo > &, double)
 
CommonMETData computeNegPUCRecoil (const CommonMETData &, const std::vector< pfCandInfo > &, const std::vector< JetInfo > &, double)
 
CommonMETData computeNegTrackRecoil (const CommonMETData &, const std::vector< pfCandInfo > &, double)
 
CommonMETData computePFCandSum (const std::vector< pfCandInfo > &, double, int)
 
CommonMETData computePUMEt (const std::vector< pfCandInfo > &, const std::vector< JetInfo > &, double)
 
CommonMETData computeSumLeptons (const std::vector< leptonInfo > &leptons, bool iCharged)
 
void finalize (CommonMETData &metData)
 
reco::Candidate::LorentzVector leadJetP4 (const std::vector< JetInfo > &)
 
 mvaMEtUtilities (const edm::ParameterSet &cfg)
 
unsigned numJetsAboveThreshold (const std::vector< JetInfo > &, double)
 
bool passesMVA (const reco::Candidate::LorentzVector &, double)
 
reco::Candidate::LorentzVector subleadJetP4 (const std::vector< JetInfo > &)
 
virtual ~mvaMEtUtilities ()
 

Protected Member Functions

reco::Candidate::LorentzVector jetP4 (const std::vector< JetInfo > &, unsigned)
 

Protected Attributes

double mvaCut_ [3][4][4]
 

Friends

bool operator< (const JetInfo &, const JetInfo &)
 

Detailed Description

Definition at line 12 of file mvaMEtUtilities.h.

Constructor & Destructor Documentation

mvaMEtUtilities::mvaMEtUtilities ( const edm::ParameterSet cfg)

Definition at line 8 of file mvaMEtUtilities.cc.

References mvaCut_.

9 {
10  // jet id
11  /* ===> this code parses an xml it uses parameter set
12  edm::ParameterSet jetConfig = iConfig.getParameter<edm::ParameterSet>("JetIdParams");
13  for(int i0 = 0; i0 < 3; i0++) {
14  std::string lCutType = "Tight";
15  if(i0 == PileupJetIdentifier::kMedium) lCutType = "Medium";
16  if(i0 == PileupJetIdentifier::kLoose) lCutType = "Loose";
17  std::vector<double> pt010 = jetConfig.getParameter<std::vector<double> >(("Pt010_" +lCutType).c_str());
18  std::vector<double> pt1020 = jetConfig.getParameter<std::vector<double> >(("Pt1020_"+lCutType).c_str());
19  std::vector<double> pt2030 = jetConfig.getParameter<std::vector<double> >(("Pt2030_"+lCutType).c_str());
20  std::vector<double> pt3050 = jetConfig.getParameter<std::vector<double> >(("Pt3050_"+lCutType).c_str());
21  for(int i1 = 0; i1 < 4; i1++) mvacut_[i0][0][i1] = pt010 [i1];
22  for(int i1 = 0; i1 < 4; i1++) mvacut_[i0][1][i1] = pt1020[i1];
23  for(int i1 = 0; i1 < 4; i1++) mvacut_[i0][2][i1] = pt2030[i1];
24  for(int i1 = 0; i1 < 4; i1++) mvacut_[i0][3][i1] = pt3050[i1];
25  }
26  */
27  //Tight Id => not used
28  mvaCut_[0][0][0] = 0.5; mvaCut_[0][0][1] = 0.6; mvaCut_[0][0][2] = 0.6; mvaCut_[0][0][3] = 0.9;
29  mvaCut_[0][1][0] = -0.2; mvaCut_[0][1][1] = 0.2; mvaCut_[0][1][2] = 0.2; mvaCut_[0][1][3] = 0.6;
30  mvaCut_[0][2][0] = 0.3; mvaCut_[0][2][1] = 0.4; mvaCut_[0][2][2] = 0.7; mvaCut_[0][2][3] = 0.8;
31  mvaCut_[0][3][0] = 0.5; mvaCut_[0][3][1] = 0.4; mvaCut_[0][3][2] = 0.8; mvaCut_[0][3][3] = 0.9;
32  //Medium id => not used
33  mvaCut_[1][0][0] = 0.2; mvaCut_[1][0][1] = 0.4; mvaCut_[1][0][2] = 0.2; mvaCut_[1][0][3] = 0.6;
34  mvaCut_[1][1][0] = -0.3; mvaCut_[1][1][1] = 0. ; mvaCut_[1][1][2] = 0. ; mvaCut_[1][1][3] = 0.5;
35  mvaCut_[1][2][0] = 0.2; mvaCut_[1][2][1] = 0.2; mvaCut_[1][2][2] = 0.5; mvaCut_[1][2][3] = 0.7;
36  mvaCut_[1][3][0] = 0.3; mvaCut_[1][3][1] = 0.2; mvaCut_[1][3][2] = 0.7; mvaCut_[1][3][3] = 0.8;
37  //Met Id => used
38  mvaCut_[2][0][0] = -0.2; mvaCut_[2][0][1] = -0.3; mvaCut_[2][0][2] = -0.5; mvaCut_[2][0][3] = -0.5;
39  mvaCut_[2][1][0] = -0.2; mvaCut_[2][1][1] = -0.2; mvaCut_[2][1][2] = -0.5; mvaCut_[2][1][3] = -0.3;
40  mvaCut_[2][2][0] = -0.2; mvaCut_[2][2][1] = -0.2; mvaCut_[2][2][2] = -0.2; mvaCut_[2][2][3] = 0.1;
41  mvaCut_[2][3][0] = -0.2; mvaCut_[2][3][1] = -0.2; mvaCut_[2][3][2] = 0. ; mvaCut_[2][3][3] = 0.2;
42 }
double mvaCut_[3][4][4]
mvaMEtUtilities::~mvaMEtUtilities ( )
virtual

Definition at line 44 of file mvaMEtUtilities.cc.

45 {
46 // nothing to be done yet...
47 }

Member Function Documentation

std::vector< mvaMEtUtilities::JetInfo > mvaMEtUtilities::cleanJets ( const std::vector< JetInfo > &  jets,
const std::vector< leptonInfo > &  leptons,
double  ptThreshold,
double  dRmatch 
)

Definition at line 98 of file mvaMEtUtilities.cc.

References reco::deltaR2(), and metsig::jet.

Referenced by PFMETAlgorithmMVA::setInput().

101 {
102 
103  double dR2match = dRmatch*dRmatch;
104  std::vector<JetInfo> retVal;
105  for ( std::vector<JetInfo>::const_iterator jet = jets.begin();
106  jet != jets.end(); ++jet ) {
107  bool isOverlap = false;
108  for ( std::vector<leptonInfo>::const_iterator lepton = leptons.begin();
109  lepton != leptons.end(); ++lepton ) {
110  if ( deltaR2(jet->p4_, lepton->p4_) < dR2match ) isOverlap = true;
111  }
112  if ( jet->p4_.pt() > ptThreshold && !isOverlap ) retVal.push_back(*jet);
113  }
114  return retVal;
115 }
vector< PseudoJet > jets
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
std::vector< mvaMEtUtilities::pfCandInfo > mvaMEtUtilities::cleanPFCands ( const std::vector< pfCandInfo > &  pfCandidates,
const std::vector< leptonInfo > &  leptons,
double  dRmatch,
bool  invert 
)

Definition at line 117 of file mvaMEtUtilities.cc.

References reco::deltaR2().

120 {
121 
122  double dR2match = dRmatch*dRmatch;
123  std::vector<pfCandInfo> retVal;
124  for ( std::vector<pfCandInfo>::const_iterator pfCandidate = pfCandidates.begin();
125  pfCandidate != pfCandidates.end(); ++pfCandidate ) {
126  bool isOverlap = false;
127  for ( std::vector<leptonInfo>::const_iterator lepton = leptons.begin();
128  lepton != leptons.end(); ++lepton ) {
129  if ( deltaR2(pfCandidate->p4_, lepton->p4_) < dR2match ) isOverlap = true;
130  }
131  if ( (!isOverlap && !invert) || (isOverlap && invert) ) retVal.push_back(*pfCandidate);
132  }
133  return retVal;
134 }
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
CommonMETData mvaMEtUtilities::computeJetSum_neutral ( const std::vector< JetInfo > &  jets,
bool  mvaPassFlag 
)

Definition at line 186 of file mvaMEtUtilities.cc.

References finalize(), metsig::jet, CommonMETData::mex, CommonMETData::mey, passesMVA(), and CommonMETData::sumet.

Referenced by computeNegNoPURecoil(), computeNegPUCRecoil(), and computePUMEt().

187 {
188  CommonMETData retVal;
189  retVal.mex = 0.;
190  retVal.mey = 0.;
191  retVal.sumet = 0.;
192  for ( std::vector<JetInfo>::const_iterator jet = jets.begin();
193  jet != jets.end(); ++jet ) {
194  bool passesMVAjetId = passesMVA(jet->p4_, jet->mva_);
195  if ( passesMVAjetId && !mvaPassFlag ) continue;
196  if ( !passesMVAjetId && mvaPassFlag ) continue;
197  retVal.mex += jet->p4_.px()*jet->neutralEnFrac_;
198  retVal.mey += jet->p4_.py()*jet->neutralEnFrac_;
199  retVal.sumet += jet->p4_.pt()*jet->neutralEnFrac_;
200  }
201  finalize(retVal);
202  return retVal;
203 }
vector< PseudoJet > jets
Structure containing data common to all types of MET.
Definition: CommonMETData.h:12
void finalize(CommonMETData &metData)
bool passesMVA(const reco::Candidate::LorentzVector &, double)
CommonMETData mvaMEtUtilities::computeNegNoPURecoil ( const CommonMETData leptons,
const std::vector< pfCandInfo > &  pfCandidates,
const std::vector< JetInfo > &  jets,
double  dZcut 
)

Definition at line 245 of file mvaMEtUtilities.cc.

References computeJetSum_neutral(), computePFCandSum(), finalize(), CommonMETData::mex, CommonMETData::mey, and CommonMETData::sumet.

Referenced by PFMETAlgorithmMVA::setInput().

248 {
249  CommonMETData retVal;
250  retVal.mex = 0.;
251  retVal.mey = 0.;
252  retVal.sumet = 0.;
253  CommonMETData trackSumNoPU = computePFCandSum(pfCandidates, dZcut, 0);
254  CommonMETData jetSumNoPU_neutral = computeJetSum_neutral(jets, true);
255  retVal.mex = -(trackSumNoPU.mex + jetSumNoPU_neutral.mex) + leptons.mex;
256  retVal.mey = -(trackSumNoPU.mey + jetSumNoPU_neutral.mey) + leptons.mey;
257  retVal.sumet = trackSumNoPU.sumet + jetSumNoPU_neutral.sumet - leptons.sumet;
258  finalize(retVal);
259  return retVal;
260 }
CommonMETData computeJetSum_neutral(const std::vector< JetInfo > &, bool)
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
vector< PseudoJet > jets
Structure containing data common to all types of MET.
Definition: CommonMETData.h:12
void finalize(CommonMETData &metData)
CommonMETData computePFCandSum(const std::vector< pfCandInfo > &, double, int)
CommonMETData mvaMEtUtilities::computeNegPFRecoil ( const CommonMETData leptons,
const std::vector< pfCandInfo > &  pfCandidates,
double  dZcut 
)

Definition at line 221 of file mvaMEtUtilities.cc.

References computePFCandSum(), finalize(), CommonMETData::mex, CommonMETData::mey, and CommonMETData::sumet.

Referenced by PFMETAlgorithmMVA::setInput().

223 {
224  CommonMETData retVal;
226  retVal.mex = -pfCandSum.mex + leptons.mex;
227  retVal.mey = -pfCandSum.mey + leptons.mey;
228  retVal.sumet = pfCandSum.sumet - leptons.sumet;
229  finalize(retVal);
230  return retVal;
231 }
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
Structure containing data common to all types of MET.
Definition: CommonMETData.h:12
void finalize(CommonMETData &metData)
CommonMETData computePFCandSum(const std::vector< pfCandInfo > &, double, int)
CommonMETData mvaMEtUtilities::computeNegPUCRecoil ( const CommonMETData leptons,
const std::vector< pfCandInfo > &  pfCandidates,
const std::vector< JetInfo > &  jets,
double  dZcut 
)

Definition at line 262 of file mvaMEtUtilities.cc.

References computeJetSum_neutral(), computePFCandSum(), finalize(), CommonMETData::mex, CommonMETData::mey, and CommonMETData::sumet.

Referenced by PFMETAlgorithmMVA::setInput().

265 {
266  CommonMETData retVal;
267  retVal.mex = 0.;
268  retVal.mey = 0.;
269  retVal.sumet = 0.;
271  CommonMETData trackSumNoPU = computePFCandSum(pfCandidates, dZcut, 1);
272  CommonMETData jetSumPU_neutral = computeJetSum_neutral(jets, false);
273  retVal.mex = -(pfCandSum.mex - (trackSumNoPU.mex + jetSumPU_neutral.mex)) + leptons.mex;
274  retVal.mey = -(pfCandSum.mey - (trackSumNoPU.mey + jetSumPU_neutral.mey)) + leptons.mey;
275  retVal.sumet = pfCandSum.sumet - (trackSumNoPU.sumet + jetSumPU_neutral.sumet) - leptons.sumet;
276  finalize(retVal);
277  return retVal;
278 }
CommonMETData computeJetSum_neutral(const std::vector< JetInfo > &, bool)
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
vector< PseudoJet > jets
Structure containing data common to all types of MET.
Definition: CommonMETData.h:12
void finalize(CommonMETData &metData)
CommonMETData computePFCandSum(const std::vector< pfCandInfo > &, double, int)
CommonMETData mvaMEtUtilities::computeNegTrackRecoil ( const CommonMETData leptons,
const std::vector< pfCandInfo > &  pfCandidates,
double  dZcut 
)

Definition at line 233 of file mvaMEtUtilities.cc.

References computePFCandSum(), finalize(), CommonMETData::mex, CommonMETData::mey, and CommonMETData::sumet.

Referenced by PFMETAlgorithmMVA::setInput().

235 {
236  CommonMETData retVal;
238  retVal.mex = -trackSum.mex + leptons.mex;
239  retVal.mey = -trackSum.mey + leptons.mey;
240  retVal.sumet = trackSum.sumet - leptons.sumet;
241  finalize(retVal);
242  return retVal;
243 }
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
Structure containing data common to all types of MET.
Definition: CommonMETData.h:12
void finalize(CommonMETData &metData)
CommonMETData computePFCandSum(const std::vector< pfCandInfo > &, double, int)
CommonMETData mvaMEtUtilities::computePFCandSum ( const std::vector< pfCandInfo > &  pfCandidates,
double  dZmax,
int  dZflag 
)

Definition at line 142 of file mvaMEtUtilities.cc.

References finalize(), CommonMETData::mex, CommonMETData::mey, and CommonMETData::sumet.

Referenced by computeNegNoPURecoil(), computeNegPFRecoil(), computeNegPUCRecoil(), computeNegTrackRecoil(), and computePUMEt().

143 {
144  // dZcut
145  // maximum distance within which tracks are considered to be associated to hard scatter vertex
146  // dZflag
147  // 0 : select charged PFCandidates originating from hard scatter vertex
148  // 1 : select charged PFCandidates originating from pile-up vertices
149  // 2 : select all PFCandidates
150  CommonMETData retVal;
151  retVal.mex = 0.;
152  retVal.mey = 0.;
153  retVal.sumet = 0.;
154  for ( std::vector<pfCandInfo>::const_iterator pfCandidate = pfCandidates.begin();
155  pfCandidate != pfCandidates.end(); ++pfCandidate ) {
156  if ( pfCandidate->dZ_ < 0. && dZflag != 2 ) continue;
157  if ( pfCandidate->dZ_ > dZmax && dZflag == 0 ) continue;
158  if ( pfCandidate->dZ_ < dZmax && dZflag == 1 ) continue;
159  retVal.mex += pfCandidate->p4_.px();
160  retVal.mey += pfCandidate->p4_.py();
161  retVal.sumet += pfCandidate->p4_.pt();
162  }
163  finalize(retVal);
164  return retVal;
165 }
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
Structure containing data common to all types of MET.
Definition: CommonMETData.h:12
void finalize(CommonMETData &metData)
CommonMETData mvaMEtUtilities::computePUMEt ( const std::vector< pfCandInfo > &  pfCandidates,
const std::vector< JetInfo > &  jets,
double  dZcut 
)

Definition at line 205 of file mvaMEtUtilities.cc.

References computeJetSum_neutral(), computePFCandSum(), finalize(), CommonMETData::mex, CommonMETData::mey, and CommonMETData::sumet.

Referenced by PFMETAlgorithmMVA::setInput().

207 {
208  CommonMETData retVal;
209  retVal.mex = 0.;
210  retVal.mey = 0.;
211  retVal.sumet = 0.;
213  CommonMETData jetSumPU_neutral = computeJetSum_neutral(jets, false);
214  retVal.mex = -(trackSumPU.mex + jetSumPU_neutral.mex);
215  retVal.mey = -(trackSumPU.mey + jetSumPU_neutral.mey);
216  retVal.sumet = trackSumPU.sumet + jetSumPU_neutral.sumet;
217  finalize(retVal);
218  return retVal;
219 }
CommonMETData computeJetSum_neutral(const std::vector< JetInfo > &, bool)
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
vector< PseudoJet > jets
Structure containing data common to all types of MET.
Definition: CommonMETData.h:12
void finalize(CommonMETData &metData)
CommonMETData computePFCandSum(const std::vector< pfCandInfo > &, double, int)
CommonMETData mvaMEtUtilities::computeSumLeptons ( const std::vector< leptonInfo > &  leptons,
bool  iCharged 
)

Definition at line 167 of file mvaMEtUtilities.cc.

References finalize(), CommonMETData::mex, CommonMETData::mey, and CommonMETData::sumet.

Referenced by PFMETAlgorithmMVA::setInput().

168 {
169  //Computes vectorial sum of charged(iCharged == true) or all(iCharged == false) components
170  CommonMETData retVal;
171  retVal.mex = 0.;
172  retVal.mey = 0.;
173  retVal.sumet = 0.;
174  for ( std::vector<leptonInfo>::const_iterator lepton = leptons.begin();
175  lepton != leptons.end(); ++lepton ) {
176  double pChargedFrac = 1;
177  if(iCharged) pChargedFrac = lepton->chargedFrac_;
178  retVal.mex += lepton->p4_.px()*pChargedFrac;
179  retVal.mey += lepton->p4_.py()*pChargedFrac;
180  retVal.sumet += lepton->p4_.pt()*pChargedFrac;
181  }
182  finalize(retVal);
183  return retVal;
184 }
Structure containing data common to all types of MET.
Definition: CommonMETData.h:12
void finalize(CommonMETData &metData)
void mvaMEtUtilities::finalize ( CommonMETData metData)

Definition at line 135 of file mvaMEtUtilities.cc.

References CommonMETData::met, CommonMETData::mex, CommonMETData::mey, CommonMETData::mez, CommonMETData::phi, and mathSSE::sqrt().

Referenced by computeJetSum_neutral(), computeNegNoPURecoil(), computeNegPFRecoil(), computeNegPUCRecoil(), computeNegTrackRecoil(), computePFCandSum(), computePUMEt(), and computeSumLeptons().

136 {
137  metData.met = sqrt(metData.mex*metData.mex + metData.mey*metData.mey);
138  metData.mez = 0.;
139  metData.phi = atan2(metData.mey, metData.mex);
140 }
T sqrt(T t)
Definition: SSEVec.h:48
reco::Candidate::LorentzVector mvaMEtUtilities::jetP4 ( const std::vector< JetInfo > &  jets,
unsigned  idx 
)
protected

Definition at line 79 of file mvaMEtUtilities.cc.

References customizeTrackingMonitorSeedNumber::idx, fwrapper::jets, and python.multivaluedict::sort().

Referenced by leadJetP4(), and subleadJetP4().

80 {
81  reco::Candidate::LorentzVector retVal(0.,0.,0.,0.);
82  if ( idx < jets.size() ) {
83  std::vector<JetInfo> jets_sorted = jets;
84  std::sort(jets_sorted.begin(), jets_sorted.end());
85  retVal = jets_sorted[idx].p4_;
86  }
87  return retVal;
88 }
vector< PseudoJet > jets
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:41
reco::Candidate::LorentzVector mvaMEtUtilities::leadJetP4 ( const std::vector< JetInfo > &  jets)

Definition at line 64 of file mvaMEtUtilities.cc.

References jetP4().

Referenced by PFMETAlgorithmMVA::setInput().

65 {
66  return jetP4(jets, 0);
67 }
reco::Candidate::LorentzVector jetP4(const std::vector< JetInfo > &, unsigned)
vector< PseudoJet > jets
unsigned mvaMEtUtilities::numJetsAboveThreshold ( const std::vector< JetInfo > &  jets,
double  ptThreshold 
)

Definition at line 89 of file mvaMEtUtilities.cc.

References metsig::jet.

Referenced by PFMETAlgorithmMVA::setInput().

90 {
91  unsigned retVal = 0;
92  for ( std::vector<JetInfo>::const_iterator jet = jets.begin();
93  jet != jets.end(); ++jet ) {
94  if ( jet->p4_.pt() > ptThreshold ) ++retVal;
95  }
96  return retVal;
97 }
vector< PseudoJet > jets
bool mvaMEtUtilities::passesMVA ( const reco::Candidate::LorentzVector jetP4,
double  mvaJetId 
)

Definition at line 49 of file mvaMEtUtilities.cc.

References funct::abs(), and mvaCut_.

Referenced by computeJetSum_neutral().

50 {
51  int ptBin = 0;
52  if ( jetP4.pt() > 10. && jetP4.pt() < 20. ) ptBin = 1;
53  if ( jetP4.pt() > 20. && jetP4.pt() < 30. ) ptBin = 2;
54  if ( jetP4.pt() > 30. ) ptBin = 3;
55 
56  int etaBin = 0;
57  if ( std::abs(jetP4.eta()) > 2.5 && std::abs(jetP4.eta()) < 2.75) etaBin = 1;
58  if ( std::abs(jetP4.eta()) > 2.75 && std::abs(jetP4.eta()) < 3.0 ) etaBin = 2;
59  if ( std::abs(jetP4.eta()) > 3.0 && std::abs(jetP4.eta()) < 5.0 ) etaBin = 3;
60 
61  return ( mvaJetId > mvaCut_[2][ptBin][etaBin] );
62 }
reco::Candidate::LorentzVector jetP4(const std::vector< JetInfo > &, unsigned)
double mvaCut_[3][4][4]
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
reco::Candidate::LorentzVector mvaMEtUtilities::subleadJetP4 ( const std::vector< JetInfo > &  jets)

Definition at line 69 of file mvaMEtUtilities.cc.

References jetP4().

Referenced by PFMETAlgorithmMVA::setInput().

70 {
71  return jetP4(jets, 1);
72 }
reco::Candidate::LorentzVector jetP4(const std::vector< JetInfo > &, unsigned)
vector< PseudoJet > jets

Friends And Related Function Documentation

bool operator< ( const JetInfo jet1,
const JetInfo jet2 
)
friend

Definition at line 74 of file mvaMEtUtilities.cc.

75 {
76  return jet1.p4_.pt() > jet2.p4_.pt();
77 }

Member Data Documentation

double mvaMEtUtilities::mvaCut_[3][4][4]
protected

Definition at line 84 of file mvaMEtUtilities.h.

Referenced by mvaMEtUtilities(), and passesMVA().