CMS 3D CMS Logo

TtSemiLepKinFitter.cc
Go to the documentation of this file.
8 
11 
13 
15 TtSemiLepKinFitter::TtSemiLepKinFitter() : TopKinFitter(), jetParam_(kEMom), lepParam_(kEMom), metParam_(kEMom) {
16  setupFitter();
17 }
18 
20  Param lepParam,
21  Param metParam,
22  int maxNrIter,
23  double maxDeltaS,
24  double maxF,
25  const std::vector<Constraint>& constraints,
26  double mW,
27  double mTop,
28  const std::vector<edm::ParameterSet>* udscResolutions,
29  const std::vector<edm::ParameterSet>* bResolutions,
30  const std::vector<edm::ParameterSet>* lepResolutions,
31  const std::vector<edm::ParameterSet>* metResolutions,
32  const std::vector<double>* jetEnergyResolutionScaleFactors,
33  const std::vector<double>* jetEnergyResolutionEtaBinning)
35  udscResolutions_(udscResolutions),
36  bResolutions_(bResolutions),
37  lepResolutions_(lepResolutions),
38  metResolutions_(metResolutions),
39  jetEnergyResolutionScaleFactors_(jetEnergyResolutionScaleFactors),
40  jetEnergyResolutionEtaBinning_(jetEnergyResolutionEtaBinning),
41  jetParam_(jetParam),
42  lepParam_(lepParam),
43  metParam_(metParam),
44  constrList_(constraints) {
45  setupFitter();
46 }
47 
49 
51  std::stringstream constr;
52  for (unsigned int i = 0; i < constrList_.size(); ++i) {
53  switch (constrList_[i]) {
54  case kWHadMass:
55  constr << " * hadronic W-mass (" << mW_ << " GeV) \n";
56  break;
57  case kWLepMass:
58  constr << " * leptonic W-mass (" << mW_ << " GeV) \n";
59  break;
60  case kTopHadMass:
61  constr << " * hadronic t-mass (" << mTop_ << " GeV) \n";
62  break;
63  case kTopLepMass:
64  constr << " * leptonic t-mass (" << mTop_ << " GeV) \n";
65  break;
66  case kNeutrinoMass:
67  constr << " * neutrino mass (0 GeV) \n";
68  break;
69  case kEqualTopMasses:
70  constr << " * equal t-masses \n";
71  break;
72  case kSumPt:
73  constr << " * summed transverse momentum \n";
74  break;
75  }
76  }
77  edm::LogVerbatim("TtSemiLepKinFitter") << "\n"
78  << "+++++++++++ TtSemiLepKinFitter Setup ++++++++++++ \n"
79  << " Parametrization: \n"
80  << " * jet : " << param(jetParam_) << "\n"
81  << " * lep : " << param(lepParam_) << "\n"
82  << " * met : " << param(metParam_) << "\n"
83  << " Constraints: \n"
84  << constr.str() << " Max(No iterations): " << maxNrIter_ << "\n"
85  << " Max(deltaS) : " << maxDeltaS_ << "\n"
86  << " Max(F) : " << maxF_ << "\n"
87  << "+++++++++++++++++++++++++++++++++++++++++++++++++ \n";
88 }
89 
91  TMatrixD empty3x3(3, 3);
92  TMatrixD empty4x4(4, 4);
93  switch (jetParam_) { // setup jets according to parameterization
94  case kEMom:
95  hadB_ = std::make_unique<TFitParticleEMomDev>("Jet1", "Jet1", nullptr, &empty4x4);
96  hadP_ = std::make_unique<TFitParticleEMomDev>("Jet2", "Jet2", nullptr, &empty4x4);
97  hadQ_ = std::make_unique<TFitParticleEMomDev>("Jet3", "Jet3", nullptr, &empty4x4);
98  lepB_ = std::make_unique<TFitParticleEMomDev>("Jet4", "Jet4", nullptr, &empty4x4);
99  break;
100  case kEtEtaPhi:
101  hadB_ = std::make_unique<TFitParticleEtEtaPhi>("Jet1", "Jet1", nullptr, &empty3x3);
102  hadP_ = std::make_unique<TFitParticleEtEtaPhi>("Jet2", "Jet2", nullptr, &empty3x3);
103  hadQ_ = std::make_unique<TFitParticleEtEtaPhi>("Jet3", "Jet3", nullptr, &empty3x3);
104  lepB_ = std::make_unique<TFitParticleEtEtaPhi>("Jet4", "Jet4", nullptr, &empty3x3);
105  break;
106  case kEtThetaPhi:
107  hadB_ = std::make_unique<TFitParticleEtThetaPhi>("Jet1", "Jet1", nullptr, &empty3x3);
108  hadP_ = std::make_unique<TFitParticleEtThetaPhi>("Jet2", "Jet2", nullptr, &empty3x3);
109  hadQ_ = std::make_unique<TFitParticleEtThetaPhi>("Jet3", "Jet3", nullptr, &empty3x3);
110  lepB_ = std::make_unique<TFitParticleEtThetaPhi>("Jet4", "Jet4", nullptr, &empty3x3);
111  break;
112  }
113 }
114 
116  TMatrixD empty3x3(3, 3);
117  switch (lepParam_) { // setup lepton according to parameterization
118  case kEMom:
119  lepton_ = std::make_unique<TFitParticleEScaledMomDev>("Lepton", "Lepton", nullptr, &empty3x3);
120  break;
121  case kEtEtaPhi:
122  lepton_ = std::make_unique<TFitParticleEtEtaPhi>("Lepton", "Lepton", nullptr, &empty3x3);
123  break;
124  case kEtThetaPhi:
125  lepton_ = std::make_unique<TFitParticleEtThetaPhi>("Lepton", "Lepton", nullptr, &empty3x3);
126  break;
127  }
128  switch (metParam_) { // setup neutrino according to parameterization
129  case kEMom:
130  neutrino_ = std::make_unique<TFitParticleEScaledMomDev>("Neutrino", "Neutrino", nullptr, &empty3x3);
131  break;
132  case kEtEtaPhi:
133  neutrino_ = std::make_unique<TFitParticleEtEtaPhi>("Neutrino", "Neutrino", nullptr, &empty3x3);
134  break;
135  case kEtThetaPhi:
136  neutrino_ = std::make_unique<TFitParticleEtThetaPhi>("Neutrino", "Neutrino", nullptr, &empty3x3);
137  break;
138  }
139 }
140 
142  massConstr_[kWHadMass] = std::make_unique<TFitConstraintM>("WMassHad", "WMassHad", nullptr, nullptr, mW_);
143  massConstr_[kWLepMass] = std::make_unique<TFitConstraintM>("WMassLep", "WMassLep", nullptr, nullptr, mW_);
144  massConstr_[kTopHadMass] = std::make_unique<TFitConstraintM>("TopMassHad", "TopMassHad", nullptr, nullptr, mTop_);
145  massConstr_[kTopLepMass] = std::make_unique<TFitConstraintM>("TopMassLep", "TopMassLep", nullptr, nullptr, mTop_);
146  massConstr_[kNeutrinoMass] = std::make_unique<TFitConstraintM>("NeutrinoMass", "NeutrinoMass", nullptr, nullptr, 0.);
148  std::make_unique<TFitConstraintM>("EqualTopMasses", "EqualTopMasses", nullptr, nullptr, 0.);
149  sumPxConstr_ = std::make_unique<TFitConstraintEp>("SumPx", "SumPx", nullptr, TFitConstraintEp::pX, 0.);
150  sumPyConstr_ = std::make_unique<TFitConstraintEp>("SumPy", "SumPy", nullptr, TFitConstraintEp::pY, 0.);
151 
152  massConstr_[kWHadMass]->addParticles1(hadP_.get(), hadQ_.get());
153  massConstr_[kWLepMass]->addParticles1(lepton_.get(), neutrino_.get());
154  massConstr_[kTopHadMass]->addParticles1(hadP_.get(), hadQ_.get(), hadB_.get());
155  massConstr_[kTopLepMass]->addParticles1(lepton_.get(), neutrino_.get(), lepB_.get());
156  massConstr_[kNeutrinoMass]->addParticle1(neutrino_.get());
157  massConstr_[kEqualTopMasses]->addParticles1(hadP_.get(), hadQ_.get(), hadB_.get());
158  massConstr_[kEqualTopMasses]->addParticles2(lepton_.get(), neutrino_.get(), lepB_.get());
159  sumPxConstr_->addParticles(lepton_.get(), neutrino_.get(), hadP_.get(), hadQ_.get(), hadB_.get(), lepB_.get());
160  sumPyConstr_->addParticles(lepton_.get(), neutrino_.get(), hadP_.get(), hadQ_.get(), hadB_.get(), lepB_.get());
161 
162  if (std::find(constrList_.begin(), constrList_.end(), kSumPt) != constrList_.end())
163  constrainSumPt_ = true;
164  constrainSumPt_ = false;
165 }
166 
168  printSetup();
169 
170  setupJets();
171  setupLeptons();
173 
174  // add measured particles
175  fitter_->addMeasParticle(hadB_.get());
176  fitter_->addMeasParticle(hadP_.get());
177  fitter_->addMeasParticle(hadQ_.get());
178  fitter_->addMeasParticle(lepB_.get());
179  fitter_->addMeasParticle(lepton_.get());
180  fitter_->addMeasParticle(neutrino_.get());
181 
182  // add constraints
183  for (unsigned int i = 0; i < constrList_.size(); i++) {
184  if (constrList_[i] != kSumPt)
185  fitter_->addConstraint(massConstr_[constrList_[i]].get());
186  }
187  if (constrainSumPt_) {
188  fitter_->addConstraint(sumPxConstr_.get());
189  fitter_->addConstraint(sumPyConstr_.get());
190  }
191 
192  // initialize helper class used to bring the resolutions into covariance matrices
193  if (!udscResolutions_->empty() && !bResolutions_->empty() && !lepResolutions_->empty() && !metResolutions_->empty())
194  covM_ = std::make_unique<CovarianceMatrix>(*udscResolutions_,
195  *bResolutions_,
200  else
201  covM_ = std::make_unique<CovarianceMatrix>();
202 }
203 
204 int TtSemiLepKinFitter::fit(const TLorentzVector& p4HadP,
205  const TLorentzVector& p4HadQ,
206  const TLorentzVector& p4HadB,
207  const TLorentzVector& p4LepB,
208  const TLorentzVector& p4Lepton,
209  const TLorentzVector& p4Neutrino,
210  const int leptonCharge,
212  // initialize covariance matrices
213  TMatrixD covHadP = covM_->setupMatrix(p4HadP, CovarianceMatrix::kUdscJet, jetParam_);
214  TMatrixD covHadQ = covM_->setupMatrix(p4HadQ, CovarianceMatrix::kUdscJet, jetParam_);
215  TMatrixD covHadB = covM_->setupMatrix(p4HadB, CovarianceMatrix::kBJet, jetParam_);
216  TMatrixD covLepB = covM_->setupMatrix(p4LepB, CovarianceMatrix::kBJet, jetParam_);
217  TMatrixD covLepton = covM_->setupMatrix(p4Lepton, leptonType, lepParam_);
218  TMatrixD covNeutrino = covM_->setupMatrix(p4Neutrino, CovarianceMatrix::kMet, metParam_);
219 
220  // now do the part that is fully independent of PAT features
221  return fit(p4HadP,
222  p4HadQ,
223  p4HadB,
224  p4LepB,
225  p4Lepton,
226  p4Neutrino,
227  covHadP,
228  covHadQ,
229  covHadB,
230  covLepB,
231  covLepton,
232  covNeutrino,
233  leptonCharge);
234 }
235 
236 int TtSemiLepKinFitter::fit(const TLorentzVector& p4HadP,
237  const TLorentzVector& p4HadQ,
238  const TLorentzVector& p4HadB,
239  const TLorentzVector& p4LepB,
240  const TLorentzVector& p4Lepton,
241  const TLorentzVector& p4Neutrino,
242  const TMatrixD& covHadP,
243  const TMatrixD& covHadQ,
244  const TMatrixD& covHadB,
245  const TMatrixD& covLepB,
246  const TMatrixD& covLepton,
247  const TMatrixD& covNeutrino,
248  const int leptonCharge) {
249  // set the kinematics of the objects to be fitted
250  hadP_->setIni4Vec(&p4HadP);
251  hadQ_->setIni4Vec(&p4HadQ);
252  hadB_->setIni4Vec(&p4HadB);
253  lepB_->setIni4Vec(&p4LepB);
254  lepton_->setIni4Vec(&p4Lepton);
255  neutrino_->setIni4Vec(&p4Neutrino);
256 
257  hadP_->setCovMatrix(&covHadP);
258  hadQ_->setCovMatrix(&covHadQ);
259  hadB_->setCovMatrix(&covHadB);
260  lepB_->setCovMatrix(&covLepB);
261  lepton_->setCovMatrix(&covLepton);
262  neutrino_->setCovMatrix(&covNeutrino);
263 
264  if (constrainSumPt_) {
265  // setup Px and Py constraint for curent event configuration so that sum Pt will be conserved
266  sumPxConstr_->setConstraint(p4HadP.Px() + p4HadQ.Px() + p4HadB.Px() + p4LepB.Px() + p4Lepton.Px() +
267  p4Neutrino.Px());
268  sumPyConstr_->setConstraint(p4HadP.Py() + p4HadQ.Py() + p4HadB.Py() + p4LepB.Py() + p4Lepton.Py() +
269  p4Neutrino.Py());
270  }
271 
272  // now do the fit
273  fitter_->fit();
274 
275  // read back the resulting particles if the fit converged
276  if (fitter_->getStatus() == 0) {
277  // read back jet kinematics
279  0,
281  hadP_->getCurr4Vec()->X(), hadP_->getCurr4Vec()->Y(), hadP_->getCurr4Vec()->Z(), hadP_->getCurr4Vec()->E()),
282  math::XYZPoint()));
284  0,
286  hadQ_->getCurr4Vec()->X(), hadQ_->getCurr4Vec()->Y(), hadQ_->getCurr4Vec()->Z(), hadQ_->getCurr4Vec()->E()),
287  math::XYZPoint()));
289  0,
291  hadB_->getCurr4Vec()->X(), hadB_->getCurr4Vec()->Y(), hadB_->getCurr4Vec()->Z(), hadB_->getCurr4Vec()->E()),
292  math::XYZPoint()));
294  0,
296  lepB_->getCurr4Vec()->X(), lepB_->getCurr4Vec()->Y(), lepB_->getCurr4Vec()->Z(), lepB_->getCurr4Vec()->E()),
297  math::XYZPoint()));
298 
299  // read back lepton kinematics
301  math::XYZTLorentzVector(lepton_->getCurr4Vec()->X(),
302  lepton_->getCurr4Vec()->Y(),
303  lepton_->getCurr4Vec()->Z(),
304  lepton_->getCurr4Vec()->E()),
305  math::XYZPoint()));
306 
307  // read back the MET kinematics
309  math::XYZTLorentzVector(neutrino_->getCurr4Vec()->X(),
310  neutrino_->getCurr4Vec()->Y(),
311  neutrino_->getCurr4Vec()->Z(),
312  neutrino_->getCurr4Vec()->E()),
313  math::XYZPoint()));
314  }
315  return fitter_->getStatus();
316 }
317 
319  TtSemiEvtSolution fitsol(*asol);
320 
321  std::vector<pat::Jet> jets;
322  jets.resize(4);
327 
328  // perform the fit, either using the electron or the muon
329  if (fitsol.getDecay() == "electron")
330  fit(jets, fitsol.getCalLepe(), fitsol.getCalLepn());
331  if (fitsol.getDecay() == "muon")
332  fit(jets, fitsol.getCalLepm(), fitsol.getCalLepn());
333 
334  // add fitted information to the solution
335  if (fitter_->getStatus() == 0) {
336  // fill the fitted particles
337  fitsol.setFitHadb(fittedHadB());
338  fitsol.setFitHadp(fittedHadP());
339  fitsol.setFitHadq(fittedHadQ());
340  fitsol.setFitLepb(fittedLepB());
341  fitsol.setFitLepl(fittedLepton());
342  fitsol.setFitLepn(fittedNeutrino());
343  // store the fit's chi2 probability
344  fitsol.setProbChi2(fitProb());
345  }
346  return fitsol;
347 }
std::string param(const Param &param) const
convert Param to human readable form
Definition: TopKinFitter.cc:18
Log< level::Info, true > LogVerbatim
void setFitHadq(const pat::Particle &aFitHadq)
void setFitHadb(const pat::Particle &aFitHadb)
std::map< Constraint, std::unique_ptr< TFitConstraintM > > massConstr_
supported constraints
std::unique_ptr< TAbsFitParticle > hadQ_
Param
supported parameterizations
Definition: TopKinFitter.h:22
pat::Jet getCalHadb() const
const std::vector< edm::ParameterSet > * udscResolutions_
resolutions
void printSetup() const
print fitter setup
pat::Jet getCalHadp() const
std::unique_ptr< TKinFitter > fitter_
kinematic fitter
Definition: TopKinFitter.h:49
std::unique_ptr< TAbsFitParticle > neutrino_
bool constrainSumPt_
internally use simple boolean for this constraint to reduce the per-event computing time ...
pat::Muon getCalLepm() const
const pat::Particle fittedHadQ() const
return hadronic light quark candidate
pat::Jet getCalHadq() const
const pat::Particle fittedNeutrino() const
return neutrino candidate
pat::Electron getCalLepe() const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
int maxNrIter_
maximal allowed number of iterations to be used for the fit
Definition: TopKinFitter.h:51
std::unique_ptr< TAbsFitParticle > lepton_
pat::Particle fittedHadP_
pat::Particle fittedLepB_
pat::MET getCalLepn() const
Param lepParam_
lepton parametrization
pat::Particle fittedHadB_
output particles
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
const std::vector< edm::ParameterSet > * bResolutions_
leptonType
LEPTON
Definition: autophobj.py:48
int fit(const std::vector< pat::Jet > &jets, const pat::Lepton< LeptonType > &leps, const pat::MET &met)
kinematic fit interface for PAT objects
void setupJets()
initialize jet inputs
const std::vector< edm::ParameterSet > * lepResolutions_
pat::Jet getCalLepb() const
const pat::Particle fittedLepB() const
return leptonic b quark candidate
double fitProb() const
return fit probability
Definition: TopKinFitter.h:39
void setFitLepb(const pat::Particle &aFitLepb)
pat::Particle fittedLepton_
void setFitHadp(const pat::Particle &aFitHadp)
const std::vector< double > * jetEnergyResolutionScaleFactors_
scale factors for the jet energy resolution
Param jetParam_
jet parametrization
std::unique_ptr< TFitConstraintEp > sumPyConstr_
double mW_
W mass value used for constraints.
Definition: TopKinFitter.h:57
TtSemiLepKinFitter()
default constructor
Param metParam_
met parametrization
std::string getDecay() const
pat::Particle fittedHadQ_
double maxDeltaS_
maximal allowed chi2 (not normalized to degrees of freedom)
Definition: TopKinFitter.h:53
TtSemiEvtSolution addKinFitInfo(TtSemiEvtSolution *asol)
add kin fit information to the old event solution (in for legacy reasons)
~TtSemiLepKinFitter()
default destructor
void setupConstraints()
initialize constraints
std::vector< Constraint > constrList_
vector of constraints to be used
std::unique_ptr< TAbsFitParticle > lepB_
void setupLeptons()
initialize lepton inputs
void setFitLepn(const pat::Particle &aFitLepn)
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
Analysis-level particle class.
Definition: Particle.h:30
std::unique_ptr< TFitConstraintEp > sumPxConstr_
std::unique_ptr< TAbsFitParticle > hadB_
input particles
const pat::Particle fittedHadP() const
return hadronic light quark candidate
void setProbChi2(double c)
void setupFitter()
setup fitter
pat::Particle fittedNeutrino_
const std::vector< double > * jetEnergyResolutionEtaBinning_
void setFitLepl(const pat::Particle &aFitLepl)
std::unique_ptr< CovarianceMatrix > covM_
object used to construct the covariance matrices for the individual particles
const pat::Particle fittedHadB() const
return hadronic b quark candidate
std::unique_ptr< TAbsFitParticle > hadP_
double mTop_
top mass value used for constraints
Definition: TopKinFitter.h:59
double maxF_
maximal allowed distance from constraints
Definition: TopKinFitter.h:55
const pat::Particle fittedLepton() const
return lepton candidate
const std::vector< edm::ParameterSet > * metResolutions_