Belle II Software  light-2212-foldex
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 #include <cstdio>
10 
11 #include <TMatrixFSym.h>
12 
13 #include <analysis/VertexFitting/KFit/MassFourCFitKFit.h>
14 #include <analysis/VertexFitting/KFit/MakeMotherKFit.h>
15 #include <analysis/utility/CLHEPToROOT.h>
16 #include <framework/gearbox/Const.h>
17 
18 using namespace std;
19 using namespace Belle2;
20 using namespace Belle2::analysis;
21 using namespace CLHEP;
22 using namespace ROOT::Math;
23 
24 MassFourCFitKFit::MassFourCFitKFit() : m_AfterVertexError(HepSymMatrix(3, 0)),
25  m_FourMomentum(PxPyPzEVector())
26 {
27  m_FlagFitted = false;
28  m_FlagTrackVertexError = false;
30  m_FlagAtDecayPoint = true;
32  m_InvariantMass = -1.0;
36 }
37 
38 
40 
42 MassFourCFitKFit::addMassConstraint(const double m, std::vector<unsigned>& childTrackId) {
43  if ((childTrackId.back() - childTrackId.front()) != childTrackId.size() - 1)
44  {
46  }
47  m_ConstraintMassList.push_back(m);
48  m_ConstraintMassChildLists.push_back(std::make_pair(childTrackId.front(), childTrackId.back()));
50 }
51 
53 MassFourCFitKFit::setVertex(const HepPoint3D& v) {
54  m_BeforeVertex = v;
55 
57 }
58 
59 
61 MassFourCFitKFit::setVertexError(const HepSymMatrix& e) {
62  if (e.num_row() != 3)
63  {
65  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
66  return m_ErrorCode;
67  }
68 
71 
73 }
74 
75 
77 MassFourCFitKFit::setInvariantMass(const double m) {
78  m_InvariantMass = m;
79 
81 }
82 
83 
85 MassFourCFitKFit::setFourMomentum(const PxPyPzEVector& m) {
86  m_FourMomentum = m;
87 
89 }
90 
91 
94  m_FlagAtDecayPoint = flag;
95 
97 }
98 
99 
100 enum KFitError::ECode
102  m_IsFixMass.push_back(true);
103 
105 }
106 
107 
108 enum KFitError::ECode
110  m_IsFixMass.push_back(false);
111 
113 }
114 
115 
116 enum KFitError::ECode
117 MassFourCFitKFit::setTrackVertexError(const HepMatrix& e) {
118  if (e.num_row() != 3 || e.num_col() != KFitConst::kNumber7)
119  {
121  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
122  return m_ErrorCode;
123  }
124 
125  m_BeforeTrackVertexError.push_back(e);
126  m_FlagTrackVertexError = true;
128 
130 }
131 
132 
133 enum KFitError::ECode
135  HepMatrix zero(3, KFitConst::kNumber7, 0);
136 
137  return this->setTrackVertexError(zero);
138 }
139 
140 
141 enum KFitError::ECode
142 MassFourCFitKFit::setCorrelation(const HepMatrix& m) {
143  return KFitBase::setCorrelation(m);
144 }
145 
146 
147 enum KFitError::ECode
150 }
151 
152 
153 const HepPoint3D
154 MassFourCFitKFit::getVertex(const int flag) const
155 {
156  if (flag == KFitConst::kAfterFit && !isFitted()) return HepPoint3D();
157 
158  switch (flag) {
160  return m_BeforeVertex;
161 
163  return m_AfterVertex;
164 
165  default:
166  KFitError::displayError(__FILE__, __LINE__, __func__, KFitError::kOutOfRange);
167  return HepPoint3D();
168  }
169 }
170 
171 
172 const HepSymMatrix
173 MassFourCFitKFit::getVertexError(const int flag) const
174 {
175  if (flag == KFitConst::kAfterFit && !isFitted()) return HepSymMatrix(3, 0);
176 
177  if (flag == KFitConst::kBeforeFit)
178  return m_BeforeVertexError;
180  return m_AfterVertexError;
181  else {
182  KFitError::displayError(__FILE__, __LINE__, __func__, KFitError::kOutOfRange);
183  return HepSymMatrix(3, 0);
184  }
185 }
186 
187 
188 double
190 {
191  return m_InvariantMass;
192 }
193 
194 
195 bool
197 {
198  return m_FlagAtDecayPoint;
199 }
200 
201 
202 bool
204 {
206 }
207 
208 
209 double
211 {
212  return m_CHIsq;
213 }
214 
215 
216 const HepMatrix
217 MassFourCFitKFit::getTrackVertexError(const int id, const int flag) const
218 {
219  if (flag == KFitConst::kAfterFit && !isFitted()) return HepMatrix(3, KFitConst::kNumber7, 0);
220  if (!isTrackIDInRange(id)) return HepMatrix(3, KFitConst::kNumber7, 0);
221 
222  if (flag == KFitConst::kBeforeFit)
223  return m_BeforeTrackVertexError[id];
225  return m_AfterTrackVertexError[id];
226  else {
227  KFitError::displayError(__FILE__, __LINE__, __func__, KFitError::kOutOfRange);
228  return HepMatrix(3, KFitConst::kNumber7, 0);
229  }
230 }
231 
232 
233 double
235 {
236  if (!isFitted()) return -1;
237  if (!isTrackIDInRange(id)) return -1;
238 
239  if (m_IsFixMass[id]) {
240 
241  HepMatrix da(m_Tracks[id].getFitParameter(KFitConst::kBeforeFit) - m_Tracks[id].getFitParameter(KFitConst::kAfterFit));
242  int err_inverse = 0;
243  const double chisq = (da.T() * (m_Tracks[id].getFitError(KFitConst::kBeforeFit).inverse(err_inverse)) * da)[0][0];
244 
245  if (err_inverse) {
246  KFitError::displayError(__FILE__, __LINE__, __func__, KFitError::kCannotGetMatrixInverse);
247  return -1;
248  }
249 
250  return chisq;
251 
252  } else {
253 
254  HepMatrix da(m_Tracks[id].getMomPos(KFitConst::kBeforeFit) - m_Tracks[id].getMomPos(KFitConst::kAfterFit));
255  int err_inverse = 0;
256  const double chisq = (da.T() * (m_Tracks[id].getError(KFitConst::kBeforeFit).inverse(err_inverse)) * da)[0][0];
257 
258  if (err_inverse) {
259  KFitError::displayError(__FILE__, __LINE__, __func__, KFitError::kCannotGetMatrixInverse);
260  return -1;
261  }
262 
263  return chisq;
264 
265  }
266 }
267 
268 
269 const HepMatrix
270 MassFourCFitKFit::getCorrelation(const int id1, const int id2, const int flag) const
271 {
272  if (flag == KFitConst::kAfterFit && !isFitted()) return HepMatrix(KFitConst::kNumber7, KFitConst::kNumber7, 0);
273  if (!isTrackIDInRange(id1)) return HepMatrix(KFitConst::kNumber7, KFitConst::kNumber7, 0);
274  if (!isTrackIDInRange(id2)) return HepMatrix(KFitConst::kNumber7, KFitConst::kNumber7, 0);
275 
276  switch (flag) {
278  return KFitBase::getCorrelation(id1, id2, flag);
279 
281  return makeError3(
282  this->getTrackMomentum(id1),
283  this->getTrackMomentum(id2),
284  m_V_al_1.sub(KFitConst::kNumber7 * id1 + 1, KFitConst::kNumber7 * (id1 + 1), KFitConst::kNumber7 * id2 + 1,
285  KFitConst::kNumber7 * (id2 + 1)),
286  m_IsFixMass[id1],
287  m_IsFixMass[id2]);
288 
289  default:
290  KFitError::displayError(__FILE__, __LINE__, __func__, KFitError::kOutOfRange);
291  return HepMatrix(KFitConst::kNumber7, KFitConst::kNumber7, 0);
292  }
293 }
294 
295 
296 enum KFitError::ECode
298  return KFitBase::doFit1();
299 }
300 
301 
302 enum KFitError::ECode
305  m_d = HepMatrix(4 + m_ConstraintMassCount, 1, 0);
306  m_V_D = HepMatrix(4 + m_ConstraintMassCount, 4 + m_ConstraintMassCount, 0);
307  m_lam = HepMatrix(4 + m_ConstraintMassCount, 1, 0);
308 
310  {
312  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
313  return m_ErrorCode;
314  }
315 
316  if (m_IsFixMass.size() == 0)
317  {
318  // If no fix_mass flag at all,
319  // all tracks are considered to be fixed at mass.
320  for (int i = 0; i < m_TrackCount; i++) this->fixMass();
321  } else if (m_IsFixMass.size() != (unsigned int)m_TrackCount)
322  {
324  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
325  return m_ErrorCode;
326  }
327 
329  {
330  int index = 0;
331  m_al_0 = HepMatrix(KFitConst::kNumber7 * m_TrackCount, 1, 0);
332  m_property = HepMatrix(m_TrackCount, 3, 0);
333  m_V_al_0 = HepSymMatrix(KFitConst::kNumber7 * m_TrackCount, 0);
334 
335  for (auto& track : m_Tracks) {
336  // momentum x,y,z and position x,y,z
337  m_al_0[index * KFitConst::kNumber7 + 0][0] = track.getMomentum(KFitConst::kBeforeFit).x();
338  m_al_0[index * KFitConst::kNumber7 + 1][0] = track.getMomentum(KFitConst::kBeforeFit).y();
339  m_al_0[index * KFitConst::kNumber7 + 2][0] = track.getMomentum(KFitConst::kBeforeFit).z();
340  m_al_0[index * KFitConst::kNumber7 + 3][0] = track.getMomentum(KFitConst::kBeforeFit).t();
341  m_al_0[index * KFitConst::kNumber7 + 4][0] = track.getPosition(KFitConst::kBeforeFit).x();
342  m_al_0[index * KFitConst::kNumber7 + 5][0] = track.getPosition(KFitConst::kBeforeFit).y();
343  m_al_0[index * KFitConst::kNumber7 + 6][0] = track.getPosition(KFitConst::kBeforeFit).z();
344  // these error
345  m_V_al_0.sub(index * KFitConst::kNumber7 + 1, track.getError(KFitConst::kBeforeFit));
346  // charge, mass, a
347  m_property[index][0] = track.getCharge();
348  m_property[index][1] = track.getMass();
349  const double c = Belle2::Const::speedOfLight * 1e-4;
350  m_property[index][2] = -c * m_MagneticField * track.getCharge();
351  index++;
352  }
353 
354  // error between track and track
355  if (m_FlagCorrelation) {
356  this->prepareCorrelation();
358  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
359  return m_ErrorCode;
360  }
361  }
362 
363  // set member matrix
364  m_al_1 = m_al_0;
365 
366  // define size of matrix
369 
370  } else
371  {
372  //TODO: Not Implemented
374  // m_FlagFitIncludingVertex == true
375  int index = 0;
376  m_al_0 = HepMatrix(KFitConst::kNumber7 * m_TrackCount + 3, 1, 0);
377  m_property = HepMatrix(m_TrackCount, 3, 0);
378  m_V_al_0 = HepSymMatrix(KFitConst::kNumber7 * m_TrackCount + 3, 0);
379 
380  for (auto& track : m_Tracks) {
381  // momentum x,y,z and position x,y,z
382  m_al_0[index * KFitConst::kNumber7 + 0][0] = track.getMomentum(KFitConst::kBeforeFit).x();
383  m_al_0[index * KFitConst::kNumber7 + 1][0] = track.getMomentum(KFitConst::kBeforeFit).y();
384  m_al_0[index * KFitConst::kNumber7 + 2][0] = track.getMomentum(KFitConst::kBeforeFit).z();
385  m_al_0[index * KFitConst::kNumber7 + 3][0] = track.getMomentum(KFitConst::kBeforeFit).t();
386  m_al_0[index * KFitConst::kNumber7 + 4][0] = track.getPosition(KFitConst::kBeforeFit).x();
387  m_al_0[index * KFitConst::kNumber7 + 5][0] = track.getPosition(KFitConst::kBeforeFit).y();
388  m_al_0[index * KFitConst::kNumber7 + 6][0] = track.getPosition(KFitConst::kBeforeFit).z();
389  // these error
390  m_V_al_0.sub(index * KFitConst::kNumber7 + 1, track.getError(KFitConst::kBeforeFit));
391  // charge, mass, a
392  m_property[index][0] = track.getCharge();
393  m_property[index][1] = track.getMass();
394  const double c = Belle2::Const::speedOfLight * 1e-4;
395  m_property[index][2] = -c * m_MagneticField * track.getCharge();
396  index++;
397  }
398 
399  // vertex
404 
405  // error between track and track
406  if (m_FlagCorrelation) {
407  this->prepareCorrelation();
409  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
410  return m_ErrorCode;
411  }
412  }
413 
414  // set member matrix
415  m_al_1 = m_al_0;
416 
417  // define size of matrix
419  m_D = m_V_al_1.sub(1, 4, 1, KFitConst::kNumber7 * m_TrackCount + 3);
420  }
421 
423 }
424 
425 
426 enum KFitError::ECode
428  char buf[1024];
429  sprintf(buf, "%s:%s(): internal error; this function should never be called", __FILE__, __func__);
430  B2FATAL(buf);
431 
432  /* NEVER REACHEd HERE */
433  return KFitError::kOutOfRange;
434 }
435 
436 
437 enum KFitError::ECode
439  if (m_BeforeCorrelation.size() != static_cast<unsigned int>(m_TrackCount * (m_TrackCount - 1) / 2))
440  {
442  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
443  return m_ErrorCode;
444  }
445 
446  int row = 0, col = 0;
447 
448  for (auto& hm : m_BeforeCorrelation)
449  {
450  // counter
451  row++;
452  if (row == m_TrackCount) {
453  col++;
454  row = col + 1;
455  }
456 
457  int ii = 0, jj = 0;
458  for (int i = KFitConst::kNumber7 * row; i < KFitConst::kNumber7 * (row + 1); i++) {
459  for (int j = KFitConst::kNumber7 * col; j < KFitConst::kNumber7 * (col + 1); j++) {
460  m_V_al_0[i][j] = hm[ii][jj];
461  jj++;
462  }
463  jj = 0;
464  ii++;
465  }
466  }
467 
469  {
470  //TODO: Not Implemented
472 
473  // ...error of vertex
475 
476  // ...error matrix between vertex and tracks
478  if (m_BeforeTrackVertexError.size() != (unsigned int)m_TrackCount) {
480  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
481  return m_ErrorCode;
482  }
483 
484  int i = 0;
485  for (auto& hm : m_BeforeTrackVertexError) {
486  for (int j = 0; j < 3; j++) for (int k = 0; k < KFitConst::kNumber7; k++) {
488  }
489  i++;
490  }
491  }
492  }
493 
495 }
496 
497 
498 enum KFitError::ECode
500  Hep3Vector h3v;
501  int index = 0;
502  for (auto& pdata : m_Tracks)
503  {
504  // tracks
505  // momentum
506  h3v.setX(m_al_1[index * KFitConst::kNumber7 + 0][0]);
507  h3v.setY(m_al_1[index * KFitConst::kNumber7 + 1][0]);
508  h3v.setZ(m_al_1[index * KFitConst::kNumber7 + 2][0]);
509  if (m_IsFixMass[index])
510  pdata.setMomentum(HepLorentzVector(h3v, sqrt(h3v.mag2() + pdata.getMass()*pdata.getMass())), KFitConst::kAfterFit);
511  else
512  pdata.setMomentum(HepLorentzVector(h3v, m_al_1[index * KFitConst::kNumber7 + 3][0]), KFitConst::kAfterFit);
513  // position
514  pdata.setPosition(HepPoint3D(
515  m_al_1[index * KFitConst::kNumber7 + 4][0],
516  m_al_1[index * KFitConst::kNumber7 + 5][0],
518  // error of the tracks
519  pdata.setError(this->makeError3(pdata.getMomentum(),
520  m_V_al_1.sub(
521  index * KFitConst::kNumber7 + 1,
522  (index + 1)*KFitConst::kNumber7,
523  index * KFitConst::kNumber7 + 1,
524  (index + 1)*KFitConst::kNumber7), m_IsFixMass[index]),
526  if (m_ErrorCode != KFitError::kNoError) break;
527  index++;
528  }
529 
531  {
532  //TODO: Not Implemented
534  // vertex
538  // error of the vertex
539  for (int i = 0; i < 3; i++) for (int j = i; j < 3; j++) {
541  }
542  // error between vertex and tracks
543  for (int i = 0; i < m_TrackCount; i++) {
544  HepMatrix hm(3, KFitConst::kNumber7, 0);
545  for (int j = 0; j < 3; j++) for (int k = 0; k < KFitConst::kNumber7; k++) {
547  }
548  if (m_IsFixMass[i])
549  m_AfterTrackVertexError.push_back(this->makeError4(m_Tracks[i].getMomentum(), hm));
550  else
551  m_AfterTrackVertexError.push_back(hm);
552  }
553  } else
554  {
555  // not fit
557  }
558 
560 }
561 
562 
563 enum KFitError::ECode
566  {
567 
568  HepMatrix al_1_prime(m_al_1);
569  HepMatrix Sum_al_1(4, 1, 0);
570  HepMatrix Sum_child_al_1(4 * m_ConstraintMassCount, 1, 0);
571  double energy[KFitConst::kMaxTrackCount2];
572  double a;
573 
574  for (int i = 0; i < m_TrackCount; i++) {
575  a = m_property[i][2];
576  if (!m_FlagAtDecayPoint) a = 0.;
577  al_1_prime[i * KFitConst::kNumber7 + 0][0] -= a * (m_BeforeVertex.y() - al_1_prime[i * KFitConst::kNumber7 + 5][0]);
578  al_1_prime[i * KFitConst::kNumber7 + 1][0] += a * (m_BeforeVertex.x() - al_1_prime[i * KFitConst::kNumber7 + 4][0]);
579  energy[i] = sqrt(al_1_prime[i * KFitConst::kNumber7 + 0][0] * al_1_prime[i * KFitConst::kNumber7 + 0][0] +
580  al_1_prime[i * KFitConst::kNumber7 + 1][0] * al_1_prime[i * KFitConst::kNumber7 + 1][0] +
581  al_1_prime[i * KFitConst::kNumber7 + 2][0] * al_1_prime[i * KFitConst::kNumber7 + 2][0] +
582  m_property[i][1] * m_property[i][1]);
583  }
584 
585  for (int i = 0; i < m_TrackCount; i++) {
586  if (m_IsFixMass[i])
587  Sum_al_1[3][0] += energy[i];
588  else
589  Sum_al_1[3][0] += al_1_prime[i * KFitConst::kNumber7 + 3][0];
590 
591  for (int j = 0; j < 3; j++) Sum_al_1[j][0] += al_1_prime[i * KFitConst::kNumber7 + j][0];
592  }
593 
594  for (int i = 0; i < m_ConstraintMassCount; i++) {
595  for (int k = m_ConstraintMassChildLists[i].first; k <= m_ConstraintMassChildLists[i].second; k++) {
596  if (m_IsFixMass[k])
597  Sum_child_al_1[i * 4 + 3][0] += energy[k];
598  else
599  Sum_child_al_1[i * 4 + 3][0] += al_1_prime[k * KFitConst::kNumber7 + 3][0];
600  for (int j = 0; j < 3; j++) Sum_child_al_1[i * 4 + j][0] += al_1_prime[k * KFitConst::kNumber7 + j][0];
601  }
602  }
603 
604  m_d[0][0] = Sum_al_1[0][0] - m_FourMomentum.Px();
605  m_d[1][0] = Sum_al_1[1][0] - m_FourMomentum.Py();
606  m_d[2][0] = Sum_al_1[2][0] - m_FourMomentum.Pz();
607  m_d[3][0] = Sum_al_1[3][0] - m_FourMomentum.E();
608 
609  for (int i = 0; i < m_ConstraintMassCount; i++) {
610  m_d[4 + i][0] =
611  + 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]
612  - 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]
614  }
615 
616  for (int i = 0; i < m_TrackCount; i++) {
617  if (energy[i] == 0) {
619  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
620  break;
621  }
622 
623  a = m_property[i][2];
624  if (!m_FlagAtDecayPoint) a = 0.;
625 
626  // four-momentum conservation constraint
627  for (int l = 0; l < 4; l++) {
628  for (int n = 0; n < 6; n++) {
629  m_D[l][i * KFitConst::kNumber7 + n] = 0;
630  }
631  }
632  if (m_IsFixMass[i]) {
633  double invE = 1. / energy[i];
634  m_D[0][i * KFitConst::kNumber7 + 0] = 1;
635  m_D[0][i * KFitConst::kNumber7 + 5] = -a;
636  m_D[1][i * KFitConst::kNumber7 + 1] = 1;
637  m_D[1][i * KFitConst::kNumber7 + 4] = a;
638  m_D[2][i * KFitConst::kNumber7 + 2] = 1;
639  m_D[3][i * KFitConst::kNumber7 + 0] = al_1_prime[i * KFitConst::kNumber7 + 0][0] * invE;
640  m_D[3][i * KFitConst::kNumber7 + 1] = al_1_prime[i * KFitConst::kNumber7 + 1][0] * invE;
641  m_D[3][i * KFitConst::kNumber7 + 2] = al_1_prime[i * KFitConst::kNumber7 + 2][0] * invE;
642  m_D[3][i * KFitConst::kNumber7 + 4] = -al_1_prime[i * KFitConst::kNumber7 + 1][0] * invE * a;
643  m_D[3][i * KFitConst::kNumber7 + 5] = al_1_prime[i * KFitConst::kNumber7 + 0][0] * invE * a;
644  } else {
645  m_D[0][i * KFitConst::kNumber7 + 0] = 1;
646  m_D[1][i * KFitConst::kNumber7 + 1] = 1;
647  m_D[2][i * KFitConst::kNumber7 + 2] = 1;
648  m_D[3][i * KFitConst::kNumber7 + 3] = 1;
649  }
650 
651  // invariant mass constraint
652  for (int l = 0; l < m_ConstraintMassCount; l++) {
653  if (i >= m_ConstraintMassChildLists[l].first && i <= m_ConstraintMassChildLists[l].second) {
654  double invE = 1. / energy[i];
655  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 -
656  Sum_child_al_1[l * 4 + 0][0]);
657  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 -
658  Sum_child_al_1[l * 4 + 1][0]);
659  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 -
660  Sum_child_al_1[l * 4 + 2][0]);
661  m_D[4 + l][i * KFitConst::kNumber7 + 3] = 0.;
662  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 -
663  Sum_child_al_1[l * 4 + 1][0]) * a;
664  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 -
665  Sum_child_al_1[l * 4 + 0][0]) * a;
666  m_D[4 + l][i * KFitConst::kNumber7 + 6] = 0.;
667  } else {
668  for (int n = 0; n < 6; n++) {
669  m_D[4 + l][i * KFitConst::kNumber7 + n] = 0;
670  }
671  }
672  }
673  }
674 
675  } else
676  {
677  //TODO: Not Implemented
679 
680  // m_FlagFitIncludingVertex == true
681  HepMatrix al_1_prime(m_al_1);
682  HepMatrix Sum_al_1(7, 1, 0);
683  double energy[KFitConst::kMaxTrackCount2];
684 
685  for (int i = 0; i < m_TrackCount; i++) {
686  const double a = m_property[i][2];
687  al_1_prime[i * KFitConst::kNumber7 + 0][0] -= a * (al_1_prime[KFitConst::kNumber7 * m_TrackCount + 1][0] - al_1_prime[i *
688  KFitConst::kNumber7 + 5][0]);
689  al_1_prime[i * KFitConst::kNumber7 + 1][0] += a * (al_1_prime[KFitConst::kNumber7 * m_TrackCount + 0][0] - al_1_prime[i *
690  KFitConst::kNumber7 + 4][0]);
691  energy[i] = sqrt(al_1_prime[i * KFitConst::kNumber7 + 0][0] * al_1_prime[i * KFitConst::kNumber7 + 0][0] +
692  al_1_prime[i * KFitConst::kNumber7 + 1][0] * al_1_prime[i * KFitConst::kNumber7 + 1][0] +
693  al_1_prime[i * KFitConst::kNumber7 + 2][0] * al_1_prime[i * KFitConst::kNumber7 + 2][0] +
694  m_property[i][1] * m_property[i][1]);
695  Sum_al_1[6][0] = + a;
696  }
697 
698  for (int i = 0; i < m_TrackCount; i++) {
699  if (energy[i] == 0) {
701  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
702  break;
703  }
704 
705  if (m_IsFixMass[i]) {
706  double invE = 1. / energy[i];
707  Sum_al_1[3][0] += energy[i];
708  Sum_al_1[4][0] += al_1_prime[i * KFitConst::kNumber7 + 1][0] * m_property[i][2] * invE;
709  Sum_al_1[5][0] += al_1_prime[i * KFitConst::kNumber7 + 0][0] * m_property[i][2] * invE;
710  } else {
711  Sum_al_1[3][0] += al_1_prime[i * KFitConst::kNumber7 + 3][0];
712  }
713 
714  for (int j = 0; j < 3; j++) Sum_al_1[j][0] += al_1_prime[i * KFitConst::kNumber7 + j][0];
715  }
716 
717  m_d[0][0] = Sum_al_1[0][0] - m_FourMomentum.Px();
718  m_d[1][0] = Sum_al_1[1][0] - m_FourMomentum.Py();
719  m_d[2][0] = Sum_al_1[2][0] - m_FourMomentum.Pz();
720  m_d[3][0] = Sum_al_1[3][0] - m_FourMomentum.E();
721 
722  for (int i = 0; i < m_TrackCount; i++) {
723  if (energy[i] == 0) {
725  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
726  break;
727  }
728 
729  for (int l = 0; l < 4; l++) {
730  for (int n = 0; n < 6; n++) {
731  if (l == n) m_D[l][i * KFitConst::kNumber7 + n] = 1;
732  else m_D[l][i * KFitConst::kNumber7 + n] = 0;
733  }
734  }
735  }
736 
737  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]);
738  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]);
739  m_D[0][KFitConst::kNumber7 * m_TrackCount + 2] = 0.;
740  }
741 
743 }
744 
745 
746 enum KFitError::ECode
749 
751 }
752 
754 {
755  MakeMotherKFit kmm;
757  unsigned n = getTrackCount();
758  for (unsigned i = 0; i < n; ++i) {
760  getTrack(i).getCharge());
761  if (getFlagFitWithVertex())
763  for (unsigned j = i + 1; j < n; ++j) {
764  kmm.setCorrelation(getCorrelation(i, j));
765  }
766  }
767  kmm.setVertex(getVertex());
768  if (getFlagFitWithVertex())
770  m_ErrorCode = kmm.doMake();
772  return m_ErrorCode;
773  double chi2 = getCHIsq();
774  int ndf = getNDF();
775  double prob = TMath::Prob(chi2, ndf);
776  //
777  bool haschi2 = mother->hasExtraInfo("chiSquared");
778  if (haschi2) {
779  mother->setExtraInfo("chiSquared", chi2);
780  mother->setExtraInfo("ndf", ndf);
781  } else {
782  mother->addExtraInfo("chiSquared", chi2);
783  mother->addExtraInfo("ndf", ndf);
784  }
785 
786  mother->updateMomentum(
787  CLHEPToROOT::getLorentzVector(kmm.getMotherMomentum()),
788  CLHEPToROOT::getXYZVector(kmm.getMotherPosition()),
789  CLHEPToROOT::getTMatrixFSym(kmm.getMotherError()),
790  prob);
792  return m_ErrorCode;
793 }
static const double speedOfLight
[cm/ns]
Definition: Const.h:685
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:1344
bool hasExtraInfo(const std::string &name) const
Return whether the extra info with the given name is set.
Definition: Particle.cc:1293
void addExtraInfo(const std::string &name, double value)
Sets the user-defined data of given name to the given value.
Definition: Particle.cc:1363
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:373
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 const int kMaxTrackCount
Maximum track size.
Definition: KFitConst.h:40
static const int kMaxTrackCount2
Maximum track size (internal use)
Definition: KFitConst.h:42
static const int kAfterFit
Input parameter to specify after-fit when setting/getting a track attribute.
Definition: KFitConst.h:37
static const int kBeforeFit
Input parameter to specify before-fit when setting/getting a track attribute.
Definition: KFitConst.h:35
static const int kNumber7
Constant 7 to check matrix size (internal use)
Definition: KFitConst.h:32