Belle II Software  light-2212-foldex
InternalParticle.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * External Contributor: Wouter Hulsbergen *
5  * *
6  * See git log for contributors and copyright holders. *
7  * This file is licensed under LGPL-3.0, see LICENSE.md. *
8  **************************************************************************/
9 
10 #include <analysis/dataobjects/Particle.h>
11 
12 #include <analysis/VertexFitting/TreeFitter/InternalParticle.h>
13 #include <analysis/VertexFitting/TreeFitter/FitParams.h>
14 #include <analysis/VertexFitting/TreeFitter/HelixUtils.h>
15 #include <framework/logging/Logger.h>
16 #include <mdst/dataobjects/Track.h>
17 
18 using std::vector;
19 
20 namespace TreeFitter {
21 
22  inline bool sortByType(const ParticleBase* lhs, const ParticleBase* rhs)
23  {
24  int lhstype = lhs->type() ;
25  int rhstype = rhs->type() ;
26  bool rc = false ;
27  if (lhstype == rhstype &&
28  lhstype == ParticleBase::TFParticleType::kRecoTrack) {
29 
30  rc = lhs->particle()->getMomentum().Rho() > rhs->particle()->getMomentum().Rho();
31  } else if (lhs->particle() && rhs->particle() && lhs->particle()->getNDaughters() > 0 &&
32  rhs->particle()->getNDaughters() > 0) {
33  rc = lhs->nFinalChargedCandidates() > rhs->nFinalChargedCandidates();
34  } else {
35  rc = lhstype < rhstype;
36  }
37  return rc;
38  }
39 
41  const ParticleBase* mother,
42  const ConstraintConfiguration& config,
43  bool forceFitAll
44  ) :
45  ParticleBase(particle, mother, &config),// config pointer here to allow final states not to have it
46  m_massconstraint(false),
47  m_beamconstraint(false),
48  m_lifetimeconstraint(false),
49  m_isconversion(false),
50  m_automatic_vertex_constraining(config.m_automatic_vertex_constraining)
51  {
52  if (particle) {
53  for (auto daughter : particle->getDaughters()) {
54  addDaughter(daughter, config, forceFitAll);
55  }
56  } else {
57  B2ERROR("Trying to create an InternalParticle from NULL. This should never happen.");
58  }
59 
60  m_massconstraint = std::find(config.m_massConstraintListPDG.begin(), config.m_massConstraintListPDG.end(),
61  std::abs(m_particle->getPDGCode())) != config.m_massConstraintListPDG.end();
62 
63  m_beamconstraint = (std::abs(m_particle->getPDGCode()) == config.m_beamConstraintPDG);
64 
66  // if this is a hadronically decaying resonance it is useful to constraint the decay vertex to its mothers decay vertex.
67  //
69  config.m_fixedToMotherVertexListPDG.end(),
70  std::abs(m_particle->getPDGCode())) != config.m_fixedToMotherVertexListPDG.end() && this->mother();
71 
72  // use geo constraint if this particle is in the list to constrain
73  m_geo_constraint = std::find(config.m_geoConstraintListPDG.begin(),
74  config.m_geoConstraintListPDG.end(),
75  std::abs(m_particle->getPDGCode())) != config.m_geoConstraintListPDG.end() && this->mother() && !m_shares_vertex_with_mother;
76  } else {
79  }
80  }
81 
83  {
84 
85  return lhs->particle()->getMomentum().Rho() > rhs->particle()->getMomentum().Rho();
86  }
87 
89  {
90  ErrCode status ;
91  for (auto daughter : m_daughters) {
92  status |= daughter->initMotherlessParticle(fitparams);
93  }
94 
95  int posindex = posIndex();
96  // FIXME update this and check if this position is already initialised
97  // lower stream vertices might be better
98  if (hasPosition()) {
99  fitparams.getStateVector().segment(posindex, 3) = Eigen::Matrix<double, 3, 1>::Zero(3);
100 
101  std::vector<ParticleBase*> alldaughters;
102  ParticleBase::collectVertexDaughters(alldaughters, posindex);
103 
104  std::vector<ParticleBase*> vtxdaughters;
105 
106  vector<RecoTrack*> trkdaughters;
107  for (auto daughter : alldaughters) {
108  if (daughter->type() == ParticleBase::TFParticleType::kRecoTrack) {
109  trkdaughters.push_back(static_cast<RecoTrack*>(daughter));
110  } else if (daughter->hasPosition()
111  && fitparams.getStateVector()(daughter->posIndex()) != 0) {
112  vtxdaughters.push_back(daughter);
113  }
114  }
115 
117 
118  if (trkdaughters.size() >= 2) {
119  std::sort(trkdaughters.begin(), trkdaughters.end(), compTrkTransverseMomentum);
120 
121  RecoTrack* dau1 = trkdaughters[0];
122  RecoTrack* dau2 = trkdaughters[1];
123 
125  dau1->particle()->getPDGCode())))->getHelix();
127  dau2->particle()->getPDGCode())))->getHelix();
128 
129  double flt1(0), flt2(0);
130  HelixUtils::helixPoca(helix1, helix2, flt1, flt2, v, m_isconversion);
131 
132  fitparams.getStateVector()(posindex) = v.X();
133  fitparams.getStateVector()(posindex + 1) = v.Y();
134  fitparams.getStateVector()(posindex + 2) = v.Z();
135 
136  dau1->setFlightLength(flt1);
137  dau2->setFlightLength(flt2);
138 
139  } else if (false && trkdaughters.size() + vtxdaughters.size() >= 2) {
140  // TODO switched off waiting for refactoring of init1 and init2 functions (does not affect performance)
141  } else if (mother() && mother()->posIndex() >= 0) {
142  const int posindexmother = mother()->posIndex();
143  const int dim = m_config->m_originDimension; //TODO access mother
144  fitparams.getStateVector().segment(posindex, dim) = fitparams.getStateVector().segment(posindexmother, dim);
145  } else {
147  fitparams.getStateVector().segment(posindex, 3) = Eigen::Matrix<double, 1, 3>::Zero(3);
148  }
149  }
150 
151  for (auto daughter : m_daughters) {
152  daughter->initParticleWithMother(fitparams);
153  }
154 
155  initMomentum(fitparams);
156  return status;
157  }
158 
160  {
161  int posindex = posIndex();
162  if (hasPosition() &&
163  mother() &&
164  fitparams.getStateVector()(posindex) == 0 &&
165  fitparams.getStateVector()(posindex + 1) == 0 && \
166  fitparams.getStateVector()(posindex + 2) == 0) {
167  const int posindexmom = mother()->posIndex();
168  const int dim = m_config->m_originDimension; //TODO access mother?
169  fitparams.getStateVector().segment(posindex, dim) = fitparams.getStateVector().segment(posindexmom, dim);
170  }
171  return initTau(fitparams);
172  }
173 
175  {
176  int momindex = momIndex();
177  fitparams.getStateVector().segment(momindex, 4) = Eigen::Matrix<double, 4, 1>::Zero(4);
178 
179  for (auto daughter : m_daughters) {
180  int daumomindex = daughter->momIndex();
181  int maxrow = daughter->hasEnergy() ? 4 : 3;
182 
183  double e2 = fitparams.getStateVector().segment(daumomindex, maxrow).squaredNorm();
184  fitparams.getStateVector().segment(momindex, maxrow) += fitparams.getStateVector().segment(daumomindex, maxrow);
185 
186  if (maxrow == 3) {
187  double mass = daughter->particle()->getPDGMass();
188  fitparams.getStateVector()(momindex + 3) += std::sqrt(e2 + mass * mass);
189  }
190  }
191  return ErrCode(ErrCode::Status::success);
192  }
193 
195  {
196  ErrCode status;
197  ParticleBase::initCovariance(fitparams);
198  for (auto daughter : m_daughters) {
199  status |= daughter->initCovariance(fitparams);
200  }
201  return status;
202  }
203 
204 
206  Projection& p) const
207  {
208  const int momindex = momIndex();
209 
210  // `this` always has an energy row
211  p.getResiduals().segment(0, 4) = fitparams.getStateVector().segment(momindex, 4);
212 
213  for (int imom = 0; imom < 4; ++imom) {
214  p.getH()(imom, momindex + imom) = 1;
215  }
216 
217  for (const auto daughter : m_daughters) {
218  const int daumomindex = daughter->momIndex();
219  const Eigen::Matrix<double, 1, 3> p3_vec = fitparams.getStateVector().segment(daumomindex, 3);
220 
221  // three momentum is easy just subtract the vectors
222  p.getResiduals().segment(0, 3) -= p3_vec;
223 
224  // energy depends on the parametrisation!
225  if (daughter->hasEnergy()) {
226  p.getResiduals()(3) -= fitparams.getStateVector()(daumomindex + 3);
227  p.getH()(3, daumomindex + 3) = -1; // d/dE -E
228  } else {
229  // m^2 + p^2 = E^2
230  // so
231  // E = sqrt(m^2 + p^2)
232  const double mass = daughter->particle()->getPDGMass();
233  const double p2 = p3_vec.squaredNorm();
234  const double energy = std::sqrt(mass * mass + p2);
235  p.getResiduals()(3) -= energy;
236 
237  for (unsigned i = 0; i < 3; ++i) {
238  // d/dpx_i sqrt(m^2 + p^2)
239  p.getH()(3, daumomindex + i) = -1 * p3_vec(i) / energy;
240  }
241  }
242 
243  // this has to be in any case
244  // d/dp_i p_i
245  for (unsigned i = 0; i < 3; ++i) {
246  p.getH()(i, daumomindex + i) = -1;
247  }
248  }
249  return ErrCode(ErrCode::Status::success);
250  }
251 
252 
254  Projection& p) const
255  {
256 
257  const int momindex = momIndex() ;
258 
259  const Eigen::Matrix<double, 4, 1> fitMomE = fitparams.getStateVector().segment(momindex, 4);
260 
261  p.getResiduals() = m_config->m_beamMomE - fitMomE;
262 
263  for (int row = 0; row < 4; ++row) {
264  p.getH()(row, momindex + row) = -1;
265  }
266 
267  p.getV() = m_config->m_beamCovariance;
268 
269  return ErrCode(ErrCode::Status::success) ;
270  }
271 
272 
274  const FitParams& fitparams,
275  Projection& p) const
276  {
277  ErrCode status;
278  switch (type) {
279  case Constraint::mass:
280  status |= projectMassConstraint(fitparams, p);
281  break;
282  case Constraint::geometric:
283  status |= projectGeoConstraint(fitparams, p);
284  break;
285  case Constraint::kinematic:
286  status |= projectKineConstraint(fitparams, p);
287  break;
288  case Constraint::beam:
289  status |= projectBeamConstraint(fitparams, p);
290  break;
291  default:
292  status |= ParticleBase::projectConstraint(type, fitparams, p);
293  }
294 
295  return status;
296  }
297 
299  {
300  // for example B0 and D* can share the same vertex
301  return m_shares_vertex_with_mother ? this->mother()->posIndex() : index();
302  }
303 
305  {
306  // { x, y, z, tau, px, py, pz, E }
307  // the last 4 always exists for composite particles
308  // tau index only exist with a vertex and geo constraint
309  //
310  if (m_shares_vertex_with_mother) { return 4; }
311  else if (!m_geo_constraint) { return 7; }
312  else { return 8; }
313  }
314 
316  {
318  return m_geo_constraint ? index() + 3 : -1;
319  }
320 
322  {
327  if (m_geo_constraint && !m_shares_vertex_with_mother) { return this->index() + 4; }
328 
329  if (m_shares_vertex_with_mother) { return this->index(); }
330 
331  if (!m_geo_constraint) {return index() + 3 ;}
332 
333  // this will crash the initialisation
334  return -1;
335  }
336 
338  {
339  // does this particle have a position index (decayvertex)
340  // true in any case
341  return true;
342  }
343 
345  int depth) const
346  {
347 
348  for (auto daughter : m_daughters) {
349  daughter->addToConstraintList(list, depth - 1);
350  }
351  if (tauIndex() >= 0 && m_lifetimeconstraint) {
352  list.push_back(Constraint(this, Constraint::lifetime, depth, 1));
353  }
354  if (momIndex() >= 0) {
355  list.push_back(Constraint(this, Constraint::kinematic, depth, 4, 3));
356  }
357  if (m_geo_constraint) {
358  assert(m_config);
359  const int dim = m_config->m_originDimension == 2 && std::abs(m_particle->getPDGCode()) == m_config->m_headOfTreePDG ? 2 : 3;
360  list.push_back(Constraint(this, Constraint::geometric, depth, dim, 3));
361  }
362  if (m_massconstraint) {
363  list.push_back(Constraint(this, Constraint::mass, depth, 1, 3));
364  }
365  if (m_beamconstraint) {
366  assert(m_config);
367  list.push_back(Constraint(this, Constraint::beam, depth, 4, 3));
368  }
369 
370  }
371 
372  std::string InternalParticle::parname(int thisindex) const
373  {
374  int id = thisindex;
375  if (!mother() && id >= 3) {++id;}
376  return ParticleBase::parname(id);
377  }
378 
380  {
381  for (const auto daughter : m_daughters) {
382  daughter->forceP4Sum(fitparams);
383  }
384  const int momindex = momIndex();
385  if (momindex > 0) {
386  const int dim = hasEnergy() ? 4 : 3;
387  Projection p(fitparams.getDimensionOfState(), dim);
388  projectKineConstraint(fitparams, p);
389  fitparams.getStateVector().segment(momindex, dim) -= p.getResiduals().segment(0, dim);
390  }
391  }
392 
393 }
394 
DataType Z() const
access variable Z (= .at(2) without boundary check)
Definition: B2Vector3.h:435
DataType X() const
access variable X (= .at(0) without boundary check)
Definition: B2Vector3.h:431
DataType Y() const
access variable Y (= .at(1) without boundary check)
Definition: B2Vector3.h:433
Provides a type-safe way to pass members of the chargedStableSet set.
Definition: Const.h:580
This class represents an ideal helix in perigee parameterization.
Definition: Helix.h:59
Class to store reconstructed particles.
Definition: Particle.h:74
const Track * getTrack() const
Returns the pointer to the Track object that was used to create this Particle (ParticleType == c_Trac...
Definition: Particle.cc:875
int getPDGCode(void) const
Returns PDG code.
Definition: Particle.h:441
std::vector< Belle2::Particle * > getDaughters() const
Returns a vector of pointers to daughter particles.
Definition: Particle.cc:669
ROOT::Math::XYZVector getMomentum() const
Returns momentum vector.
Definition: Particle.h:541
Helix getHelix() const
Conversion to framework Helix (without covariance).
const TrackFitResult * getTrackFitResultWithClosestMass(const Const::ChargedStable &requestedType) const
Return the track fit for a fit hypothesis with the closest mass.
Definition: Track.cc:104
const int m_originDimension
dimension of the origin constraint and ALL geometric gcosntraints
Eigen::Matrix< double, 4, 1 > m_beamMomE
Beam four-momentum.
const std::vector< int > m_geoConstraintListPDG
list of pdg codes to mass constrain
const std::vector< int > m_fixedToMotherVertexListPDG
list of pdg codes to mass constrain
int m_headOfTreePDG
PDG code of the head particle.
const std::vector< int > m_massConstraintListPDG
list of pdg codes to mass constrain
const int m_beamConstraintPDG
PDG code to beam constraint.
Eigen::Matrix< double, 4, 4 > m_beamCovariance
Beam Covariance.
class to manage the order of constraints and their filtering
Definition: Constraint.h:20
Type
type of constraints the order of these constraints is important: it is the order in which they are ap...
Definition: Constraint.h:27
abstract errorocode be aware that the default is success
Definition: ErrCode.h:14
Class to store and manage fitparams (statevector)
Definition: FitParams.h:20
Eigen::Matrix< double, -1, 1, 0, MAX_MATRIX_SIZE, 1 > & getStateVector()
getter for the fit parameters/statevector
Definition: FitParams.h:65
int getDimensionOfState() const
get the states dimension
Definition: FitParams.h:89
static double helixPoca(const Belle2::Helix &helix1, const Belle2::Helix &helix2, double &flt1, double &flt2, Belle2::B2Vector3D &vertex, bool parallel=false)
POCA between two tracks.
Definition: HelixUtils.cc:210
virtual void forceP4Sum(FitParams &) const override
enforce conservation of momentum sum
virtual int tauIndex() const override
tau index in fit params only if it has a mother
bool m_isconversion
is conversion
virtual std::string parname(int index) const override
name
ErrCode projectKineConstraint(const FitParams &, Projection &) const
project kinematical constraint
virtual int dim() const override
space reserved in fit params, if has mother then it has tau
bool m_beamconstraint
has beam constraint
bool m_massconstraint
has mass constraint
ErrCode projectBeamConstraint(const FitParams &, Projection &) const
project beam four momentum constraint
virtual ErrCode initCovariance(FitParams &) const override
init covariance
virtual int type() const override
type
bool m_automatic_vertex_constraining
automatically figure out if mother and particle vertex should be the same and also add geometric cons...
virtual bool hasEnergy() const override
has energy in fitparams
virtual bool hasPosition() const override
has position index
ErrCode projectConstraint(const Constraint::Type type, const FitParams &fitparams, Projection &p) const override
find out which constraint it is and project
virtual ErrCode initParticleWithMother(FitParams &fitparams) override
init particle in case it has a mother
virtual int momIndex() const override
momentum index in fit params depending on whether it has a mother
ErrCode initMomentum(FitParams &fitparams) const
init momentum of *this and daughters
virtual void addToConstraintList(constraintlist &list, int depth) const override
add to constraint list
virtual ErrCode initMotherlessParticle(FitParams &fitparams) override
init particle in case it has no mother
static bool compTrkTransverseMomentum(const RecoTrack *lhs, const RecoTrack *rhs)
compare transverse track momentum
virtual int posIndex() const override
position index in fit params
InternalParticle(Belle2::Particle *particle, const ParticleBase *mother, const ConstraintConfiguration &config, bool forceFitAll)
constructor
bool m_lifetimeconstraint
has lifetime constraint
bool m_geo_constraint
use a geo metric constraint
bool m_shares_vertex_with_mother
shares vertex with mother, that means decay vertex = productionvertex
base class for all particles
Definition: ParticleBase.h:25
Belle2::Particle * particle() const
get basf2 particle
Definition: ParticleBase.h:96
ErrCode initTau(FitParams &par) const
initialises tau as a length
virtual ParticleBase * addDaughter(Belle2::Particle *, const ConstraintConfiguration &config, bool forceFitAll=false)
add daughter
Definition: ParticleBase.cc:82
virtual ErrCode projectMassConstraint(const FitParams &, Projection &) const
project mass constraint abstract
const ConstraintConfiguration * m_config
has all the constraint config
Definition: ParticleBase.h:208
virtual std::string parname(int index) const
get name of parameter i
virtual ErrCode projectConstraint(Constraint::Type, const FitParams &, Projection &) const
project constraint.
bool m_isStronglyDecayingResonance
decay length less than 1 micron
Definition: ParticleBase.h:205
virtual ErrCode initCovariance(FitParams &) const
init covariance matrix
Belle2::Particle * m_particle
pointer to framework type
Definition: ParticleBase.h:196
virtual ErrCode projectGeoConstraint(const FitParams &, Projection &) const
project geometrical constraint
virtual int posIndex() const
get vertex index (in statevector!)
Definition: ParticleBase.h:126
int index() const
get index
Definition: ParticleBase.h:99
std::vector< ParticleBase * > m_daughters
daughter container
Definition: ParticleBase.h:202
void collectVertexDaughters(std::vector< ParticleBase * > &particles, int posindex)
get vertex daughters
std::vector< Constraint > constraintlist
alias
Definition: ParticleBase.h:56
const ParticleBase * mother() const
getMother() / hasMother()
Definition: ParticleBase.h:102
class to store the projected residuals and the corresponding jacobian as well as the covariance matri...
Definition: Projection.h:18
representation of all charged final states FIXME rename since this name is taken in tracking
Definition: RecoTrack.h:18
void setFlightLength(double flt)
setter for the flight length
Definition: RecoTrack.h:64