CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CovarianceMatrix.h
Go to the documentation of this file.
1 #ifndef CovarianceMatrix_h
2 #define CovarianceMatrix_h
3 
5 
8 
14 
20 
22 
23  private:
24 
26  std::vector<std::string> binsUdsc_, binsB_, binsLep_, binsMet_;
28  std::vector<std::string> funcEtUdsc_ , funcEtB_ , funcEtLep_ , funcEtMet_;
29  std::vector<std::string> funcEtaUdsc_, funcEtaB_, funcEtaLep_, funcEtaMet_;
30  std::vector<std::string> funcPhiUdsc_, funcPhiB_, funcPhiLep_, funcPhiMet_;
31 
32  public:
33 
35  CovarianceMatrix(const std::vector<edm::ParameterSet> udscResolutions, const std::vector<edm::ParameterSet> bResolutions){
36  for(std::vector<edm::ParameterSet>::const_iterator iSet = udscResolutions.begin(); iSet != udscResolutions.end(); ++iSet){
37  if(iSet->exists("bin")) binsUdsc_.push_back(iSet->getParameter<std::string>("bin"));
38  else if(udscResolutions.size()==1) binsUdsc_.push_back("");
39  else throw cms::Exception("WrongConfig") << "Parameter 'bin' is needed if more than one PSet is specified!\n";
40 
41  funcEtUdsc_.push_back(iSet->getParameter<std::string>("et"));
42  funcEtaUdsc_.push_back(iSet->getParameter<std::string>("eta"));
43  funcPhiUdsc_.push_back(iSet->getParameter<std::string>("phi"));
44  }
45  for(std::vector<edm::ParameterSet>::const_iterator iSet = bResolutions.begin(); iSet != bResolutions.end(); ++iSet){
46  if(iSet->exists("bin")) binsB_.push_back(iSet->getParameter<std::string>("bin"));
47  else if(bResolutions.size()==1) binsB_.push_back("");
48  else throw cms::Exception("WrongConfig") << "Parameter 'bin' is needed if more than one PSet is specified!\n";
49 
50  funcEtB_.push_back(iSet->getParameter<std::string>("et"));
51  funcEtaB_.push_back(iSet->getParameter<std::string>("eta"));
52  funcPhiB_.push_back(iSet->getParameter<std::string>("phi"));
53  }
54  };
55  CovarianceMatrix(const std::vector<edm::ParameterSet> udscResolutions, const std::vector<edm::ParameterSet> bResolutions, const std::vector<edm::ParameterSet> lepResolutions, const std::vector<edm::ParameterSet> metResolutions){
56  for(std::vector<edm::ParameterSet>::const_iterator iSet = udscResolutions.begin(); iSet != udscResolutions.end(); ++iSet){
57  if(iSet->exists("bin")) binsUdsc_.push_back(iSet->getParameter<std::string>("bin"));
58  else if(udscResolutions.size()==1) binsUdsc_.push_back("");
59  else throw cms::Exception("WrongConfig") << "Parameter 'bin' is needed if more than one PSet is specified!\n";
60 
61  funcEtUdsc_.push_back(iSet->getParameter<std::string>("et"));
62  funcEtaUdsc_.push_back(iSet->getParameter<std::string>("eta"));
63  funcPhiUdsc_.push_back(iSet->getParameter<std::string>("phi"));
64  }
65  for(std::vector<edm::ParameterSet>::const_iterator iSet = bResolutions.begin(); iSet != bResolutions.end(); ++iSet){
66  if(iSet->exists("bin")) binsB_.push_back(iSet->getParameter<std::string>("bin"));
67  else if(bResolutions.size()==1) binsB_.push_back("");
68  else throw cms::Exception("WrongConfig") << "Parameter 'bin' is needed if more than one PSet is specified!\n";
69 
70  funcEtB_.push_back(iSet->getParameter<std::string>("et"));
71  funcEtaB_.push_back(iSet->getParameter<std::string>("eta"));
72  funcPhiB_.push_back(iSet->getParameter<std::string>("phi"));
73  }
74  for(std::vector<edm::ParameterSet>::const_iterator iSet = lepResolutions.begin(); iSet != lepResolutions.end(); ++iSet){
75  if(iSet->exists("bin")) binsLep_.push_back(iSet->getParameter<std::string>("bin"));
76  else if(lepResolutions.size()==1) binsLep_.push_back("");
77  else throw cms::Exception("WrongConfig") << "Parameter 'bin' is needed if more than one PSet is specified!\n";
78 
79  funcEtLep_.push_back(iSet->getParameter<std::string>("et"));
80  funcEtaLep_.push_back(iSet->getParameter<std::string>("eta"));
81  funcPhiLep_.push_back(iSet->getParameter<std::string>("phi"));
82  }
83  for(std::vector<edm::ParameterSet>::const_iterator iSet = metResolutions.begin(); iSet != metResolutions.end(); ++iSet){
84  if(iSet->exists("bin")) binsMet_.push_back(iSet->getParameter<std::string>("bin"));
85  else if(metResolutions.size()==1) binsMet_.push_back("");
86  else throw cms::Exception("WrongConfig") << "Parameter 'bin' is needed if more than one PSet is specified!\n";
87 
88  funcEtMet_.push_back(iSet->getParameter<std::string>("et"));
89  funcEtaMet_.push_back(iSet->getParameter<std::string>("eta"));
90  funcPhiMet_.push_back(iSet->getParameter<std::string>("phi"));
91  }
92  };
94 
95  template <class ObjectType>
96  TMatrixD setupMatrix(const pat::PATObject<ObjectType>& object, TopKinFitter::Param param, std::string resolutionProvider);
97 
98  template <class ObjectType>
99  double getResolution(const pat::PATObject<ObjectType>& object, const std::string whichResolution, bool isBJet);
100 };
101 
102 template <class ObjectType>
103 double CovarianceMatrix::getResolution(const pat::PATObject<ObjectType>& object, const std::string whichResolution = "", bool isBJet = false)
104 {
105  std::vector<std::string> * bins_, * funcEt_, * funcEta_, * funcPhi_;
106 
107  if( dynamic_cast<const reco::Jet*>(&object) && !isBJet ) {
108  bins_ = &binsUdsc_;
109  funcEt_ = &funcEtUdsc_;
110  funcEta_ = &funcEtaUdsc_;
111  funcPhi_ = &funcPhiUdsc_;
112  }
113  else if( dynamic_cast<const reco::Jet*>(&object) && isBJet ) {
114  bins_ = &binsB_;
115  funcEt_ = &funcEtB_;
116  funcEta_ = &funcEtaB_;
117  funcPhi_ = &funcPhiB_;
118  }
119  else if( dynamic_cast<const reco::Muon*>(&object) || dynamic_cast<const reco::GsfElectron*>(&object) ) {
120  bins_ = &binsLep_;
121  funcEt_ = &funcEtLep_;
122  funcEta_ = &funcEtaLep_;
123  funcPhi_ = &funcPhiLep_;
124  }
125  else if( dynamic_cast<const reco::MET*>(&object) ) {
126  bins_ = &binsMet_;
127  funcEt_ = &funcEtMet_;
128  funcEta_ = &funcEtaMet_;
129  funcPhi_ = &funcPhiMet_;
130  }
131  else{
132  throw cms::Exception("UnsupportedObject") << "The object given is not supported!\n";
133  }
134 
135  int selectedBin=-1;
136  reco::LeafCandidate candidate;
137  for(unsigned int i=0; i<bins_->size(); ++i){
139  candidate = reco::LeafCandidate( 0, reco::LeafCandidate::LorentzVector(object.px(), object.py(), object.pz(), object.energy()) );
140  if(select_(candidate)){
141  selectedBin = i;
142  break;
143  }
144  }
145  if(selectedBin>=0){
146  if(whichResolution == "et") return StringObjectFunction<reco::LeafCandidate>(funcEt_ ->at(selectedBin)).operator()(candidate);
147  else if(whichResolution == "eta") return StringObjectFunction<reco::LeafCandidate>(funcEta_->at(selectedBin)).operator()(candidate);
148  else if(whichResolution == "phi") return StringObjectFunction<reco::LeafCandidate>(funcPhi_->at(selectedBin)).operator()(candidate);
149  else throw cms::Exception("ProgrammingError") << "Only 'et', 'eta' and 'phi' resolutions supported!\n";
150  }
151  return 0;
152 }
153 
154 
155 template <class ObjectType>
156 TMatrixD CovarianceMatrix::setupMatrix(const pat::PATObject<ObjectType>& object, TopKinFitter::Param param, std::string resolutionProvider = "")
157 {
158  TMatrixD CovM3 (3,3); CovM3.Zero();
159  TMatrixD CovM4 (4,4); CovM4.Zero();
160  TMatrixD* CovM = &CovM3;
161  // This part is for pat objects with resolutions embedded
162  if(object.hasKinResolution())
163  {
164  switch(param){
166  CovM3(0,0) = pow(object.resolEt(resolutionProvider) , 2);
167  if( dynamic_cast<const reco::MET*>(&object) ) CovM3(1,1) = pow(9999., 2);
168  else CovM3(1,1) = pow(object.resolEta(resolutionProvider), 2);
169  CovM3(2,2) = pow(object.resolPhi(resolutionProvider), 2);
170  CovM = &CovM3;
171  break;
173  CovM3(0,0) = pow(object.resolEt(resolutionProvider) , 2);
174  CovM3(1,1) = pow(object.resolTheta(resolutionProvider), 2);
175  CovM3(2,2) = pow(object.resolPhi(resolutionProvider) , 2);
176  CovM = &CovM3;
177  break;
178  case TopKinFitter::kEMom :
179  CovM4(0,0) = pow(1, 2);
180  CovM4(1,1) = pow(1, 2);
181  CovM4(2,2) = pow(1, 2);
182  CovM4(3,3) = pow(1, 2);
183  CovM = &CovM4;
184  break;
185  }
186  }
187  // This part is for objects without resolutions embedded
188  else
189  {
190  double pt = object.pt(), eta = object.eta();
191  // if object is a jet
192  if( dynamic_cast<const reco::Jet*>(&object) ) {
193  res::HelperJet jetRes;
194  switch(param){
195  case TopKinFitter::kEMom :
196  if(resolutionProvider == "bjets") {
197  CovM4(0,0) = pow(jetRes.a (pt, eta, res::HelperJet::kB ), 2);
198  CovM4(1,1) = pow(jetRes.b (pt, eta, res::HelperJet::kB ), 2);
199  CovM4(2,2) = pow(jetRes.c (pt, eta, res::HelperJet::kB ), 2);
200  CovM4(3,3) = pow(jetRes.d (pt, eta, res::HelperJet::kB ), 2);
201  }
202  else {
203  CovM4(0,0) = pow(jetRes.a (pt, eta, res::HelperJet::kUds), 2);
204  CovM4(1,1) = pow(jetRes.b (pt, eta, res::HelperJet::kUds), 2);
205  CovM4(2,2) = pow(jetRes.c (pt, eta, res::HelperJet::kUds), 2);
206  CovM4(3,3) = pow(jetRes.d (pt, eta, res::HelperJet::kUds), 2);
207  }
208  CovM = &CovM4;
209  break;
211  if(resolutionProvider == "bjets") {
212  if(!binsB_.size()){
213  CovM3(0,0) = pow(jetRes.et (pt, eta, res::HelperJet::kB ), 2);
214  CovM3(1,1) = pow(jetRes.eta(pt, eta, res::HelperJet::kB ), 2);
215  CovM3(2,2) = pow(jetRes.phi(pt, eta, res::HelperJet::kB ), 2);
216  }
217  else{
218  CovM3(0,0) = pow(getResolution(object, "et" , true), 2);
219  CovM3(1,1) = pow(getResolution(object, "eta", true), 2);
220  CovM3(2,2) = pow(getResolution(object, "phi", true), 2);
221  }
222  }
223  else {
224  if(!binsUdsc_.size()){
225  CovM3(0,0) = pow(jetRes.et (pt, eta, res::HelperJet::kUds), 2);
226  CovM3(1,1) = pow(jetRes.eta(pt, eta, res::HelperJet::kUds), 2);
227  CovM3(2,2) = pow(jetRes.phi(pt, eta, res::HelperJet::kUds), 2);
228  }
229  else{
230  CovM3(0,0) = pow(getResolution(object, "et") , 2);
231  CovM3(1,1) = pow(getResolution(object, "eta"), 2);
232  CovM3(2,2) = pow(getResolution(object, "phi"), 2);
233  }
234  }
235  CovM = &CovM3;
236  break;
238  if(resolutionProvider == "bjets") {
239  CovM3(0,0) = pow(jetRes.et (pt, eta, res::HelperJet::kB ), 2);
240  CovM3(1,1) = pow(jetRes.theta(pt, eta, res::HelperJet::kB ), 2);
241  CovM3(2,2) = pow(jetRes.phi (pt, eta, res::HelperJet::kB ), 2);
242  }
243  else {
244  CovM3(0,0) = pow(jetRes.et (pt, eta, res::HelperJet::kUds), 2);
245  CovM3(1,1) = pow(jetRes.theta(pt, eta, res::HelperJet::kUds), 2);
246  CovM3(2,2) = pow(jetRes.phi (pt, eta, res::HelperJet::kUds), 2);
247  }
248  CovM = &CovM3;
249  break;
250  }
251  }
252  // if object is an electron
253  else if( dynamic_cast<const reco::GsfElectron*>(&object) ) {
254  res::HelperElectron elecRes;
255  switch(param){
256  case TopKinFitter::kEMom :
257  CovM3(0,0) = pow(elecRes.a (pt, eta), 2);
258  CovM3(1,1) = pow(elecRes.b (pt, eta), 2);
259  CovM3(2,2) = pow(elecRes.c (pt, eta), 2);
260  CovM = &CovM3;
261  break;
263  if(!binsLep_.size()){
264  CovM3(0,0) = pow(elecRes.et (pt, eta), 2);
265  CovM3(1,1) = pow(elecRes.eta(pt, eta), 2);
266  CovM3(2,2) = pow(elecRes.phi(pt, eta), 2);
267  }
268  else{
269  CovM3(0,0) = pow(getResolution(object, "et") , 2);
270  CovM3(1,1) = pow(getResolution(object, "eta"), 2);
271  CovM3(2,2) = pow(getResolution(object, "phi"), 2);
272  }
273  CovM = &CovM3;
274  break;
276  CovM3(0,0) = pow(elecRes.et (pt, eta), 2);
277  CovM3(1,1) = pow(elecRes.theta(pt, eta), 2);
278  CovM3(2,2) = pow(elecRes.phi (pt, eta), 2);
279  CovM = &CovM3;
280  break;
281  }
282  }
283  // if object is a muon
284  else if( dynamic_cast<const reco::Muon*>(&object) ) {
285  res::HelperMuon muonRes;
286  switch(param){
287  case TopKinFitter::kEMom :
288  CovM3(0,0) = pow(muonRes.a (pt, eta), 2);
289  CovM3(1,1) = pow(muonRes.b (pt, eta), 2);
290  CovM3(2,2) = pow(muonRes.c (pt, eta), 2);
291  CovM = &CovM3;
292  break;
294  if(!binsLep_.size()){
295  CovM3(0,0) = pow(muonRes.et (pt, eta), 2);
296  CovM3(1,1) = pow(muonRes.eta(pt, eta), 2);
297  CovM3(2,2) = pow(muonRes.phi(pt, eta), 2);
298  }
299  else{
300  CovM3(0,0) = pow(getResolution(object, "et") , 2);
301  CovM3(1,1) = pow(getResolution(object, "eta"), 2);
302  CovM3(2,2) = pow(getResolution(object, "phi"), 2);
303  }
304  CovM = &CovM3;
305  break;
307  CovM3(0,0) = pow(muonRes.et (pt, eta), 2);
308  CovM3(1,1) = pow(muonRes.theta(pt, eta), 2);
309  CovM3(2,2) = pow(muonRes.phi (pt, eta), 2);
310  CovM = &CovM3;
311  break;
312  }
313  }
314  // if object is met
315  else if( dynamic_cast<const reco::MET*>(&object) ) {
316  res::HelperMET metRes;
317  switch(param){
318  case TopKinFitter::kEMom :
319  CovM3(0,0) = pow(metRes.a(pt), 2);
320  CovM3(1,1) = pow(metRes.b(pt), 2);
321  CovM3(2,2) = pow(metRes.c(pt), 2);
322  CovM = &CovM3;
323  break;
325  if(!binsMet_.size()){
326  CovM3(0,0) = pow(metRes.et(pt) , 2);
327  CovM3(1,1) = pow( 9999. , 2);
328  CovM3(2,2) = pow(metRes.phi(pt), 2);
329  }
330  else{
331  CovM3(0,0) = pow(getResolution(object, "et") , 2);
332  CovM3(1,1) = pow(getResolution(object, "eta"), 2);
333  CovM3(2,2) = pow(getResolution(object, "phi"), 2);
334  }
335  CovM = &CovM3;
336  break;
338  CovM3(0,0) = pow(metRes.et(pt) , 2);
339  CovM3(1,1) = pow( 9999. , 2);
340  CovM3(2,2) = pow(metRes.phi(pt), 2);
341  CovM = &CovM3;
342  break;
343  }
344  }
345  }
346  return *CovM;
347 }
348 
349 #endif
int i
Definition: DBlmapReader.cc:9
Param
supported parameterizations
Definition: TopKinFitter.h:22
double c(double pt, double eta)
Definition: Electron.h:67
std::vector< std::string > binsB_
double b(double pt)
Definition: MET.h:34
double et(double pt, double eta)
Definition: Electron.h:131
double phi(double pt, double eta, Flavor flav)
Definition: Jet.h:183
double a(double pt)
Definition: MET.h:28
double phi(double pt, double eta)
Definition: Electron.h:115
std::vector< std::string > funcEtUdsc_
vectors for the resolution functions
double phi(double pt)
Definition: MET.h:58
double b(double pt, double eta)
Definition: Muon.h:45
std::vector< std::string > binsMet_
double theta(double pt, double eta, Flavor flav)
Definition: Jet.h:154
double c(double pt, double eta)
Definition: Muon.h:61
TMatrixD setupMatrix(const pat::PATObject< ObjectType > &object, TopKinFitter::Param param, std::string resolutionProvider)
double eta(double pt, double eta)
Definition: Electron.h:147
CovarianceMatrix(const std::vector< edm::ParameterSet > udscResolutions, const std::vector< edm::ParameterSet > bResolutions)
double et(double pt, double eta, Flavor flav)
Definition: Jet.h:212
double getResolution(const pat::PATObject< ObjectType > &object, const std::string whichResolution, bool isBJet)
std::vector< std::string > funcEtB_
T eta() const
std::vector< std::string > binsUdsc_
vector of strings for the binning of the resolutions
std::vector< std::string > funcPhiB_
double b(double pt, double eta, Flavor flav)
Definition: Jet.h:67
CovarianceMatrix(const std::vector< edm::ParameterSet > udscResolutions, const std::vector< edm::ParameterSet > bResolutions, const std::vector< edm::ParameterSet > lepResolutions, const std::vector< edm::ParameterSet > metResolutions)
double d(double pt, double eta, Flavor flav)
Definition: Jet.h:125
std::vector< std::string > funcEtaB_
double a(double pt, double eta)
Definition: Electron.h:35
double theta(double pt, double eta)
Definition: Muon.h:93
double et(double pt, double eta)
Definition: Muon.h:125
std::vector< std::string > funcPhiUdsc_
std::vector< std::string > funcPhiLep_
double c(double pt)
Definition: MET.h:40
double eta(double pt, double eta, Flavor flav)
Definition: Jet.h:241
std::vector< std::string > funcEtaUdsc_
std::vector< std::string > funcEtLep_
double theta(double pt, double eta)
Definition: Electron.h:99
std::vector< std::string > funcEtMet_
double eta(double pt, double eta)
Definition: Muon.h:141
double a(double pt, double eta)
Definition: Muon.h:29
std::vector< std::string > funcEtaLep_
std::vector< std::string > binsLep_
Templated PAT object container.
Definition: PATObject.h:42
double c(double pt, double eta, Flavor flav)
Definition: Jet.h:96
double et(double pt)
Definition: MET.h:64
std::vector< std::string > funcEtaMet_
double phi(double pt, double eta)
Definition: Muon.h:109
std::vector< std::string > funcPhiMet_
double a(double pt, double eta, Flavor flav)
Definition: Jet.h:38
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: LeafCandidate.h:25
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
double b(double pt, double eta)
Definition: Electron.h:51