CMS 3D CMS Logo

CovarianceMatrix.cc
Go to the documentation of this file.
2 
5 
6 CovarianceMatrix::CovarianceMatrix(const std::vector<edm::ParameterSet>& udscResolutions,
7  const std::vector<edm::ParameterSet>& bResolutions,
8  const std::vector<double>& jetEnergyResolutionScaleFactors,
9  const std::vector<double>& jetEnergyResolutionEtaBinning)
10  : jetEnergyResolutionScaleFactors_(jetEnergyResolutionScaleFactors),
11  jetEnergyResolutionEtaBinning_(jetEnergyResolutionEtaBinning) {
12  for (std::vector<edm::ParameterSet>::const_iterator iSet = udscResolutions.begin(); iSet != udscResolutions.end();
13  ++iSet) {
14  if (iSet->exists("bin"))
15  binsUdsc_.push_back(iSet->getParameter<std::string>("bin"));
16  else if (udscResolutions.size() == 1)
17  binsUdsc_.push_back("");
18  else
19  throw cms::Exception("Configuration") << "Parameter 'bin' is needed if more than one PSet is specified!\n";
20 
21  funcEtUdsc_.push_back(iSet->getParameter<std::string>("et"));
22  funcEtaUdsc_.push_back(iSet->getParameter<std::string>("eta"));
23  funcPhiUdsc_.push_back(iSet->getParameter<std::string>("phi"));
24  }
25  for (std::vector<edm::ParameterSet>::const_iterator iSet = bResolutions.begin(); iSet != bResolutions.end(); ++iSet) {
26  if (iSet->exists("bin"))
27  binsB_.push_back(iSet->getParameter<std::string>("bin"));
28  else if (bResolutions.size() == 1)
29  binsB_.push_back("");
30  else
31  throw cms::Exception("Configuration") << "Parameter 'bin' is needed if more than one PSet is specified!\n";
32 
33  funcEtB_.push_back(iSet->getParameter<std::string>("et"));
34  funcEtaB_.push_back(iSet->getParameter<std::string>("eta"));
35  funcPhiB_.push_back(iSet->getParameter<std::string>("phi"));
36  }
37 }
38 
39 CovarianceMatrix::CovarianceMatrix(const std::vector<edm::ParameterSet>& udscResolutions,
40  const std::vector<edm::ParameterSet>& bResolutions,
41  const std::vector<edm::ParameterSet>& lepResolutions,
42  const std::vector<edm::ParameterSet>& metResolutions,
43  const std::vector<double>& jetEnergyResolutionScaleFactors,
44  const std::vector<double>& jetEnergyResolutionEtaBinning)
45  : jetEnergyResolutionScaleFactors_(jetEnergyResolutionScaleFactors),
46  jetEnergyResolutionEtaBinning_(jetEnergyResolutionEtaBinning) {
48  throw cms::Exception("Configuration") << "The number of scale factors does not fit to the number of eta bins!\n";
49  for (unsigned int i = 0; i < jetEnergyResolutionEtaBinning_.size(); i++)
51  throw cms::Exception("Configuration") << "eta binning in absolut values required!\n";
52 
53  for (std::vector<edm::ParameterSet>::const_iterator iSet = udscResolutions.begin(); iSet != udscResolutions.end();
54  ++iSet) {
55  if (iSet->exists("bin"))
56  binsUdsc_.push_back(iSet->getParameter<std::string>("bin"));
57  else if (udscResolutions.size() == 1)
58  binsUdsc_.push_back("");
59  else
60  throw cms::Exception("Configuration") << "Parameter 'bin' is needed if more than one PSet is specified!\n";
61 
62  funcEtUdsc_.push_back(iSet->getParameter<std::string>("et"));
63  funcEtaUdsc_.push_back(iSet->getParameter<std::string>("eta"));
64  funcPhiUdsc_.push_back(iSet->getParameter<std::string>("phi"));
65  }
66  for (std::vector<edm::ParameterSet>::const_iterator iSet = bResolutions.begin(); iSet != bResolutions.end(); ++iSet) {
67  if (iSet->exists("bin"))
68  binsB_.push_back(iSet->getParameter<std::string>("bin"));
69  else if (bResolutions.size() == 1)
70  binsB_.push_back("");
71  else
72  throw cms::Exception("Configuration") << "Parameter 'bin' is needed if more than one PSet is specified!\n";
73 
74  funcEtB_.push_back(iSet->getParameter<std::string>("et"));
75  funcEtaB_.push_back(iSet->getParameter<std::string>("eta"));
76  funcPhiB_.push_back(iSet->getParameter<std::string>("phi"));
77  }
78  for (std::vector<edm::ParameterSet>::const_iterator iSet = lepResolutions.begin(); iSet != lepResolutions.end();
79  ++iSet) {
80  if (iSet->exists("bin"))
81  binsLep_.push_back(iSet->getParameter<std::string>("bin"));
82  else if (lepResolutions.size() == 1)
83  binsLep_.push_back("");
84  else
85  throw cms::Exception("Configuration") << "Parameter 'bin' is needed if more than one PSet is specified!\n";
86 
87  funcEtLep_.push_back(iSet->getParameter<std::string>("et"));
88  funcEtaLep_.push_back(iSet->getParameter<std::string>("eta"));
89  funcPhiLep_.push_back(iSet->getParameter<std::string>("phi"));
90  }
91  for (std::vector<edm::ParameterSet>::const_iterator iSet = metResolutions.begin(); iSet != metResolutions.end();
92  ++iSet) {
93  if (iSet->exists("bin"))
94  binsMet_.push_back(iSet->getParameter<std::string>("bin"));
95  else if (metResolutions.size() == 1)
96  binsMet_.push_back("");
97  else
98  throw cms::Exception("Configuration") << "Parameter 'bin' is needed if more than one PSet is specified!\n";
99 
100  funcEtMet_.push_back(iSet->getParameter<std::string>("et"));
101  funcEtaMet_.push_back(iSet->getParameter<std::string>("eta"));
102  funcPhiMet_.push_back(iSet->getParameter<std::string>("phi"));
103  }
104 }
105 
106 double CovarianceMatrix::getResolution(const TLorentzVector& object,
107  const ObjectType objType,
108  const std::string& whichResolution) {
109  std::vector<std::string>*bins_, *funcEt_, *funcEta_, *funcPhi_;
110 
111  switch (objType) {
112  case kUdscJet:
113  bins_ = &binsUdsc_;
114  funcEt_ = &funcEtUdsc_;
115  funcEta_ = &funcEtaUdsc_;
116  funcPhi_ = &funcPhiUdsc_;
117  break;
118  case kBJet:
119  bins_ = &binsB_;
120  funcEt_ = &funcEtB_;
121  funcEta_ = &funcEtaB_;
122  funcPhi_ = &funcPhiB_;
123  break;
124  case kMuon:
125  bins_ = &binsLep_;
126  funcEt_ = &funcEtLep_;
127  funcEta_ = &funcEtaLep_;
128  funcPhi_ = &funcPhiLep_;
129  break;
130  case kElectron:
131  bins_ = &binsLep_;
132  funcEt_ = &funcEtLep_;
133  funcEta_ = &funcEtaLep_;
134  funcPhi_ = &funcPhiLep_;
135  break;
136  case kMet:
137  bins_ = &binsMet_;
138  funcEt_ = &funcEtMet_;
139  funcEta_ = &funcEtaMet_;
140  funcPhi_ = &funcPhiMet_;
141  break;
142  default:
143  throw cms::Exception("UnsupportedObject") << "The object given is not supported!\n";
144  }
145 
146  int selectedBin = -1;
147  reco::LeafCandidate candidate;
148  for (unsigned int i = 0; i < bins_->size(); ++i) {
150  candidate = reco::LeafCandidate(
151  0, reco::LeafCandidate::LorentzVector(object.Px(), object.Py(), object.Pz(), object.Energy()));
152  if (select_(candidate)) {
153  selectedBin = i;
154  break;
155  }
156  }
157  double res = 0;
158  if (selectedBin >= 0) {
159  if (whichResolution == "et")
160  res = StringObjectFunction<reco::LeafCandidate>(funcEt_->at(selectedBin)).operator()(candidate);
161  else if (whichResolution == "eta")
162  res = StringObjectFunction<reco::LeafCandidate>(funcEta_->at(selectedBin)).operator()(candidate);
163  else if (whichResolution == "phi")
164  res = StringObjectFunction<reco::LeafCandidate>(funcPhi_->at(selectedBin)).operator()(candidate);
165  else
166  throw cms::Exception("ProgrammingError") << "Only 'et', 'eta' and 'phi' resolutions supported!\n";
167  }
168  return res;
169 }
170 
171 TMatrixD CovarianceMatrix::setupMatrix(const TLorentzVector& object,
172  const ObjectType objType,
173  const TopKinFitter::Param param) {
174  TMatrixD CovM3(3, 3);
175  CovM3.Zero();
176  TMatrixD CovM4(4, 4);
177  CovM4.Zero();
178  const double pt = object.Pt();
179  const double eta = object.Eta();
180  switch (objType) {
181  case kUdscJet:
182  // if object is a non-b jet
183  {
184  res::HelperJet jetRes;
185  switch (param) {
186  case TopKinFitter::kEMom:
187  CovM4(0, 0) = pow(jetRes.a(pt, eta, res::HelperJet::kUds), 2);
188  CovM4(1, 1) = pow(jetRes.b(pt, eta, res::HelperJet::kUds), 2);
189  CovM4(2, 2) = pow(jetRes.c(pt, eta, res::HelperJet::kUds), 2);
190  CovM4(3, 3) = pow(jetRes.d(pt, eta, res::HelperJet::kUds), 2);
191  return CovM4;
193  if (binsUdsc_.empty()) {
194  CovM3(0, 0) = pow(jetRes.et(pt, eta, res::HelperJet::kUds), 2);
195  CovM3(0, 0) *= pow(getEtaDependentScaleFactor(object), 2);
196  CovM3(1, 1) = pow(jetRes.eta(pt, eta, res::HelperJet::kUds), 2);
197  CovM3(2, 2) = pow(jetRes.phi(pt, eta, res::HelperJet::kUds), 2);
198  } else {
199  CovM3(0, 0) = pow(getResolution(object, objType, "et"), 2);
200  CovM3(0, 0) *= pow(getEtaDependentScaleFactor(object), 2);
201  CovM3(1, 1) = pow(getResolution(object, objType, "eta"), 2);
202  CovM3(2, 2) = pow(getResolution(object, objType, "phi"), 2);
203  }
204  return CovM3;
206  CovM3(0, 0) = pow(jetRes.et(pt, eta, res::HelperJet::kUds), 2);
207  CovM3(0, 0) *= pow(getEtaDependentScaleFactor(object), 2);
208  CovM3(1, 1) = pow(jetRes.theta(pt, eta, res::HelperJet::kUds), 2);
209  CovM3(2, 2) = pow(jetRes.phi(pt, eta, res::HelperJet::kUds), 2);
210  return CovM3;
211  }
212  }
213  break;
214  case kBJet:
215  // if object is a b jet
216  {
217  res::HelperJet jetRes;
218  switch (param) {
219  case TopKinFitter::kEMom:
220  CovM4(0, 0) = pow(jetRes.a(pt, eta, res::HelperJet::kB), 2);
221  CovM4(1, 1) = pow(jetRes.b(pt, eta, res::HelperJet::kB), 2);
222  CovM4(2, 2) = pow(jetRes.c(pt, eta, res::HelperJet::kB), 2);
223  CovM4(3, 3) = pow(jetRes.d(pt, eta, res::HelperJet::kB), 2);
224  return CovM4;
226  if (binsUdsc_.empty()) {
227  CovM3(0, 0) = pow(jetRes.et(pt, eta, res::HelperJet::kB), 2);
228  CovM3(0, 0) *= pow(getEtaDependentScaleFactor(object), 2);
229  CovM3(1, 1) = pow(jetRes.eta(pt, eta, res::HelperJet::kB), 2);
230  CovM3(2, 2) = pow(jetRes.phi(pt, eta, res::HelperJet::kB), 2);
231  } else {
232  CovM3(0, 0) = pow(getResolution(object, objType, "et"), 2);
233  CovM3(0, 0) *= pow(getEtaDependentScaleFactor(object), 2);
234  CovM3(1, 1) = pow(getResolution(object, objType, "eta"), 2);
235  CovM3(2, 2) = pow(getResolution(object, objType, "phi"), 2);
236  }
237  return CovM3;
239  CovM3(0, 0) = pow(jetRes.et(pt, eta, res::HelperJet::kB), 2);
240  CovM3(0, 0) *= pow(getEtaDependentScaleFactor(object), 2);
241  CovM3(1, 1) = pow(jetRes.theta(pt, eta, res::HelperJet::kB), 2);
242  CovM3(2, 2) = pow(jetRes.phi(pt, eta, res::HelperJet::kB), 2);
243  return CovM3;
244  }
245  }
246  break;
247  case kMuon:
248  // if object is a muon
249  {
250  res::HelperMuon muonRes;
251  switch (param) {
252  case TopKinFitter::kEMom:
253  CovM3(0, 0) = pow(muonRes.a(pt, eta), 2);
254  CovM3(1, 1) = pow(muonRes.b(pt, eta), 2);
255  CovM3(2, 2) = pow(muonRes.c(pt, eta), 2);
256  return CovM3;
258  if (binsLep_.empty()) {
259  CovM3(0, 0) = pow(muonRes.et(pt, eta), 2);
260  CovM3(1, 1) = pow(muonRes.eta(pt, eta), 2);
261  CovM3(2, 2) = pow(muonRes.phi(pt, eta), 2);
262  } else {
263  CovM3(0, 0) = pow(getResolution(object, objType, "et"), 2);
264  CovM3(1, 1) = pow(getResolution(object, objType, "eta"), 2);
265  CovM3(2, 2) = pow(getResolution(object, objType, "phi"), 2);
266  }
267  return CovM3;
269  CovM3(0, 0) = pow(muonRes.et(pt, eta), 2);
270  CovM3(1, 1) = pow(muonRes.theta(pt, eta), 2);
271  CovM3(2, 2) = pow(muonRes.phi(pt, eta), 2);
272  return CovM3;
273  }
274  }
275  break;
276  case kElectron: {
277  // if object is an electron
278  res::HelperElectron elecRes;
279  switch (param) {
280  case TopKinFitter::kEMom:
281  CovM3(0, 0) = pow(elecRes.a(pt, eta), 2);
282  CovM3(1, 1) = pow(elecRes.b(pt, eta), 2);
283  CovM3(2, 2) = pow(elecRes.c(pt, eta), 2);
284  return CovM3;
286  if (binsLep_.empty()) {
287  CovM3(0, 0) = pow(elecRes.et(pt, eta), 2);
288  CovM3(1, 1) = pow(elecRes.eta(pt, eta), 2);
289  CovM3(2, 2) = pow(elecRes.phi(pt, eta), 2);
290  } else {
291  CovM3(0, 0) = pow(getResolution(object, objType, "et"), 2);
292  CovM3(1, 1) = pow(getResolution(object, objType, "eta"), 2);
293  CovM3(2, 2) = pow(getResolution(object, objType, "phi"), 2);
294  }
295  return CovM3;
297  CovM3(0, 0) = pow(elecRes.et(pt, eta), 2);
298  CovM3(1, 1) = pow(elecRes.theta(pt, eta), 2);
299  CovM3(2, 2) = pow(elecRes.phi(pt, eta), 2);
300  return CovM3;
301  }
302  } break;
303  case kMet:
304  // if object is met
305  {
306  res::HelperMET metRes;
307  switch (param) {
308  case TopKinFitter::kEMom:
309  CovM3(0, 0) = pow(metRes.a(pt), 2);
310  CovM3(1, 1) = pow(metRes.b(pt), 2);
311  CovM3(2, 2) = pow(metRes.c(pt), 2);
312  return CovM3;
314  if (binsMet_.empty()) {
315  CovM3(0, 0) = pow(metRes.et(pt), 2);
316  CovM3(1, 1) = pow(9999., 2);
317  CovM3(2, 2) = pow(metRes.phi(pt), 2);
318  } else {
319  CovM3(0, 0) = pow(getResolution(object, objType, "et"), 2);
320  CovM3(1, 1) = pow(getResolution(object, objType, "eta"), 2);
321  CovM3(2, 2) = pow(getResolution(object, objType, "phi"), 2);
322  }
323  return CovM3;
325  CovM3(0, 0) = pow(metRes.et(pt), 2);
326  CovM3(1, 1) = pow(9999., 2);
327  CovM3(2, 2) = pow(metRes.phi(pt), 2);
328  return CovM3;
329  }
330  }
331  break;
332  }
333  cms::Exception("Logic") << "Something went wrong while trying to setup a covariance matrix!\n";
334  return CovM4; //should never get here
335 }
336 
337 double CovarianceMatrix::getEtaDependentScaleFactor(const TLorentzVector& object) {
338  double etaDependentScaleFactor = 1.;
339  for (unsigned int i = 0; i < jetEnergyResolutionEtaBinning_.size(); i++) {
341  if (i == jetEnergyResolutionEtaBinning_.size() - 1) {
342  edm::LogWarning("CovarianceMatrix") << "object eta (" << std::abs(object.Eta()) << ") beyond last eta bin ("
343  << jetEnergyResolutionEtaBinning_[i] << ") using scale factor 1.0!";
344  etaDependentScaleFactor = 1.;
345  break;
346  }
347  etaDependentScaleFactor = jetEnergyResolutionScaleFactors_[i];
348  } else
349  break;
350  }
351  return etaDependentScaleFactor;
352 }
TopKinFitter::kEtEtaPhi
Definition: TopKinFitter.h:20
res::HelperElectron::et
double et(double pt, double eta)
Definition: Electron.h:185
CovarianceMatrix::jetEnergyResolutionScaleFactors_
const std::vector< double > jetEnergyResolutionScaleFactors_
scale factors for the jet energy resolution
Definition: CovarianceMatrix.h:73
CovarianceMatrix::funcPhiUdsc_
std::vector< std::string > funcPhiUdsc_
Definition: CovarianceMatrix.h:71
res::HelperMuon
Definition: Muon.h:7
res::HelperMET
Definition: MET.h:7
mps_fire.i
i
Definition: mps_fire.py:355
StringObjectFunction
Definition: StringObjectFunction.h:16
res::HelperMuon::c
double c(double pt, double eta)
Definition: Muon.h:79
EcalCondDBWriter_cfi.Energy
Energy
Definition: EcalCondDBWriter_cfi.py:152
res::HelperJet::c
double c(double pt, double eta, Flavor flav)
Definition: Jet.h:138
TtFullHadKinFitProducer_cfi.udscResolutions
udscResolutions
Definition: TtFullHadKinFitProducer_cfi.py:83
res::HelperMET::b
double b(double pt)
Definition: MET.h:31
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
res::HelperJet
Definition: Jet.h:7
CovarianceMatrix::funcEtMet_
std::vector< std::string > funcEtMet_
Definition: CovarianceMatrix.h:69
CovarianceMatrix::CovarianceMatrix
CovarianceMatrix()
default constructor
Definition: CovarianceMatrix.h:32
CovarianceMatrix::binsLep_
std::vector< std::string > binsLep_
Definition: CovarianceMatrix.h:67
CovarianceMatrix::funcPhiLep_
std::vector< std::string > funcPhiLep_
Definition: CovarianceMatrix.h:71
res::HelperMET::phi
double phi(double pt)
Definition: MET.h:51
res::HelperMET::et
double et(double pt)
Definition: MET.h:56
res::HelperMuon::b
double b(double pt, double eta)
Definition: Muon.h:54
res::HelperJet::phi
double phi(double pt, double eta, Flavor flav)
Definition: Jet.h:282
TopKinFitter::kEMom
Definition: TopKinFitter.h:20
res::HelperElectron::phi
double phi(double pt, double eta)
Definition: Electron.h:160
CovarianceMatrix.h
res::HelperMuon::theta
double theta(double pt, double eta)
Definition: Muon.h:129
res::HelperElectron::eta
double eta(double pt, double eta)
Definition: Electron.h:210
CovarianceMatrix::kElectron
Definition: CovarianceMatrix.h:29
CovarianceMatrix::funcPhiMet_
std::vector< std::string > funcPhiMet_
Definition: CovarianceMatrix.h:71
res::HelperJet::kB
Definition: Jet.h:9
CovarianceMatrix::funcEtaLep_
std::vector< std::string > funcEtaLep_
Definition: CovarianceMatrix.h:70
PVValHelper::eta
Definition: PVValidationHelpers.h:69
CovarianceMatrix::getResolution
double getResolution(const TLorentzVector &object, const ObjectType objType, const std::string &whichResolution="")
get resolution for a given component of an object
Definition: CovarianceMatrix.cc:106
res::HelperMuon::et
double et(double pt, double eta)
Definition: Muon.h:179
CovarianceMatrix::kBJet
Definition: CovarianceMatrix.h:29
res::HelperJet::d
double d(double pt, double eta, Flavor flav)
Definition: Jet.h:186
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
edm::LogWarning
Definition: MessageLogger.h:141
CovarianceMatrix::funcEtaMet_
std::vector< std::string > funcEtaMet_
Definition: CovarianceMatrix.h:70
CovarianceMatrix::kMet
Definition: CovarianceMatrix.h:29
CovarianceMatrix::binsB_
std::vector< std::string > binsB_
Definition: CovarianceMatrix.h:67
electrons_cff.objType
objType
Definition: electrons_cff.py:491
res::HelperMET::c
double c(double pt)
Definition: MET.h:36
CovarianceMatrix::jetEnergyResolutionEtaBinning_
const std::vector< double > jetEnergyResolutionEtaBinning_
Definition: CovarianceMatrix.h:74
res::HelperJet::b
double b(double pt, double eta, Flavor flav)
Definition: Jet.h:90
CovarianceMatrix::kMuon
Definition: CovarianceMatrix.h:29
res::HelperElectron
Definition: Electron.h:7
res::HelperJet::theta
double theta(double pt, double eta, Flavor flav)
Definition: Jet.h:234
CovarianceMatrix::getEtaDependentScaleFactor
double getEtaDependentScaleFactor(const pat::PATObject< T > &object)
get eta dependent smear factor for a PAT object
Definition: CovarianceMatrix.h:161
CovarianceMatrix::funcEtUdsc_
std::vector< std::string > funcEtUdsc_
vectors for the resolution functions
Definition: CovarianceMatrix.h:69
res
Definition: Electron.h:6
TtFullHadKinFitProducer_cfi.jetEnergyResolutionEtaBinning
jetEnergyResolutionEtaBinning
Definition: TtFullHadKinFitProducer_cfi.py:95
TopKinFitter::kEtThetaPhi
Definition: TopKinFitter.h:20
TtFullHadKinFitProducer_cfi.jetEnergyResolutionScaleFactors
jetEnergyResolutionScaleFactors
Definition: TtFullHadKinFitProducer_cfi.py:94
CovarianceMatrix::funcEtB_
std::vector< std::string > funcEtB_
Definition: CovarianceMatrix.h:69
StringCutObjectSelector.h
res::HelperElectron::a
double a(double pt, double eta)
Definition: Electron.h:35
res::HelperMuon::eta
double eta(double pt, double eta)
Definition: Muon.h:204
res::HelperElectron::theta
double theta(double pt, double eta)
Definition: Electron.h:135
CovarianceMatrix::funcPhiB_
std::vector< std::string > funcPhiB_
Definition: CovarianceMatrix.h:71
StringCutObjectSelector
Definition: StringCutObjectSelector.h:16
Exception
Definition: hltDiff.cc:246
TtFullHadKinFitProducer_cfi.bResolutions
bResolutions
Definition: TtFullHadKinFitProducer_cfi.py:84
res::HelperJet::eta
double eta(double pt, double eta, Flavor flav)
Definition: Jet.h:378
CovarianceMatrix::kUdscJet
Definition: CovarianceMatrix.h:29
CovarianceMatrix::binsMet_
std::vector< std::string > binsMet_
Definition: CovarianceMatrix.h:67
res::HelperMuon::phi
double phi(double pt, double eta)
Definition: Muon.h:154
res::HelperElectron::b
double b(double pt, double eta)
Definition: Electron.h:60
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:30
CovarianceMatrix::setupMatrix
TMatrixD setupMatrix(const pat::PATObject< T > &object, const TopKinFitter::Param param, const std::string &resolutionProvider="")
return covariance matrix for a PAT object
Definition: CovarianceMatrix.h:87
CovarianceMatrix::binsUdsc_
std::vector< std::string > binsUdsc_
vector of strings for the binning of the resolutions
Definition: CovarianceMatrix.h:67
reco::LeafCandidate
Definition: LeafCandidate.h:16
cms::Exception
Definition: Exception.h:70
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
res::HelperMuon::a
double a(double pt, double eta)
Definition: Muon.h:29
CovarianceMatrix::ObjectType
ObjectType
Definition: CovarianceMatrix.h:29
res::HelperMET::a
double a(double pt)
Definition: MET.h:26
res::HelperJet::a
double a(double pt, double eta, Flavor flav)
Definition: Jet.h:42
TopKinFitter::Param
Param
supported parameterizations
Definition: TopKinFitter.h:20
res::HelperJet::et
double et(double pt, double eta, Flavor flav)
Definition: Jet.h:330
StringObjectFunction.h
res::HelperJet::kUds
Definition: Jet.h:9
res::HelperElectron::c
double c(double pt, double eta)
Definition: Electron.h:85
CovarianceMatrix::funcEtaUdsc_
std::vector< std::string > funcEtaUdsc_
Definition: CovarianceMatrix.h:70
reco::LeafCandidate::LorentzVector
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: LeafCandidate.h:23
CovarianceMatrix::funcEtLep_
std::vector< std::string > funcEtLep_
Definition: CovarianceMatrix.h:69
CovarianceMatrix::funcEtaB_
std::vector< std::string > funcEtaB_
Definition: CovarianceMatrix.h:70