Belle II Software  light-2205-abys
MassFourCFitKFit.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 
10 #include <cstdio>
11 
12 #include <TMatrixFSym.h>
13 
14 #include <analysis/VertexFitting/KFit/MassFourCFitKFit.h>
15 #include <analysis/VertexFitting/KFit/MakeMotherKFit.h>
16 #include <analysis/utility/CLHEPToROOT.h>
17 
18 
19 using namespace std;
20 using namespace Belle2;
21 using namespace Belle2::analysis;
22 using namespace CLHEP;
23 using namespace ROOT::Math;
24 
25 MassFourCFitKFit::MassFourCFitKFit() : m_AfterVertexError(HepSymMatrix(3, 0)),
26  m_FourMomentum(PxPyPzEVector())
27 {
28  m_FlagFitted = false;
29  m_FlagTrackVertexError = false;
31  m_FlagAtDecayPoint = true;
33  m_InvariantMass = -1.0;
37 }
38 
39 
41 
43 MassFourCFitKFit::addMassConstraint(const double m, std::vector<unsigned>& childTrackId) {
44  if ((childTrackId.back() - childTrackId.front()) != childTrackId.size() - 1)
45  {
47  }
48  m_ConstraintMassList.push_back(m);
49  m_ConstraintMassChildLists.push_back(std::make_pair(childTrackId.front(), childTrackId.back()));
51 }
52 
54 MassFourCFitKFit::setVertex(const HepPoint3D& v) {
55  m_BeforeVertex = v;
56 
58 }
59 
60 
62 MassFourCFitKFit::setVertexError(const HepSymMatrix& e) {
63  if (e.num_row() != 3)
64  {
66  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
67  return m_ErrorCode;
68  }
69 
72 
74 }
75 
76 
78 MassFourCFitKFit::setInvariantMass(const double m) {
79  m_InvariantMass = m;
80 
82 }
83 
84 
86 MassFourCFitKFit::setFourMomentum(const PxPyPzEVector& m) {
87  m_FourMomentum = m;
88 
90 }
91 
92 
95  m_FlagAtDecayPoint = flag;
96 
98 }
99 
100 
101 enum KFitError::ECode
103  m_IsFixMass.push_back(true);
104 
106 }
107 
108 
109 enum KFitError::ECode
111  m_IsFixMass.push_back(false);
112 
114 }
115 
116 
117 enum KFitError::ECode
118 MassFourCFitKFit::setTrackVertexError(const HepMatrix& e) {
119  if (e.num_row() != 3 || e.num_col() != KFitConst::kNumber7)
120  {
122  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
123  return m_ErrorCode;
124  }
125 
126  m_BeforeTrackVertexError.push_back(e);
127  m_FlagTrackVertexError = true;
129 
131 }
132 
133 
134 enum KFitError::ECode
136  HepMatrix zero(3, KFitConst::kNumber7, 0);
137 
138  return this->setTrackVertexError(zero);
139 }
140 
141 
142 enum KFitError::ECode
143 MassFourCFitKFit::setCorrelation(const HepMatrix& m) {
144  return KFitBase::setCorrelation(m);
145 }
146 
147 
148 enum KFitError::ECode
151 }
152 
153 
154 const HepPoint3D
155 MassFourCFitKFit::getVertex(const int flag) const
156 {
157  if (flag == KFitConst::kAfterFit && !isFitted()) return HepPoint3D();
158 
159  switch (flag) {
161  return m_BeforeVertex;
162 
164  return m_AfterVertex;
165 
166  default:
167  KFitError::displayError(__FILE__, __LINE__, __func__, KFitError::kOutOfRange);
168  return HepPoint3D();
169  }
170 }
171 
172 
173 const HepSymMatrix
174 MassFourCFitKFit::getVertexError(const int flag) const
175 {
176  if (flag == KFitConst::kAfterFit && !isFitted()) return HepSymMatrix(3, 0);
177 
178  if (flag == KFitConst::kBeforeFit)
179  return m_BeforeVertexError;
181  return m_AfterVertexError;
182  else {
183  KFitError::displayError(__FILE__, __LINE__, __func__, KFitError::kOutOfRange);
184  return HepSymMatrix(3, 0);
185  }
186 }
187 
188 
189 double
191 {
192  return m_InvariantMass;
193 }
194 
195 
196 bool
198 {
199  return m_FlagAtDecayPoint;
200 }
201 
202 
203 bool
205 {
207 }
208 
209 
210 double
212 {
213  return m_CHIsq;
214 }
215 
216 
217 const HepMatrix
218 MassFourCFitKFit::getTrackVertexError(const int id, const int flag) const
219 {
220  if (flag == KFitConst::kAfterFit && !isFitted()) return HepMatrix(3, KFitConst::kNumber7, 0);
221  if (!isTrackIDInRange(id)) return HepMatrix(3, KFitConst::kNumber7, 0);
222 
223  if (flag == KFitConst::kBeforeFit)
224  return m_BeforeTrackVertexError[id];
226  return m_AfterTrackVertexError[id];
227  else {
228  KFitError::displayError(__FILE__, __LINE__, __func__, KFitError::kOutOfRange);
229  return HepMatrix(3, KFitConst::kNumber7, 0);
230  }
231 }
232 
233 
234 double
236 {
237  if (!isFitted()) return -1;
238  if (!isTrackIDInRange(id)) return -1;
239 
240  if (m_IsFixMass[id]) {
241 
242  HepMatrix da(m_Tracks[id].getFitParameter(KFitConst::kBeforeFit) - m_Tracks[id].getFitParameter(KFitConst::kAfterFit));
243  int err_inverse = 0;
244  const double chisq = (da.T() * (m_Tracks[id].getFitError(KFitConst::kBeforeFit).inverse(err_inverse)) * da)[0][0];
245 
246  if (err_inverse) {
247  KFitError::displayError(__FILE__, __LINE__, __func__, KFitError::kCannotGetMatrixInverse);
248  return -1;
249  }
250 
251  return chisq;
252 
253  } else {
254 
255  HepMatrix da(m_Tracks[id].getMomPos(KFitConst::kBeforeFit) - m_Tracks[id].getMomPos(KFitConst::kAfterFit));
256  int err_inverse = 0;
257  const double chisq = (da.T() * (m_Tracks[id].getError(KFitConst::kBeforeFit).inverse(err_inverse)) * da)[0][0];
258 
259  if (err_inverse) {
260  KFitError::displayError(__FILE__, __LINE__, __func__, KFitError::kCannotGetMatrixInverse);
261  return -1;
262  }
263 
264  return chisq;
265 
266  }
267 }
268 
269 
270 const HepMatrix
271 MassFourCFitKFit::getCorrelation(const int id1, const int id2, const int flag) const
272 {
273  if (flag == KFitConst::kAfterFit && !isFitted()) return HepMatrix(KFitConst::kNumber7, KFitConst::kNumber7, 0);
274  if (!isTrackIDInRange(id1)) return HepMatrix(KFitConst::kNumber7, KFitConst::kNumber7, 0);
275  if (!isTrackIDInRange(id2)) return HepMatrix(KFitConst::kNumber7, KFitConst::kNumber7, 0);
276 
277  switch (flag) {
279  return KFitBase::getCorrelation(id1, id2, flag);
280 
282  return makeError3(
283  this->getTrackMomentum(id1),
284  this->getTrackMomentum(id2),
285  m_V_al_1.sub(KFitConst::kNumber7 * id1 + 1, KFitConst::kNumber7 * (id1 + 1), KFitConst::kNumber7 * id2 + 1,
286  KFitConst::kNumber7 * (id2 + 1)),
287  m_IsFixMass[id1],
288  m_IsFixMass[id2]);
289 
290  default:
291  KFitError::displayError(__FILE__, __LINE__, __func__, KFitError::kOutOfRange);
292  return HepMatrix(KFitConst::kNumber7, KFitConst::kNumber7, 0);
293  }
294 }
295 
296 
297 enum KFitError::ECode
299  return KFitBase::doFit1();
300 }
301 
302 
303 enum KFitError::ECode
306  m_d = HepMatrix(4 + m_ConstraintMassCount, 1, 0);
307  m_V_D = HepMatrix(4 + m_ConstraintMassCount, 4 + m_ConstraintMassCount, 0);
308  m_lam = HepMatrix(4 + m_ConstraintMassCount, 1, 0);
309 
311  {
313  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
314  return m_ErrorCode;
315  }
316 
317  if (m_IsFixMass.size() == 0)
318  {
319  // If no fix_mass flag at all,
320  // all tracks are considered to be fixed at mass.
321  for (int i = 0; i < m_TrackCount; i++) this->fixMass();
322  } else if (m_IsFixMass.size() != (unsigned int)m_TrackCount)
323  {
325  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
326  return m_ErrorCode;
327  }
328 
330  {
331  int index = 0;
332  m_al_0 = HepMatrix(KFitConst::kNumber7 * m_TrackCount, 1, 0);
333  m_property = HepMatrix(m_TrackCount, 3, 0);
334  m_V_al_0 = HepSymMatrix(KFitConst::kNumber7 * m_TrackCount, 0);
335 
336  for (auto& track : m_Tracks) {
337  // momentum x,y,z and position x,y,z
338  m_al_0[index * KFitConst::kNumber7 + 0][0] = track.getMomentum(KFitConst::kBeforeFit).x();
339  m_al_0[index * KFitConst::kNumber7 + 1][0] = track.getMomentum(KFitConst::kBeforeFit).y();
340  m_al_0[index * KFitConst::kNumber7 + 2][0] = track.getMomentum(KFitConst::kBeforeFit).z();
341  m_al_0[index * KFitConst::kNumber7 + 3][0] = track.getMomentum(KFitConst::kBeforeFit).t();
342  m_al_0[index * KFitConst::kNumber7 + 4][0] = track.getPosition(KFitConst::kBeforeFit).x();
343  m_al_0[index * KFitConst::kNumber7 + 5][0] = track.getPosition(KFitConst::kBeforeFit).y();
344  m_al_0[index * KFitConst::kNumber7 + 6][0] = track.getPosition(KFitConst::kBeforeFit).z();
345  // these error
346  m_V_al_0.sub(index * KFitConst::kNumber7 + 1, track.getError(KFitConst::kBeforeFit));
347  // charge, mass, a
348  m_property[index][0] = track.getCharge();
349  m_property[index][1] = track.getMass();
350  const double c = KFitConst::kLightSpeed; // C++ bug?
351  // m_property[index][2] = -KFitConst::kLightSpeed * m_MagneticField * it->getCharge();
352  m_property[index][2] = -c * m_MagneticField * track.getCharge();
353  index++;
354  }
355 
356  // error between track and track
357  if (m_FlagCorrelation) {
358  this->prepareCorrelation();
360  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
361  return m_ErrorCode;
362  }
363  }
364 
365  // set member matrix
366  m_al_1 = m_al_0;
367 
368  // define size of matrix
371 
372  } else
373  {
374  //TODO: Not Implemented
376  // m_FlagFitIncludingVertex == true
377  int index = 0;
378  m_al_0 = HepMatrix(KFitConst::kNumber7 * m_TrackCount + 3, 1, 0);
379  m_property = HepMatrix(m_TrackCount, 3, 0);
380  m_V_al_0 = HepSymMatrix(KFitConst::kNumber7 * m_TrackCount + 3, 0);
381 
382  for (auto& track : m_Tracks) {
383  // momentum x,y,z and position x,y,z
384  m_al_0[index * KFitConst::kNumber7 + 0][0] = track.getMomentum(KFitConst::kBeforeFit).x();
385  m_al_0[index * KFitConst::kNumber7 + 1][0] = track.getMomentum(KFitConst::kBeforeFit).y();
386  m_al_0[index * KFitConst::kNumber7 + 2][0] = track.getMomentum(KFitConst::kBeforeFit).z();
387  m_al_0[index * KFitConst::kNumber7 + 3][0] = track.getMomentum(KFitConst::kBeforeFit).t();
388  m_al_0[index * KFitConst::kNumber7 + 4][0] = track.getPosition(KFitConst::kBeforeFit).x();
389  m_al_0[index * KFitConst::kNumber7 + 5][0] = track.getPosition(KFitConst::kBeforeFit).y();
390  m_al_0[index * KFitConst::kNumber7 + 6][0] = track.getPosition(KFitConst::kBeforeFit).z();
391  // these error
392  m_V_al_0.sub(index * KFitConst::kNumber7 + 1, track.getError(KFitConst::kBeforeFit));
393  // charge, mass, a
394  m_property[index][0] = track.getCharge();
395  m_property[index][1] = track.getMass();
396  const double c = KFitConst::kLightSpeed; // C++ bug?
397  // m_property[index][2] = -KFitConst::kLightSpeed * m_MagneticField * it->getCharge();
398  m_property[index][2] = -c * m_MagneticField * track.getCharge();
399  index++;
400  }
401 
402  // vertex
407 
408  // error between track and track
409  if (m_FlagCorrelation) {
410  this->prepareCorrelation();
412  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
413  return m_ErrorCode;
414  }
415  }
416 
417  // set member matrix
418  m_al_1 = m_al_0;
419 
420  // define size of matrix
422  m_D = m_V_al_1.sub(1, 4, 1, KFitConst::kNumber7 * m_TrackCount + 3);
423  }
424 
426 }
427 
428 
429 enum KFitError::ECode
431  char buf[1024];
432  sprintf(buf, "%s:%s(): internal error; this function should never be called", __FILE__, __func__);
433  B2FATAL(buf);
434 
435  /* NEVER REACHEd HERE */
436  return KFitError::kOutOfRange;
437 }
438 
439 
440 enum KFitError::ECode
442  if (m_BeforeCorrelation.size() != static_cast<unsigned int>(m_TrackCount * (m_TrackCount - 1) / 2))
443  {
445  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
446  return m_ErrorCode;
447  }
448 
449  int row = 0, col = 0;
450 
451  for (auto& hm : m_BeforeCorrelation)
452  {
453  // counter
454  row++;
455  if (row == m_TrackCount) {
456  col++;
457  row = col + 1;
458  }
459 
460  int ii = 0, jj = 0;
461  for (int i = KFitConst::kNumber7 * row; i < KFitConst::kNumber7 * (row + 1); i++) {
462  for (int j = KFitConst::kNumber7 * col; j < KFitConst::kNumber7 * (col + 1); j++) {
463  m_V_al_0[i][j] = hm[ii][jj];
464  jj++;
465  }
466  jj = 0;
467  ii++;
468  }
469  }
470 
472  {
473  //TODO: Not Implemented
475 
476  // ...error of vertex
478 
479  // ...error matrix between vertex and tracks
481  if (m_BeforeTrackVertexError.size() != (unsigned int)m_TrackCount) {
483  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
484  return m_ErrorCode;
485  }
486 
487  int i = 0;
488  for (auto& hm : m_BeforeTrackVertexError) {
489  for (int j = 0; j < 3; j++) for (int k = 0; k < KFitConst::kNumber7; k++) {
491  }
492  i++;
493  }
494  }
495  }
496 
498 }
499 
500 
501 enum KFitError::ECode
503  Hep3Vector h3v;
504  int index = 0;
505  for (auto& pdata : m_Tracks)
506  {
507  // tracks
508  // momentum
509  h3v.setX(m_al_1[index * KFitConst::kNumber7 + 0][0]);
510  h3v.setY(m_al_1[index * KFitConst::kNumber7 + 1][0]);
511  h3v.setZ(m_al_1[index * KFitConst::kNumber7 + 2][0]);
512  if (m_IsFixMass[index])
513  pdata.setMomentum(HepLorentzVector(h3v, sqrt(h3v.mag2() + pdata.getMass()*pdata.getMass())), KFitConst::kAfterFit);
514  else
515  pdata.setMomentum(HepLorentzVector(h3v, m_al_1[index * KFitConst::kNumber7 + 3][0]), KFitConst::kAfterFit);
516  // position
517  pdata.setPosition(HepPoint3D(
518  m_al_1[index * KFitConst::kNumber7 + 4][0],
519  m_al_1[index * KFitConst::kNumber7 + 5][0],
521  // error of the tracks
522  pdata.setError(this->makeError3(pdata.getMomentum(),
523  m_V_al_1.sub(
524  index * KFitConst::kNumber7 + 1,
525  (index + 1)*KFitConst::kNumber7,
526  index * KFitConst::kNumber7 + 1,
527  (index + 1)*KFitConst::kNumber7), m_IsFixMass[index]),
529  if (m_ErrorCode != KFitError::kNoError) break;
530  index++;
531  }
532 
534  {
535  //TODO: Not Implemented
537  // vertex
541  // error of the vertex
542  for (int i = 0; i < 3; i++) for (int j = i; j < 3; j++) {
544  }
545  // error between vertex and tracks
546  for (int i = 0; i < m_TrackCount; i++) {
547  HepMatrix hm(3, KFitConst::kNumber7, 0);
548  for (int j = 0; j < 3; j++) for (int k = 0; k < KFitConst::kNumber7; k++) {
550  }
551  if (m_IsFixMass[i])
552  m_AfterTrackVertexError.push_back(this->makeError4(m_Tracks[i].getMomentum(), hm));
553  else
554  m_AfterTrackVertexError.push_back(hm);
555  }
556  } else
557  {
558  // not fit
560  }
561 
563 }
564 
565 
566 enum KFitError::ECode
569  {
570 
571  HepMatrix al_1_prime(m_al_1);
572  HepMatrix Sum_al_1(4, 1, 0);
573  HepMatrix Sum_child_al_1(4 * m_ConstraintMassCount, 1, 0);
574  double energy[KFitConst::kMaxTrackCount2];
575  double a;
576 
577  for (int i = 0; i < m_TrackCount; i++) {
578  a = m_property[i][2];
579  if (!m_FlagAtDecayPoint) a = 0.;
580  al_1_prime[i * KFitConst::kNumber7 + 0][0] -= a * (m_BeforeVertex.y() - al_1_prime[i * KFitConst::kNumber7 + 5][0]);
581  al_1_prime[i * KFitConst::kNumber7 + 1][0] += a * (m_BeforeVertex.x() - al_1_prime[i * KFitConst::kNumber7 + 4][0]);
582  energy[i] = sqrt(al_1_prime[i * KFitConst::kNumber7 + 0][0] * al_1_prime[i * KFitConst::kNumber7 + 0][0] +
583  al_1_prime[i * KFitConst::kNumber7 + 1][0] * al_1_prime[i * KFitConst::kNumber7 + 1][0] +
584  al_1_prime[i * KFitConst::kNumber7 + 2][0] * al_1_prime[i * KFitConst::kNumber7 + 2][0] +
585  m_property[i][1] * m_property[i][1]);
586  }
587 
588  for (int i = 0; i < m_TrackCount; i++) {
589  if (m_IsFixMass[i])
590  Sum_al_1[3][0] += energy[i];
591  else
592  Sum_al_1[3][0] += al_1_prime[i * KFitConst::kNumber7 + 3][0];
593 
594  for (int j = 0; j < 3; j++) Sum_al_1[j][0] += al_1_prime[i * KFitConst::kNumber7 + j][0];
595  }
596 
597  for (int i = 0; i < m_ConstraintMassCount; i++) {
598  for (int k = m_ConstraintMassChildLists[i].first; k <= m_ConstraintMassChildLists[i].second; k++) {
599  if (m_IsFixMass[k])
600  Sum_child_al_1[i * 4 + 3][0] += energy[k];
601  else
602  Sum_child_al_1[i * 4 + 3][0] += al_1_prime[k * KFitConst::kNumber7 + 3][0];
603  for (int j = 0; j < 3; j++) Sum_child_al_1[i * 4 + j][0] += al_1_prime[k * KFitConst::kNumber7 + j][0];
604  }
605  }
606 
607  m_d[0][0] = Sum_al_1[0][0] - m_FourMomentum.Px();
608  m_d[1][0] = Sum_al_1[1][0] - m_FourMomentum.Py();
609  m_d[2][0] = Sum_al_1[2][0] - m_FourMomentum.Pz();
610  m_d[3][0] = Sum_al_1[3][0] - m_FourMomentum.E();
611 
612  for (int i = 0; i < m_ConstraintMassCount; i++) {
613  m_d[4 + i][0] =
614  + Sum_child_al_1[i * 4 + 3][0] * Sum_child_al_1[i * 4 + 3][0] - Sum_child_al_1[i * 4 + 0][0] * Sum_child_al_1[i * 4 + 0][0]
615  - Sum_child_al_1[i * 4 + 1][0] * Sum_child_al_1[i * 4 + 1][0] - Sum_child_al_1[i * 4 + 2][0] * Sum_child_al_1[i * 4 + 2][0]
617  }
618 
619  for (int i = 0; i < m_TrackCount; i++) {
620  if (energy[i] == 0) {
622  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
623  break;
624  }
625 
626  a = m_property[i][2];
627  if (!m_FlagAtDecayPoint) a = 0.;
628 
629  // four-momentum conservation constraint
630  for (int l = 0; l < 4; l++) {
631  for (int n = 0; n < 6; n++) {
632  m_D[l][i * KFitConst::kNumber7 + n] = 0;
633  }
634  }
635  if (m_IsFixMass[i]) {
636  double invE = 1. / energy[i];
637  m_D[0][i * KFitConst::kNumber7 + 0] = 1;
638  m_D[0][i * KFitConst::kNumber7 + 5] = -a;
639  m_D[1][i * KFitConst::kNumber7 + 1] = 1;
640  m_D[1][i * KFitConst::kNumber7 + 4] = a;
641  m_D[2][i * KFitConst::kNumber7 + 2] = 1;
642  m_D[3][i * KFitConst::kNumber7 + 0] = al_1_prime[i * KFitConst::kNumber7 + 0][0] * invE;
643  m_D[3][i * KFitConst::kNumber7 + 1] = al_1_prime[i * KFitConst::kNumber7 + 1][0] * invE;
644  m_D[3][i * KFitConst::kNumber7 + 2] = al_1_prime[i * KFitConst::kNumber7 + 2][0] * invE;
645  m_D[3][i * KFitConst::kNumber7 + 4] = -al_1_prime[i * KFitConst::kNumber7 + 1][0] * invE * a;
646  m_D[3][i * KFitConst::kNumber7 + 5] = al_1_prime[i * KFitConst::kNumber7 + 0][0] * invE * a;
647  } else {
648  m_D[0][i * KFitConst::kNumber7 + 0] = 1;
649  m_D[1][i * KFitConst::kNumber7 + 1] = 1;
650  m_D[2][i * KFitConst::kNumber7 + 2] = 1;
651  m_D[3][i * KFitConst::kNumber7 + 3] = 1;
652  }
653 
654  // invariant mass constraint
655  for (int l = 0; l < m_ConstraintMassCount; l++) {
656  if (i >= m_ConstraintMassChildLists[l].first && i <= m_ConstraintMassChildLists[l].second) {
657  double invE = 1. / energy[i];
658  m_D[4 + l][i * KFitConst::kNumber7 + 0] = 2.*(Sum_child_al_1[l * 4 + 3][0] * al_1_prime[i * KFitConst::kNumber7 + 0][0] * invE -
659  Sum_child_al_1[l * 4 + 0][0]);
660  m_D[4 + l][i * KFitConst::kNumber7 + 1] = 2.*(Sum_child_al_1[l * 4 + 3][0] * al_1_prime[i * KFitConst::kNumber7 + 1][0] * invE -
661  Sum_child_al_1[l * 4 + 1][0]);
662  m_D[4 + l][i * KFitConst::kNumber7 + 2] = 2.*(Sum_child_al_1[l * 4 + 3][0] * al_1_prime[i * KFitConst::kNumber7 + 2][0] * invE -
663  Sum_child_al_1[l * 4 + 2][0]);
664  m_D[4 + l][i * KFitConst::kNumber7 + 3] = 0.;
665  m_D[4 + l][i * KFitConst::kNumber7 + 4] = -2.*(Sum_child_al_1[l * 4 + 3][0] * al_1_prime[i * KFitConst::kNumber7 + 1][0] * invE -
666  Sum_child_al_1[l * 4 + 1][0]) * a;
667  m_D[4 + l][i * KFitConst::kNumber7 + 5] = 2.*(Sum_child_al_1[l * 4 + 3][0] * al_1_prime[i * KFitConst::kNumber7 + 0][0] * invE -
668  Sum_child_al_1[l * 4 + 0][0]) * a;
669  m_D[4 + l][i * KFitConst::kNumber7 + 6] = 0.;
670  } else {
671  for (int n = 0; n < 6; n++) {
672  m_D[4 + l][i * KFitConst::kNumber7 + n] = 0;
673  }
674  }
675  }
676  }
677 
678  } else
679  {
680  //TODO: Not Implemented
682 
683  // m_FlagFitIncludingVertex == true
684  HepMatrix al_1_prime(m_al_1);
685  HepMatrix Sum_al_1(7, 1, 0);
686  double energy[KFitConst::kMaxTrackCount2];
687 
688  for (int i = 0; i < m_TrackCount; i++) {
689  const double a = m_property[i][2];
690  al_1_prime[i * KFitConst::kNumber7 + 0][0] -= a * (al_1_prime[KFitConst::kNumber7 * m_TrackCount + 1][0] - al_1_prime[i *
691  KFitConst::kNumber7 + 5][0]);
692  al_1_prime[i * KFitConst::kNumber7 + 1][0] += a * (al_1_prime[KFitConst::kNumber7 * m_TrackCount + 0][0] - al_1_prime[i *
693  KFitConst::kNumber7 + 4][0]);
694  energy[i] = sqrt(al_1_prime[i * KFitConst::kNumber7 + 0][0] * al_1_prime[i * KFitConst::kNumber7 + 0][0] +
695  al_1_prime[i * KFitConst::kNumber7 + 1][0] * al_1_prime[i * KFitConst::kNumber7 + 1][0] +
696  al_1_prime[i * KFitConst::kNumber7 + 2][0] * al_1_prime[i * KFitConst::kNumber7 + 2][0] +
697  m_property[i][1] * m_property[i][1]);
698  Sum_al_1[6][0] = + a;
699  }
700 
701  for (int i = 0; i < m_TrackCount; i++) {
702  if (energy[i] == 0) {
704  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
705  break;
706  }
707 
708  if (m_IsFixMass[i]) {
709  double invE = 1. / energy[i];
710  Sum_al_1[3][0] += energy[i];
711  Sum_al_1[4][0] += al_1_prime[i * KFitConst::kNumber7 + 1][0] * m_property[i][2] * invE;
712  Sum_al_1[5][0] += al_1_prime[i * KFitConst::kNumber7 + 0][0] * m_property[i][2] * invE;
713  } else {
714  Sum_al_1[3][0] += al_1_prime[i * KFitConst::kNumber7 + 3][0];
715  }
716 
717  for (int j = 0; j < 3; j++) Sum_al_1[j][0] += al_1_prime[i * KFitConst::kNumber7 + j][0];
718  }
719 
720  m_d[0][0] = Sum_al_1[0][0] - m_FourMomentum.Px();
721  m_d[1][0] = Sum_al_1[1][0] - m_FourMomentum.Py();
722  m_d[2][0] = Sum_al_1[2][0] - m_FourMomentum.Pz();
723  m_d[3][0] = Sum_al_1[3][0] - m_FourMomentum.E();
724 
725  for (int i = 0; i < m_TrackCount; i++) {
726  if (energy[i] == 0) {
728  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
729  break;
730  }
731 
732  for (int l = 0; l < 4; l++) {
733  for (int n = 0; n < 6; n++) {
734  if (l == n) m_D[l][i * KFitConst::kNumber7 + n] = 1;
735  else m_D[l][i * KFitConst::kNumber7 + n] = 0;
736  }
737  }
738  }
739 
740  m_D[0][KFitConst::kNumber7 * m_TrackCount + 0] = 2.*(Sum_al_1[3][0] * Sum_al_1[4][0] - Sum_al_1[1][0] * Sum_al_1[6][0]);
741  m_D[0][KFitConst::kNumber7 * m_TrackCount + 1] = -2.*(Sum_al_1[3][0] * Sum_al_1[5][0] - Sum_al_1[0][0] * Sum_al_1[6][0]);
742  m_D[0][KFitConst::kNumber7 * m_TrackCount + 2] = 0.;
743  }
744 
746 }
747 
748 
749 enum KFitError::ECode
752 
754 }
755 
757 {
758  MakeMotherKFit kmm;
760  unsigned n = getTrackCount();
761  for (unsigned i = 0; i < n; ++i) {
763  getTrack(i).getCharge());
764  if (getFlagFitWithVertex())
766  for (unsigned j = i + 1; j < n; ++j) {
767  kmm.setCorrelation(getCorrelation(i, j));
768  }
769  }
770  kmm.setVertex(getVertex());
771  if (getFlagFitWithVertex())
773  m_ErrorCode = kmm.doMake();
775  return m_ErrorCode;
776  double chi2 = getCHIsq();
777  int ndf = getNDF();
778  double prob = TMath::Prob(chi2, ndf);
779  //
780  bool haschi2 = mother->hasExtraInfo("chiSquared");
781  if (haschi2) {
782  mother->setExtraInfo("chiSquared", chi2);
783  mother->setExtraInfo("ndf", ndf);
784  } else {
785  mother->addExtraInfo("chiSquared", chi2);
786  mother->addExtraInfo("ndf", ndf);
787  }
788 
789  mother->updateMomentum(
790  CLHEPToROOT::getLorentzVector(kmm.getMotherMomentum()),
791  CLHEPToROOT::getXYZVector(kmm.getMotherPosition()),
792  CLHEPToROOT::getTMatrixFSym(kmm.getMotherError()),
793  prob);
795  return m_ErrorCode;
796 }
Class to store reconstructed particles.
Definition: Particle.h:74
void setExtraInfo(const std::string &name, double value)
Sets the user-defined data of given name to the given value.
Definition: Particle.cc:1306
bool hasExtraInfo(const std::string &name) const
Return whether the extra info with the given name is set.
Definition: Particle.cc:1255
void addExtraInfo(const std::string &name, double value)
Sets the user-defined data of given name to the given value.
Definition: Particle.cc:1325
void updateMomentum(const ROOT::Math::PxPyPzEVector &p4, const ROOT::Math::XYZVector &vertex, const TMatrixFSym &errMatrix, double pValue)
Sets Lorentz vector, position, 7x7 error matrix and p-value.
Definition: Particle.h:351
int m_NecessaryTrackCount
Number needed tracks to perform fit.
Definition: KFitBase.h:303
double m_MagneticField
Magnetic field.
Definition: KFitBase.h:311
CLHEP::HepMatrix m_al_1
See J.Tanaka Ph.D (2001) p136 for definition.
Definition: KFitBase.h:259
virtual enum KFitError::ECode setCorrelation(const CLHEP::HepMatrix &c)
Set a correlation matrix.
Definition: KFitBase.cc:70
const CLHEP::HepSymMatrix getTrackError(const int id) const
Get an error matrix of the track.
Definition: KFitBase.cc:168
const CLHEP::HepLorentzVector getTrackMomentum(const int id) const
Get a Lorentz vector of the track.
Definition: KFitBase.cc:154
CLHEP::HepMatrix m_lam
See J.Tanaka Ph.D (2001) p137 for definition.
Definition: KFitBase.h:276
const HepPoint3D getTrackPosition(const int id) const
Get a position of the track.
Definition: KFitBase.cc:161
CLHEP::HepMatrix m_property
Container of charges and masses.
Definition: KFitBase.h:263
enum KFitError::ECode m_ErrorCode
Error code.
Definition: KFitBase.h:243
virtual enum KFitError::ECode setZeroCorrelation(void)
Indicate no correlation between tracks.
Definition: KFitBase.cc:85
const CLHEP::HepSymMatrix makeError3(const CLHEP::HepLorentzVector &p, const CLHEP::HepMatrix &e, const bool is_fix_mass) const
Rebuild an error matrix from a Lorentz vector and an error matrix.
CLHEP::HepMatrix m_V_al_1
See J.Tanaka Ph.D (2001) p138 for definition.
Definition: KFitBase.h:274
virtual int getNDF(void) const
Get an NDF of the fit.
Definition: KFitBase.cc:114
CLHEP::HepMatrix m_d
See J.Tanaka Ph.D (2001) p137 for definition.
Definition: KFitBase.h:268
bool isFitted(void) const
Return false if fit is not performed yet or performed fit is failed; otherwise true.
Definition: KFitBase.cc:727
CLHEP::HepMatrix m_D
See J.Tanaka Ph.D (2001) p137 for definition.
Definition: KFitBase.h:266
CLHEP::HepMatrix m_V_D
See J.Tanaka Ph.D (2001) p138 for definition.
Definition: KFitBase.h:271
bool isTrackIDInRange(const int id) const
Check if the id is in the range.
Definition: KFitBase.cc:738
virtual const CLHEP::HepMatrix getCorrelation(const int id1, const int id2, const int flag=KFitConst::kAfterFit) const
Get a correlation matrix between two tracks.
Definition: KFitBase.cc:183
bool m_FlagCorrelation
Flag whether a correlation among tracks exists.
Definition: KFitBase.h:306
CLHEP::HepSymMatrix m_V_al_0
See J.Tanaka Ph.D (2001) p137 for definition.
Definition: KFitBase.h:255
const KFitTrack getTrack(const int id) const
Get a specified track object.
Definition: KFitBase.cc:175
std::vector< CLHEP::HepMatrix > m_BeforeCorrelation
Container of input correlation matrices.
Definition: KFitBase.h:251
bool m_FlagFitted
Flag to indicate if the fit is performed and succeeded.
Definition: KFitBase.h:245
double m_CHIsq
chi-square of the fit.
Definition: KFitBase.h:297
int getTrackCount(void) const
Get the number of added tracks.
Definition: KFitBase.cc:107
int m_NDF
NDF of the fit.
Definition: KFitBase.h:295
std::vector< KFitTrack > m_Tracks
Container of input tracks.
Definition: KFitBase.h:249
const CLHEP::HepMatrix makeError4(const CLHEP::HepLorentzVector &p, const CLHEP::HepMatrix &e) const
Rebuild an error matrix from a Lorentz vector and an error matrix.
Definition: KFitBase.cc:438
int m_TrackCount
Number of tracks.
Definition: KFitBase.h:301
CLHEP::HepMatrix m_al_0
See J.Tanaka Ph.D (2001) p136 for definition.
Definition: KFitBase.h:257
enum KFitError::ECode doFit1(void)
Perform a fit (used in MassFitKFit::doFit()).
Definition: KFitBase.cc:501
static void displayError(const char *file, const int line, const char *func, const enum ECode code)
Display a description of error and its location.
Definition: KFitError.h:72
ECode
ECode is a error code enumerate.
Definition: KFitError.h:34
@ kUnimplemented
Unprepared.
Definition: KFitError.h:44
@ kCannotGetMatrixInverse
Cannot calculate matrix inverse (bad track property or internal error)
Definition: KFitError.h:58
@ kOutOfRange
Specified track-id out of range.
Definition: KFitError.h:42
@ kDivisionByZero
Division by zero (bad track property or internal error)
Definition: KFitError.h:56
@ kBadTrackSize
Track count too small to perform fit.
Definition: KFitError.h:47
@ kBadMatrixSize
Wrong correlation matrix size.
Definition: KFitError.h:49
@ kBadCorrelationSize
Wrong correlation matrix size (internal error)
Definition: KFitError.h:51
MakeMotherKFit is a class to build mother particle from kinematically fitted daughters.
enum KFitError::ECode setVertex(const HepPoint3D &v)
Set a vertex position of the mother particle.
enum KFitError::ECode addTrack(const KFitTrack &kp)
Add a track to the make-mother object.
enum KFitError::ECode doMake(void)
Perform a reconstruction of mother particle.
const CLHEP::HepSymMatrix getMotherError(void) const
Get an error matrix of the mother particle.
enum KFitError::ECode setCorrelation(const CLHEP::HepMatrix &e)
Set a correlation matrix.
const HepPoint3D getMotherPosition(void) const
Get a position of the mother particle.
enum KFitError::ECode setVertexError(const CLHEP::HepSymMatrix &e)
Set a vertex error matrix of the mother particle.
enum KFitError::ECode setTrackVertexError(const CLHEP::HepMatrix &e)
Set a vertex error matrix of the child particle in the addTrack'ed order.
const CLHEP::HepLorentzVector getMotherMomentum(void) const
Get a Lorentz vector of the mother particle.
enum KFitError::ECode setMagneticField(const double mf)
Change a magnetic field from the default value KFitConst::kDefaultMagneticField.
enum KFitError::ECode setZeroCorrelation(void) override
Indicate no correlation between tracks.
enum KFitError::ECode setVertex(const HepPoint3D &v)
Set an initial vertex position for the mass-four-momentum-constraint fit.
bool m_FlagTrackVertexError
Flag to indicate if the vertex error matrix of the child particle is preset.
enum KFitError::ECode setFourMomentum(const ROOT::Math::PxPyPzEVector &m)
Set an 4 Momentum for the mass-four-constraint fit.
bool m_FlagAtDecayPoint
Flag controlled by setFlagAtDecayPoint().
bool getFlagAtDecayPoint(void) const
Get a flag if to constraint at the decay point in the mass-four-momentum-constraint fit.
enum KFitError::ECode addMassConstraint(const double m, std::vector< unsigned > &childTrackId)
Set an invariant mass of daughter particle for the mass-four-momentum-constraint fit.
enum KFitError::ECode prepareInputMatrix(void) override
Build grand matrices for minimum search from input-track properties.
bool getFlagFitWithVertex(void) const
Get a flag if the fit is allowed with moving the vertex position.
enum KFitError::ECode calculateNDF(void) override
Calculate an NDF of the fit.
double getCHIsq(void) const override
Get a chi-square of the fit.
enum KFitError::ECode setFlagAtDecayPoint(const bool flag)
Set a flag if to constraint at the decay point in the mass-four-momentum-constraint fit.
enum KFitError::ECode setTrackZeroVertexError(void)
Indicate no vertex uncertainty in the child particle in the addTrack'ed order.
std::vector< int > m_IsFixMass
Array of flags whether the track property is fixed at the mass.
enum KFitError::ECode updateMother(Particle *mother)
Update mother particle.
enum KFitError::ECode unfixMass(void)
Tell the object to unfix the last added track property at the invariant mass.
enum KFitError::ECode prepareInputSubMatrix(void) override
Build sub-matrices for minimum search from input-track properties.
bool m_FlagFitIncludingVertex
Flag to indicate if the fit is allowed with moving the vertex position.
~MassFourCFitKFit(void)
Destruct the object.
enum KFitError::ECode doFit(void)
Perform a mass-four-momentum-constraint fit.
double m_InvariantMass
Invariant mass.
enum KFitError::ECode setInvariantMass(const double m)
Set an invariant mass for the mass-four-momentum-constraint fit.
std::vector< CLHEP::HepMatrix > m_BeforeTrackVertexError
array of vertex error matrices before the fit.
enum KFitError::ECode makeCoreMatrix(void) override
Build matrices using the kinematical constraint.
enum KFitError::ECode setVertexError(const CLHEP::HepSymMatrix &e)
Set an initial vertex error matrix for the mass-four-momentum-constraint fit.
enum KFitError::ECode setTrackVertexError(const CLHEP::HepMatrix &e)
Set a vertex error matrix of the child particle in the addTrack'ed order.
enum KFitError::ECode prepareCorrelation(void) override
Build a grand correlation matrix from input-track properties.
HepPoint3D m_AfterVertex
Vertex position after the fit.
enum KFitError::ECode setCorrelation(const CLHEP::HepMatrix &m) override
Set a correlation matrix.
std::vector< std::pair< int, int > > m_ConstraintMassChildLists
Daughter track id of constrainted particle.
std::vector< CLHEP::HepMatrix > m_AfterTrackVertexError
array of vertex error matrices after the fit.
const CLHEP::HepMatrix getCorrelation(const int id1, const int id2, const int flag=KFitConst::kAfterFit) const override
Get a correlation matrix between two tracks.
double getTrackCHIsq(const int id) const override
Get a chi-square of the track.
const CLHEP::HepSymMatrix getVertexError(const int flag=KFitConst::kAfterFit) const
Get a vertex error matrix.
const HepPoint3D getVertex(const int flag=KFitConst::kAfterFit) const
Get a vertex position.
std::vector< double > m_ConstraintMassList
constrainted mass
const CLHEP::HepMatrix getTrackVertexError(const int id, const int flag=KFitConst::kAfterFit) const
Get a vertex error matrix of the track.
CLHEP::HepSymMatrix m_AfterVertexError
Vertex error matrix after the fit.
enum KFitError::ECode prepareOutputMatrix(void) override
Build an output error matrix.
CLHEP::HepSymMatrix m_BeforeVertexError
Vertex error matrix before the fit.
HepPoint3D m_BeforeVertex
Vertex position before the fit.
enum KFitError::ECode fixMass(void)
Tell the object to fix the last added track property at the invariant mass.
ROOT::Math::PxPyPzEVector m_FourMomentum
Four Momentum.
double getInvariantMass(void) const
Get an invariant mass.
Abstract base class for different kinds of events.
Definition: ClusterUtils.h:23
static constexpr double kLightSpeed
Speed of light.
Definition: KFitConst.h:53
static const int kMaxTrackCount
Maximum track size.
Definition: KFitConst.h:39
static const int kMaxTrackCount2
Maximum track size (internal use)
Definition: KFitConst.h:41
static const int kAfterFit
Input parameter to specify after-fit when setting/getting a track attribute.
Definition: KFitConst.h:36
static const int kBeforeFit
Input parameter to specify before-fit when setting/getting a track attribute.
Definition: KFitConst.h:34
static const int kNumber7
Constant 7 to check matrix size (internal use)
Definition: KFitConst.h:31