Belle II Software  release-06-02-00
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 #include <TLorentzVector.h>
18 
19 
20 using namespace std;
21 using namespace Belle2;
22 using namespace Belle2::analysis;
23 using namespace CLHEP;
24 
25 MassFourCFitKFit::MassFourCFitKFit() : m_AfterVertexError(HepSymMatrix(3, 0)),
26  m_FourMomentum(TLorentzVector())
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 
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 TLorentzVector& 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  //TODO: Not Implemented
375  // m_FlagFitIncludingVertex == true
376  int index = 0;
377  m_al_0 = HepMatrix(KFitConst::kNumber7 * m_TrackCount + 3, 1, 0);
378  m_property = HepMatrix(m_TrackCount, 3, 0);
379  m_V_al_0 = HepSymMatrix(KFitConst::kNumber7 * m_TrackCount + 3, 0);
380 
381  for (auto& track : m_Tracks)
382  {
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  {
411  this->prepareCorrelation();
413  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
414  return m_ErrorCode;
415  }
416  }
417 
418  // set member matrix
419  m_al_1 = m_al_0;
420 
421  // define size of matrix
423  m_D = m_V_al_1.sub(1, 4, 1, KFitConst::kNumber7 * m_TrackCount + 3);
424  }
425 
427 }
428 
429 
430 enum KFitError::ECode
432  char buf[1024];
433  sprintf(buf, "%s:%s(): internal error; this function should never be called", __FILE__, __func__);
434  B2FATAL(buf);
435 
436  /* NEVER REACHEd HERE */
437  return KFitError::kOutOfRange;
438 }
439 
440 
441 enum KFitError::ECode
443  if (m_BeforeCorrelation.size() != static_cast<unsigned int>(m_TrackCount * (m_TrackCount - 1) / 2))
444  {
446  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
447  return m_ErrorCode;
448  }
449 
450  int row = 0, col = 0;
451 
452  for (auto& hm : m_BeforeCorrelation)
453  {
454  // counter
455  row++;
456  if (row == m_TrackCount) {
457  col++;
458  row = col + 1;
459  }
460 
461  int ii = 0, jj = 0;
462  for (int i = KFitConst::kNumber7 * row; i < KFitConst::kNumber7 * (row + 1); i++) {
463  for (int j = KFitConst::kNumber7 * col; j < KFitConst::kNumber7 * (col + 1); j++) {
464  m_V_al_0[i][j] = hm[ii][jj];
465  jj++;
466  }
467  jj = 0;
468  ii++;
469  }
470  }
471 
473  {
474  //TODO: Not Implemented
476 
477  // ...error of vertex
479 
480  // ...error matrix between vertex and tracks
482  if (m_BeforeTrackVertexError.size() != (unsigned int)m_TrackCount) {
484  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
485  return m_ErrorCode;
486  }
487 
488  int i = 0;
489  for (auto& hm : m_BeforeTrackVertexError) {
490  for (int j = 0; j < 3; j++) for (int k = 0; k < KFitConst::kNumber7; k++) {
492  }
493  i++;
494  }
495  }
496  }
497 
499 }
500 
501 
502 enum KFitError::ECode
504  Hep3Vector h3v;
505  int index = 0;
506  for (auto& pdata : m_Tracks)
507  {
508  // tracks
509  // momentum
510  h3v.setX(m_al_1[index * KFitConst::kNumber7 + 0][0]);
511  h3v.setY(m_al_1[index * KFitConst::kNumber7 + 1][0]);
512  h3v.setZ(m_al_1[index * KFitConst::kNumber7 + 2][0]);
513  pdata.setMomentum(HepLorentzVector(h3v, m_al_1[index * KFitConst::kNumber7 + 3][0]), KFitConst::kAfterFit);
514  // position
515  pdata.setPosition(HepPoint3D(
516  m_al_1[index * KFitConst::kNumber7 + 4][0],
517  m_al_1[index * KFitConst::kNumber7 + 5][0],
519  // error of the tracks
520  pdata.setError(this->makeError3(pdata.getMomentum(),
521  m_V_al_1.sub(
522  index * KFitConst::kNumber7 + 1,
523  (index + 1)*KFitConst::kNumber7,
524  index * KFitConst::kNumber7 + 1,
525  (index + 1)*KFitConst::kNumber7), m_IsFixMass[index]),
527  if (m_ErrorCode != KFitError::kNoError) break;
528  index++;
529  }
530 
532  {
533  //TODO: Not Implemented
535  // vertex
539  // error of the vertex
540  for (int i = 0; i < 3; i++) for (int j = i; j < 3; j++) {
542  }
543  // error between vertex and tracks
544  for (int i = 0; i < m_TrackCount; i++) {
545  HepMatrix hm(3, KFitConst::kNumber7, 0);
546  for (int j = 0; j < 3; j++) for (int k = 0; k < KFitConst::kNumber7; k++) {
548  }
549  if (m_IsFixMass[i])
550  m_AfterTrackVertexError.push_back(this->makeError4(m_Tracks[i].getMomentum(), hm));
551  else
552  m_AfterTrackVertexError.push_back(hm);
553  }
554  } else {
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  energy[i] = sqrt(al_1_prime[i * KFitConst::kNumber7 + 0][0] * al_1_prime[i * KFitConst::kNumber7 + 0][0] +
576  al_1_prime[i * KFitConst::kNumber7 + 1][0] * al_1_prime[i * KFitConst::kNumber7 + 1][0] +
577  al_1_prime[i * KFitConst::kNumber7 + 2][0] * al_1_prime[i * KFitConst::kNumber7 + 2][0] +
578  m_property[i][1] * m_property[i][1]);
579  }
580 
581  for (int i = 0; i < m_TrackCount; i++) {
582  // 3->4
583  for (int j = 0; j < 4; j++) Sum_al_1[j][0] += al_1_prime[i * KFitConst::kNumber7 + j][0];
584  }
585 
586  for (int i = 0; i < m_ConstraintMassCount; i++) {
587  for (int k = m_ConstraintMassChildLists[i].first; k <= m_ConstraintMassChildLists[i].second; k++) {
588  for (int j = 0; j < 4; j++) Sum_child_al_1[i * 4 + j][0] += al_1_prime[k * KFitConst::kNumber7 + j][0];
589  }
590  }
591 
592  m_d[0][0] = Sum_al_1[0][0] - m_FourMomentum.Px();
593  m_d[1][0] = Sum_al_1[1][0] - m_FourMomentum.Py();
594  m_d[2][0] = Sum_al_1[2][0] - m_FourMomentum.Pz();
595  m_d[3][0] = Sum_al_1[3][0] - m_FourMomentum.E();
596 
597  for (int i = 0; i < m_ConstraintMassCount; i++) {
598  m_d[4 + i][0] =
599  + 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]
600  - 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]
602  }
603 
604  for (int i = 0; i < m_TrackCount; i++) {
605  if (energy[i] == 0) {
607  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
608  break;
609  }
610 
611  // four-momentum conservation constraint
612  for (int l = 0; l < 4; l++) {
613  for (int n = 0; n < 6; n++) {
614  if (l == n) m_D[l][i * KFitConst::kNumber7 + n] = 1;
615  else m_D[l][i * KFitConst::kNumber7 + n] = 0;
616  }
617  }
618 
619  // invariant mass constraint
620  for (int l = 0; l < m_ConstraintMassCount; l++) {
621  if (i >= m_ConstraintMassChildLists[l].first && i <= m_ConstraintMassChildLists[l].second) {
622  a = m_property[i][2];
623  if (!m_FlagAtDecayPoint) a = 0.;
624  double invE = 1. / energy[i];
625  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 -
626  Sum_child_al_1[l * 4 + 0][0]);
627  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 -
628  Sum_child_al_1[l * 4 + 1][0]);
629  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 -
630  Sum_child_al_1[l * 4 + 2][0]);
631  m_D[4 + l][i * KFitConst::kNumber7 + 3] = 0.;
632  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 -
633  Sum_child_al_1[l * 4 + 1][0]) * a;
634  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 -
635  Sum_child_al_1[l * 4 + 0][0]) * a;
636  m_D[4 + l][i * KFitConst::kNumber7 + 6] = 0.;
637  } else {
638  for (int n = 0; n < 6; n++) {
639  m_D[4 + l][i * KFitConst::kNumber7 + n] = 0;
640  }
641  }
642  }
643  }
644 
645  } else {
646  //TODO: Not Implemented
648 
649  // m_FlagFitIncludingVertex == true
650  HepMatrix al_1_prime(m_al_1);
651  HepMatrix Sum_al_1(7, 1, 0);
652  double energy[KFitConst::kMaxTrackCount2];
653 
654  for (int i = 0; i < m_TrackCount; i++)
655  {
656  const double a = m_property[i][2];
657  al_1_prime[i * KFitConst::kNumber7 + 0][0] -= a * (al_1_prime[KFitConst::kNumber7 * m_TrackCount + 1][0] - al_1_prime[i *
658  KFitConst::kNumber7 + 5][0]);
659  al_1_prime[i * KFitConst::kNumber7 + 1][0] += a * (al_1_prime[KFitConst::kNumber7 * m_TrackCount + 0][0] - al_1_prime[i *
660  KFitConst::kNumber7 + 4][0]);
661  energy[i] = sqrt(al_1_prime[i * KFitConst::kNumber7 + 0][0] * al_1_prime[i * KFitConst::kNumber7 + 0][0] +
662  al_1_prime[i * KFitConst::kNumber7 + 1][0] * al_1_prime[i * KFitConst::kNumber7 + 1][0] +
663  al_1_prime[i * KFitConst::kNumber7 + 2][0] * al_1_prime[i * KFitConst::kNumber7 + 2][0] +
664  m_property[i][1] * m_property[i][1]);
665  Sum_al_1[6][0] = + a;
666  }
667 
668  for (int i = 0; i < m_TrackCount; i++)
669  {
670  if (energy[i] == 0) {
672  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
673  break;
674  }
675 
676  if (m_IsFixMass[i]) {
677  double invE = 1. / energy[i];
678  Sum_al_1[3][0] += energy[i];
679  Sum_al_1[4][0] += al_1_prime[i * KFitConst::kNumber7 + 1][0] * m_property[i][2] * invE;
680  Sum_al_1[5][0] += al_1_prime[i * KFitConst::kNumber7 + 0][0] * m_property[i][2] * invE;
681  } else {
682  Sum_al_1[3][0] += al_1_prime[i * KFitConst::kNumber7 + 3][0];
683  }
684 
685  for (int j = 0; j < 3; j++) Sum_al_1[j][0] += al_1_prime[i * KFitConst::kNumber7 + j][0];
686  }
687 
688  m_d[0][0] = Sum_al_1[0][0] - m_FourMomentum.Px();
689  m_d[1][0] = Sum_al_1[1][0] - m_FourMomentum.Py();
690  m_d[2][0] = Sum_al_1[2][0] - m_FourMomentum.Pz();
691  m_d[3][0] = Sum_al_1[3][0] - m_FourMomentum.E();
692 
693  for (int i = 0; i < m_TrackCount; i++)
694  {
695  if (energy[i] == 0) {
697  KFitError::displayError(__FILE__, __LINE__, __func__, m_ErrorCode);
698  break;
699  }
700 
701  for (int l = 0; l < 4; l++) {
702  for (int n = 0; n < 6; n++) {
703  if (l == n) m_D[l][i * KFitConst::kNumber7 + n] = 1;
704  else m_D[l][i * KFitConst::kNumber7 + n] = 0;
705  }
706  }
707  }
708 
709  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]);
710  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]);
711  m_D[0][KFitConst::kNumber7 * m_TrackCount + 2] = 0.;
712  }
713 
715 }
716 
717 
718 enum KFitError::ECode
721 
723 }
724 
726 {
727  MakeMotherKFit kmm;
729  unsigned n = getTrackCount();
730  for (unsigned i = 0; i < n; ++i) {
732  getTrack(i).getCharge());
733  if (getFlagFitWithVertex())
735  for (unsigned j = i + 1; j < n; ++j) {
736  kmm.setCorrelation(getCorrelation(i, j));
737  }
738  }
739  kmm.setVertex(getVertex());
740  if (getFlagFitWithVertex())
742  m_ErrorCode = kmm.doMake();
744  return m_ErrorCode;
745  double chi2 = getCHIsq();
746  int ndf = getNDF();
747  double prob = TMath::Prob(chi2, ndf);
748  //
749  bool haschi2 = mother->hasExtraInfo("chiSquared");
750  if (haschi2) {
751  mother->setExtraInfo("chiSquared", chi2);
752  mother->setExtraInfo("ndf", ndf);
753  } else {
754  mother->addExtraInfo("chiSquared", chi2);
755  mother->addExtraInfo("ndf", ndf);
756  }
757 
758  mother->updateMomentum(
759  CLHEPToROOT::getTLorentzVector(kmm.getMotherMomentum()),
760  CLHEPToROOT::getTVector3(kmm.getMotherPosition()),
761  CLHEPToROOT::getTMatrixFSym(kmm.getMotherError()),
762  prob);
764  return m_ErrorCode;
765 }
Class to store reconstructed particles.
Definition: Particle.h:74
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.
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 setFourMomentum(const TLorentzVector &m)
Set an 4 Momentum for the mass-four-constraint fit.
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.
TLorentzVector m_FourMomentum
Four Momentum.
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.
double getInvariantMass(void) const
Get an invariant mass.
Abstract base class for different kinds of events.
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