CMS 3D CMS Logo

CandKinematicVertexFitter.cc
Go to the documentation of this file.
9 #include <sstream>
10 #include <iostream>
11 using namespace reco;
12 using namespace std;
13 
14 // perform the kinematic fit
15 bool CandKinematicVertexFitter::fit(const vector<RefCountedKinematicParticle> &particles) const {
16  try {
17  tree_ = fitter_.fit(particles);
18  } catch (std::exception &err) {
19  std::cerr << ">>> exception thrown by KinematicParticleVertexFitter:\n"
20  << err.what() << "\n"
21  << ">>> candidate not fitted to common vertex" << std::endl;
22  return false;
23  }
24  //check tree_ is valid here!
25  if (tree_->isValid())
26  return true;
27  else
28  return false;
29 }
30 
31 // main method called by CandProducer sets the VertexCompositeCandidate
33  if (bField_ == nullptr)
35  << "B-Field was not set up CandKinematicVertexFitter.\n"
36  << "the following method must be called before fitting a candidate:\n"
37  << " CandKinematicVertexFitter:.set( const MagneticField * )" << endl;
38  vector<RefCountedKinematicParticle> particles;
39  vector<Candidate *> daughters;
40  vector<RecoCandidate::TrackType> trackTypes;
41  // fill particles with KinematicParticles and daughters with Candidates of the daughters of c
42  fill(particles, daughters, trackTypes, c);
43  assert(particles.size() == daughters.size());
44 
45  // attempt to fit the KinematicParticles, particles
46  if (fit(particles)) {
47  // after the fit, tree_ contains the KinematicTree from the fit
48  tree_->movePointerToTheTop();
49  // set the kinematic properties of the daughters from the fit
50  RefCountedKinematicVertex vertex = tree_->currentDecayVertex();
51  if (vertex->vertexIsValid()) {
52  Candidate::Point vtx(vertex->position());
53  c.setVertex(vtx);
54  vector<RefCountedKinematicParticle> treeParticles = tree_->daughterParticles();
55  vector<RefCountedKinematicParticle>::const_iterator particleIt = treeParticles.begin();
56  vector<Candidate *>::const_iterator daughterIt = daughters.begin(), daughtersEnd = daughters.end();
57  vector<RecoCandidate::TrackType>::const_iterator trackTypeIt = trackTypes.begin();
58  Candidate::LorentzVector mp4(0, 0, 0, 0);
59  for (; daughterIt != daughtersEnd; ++particleIt, ++daughterIt, ++trackTypeIt) {
60  Candidate &daughter = **daughterIt;
61  GlobalVector p3 = (*particleIt)->currentState().globalMomentum();
62  double px = p3.x(), py = p3.y(), pz = p3.z(), p = p3.mag();
63  double energy;
64 
65  if (!daughter.longLived())
66  daughter.setVertex(vtx);
67  double scale;
68  switch (*trackTypeIt) {
70  //gsf used for electron tracks
71  energy = daughter.energy();
72  scale = energy / p;
73  px *= scale;
74  py *= scale;
75  pz *= scale;
76  [[fallthrough]];
77  default:
78  double mass = (*particleIt)->currentState().mass();
79  energy = sqrt(p * p + mass * mass);
80  };
82  daughter.setP4(dp4);
83  mp4 += dp4;
84  }
85  c.setP4(mp4);
86  c.setChi2AndNdof(chi2_ = vertex->chiSquared(), ndof_ = vertex->degreesOfFreedom());
87  GlobalError err = vertex->error();
88  cov_(0, 0) = err.cxx();
89  cov_(0, 1) = err.cyx();
90  cov_(0, 2) = err.czx();
91  cov_(1, 2) = err.czy();
92  cov_(1, 1) = err.cyy();
93  cov_(2, 2) = err.czz();
94  c.setCovariance(cov_);
95  }
96  } else {
97  c.setChi2AndNdof(chi2_ = -1, ndof_ = 0);
98  c.setCovariance(cov_ = CovarianceMatrix(ROOT::Math::SMatrixIdentity()));
99  }
100 }
101 
102 // methond to fill the properties of a CompositeCandidate's daughters
103 void CandKinematicVertexFitter::fill(vector<RefCountedKinematicParticle> &particles,
104  vector<Candidate *> &daughters,
105  vector<RecoCandidate::TrackType> &trackTypes,
106  Candidate &c) const {
107  size_t nDau = c.numberOfDaughters();
108  // loop through CompositeCandidate daughters
109  for (unsigned int j = 0; j < nDau; ++j) {
110  Candidate *d = c.daughter(j);
111  if (d == nullptr) {
112  ostringstream message;
113  message << "Can't access in write mode candidate daughters. "
114  << "pdgId = " << c.pdgId() << ".\n";
115  const Candidate *d1 = c.daughter(j);
116  if (d1 == nullptr)
117  message << "Null daughter also found in read-only mode\n";
118  else
119  message << "Daughter found in read-only mode with id: " << d1->pdgId() << "\n";
120  throw edm::Exception(edm::errors::InvalidReference) << message.str();
121  }
122  //check for a daughter which itself is a composite
123  if (d->numberOfDaughters() > 0) {
124  //try to cast to VertexCompositeCandiate
125  VertexCompositeCandidate *vtxDau = dynamic_cast<VertexCompositeCandidate *>(d);
126  if (vtxDau != nullptr && vtxDau->vertexChi2() > 0) {
127  // if VertexCompositeCandidate refit vtxDau via the set method
128  (*this).set(*vtxDau);
129  // if mass constraint is desired, do it here
130  if (vtxDau->massConstraint()) {
131  KinematicParticleFitter csFitter;
132  //get particle mass from pdg table via pdgid number
133  const ParticleData *data = pdt_->particle(vtxDau->pdgId());
134  ParticleMass mass = data->mass();
135  float mass_sigma = mass * 0.000001; //needs a sigma for the fit
136  // create a KinematicConstraint and refit the tree with it
137  //KinematicConstraint * mass_c = new MassKinematicConstraint(mass,mass_sigma);
138  MassKinematicConstraint mkc(mass, mass_sigma);
139  KinematicConstraint *mass_c(&mkc);
140  tree_ = csFitter.fit(mass_c, tree_);
141  //CHECK THIS! the following works, but might not be safe
142  //tree_ = csFitter.fit(&(MassKinematicConstraint(mass,mass_sigma)),tree_);
143  }
144  // add the kinematic particle from the fit to particles
145  RefCountedKinematicParticle current = (*this).currentParticle();
146  particles.push_back(current);
147  daughters.push_back(d);
148  trackTypes.push_back(RecoCandidate::noTrackType);
149  } else {
150  fill(particles, daughters, trackTypes, *d);
151  }
152  } else {
153  //get track, make KinematicParticle and add to particles so it can be fit
154  TrackRef trk = d->get<TrackRef>();
156  if (!trk.isNull()) {
157  TransientTrack trTrk(trk, bField_);
158  float chi2 = 0, ndof = 0;
159  ParticleMass mass = d->mass();
160  float sigma = mass * 1.e-6;
161  particles.push_back(factory_.particle(trTrk, mass, chi2, ndof, sigma));
162  daughters.push_back(d);
163  trackTypes.push_back(type);
164  } else {
165  cerr << ">>> warning: candidate of type " << d->pdgId() << " has no track reference." << endl;
166  }
167  }
168  }
169 }
double vertexChi2() const override
chi-squares
virtual double energy() const =0
energy
virtual void setP4(const LorentzVector &p4)=0
set 4-momentum
TrackType
track type
Definition: RecoCandidate.h:56
double ParticleMass
Definition: ParticleMass.h:4
math::Error< 5 >::type CovarianceMatrix
assert(be >=bs)
int pdgId() const final
PDG identifier.
bool fit(const std::vector< RefCountedKinematicParticle > &tracks) const
bool massConstraint() const final
do mass constraint?
T sqrt(T t)
Definition: SSEVec.h:19
virtual void setVertex(const Point &vertex)=0
set vertex
bool isNull() const
Checks for null.
Definition: Ref.h:235
HepPDT::ParticleData ParticleData
d
Definition: ztail.py:151
void fill(std::vector< RefCountedKinematicParticle > &, std::vector< reco::Candidate *> &, std::vector< reco::RecoCandidate::TrackType > &, reco::Candidate &) const
std::vector< RefCountedKinematicTree > fit(KinematicConstraint *cs, const std::vector< RefCountedKinematicTree > &trees) const
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
void set(const MagneticField *bField)
fixed size matrix
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:80
math::XYZPoint Point
point in the space
Definition: Candidate.h:40
virtual bool longLived() const =0
is long lived?
static constexpr float d1