CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
CandMassKinFitter Struct Reference

#include <CandMassKinFitter.h>

Inheritance diagram for CandMassKinFitter:
CustomKinFitter

Public Member Functions

 CandMassKinFitter (double m)
 
reco::FitQuality set (reco::Candidate &) const
 

Private Attributes

double mass_
 

Detailed Description

Definition at line 6 of file CandMassKinFitter.h.

Constructor & Destructor Documentation

CandMassKinFitter::CandMassKinFitter ( double  m)
inlineexplicit

Member Function Documentation

FitQuality CandMassKinFitter::set ( reco::Candidate c) const

Definition at line 11 of file CandMassKinFitter.cc.

References TKinFitter::addConstraint(), TKinFitter::addMeasParticle(), TFitConstraintM::addParticle1(), HLT_2018_cff::constraint, reco::Candidate::daughter(), MillePedeFileConverter_cfg::e, runTheMatrix::err, MessageLogger_cfi::errors, TKinFitter::fit(), reco::Candidate::get(), TAbsFitParticle::getCurr4Vec(), TKinFitter::getNDF(), TKinFitter::getS(), mps_fire::i, reco::Candidate::mass(), Skims_PA_cff::name, reco::Candidate::numberOfDaughters(), reco::Candidate::p4(), p4, ecalTrigSettings_cff::particles, TKinFitter::setMaxDeltaS(), TKinFitter::setMaxF(), TKinFitter::setMaxNbIter(), reco::Candidate::setP4(), and TKinFitter::setVerbosity().

11  {
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);
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 }
virtual void setP4(const LorentzVector &p4)=0
set 4-momentum
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector
double p4[4]
Definition: TauolaWrapper.h:92
part
Definition: HCALResponse.h:20
virtual double mass() const =0
mass
const TLorentzVector * getCurr4Vec()
T get() const
get a component
Definition: Candidate.h:222
Definition: errors.py:1
FitQuality
Definition: Fit.h:30
virtual size_type numberOfDaughters() const =0
number of daughters
math::PtEtaPhiELorentzVectorF LorentzVector

Member Data Documentation

double CandMassKinFitter::mass_
private

Definition at line 11 of file CandMassKinFitter.h.