CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
CandMassKinFitter.cc
Go to the documentation of this file.
6 #include "TMatrixD.h"
7 #include <iostream>
8 using namespace reco;
9 using namespace std;
10 
12  TKinFitter fitter("CandMassFit", "CandMassFit");
13  TString name("dau0");
14  size_t daus = c.numberOfDaughters();
15  vector<TMatrixD> errors(daus, TMatrix(3, 3));
16  vector<TVector3> momenta(daus);
17  vector<TFitParticleMCCart*> particles(daus, nullptr);
18  TFitConstraintM constraint("MassConstraint", "MassConstraint", nullptr, nullptr, mass_);
19  for (size_t i = 0; i < daus; ++i) {
20  const Candidate& dau = *c.daughter(i);
21  const Particle::LorentzVector& p4 = dau.p4();
22  TMatrixD& err = errors[i];
23  TVector3& mom = momenta[i];
24  mom = TVector3(p4.px(), p4.py(), p4.pz());
25  TrackRef trk = dau.get<TrackRef>();
26  // dummy errors for now...
27  // should play with track parametrization...
28  err.Zero();
29  err(0, 0) = 0.1;
30  err(1, 1) = 0.1;
31  err(2, 2) = 0.1;
32  fitter.addMeasParticle(particles[i] = new TFitParticleMCCart(name, name, &mom, dau.mass(), &err));
33  name[3]++;
34  constraint.addParticle1(particles[i]);
35  }
36  fitter.addConstraint(&constraint);
37  fitter.setMaxNbIter(30);
38  fitter.setMaxDeltaS(1e-2);
39  fitter.setMaxF(1e-1);
40  fitter.setVerbosity(0);
41  fitter.fit();
42  // if ( fitter->getStatus() != 0 ) throw ...
43  TLorentzVector sum(0, 0, 0, 0);
44  for (size_t i = 0; i < daus; ++i) {
45  Candidate& dau = *c.daughter(i);
46  TFitParticleMCCart* part = particles[i];
47  const TLorentzVector* p4 = part->getCurr4Vec();
48  dau.setP4(Particle::LorentzVector(p4->X(), p4->Y(), p4->Z(), p4->T()));
49  sum += *p4;
50  delete particles[i];
51  }
52  c.setP4(Particle::LorentzVector(sum.X(), sum.Y(), sum.Z(), sum.T()));
53  return FitQuality(fitter.getS(), fitter.getNDF());
54 }
void setMaxF(Double_t maxF)
Definition: TKinFitter.h:54
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
const edm::EventSetup & c
Int_t fit()
Definition: TKinFitter.cc:318
virtual void setP4(const LorentzVector &p4)=0
set 4-momentum
void setMaxDeltaS(Double_t maxDeltaS)
Definition: TKinFitter.h:52
virtual double mass() const =0
mass
reco::FitQuality set(reco::Candidate &) const
Double_t getS()
Definition: TKinFitter.cc:1066
virtual size_type numberOfDaughters() const =0
number of daughters
void addConstraint(TAbsFitConstraint *constraint)
Definition: TKinFitter.cc:292
Int_t getNDF()
Definition: TKinFitter.h:45
void addMeasParticle(TAbsFitParticle *particle)
Definition: TKinFitter.cc:194
void setVerbosity(Int_t verbosity=1)
Definition: TKinFitter.cc:304
part
Definition: HCALResponse.h:20
const TLorentzVector * getCurr4Vec()
void setMaxNbIter(Int_t maxNbIter)
Definition: TKinFitter.h:48
void addParticle1(TAbsFitParticle *particle)
T get() const
get a component
Definition: Candidate.h:221
FitQuality
Definition: Fit.h:30
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:21
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector