CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
StKinFitter.cc
Go to the documentation of this file.
1 //
2 // $Id: StKinFitter.cc,v 1.8 2010/09/06 13:46:16 snaumann Exp $
3 //
4 
12 
15 
16 //introduced to repair kinFit w/o resolutions from pat
21 
22 /* other parametrizations and constraints
23 #include "PhysicsTools/KinFitter/interface/TFitParticleESpher.h"
24 #include "PhysicsTools/KinFitter/interface/TFitParticleMCPInvSpher.h"
25 #include "PhysicsTools/KinFitter/interface/TFitConstraintMGaus.h"
26 #include "PhysicsTools/KinFitter/interface/TFitConstraintEp.h"*/
27 
29  TopKinFitter(),
30  jetParam_(kEMom),
31  lepParam_(kEMom),
32  metParam_(kEMom)
33 {
34  setupFitter();
35 }
36 
37 StKinFitter::StKinFitter(int jetParam, int lepParam, int metParam,
38  int maxNrIter, double maxDeltaS, double maxF, std::vector<int> constraints) :
39  TopKinFitter(maxNrIter, maxDeltaS, maxF),
40  jetParam_((Param) jetParam),
41  lepParam_((Param) lepParam),
42  metParam_((Param) metParam),
43  constraints_(constraints)
44 {
45  setupFitter();
46 }
47 
48 StKinFitter::StKinFitter(Param jetParam, Param lepParam, Param metParam,
49  int maxNrIter, double maxDeltaS, double maxF, std::vector<int> constraints) :
50  TopKinFitter(maxNrIter, maxDeltaS, maxF),
51  jetParam_(jetParam),
52  lepParam_(lepParam),
53  metParam_(metParam),
54  constraints_(constraints)
55 {
56  setupFitter();
57 }
58 
60 {
61  delete cons1_; delete cons2_; delete cons3_;
62  delete fitBottom_; delete fitLight_; delete fitLepton_; delete fitNeutrino_;
63 }
64 
66 {
67  StEvtSolution fitsol(*asol);
68 
69  TMatrixD m1(3,3), m2(3,3);
70  TMatrixD m1b(4,4), m2b(4,4);
71  TMatrixD m3(3,3), m4(3,3);
72  m1.Zero(); m2.Zero();
73  m1b.Zero(); m2b.Zero();
74  m3.Zero(); m4.Zero();
75 
76  TLorentzVector bottomVec(fitsol.getBottom().px(),fitsol.getBottom().py(),
77  fitsol.getBottom().pz(),fitsol.getBottom().energy());
78  TLorentzVector lightVec(fitsol.getLight().px(),fitsol.getLight().py(),
79  fitsol.getLight().pz(),fitsol.getLight().energy());
80  TLorentzVector leplVec;
81  if(fitsol.getDecay()== "electron") leplVec = TLorentzVector(fitsol.getElectron().px(), fitsol.getElectron().py(),
82  fitsol.getElectron().pz(), fitsol.getElectron().energy());
83  if(fitsol.getDecay()== "muon") leplVec = TLorentzVector(fitsol.getMuon().px(), fitsol.getMuon().py(),
84  fitsol.getMuon().pz(), fitsol.getMuon().energy());
85  TLorentzVector lepnVec(fitsol.getNeutrino().px(), fitsol.getNeutrino().py(),
86  0, fitsol.getNeutrino().et());
87 
88  // jet resolutions
89  {
90  //FIXME this dirty hack needs a clean solution soon!
91  double pt = fitsol.getBottom().pt ();
92  double eta = fitsol.getBottom().eta();
93  res::HelperJet jetRes;
94  if (jetParam_ == kEMom) {
95  m1b(0,0) = pow(jetRes.pt (pt, eta, res::HelperJet::kB ), 2);
96  m1b(1,1) = pow(jetRes.pt (pt, eta, res::HelperJet::kB ), 2);
97  m1b(2,2) = pow(jetRes.pt (pt, eta, res::HelperJet::kB ), 2);
98  m1b(3,3) = pow(jetRes.pt (pt, eta, res::HelperJet::kB ), 2);
99  m2b(0,0) = pow(jetRes.pt (pt, eta, res::HelperJet::kUds), 2);
100  m2b(1,1) = pow(jetRes.pt (pt, eta, res::HelperJet::kUds), 2);
101  m2b(2,2) = pow(jetRes.pt (pt, eta, res::HelperJet::kUds), 2);
102  m2b(3,3) = pow(jetRes.pt (pt, eta, res::HelperJet::kUds), 2);
103  } else if (jetParam_ == kEtEtaPhi) {
104  m1 (0,0) = pow(jetRes.pt (pt, eta, res::HelperJet::kB ), 2);
105  m1 (1,1) = pow(jetRes.eta(pt, eta, res::HelperJet::kB ), 2);
106  m1 (2,2) = pow(jetRes.phi(pt, eta, res::HelperJet::kB ), 2);
107  m2 (0,0) = pow(jetRes.pt (pt, eta, res::HelperJet::kUds), 2);
108  m2 (1,1) = pow(jetRes.eta(pt, eta, res::HelperJet::kUds), 2);
109  m2 (2,2) = pow(jetRes.phi(pt, eta, res::HelperJet::kUds), 2);
110  } else if (jetParam_ == kEtThetaPhi) {
111  m1 (0,0) = pow(jetRes.pt (pt, eta, res::HelperJet::kB ), 2);
112  m1 (1,1) = pow(jetRes.eta(pt, eta, res::HelperJet::kB ), 2);
113  m1 (2,2) = pow(jetRes.phi(pt, eta, res::HelperJet::kB ), 2);
114  m2 (0,0) = pow(jetRes.pt (pt, eta, res::HelperJet::kUds), 2);
115  m2 (1,1) = pow(jetRes.eta(pt, eta, res::HelperJet::kUds), 2);
116  m2 (2,2) = pow(jetRes.phi(pt, eta, res::HelperJet::kUds), 2);
117  }
118  }
119  // lepton resolutions
120  {
121  //FIXME this dirty hack needs a clean solution soon!
122  double pt = fitsol.getElectron().pt ();
123  double eta = fitsol.getElectron().eta();
124  res::HelperMuon muonRes;
125  res::HelperElectron elecRes;
126  if (lepParam_ == kEMom) {
127  if(fitsol.getDecay()== "electron"){
128  m3(0,0) = pow(elecRes.pt (pt, eta), 2);
129  m3(1,1) = pow(elecRes.pt (pt, eta), 2);
130  m3(2,2) = pow(elecRes.pt (pt, eta), 2);
131  }
132  if(fitsol.getDecay()== "muon"){
133  m3(0,0) = pow(muonRes.pt (pt, eta), 2);
134  m3(1,1) = pow(muonRes.pt (pt, eta), 2);
135  m3(2,2) = pow(muonRes.pt (pt, eta), 2);
136  }
137  } else if (lepParam_ == kEtEtaPhi) {
138  if(fitsol.getDecay()== "electron"){
139  m3(0,0) = pow(elecRes.pt (pt, eta), 2);
140  m3(1,1) = pow(elecRes.eta(pt, eta), 2);
141  m3(2,2) = pow(elecRes.phi(pt, eta), 2);
142  }
143  if(fitsol.getDecay()== "muon"){
144  m3(0,0) = pow(muonRes.pt (pt, eta), 2);
145  m3(1,1) = pow(muonRes.eta(pt, eta), 2);
146  m3(2,2) = pow(muonRes.phi(pt, eta), 2);
147  }
148  } else if (lepParam_ == kEtThetaPhi) {
149  if(fitsol.getDecay()== "electron"){
150  m3(0,0) = pow(elecRes.pt (pt, eta), 2);
151  m3(1,1) = pow(elecRes.eta(pt, eta), 2);
152  m3(2,2) = pow(elecRes.phi(pt, eta), 2);
153  }
154  if(fitsol.getDecay()== "muon"){
155  m3(0,0) = pow(muonRes.pt (pt, eta), 2);
156  m3(1,1) = pow(muonRes.eta(pt, eta), 2);
157  m3(2,2) = pow(muonRes.phi(pt, eta), 2);
158  }
159  }
160  }
161  // neutrino resolutions
162  {
163  //FIXME this dirty hack needs a clean solution soon!
164  double met = fitsol.getNeutrino().pt();
165  res::HelperMET metRes;
166  if (metParam_ == kEMom) {
167  m4(0,0) = pow(metRes.met(met), 2);
168  m4(1,1) = pow( 9999., 2);
169  m4(2,2) = pow(metRes.met(met), 2);
170  } else if (metParam_ == kEtEtaPhi) {
171  m4(0,0) = pow(metRes.met(met), 2);
172  m4(1,1) = pow( 9999., 2);
173  m4(2,2) = pow(metRes.phi(met), 2);
174  } else if (metParam_ == kEtThetaPhi) {
175  m4(0,0) = pow(metRes.met(met), 2);
176  m4(1,1) = pow( 9999., 2);
177  m4(2,2) = pow(metRes.phi(met), 2);
178  }
179  }
180  // set the kinematics of the objects to be fitted
181  fitBottom_->setIni4Vec(&bottomVec);
182  fitLight_->setIni4Vec(&lightVec);
183  fitLepton_->setIni4Vec(&leplVec);
184  fitNeutrino_->setIni4Vec(&lepnVec);
185  if (jetParam_ == kEMom) {
186  fitBottom_->setCovMatrix(&m1b);
187  fitLight_->setCovMatrix(&m2b);
188  } else {
189  fitBottom_->setCovMatrix(&m1);
190  fitLight_->setCovMatrix(&m2);
191  }
192  fitLepton_->setCovMatrix(&m3);
194 
195  // perform the fit!
196  fitter_->fit();
197 
198  // add fitted information to the solution
199  if (fitter_->getStatus() == 0) {
200  // read back the jet kinematics and resolutions
203 
204  // read back the lepton kinematics and resolutions
206 
207  // read back the MET kinematics and resolutions
209 
210  // finally fill the fitted particles
211  fitsol.setFitBottom(aFitBottom);
212  fitsol.setFitLight(aFitLight);
213  fitsol.setFitLepton(aFitLepton);
214  fitsol.setFitNeutrino(aFitNeutrino);
215 
216  // store the fit's chi2 probability
217  fitsol.setChi2Prob( fitProb() );
218  }
219 
220  return fitsol;
221 
222 }
223 
224 //
225 // Setup the fitter
226 //
228 
229  // FIXME: replace by messagelogger!!!
230 
231  std::cout<<std::endl<<std::endl<<"+++++++++++ KINFIT SETUP ++++++++++++"<<std::endl;
232  std::cout<<" jet parametrisation: " << param(jetParam_) << std::endl;
233  std::cout<<" lepton parametrisation: " << param(lepParam_) << std::endl;
234  std::cout<<" met parametrisation: " << param(metParam_) << std::endl;
235  std::cout<<" constraints: "<<std::endl;
236  for(unsigned int i=0; i<constraints_.size(); i++){
237  if(constraints_[i] == 1) std::cout<<" - hadronic W-mass"<<std::endl;
238  if(constraints_[i] == 2) std::cout<<" - leptonic W-mass"<<std::endl;
239  if(constraints_[i] == 3) std::cout<<" - hadronic top mass"<<std::endl;
240  if(constraints_[i] == 4) std::cout<<" - leptonic top mass"<<std::endl;
241  if(constraints_[i] == 5) std::cout<<" - neutrino mass"<<std::endl;
242  }
243  std::cout<<"Max. number of iterations: "<<maxNrIter_<<std::endl;
244  std::cout<<"Max. deltaS: "<<maxDeltaS_<<std::endl;
245  std::cout<<"Max. F: "<<maxF_<<std::endl;
246  std::cout<<"++++++++++++++++++++++++++++++++++++++++++++"<<std::endl<<std::endl<<std::endl;
247 
248  TMatrixD empty3(3,3); TMatrixD empty4(4,4);
249  if (jetParam_ == kEMom) {
250  fitBottom_ = new TFitParticleEMomDev("Jet1", "Jet1", 0, &empty4);
251  fitLight_ = new TFitParticleEMomDev("Jet2", "Jet2", 0, &empty4);
252  } else if (jetParam_ == kEtEtaPhi) {
253  fitBottom_ = new TFitParticleEtEtaPhi("Jet1", "Jet1", 0, &empty3);
254  fitLight_ = new TFitParticleEtEtaPhi("Jet2", "Jet2", 0, &empty3);
255  } else if (jetParam_ == kEtThetaPhi) {
256  fitBottom_ = new TFitParticleEtThetaPhi("Jet1", "Jet1", 0, &empty3);
257  fitLight_ = new TFitParticleEtThetaPhi("Jet2", "Jet2", 0, &empty3);
258  }
259  if (lepParam_ == kEMom) {
260  fitLepton_ = new TFitParticleEScaledMomDev("Lepton", "Lepton", 0, &empty3);
261  } else if (lepParam_ == kEtEtaPhi) {
262  fitLepton_ = new TFitParticleEtEtaPhi("Lepton", "Lepton", 0, &empty3);
263  } else if (lepParam_ == kEtThetaPhi) {
264  fitLepton_ = new TFitParticleEtThetaPhi("Lepton", "Lepton", 0, &empty3);
265  }
266  if (metParam_ == kEMom) {
267  fitNeutrino_ = new TFitParticleEScaledMomDev("Neutrino", "Neutrino", 0, &empty3);
268  } else if (metParam_ == kEtEtaPhi) {
269  fitNeutrino_ = new TFitParticleEtEtaPhi("Neutrino", "Neutrino", 0, &empty3);
270  } else if (metParam_ == kEtThetaPhi) {
271  fitNeutrino_ = new TFitParticleEtThetaPhi("Neutrino", "Neutrino", 0, &empty3);
272  }
273 
274  cons1_ = new TFitConstraintM("MassConstraint", "Mass-Constraint", 0, 0 , mW_);
276  cons2_ = new TFitConstraintM("MassConstraint", "Mass-Constraint", 0, 0, mTop_);
278  cons3_ = new TFitConstraintM("MassConstraint", "Mass-Constraint", 0, 0, 0.);
280 
281  for (unsigned int i=0; i<constraints_.size(); i++) {
285  }
290 
291 }
double pt(double pt, double eta, Flavor flav)
Definition: Jet.h:25
int i
Definition: DBlmapReader.cc:9
TAbsFitParticle * fitLight_
Definition: StKinFitter.h:39
pat::Muon getMuon() const
Definition: StEvtSolution.h:37
TFitConstraintM * cons2_
Definition: StKinFitter.h:44
void addParticles1(TAbsFitParticle *p1, TAbsFitParticle *p2=0, TAbsFitParticle *p3=0, TAbsFitParticle *p4=0, TAbsFitParticle *p5=0, TAbsFitParticle *p6=0, TAbsFitParticle *p7=0, TAbsFitParticle *p8=0, TAbsFitParticle *p9=0, TAbsFitParticle *p10=0)
Int_t fit()
Definition: TKinFitter.cc:309
virtual double et() const
transverse energy
void setupFitter()
Definition: StKinFitter.cc:227
double phi(double pt, double eta, Flavor flav)
Definition: Jet.h:183
virtual void setIni4Vec(const TLorentzVector *pini)=0
double phi(double pt, double eta)
Definition: Electron.h:115
double phi(double pt)
Definition: MET.h:58
std::string param(const Param &param) const
convert Param to human readable form
Definition: TopKinFitter.cc:23
double pt(double pt, double eta)
Definition: Electron.h:23
double eta(double pt, double eta)
Definition: Electron.h:147
T eta() const
int maxNrIter_
maximal allowed number of iterations to be used for the fit
Definition: TopKinFitter.h:48
void setFitBottom(const pat::Particle &part)
virtual double eta() const
momentum pseudorapidity
void setChi2Prob(double prob)
void setFitNeutrino(const pat::Particle &part)
double fitProb() const
return fit probability
Definition: TopKinFitter.h:36
TFitConstraintM * cons3_
Definition: StKinFitter.h:45
TAbsFitParticle * fitBottom_
Definition: StKinFitter.h:38
pat::Jet getBottom() const
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
TAbsFitParticle * fitNeutrino_
Definition: StKinFitter.h:41
virtual double energy() const
energy
Param metParam_
Definition: StKinFitter.h:47
pat::Jet getLight() const
void addConstraint(TAbsFitConstraint *constraint)
Definition: TKinFitter.cc:281
Int_t getStatus()
Definition: TKinFitter.h:41
virtual void setCovMatrix(const TMatrixD *theCovMatrix)
std::vector< int > constraints_
Definition: StKinFitter.h:48
void setFitLepton(const pat::Particle &part)
double mW_
W mass value used for constraints.
Definition: TopKinFitter.h:54
void addMeasParticle(TAbsFitParticle *particle)
Definition: TKinFitter.cc:209
double eta(double pt, double eta, Flavor flav)
Definition: Jet.h:241
void setFitLight(const pat::Particle &part)
std::string getDecay() const
Definition: StEvtSolution.h:78
double maxDeltaS_
maximal allowed chi2 (not normalized to degrees of freedom)
Definition: TopKinFitter.h:50
pat::MET getNeutrino() const
Definition: StEvtSolution.h:39
virtual double px() const
x coordinate of momentum vector
virtual double pt() const
transverse momentum
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
Analysis-level particle class.
Definition: Particle.h:34
Param lepParam_
Definition: StKinFitter.h:47
double eta(double pt, double eta)
Definition: Muon.h:141
virtual double pz() const
z coordinate of momentum vector
TFitConstraintM * cons1_
Definition: StKinFitter.h:43
const TLorentzVector * getCurr4Vec()
void addParticle1(TAbsFitParticle *particle)
double pt(double pt, double eta)
Definition: Muon.h:23
StEvtSolution addKinFitInfo(StEvtSolution *asol)
Definition: StKinFitter.cc:65
tuple cout
Definition: gather_cfg.py:121
TAbsFitParticle * fitLepton_
Definition: StKinFitter.h:40
TKinFitter * fitter_
kinematic fitter
Definition: TopKinFitter.h:46
Param jetParam_
Definition: StKinFitter.h:47
double phi(double pt, double eta)
Definition: Muon.h:109
double mTop_
top mass value used for constraints
Definition: TopKinFitter.h:56
double met(double met)
Definition: MET.h:23
virtual double py() const
y coordinate of momentum vector
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
pat::Electron getElectron() const
Definition: StEvtSolution.h:38
double maxF_
maximal allowed distance from constraints
Definition: TopKinFitter.h:52