CMS 3D CMS Logo

GblTrajectory.cc
Go to the documentation of this file.
1 /*
2  * GblTrajectory.cpp
3  *
4  * Created on: Aug 18, 2011
5  * Author: kleinwrt
6  */
7 
139 using namespace Eigen;
140 
142 namespace gbl {
143 
145 
153  GblTrajectory::GblTrajectory(const std::vector<GblPoint> &aPointList,
154  bool flagCurv, bool flagU1dir, bool flagU2dir) :
155  numAllPoints(aPointList.size()), numPoints(), numOffsets(0),
156  numInnerTrans(0), numCurvature(flagCurv ? 1 : 0), numParameters(0),
157  numLocals(0), numMeasurements(0), externalPoint(0), skippedMeasLabel(0),
158  maxNumGlobals(0), theDimension(0), thePoints(), theData(), measDataIndex(),
159  scatDataIndex(), externalSeed(), innerTransformations(),
160  externalDerivatives(), externalMeasurements(), externalPrecisions()
161  {
162 
163  if (flagU1dir)
164  theDimension.push_back(0);
165  if (flagU2dir)
166  theDimension.push_back(1);
167  // simple (single) trajectory
168  thePoints.emplace_back(std::move(aPointList));
169  numPoints.push_back(numAllPoints);
170  construct(); // construct trajectory
171  }
172 
174 
178  GblTrajectory::GblTrajectory(const std::vector<std::pair<std::vector<GblPoint>, MatrixXd> > &aPointsAndTransList) :
180  numInnerTrans(aPointsAndTransList.size()), numParameters(0), numLocals(0),
185  {
186 
187  for (unsigned int iTraj = 0; iTraj < aPointsAndTransList.size(); ++iTraj) {
188  thePoints.push_back(aPointsAndTransList[iTraj].first);
189  numPoints.push_back(thePoints.back().size());
190  numAllPoints += numPoints.back();
191  innerTransformations.push_back(aPointsAndTransList[iTraj].second);
192  }
193  theDimension.push_back(0);
194  theDimension.push_back(1);
196  construct(); // construct (composed) trajectory
197  }
198 
199 #ifdef GBL_EIGEN_SUPPORT_ROOT
200 
212  GblTrajectory::GblTrajectory(const std::vector<GblPoint> &aPointList,
213  unsigned int aLabel, const TMatrixDSym &aSeed,
214  bool flagCurv, bool flagU1dir, bool flagU2dir) :
215  numAllPoints(aPointList.size()), numPoints(), numOffsets(0),
216  numInnerTrans(0), numCurvature(flagCurv ? 1 : 0), numParameters(0),
217  numLocals(0), numMeasurements(0), externalPoint(aLabel),
222  {
223 
224  if (flagU1dir)
225  theDimension.push_back(0);
226  if (flagU2dir)
227  theDimension.push_back(1);
228  // convert from ROOT
229  unsigned int nParSeed = aSeed.GetNrows();
230  externalSeed.resize(nParSeed, nParSeed);
231  for (unsigned int i = 0; i < nParSeed; ++i) {
232  for (unsigned int j = 0; j < nParSeed; ++j) {
233  externalSeed(i, j) = aSeed(i, j);
234  }
235  }
236  // simple (single) trajectory
237  thePoints.emplace_back(std::move(aPointList));
238  numPoints.push_back(numAllPoints);
239  construct(); // construct trajectory
240  }
241 
243 
247  GblTrajectory::GblTrajectory(const std::vector<std::pair<std::vector<GblPoint>, TMatrixD> > &aPointsAndTransList) :
249  numInnerTrans(aPointsAndTransList.size()), numParameters(0), numLocals(0),
254  {
255 
256  for (unsigned int iTraj = 0; iTraj < aPointsAndTransList.size(); ++iTraj) {
257  thePoints.emplace_back(std::move(aPointsAndTransList[iTraj].first));
258  numPoints.push_back(thePoints.back().size());
259  numAllPoints += numPoints.back();
260  // convert from ROOT
261  unsigned int nRowTrans = aPointsAndTransList[iTraj].second.GetNrows();
262  unsigned int nColTrans = aPointsAndTransList[iTraj].second.GetNcols();
263  MatrixXd aTrans(nRowTrans, nColTrans);
264  for (unsigned int i = 0; i < nRowTrans; ++i) {
265  for (unsigned int j = 0; j < nColTrans; ++j) {
266  aTrans(i, j) = aPointsAndTransList[iTraj].second(i, j);
267  }
268  }
269  innerTransformations.emplace_back(std::move(aTrans));
270  }
271  theDimension.push_back(0);
272  theDimension.push_back(1);
274  construct(); // construct (composed) trajectory
275  }
276 
278 
285  GblTrajectory::GblTrajectory(const std::vector<std::pair<std::vector<GblPoint>, TMatrixD> > &aPointsAndTransList,
286  const TMatrixD &extDerivatives,
287  const TVectorD &extMeasurements,
288  const TVectorD &extPrecisions) :
290  numInnerTrans(aPointsAndTransList.size()), numParameters(0), numLocals(0),
295  {
296 
297  // convert from ROOT
298  unsigned int nExtMeas = extDerivatives.GetNrows();
299  unsigned int nExtPar = extDerivatives.GetNcols();
300  externalDerivatives.resize(nExtMeas, nExtPar);
301  externalMeasurements.resize(nExtMeas);
302  externalPrecisions.resize(nExtMeas);
303  for (unsigned int i = 0; i < nExtMeas; ++i) {
304  externalMeasurements(i) = extMeasurements[i];
305  externalPrecisions(i) = extPrecisions[i];
306  for (unsigned int j = 0; j < nExtPar; ++j) {
307  externalDerivatives(i, j) = extDerivatives(i, j);
308  }
309  }
310  for (unsigned int iTraj = 0; iTraj < aPointsAndTransList.size(); ++iTraj) {
311  thePoints.emplace_back(std::move(aPointsAndTransList[iTraj].first));
312  numPoints.push_back(thePoints.back().size());
313  numAllPoints += numPoints.back();
314  // convert from ROOT
315  unsigned int nRowTrans = aPointsAndTransList[iTraj].second.GetNrows();
316  unsigned int nColTrans = aPointsAndTransList[iTraj].second.GetNcols();
317  MatrixXd aTrans(nRowTrans, nColTrans);
318  for (unsigned int i = 0; i < nRowTrans; ++i) {
319  for (unsigned int j = 0; j < nColTrans; ++j) {
320  aTrans(i, j) = aPointsAndTransList[iTraj].second(i, j);
321  }
322  }
323  innerTransformations.emplace_back(std::move(aTrans));
324  }
325  theDimension.push_back(0);
326  theDimension.push_back(1);
328  construct(); // construct (composed) trajectory
329  }
330 
332 
339  GblTrajectory::GblTrajectory(const std::vector<std::pair<std::vector<GblPoint>, TMatrixD> > &aPointsAndTransList,
340  const TMatrixD &extDerivatives,
341  const TVectorD &extMeasurements,
342  const TMatrixDSym &extPrecisions) :
344  numInnerTrans(aPointsAndTransList.size()), numParameters(0), numLocals(0),
348  {
349 
350  // diagonalize external measurement
351  TMatrixDSymEigen extEigen(extPrecisions);
352  TMatrixD extTransformation = extEigen.GetEigenVectors();
353  extTransformation.T();
354  TMatrixD aDerivatives = extTransformation * extDerivatives;
355  TVectorD aMeasurements = extTransformation * extMeasurements;
356  TVectorD aPrecisions = extEigen.GetEigenValues();
357  // convert from ROOT
358  unsigned int nExtMeas = aDerivatives.GetNrows();
359  unsigned int nExtPar = aDerivatives.GetNcols();
360  externalDerivatives.resize(nExtMeas, nExtPar);
361  externalMeasurements.resize(nExtMeas);
362  externalPrecisions.resize(nExtMeas);
363  for (unsigned int i = 0; i < nExtMeas; ++i) {
364  externalMeasurements(i) = aMeasurements[i];
365  externalPrecisions(i) = aPrecisions[i];
366  for (unsigned int j = 0; j < nExtPar; ++j) {
367  externalDerivatives(i, j) = aDerivatives(i, j);
368  }
369  }
370  for (unsigned int iTraj = 0; iTraj < aPointsAndTransList.size(); ++iTraj) {
371  thePoints.push_back(aPointsAndTransList[iTraj].first);
372  numPoints.push_back(thePoints.back().size());
373  numAllPoints += numPoints.back();
374  // convert from ROOT
375  unsigned int nRowTrans = aPointsAndTransList[iTraj].second.GetNrows();
376  unsigned int nColTrans = aPointsAndTransList[iTraj].second.GetNcols();
377  MatrixXd aTrans(nRowTrans, nColTrans);
378  for (unsigned int i = 0; i < nRowTrans; ++i) {
379  for (unsigned int j = 0; j < nColTrans; ++j) {
380  aTrans(i, j) = aPointsAndTransList[iTraj].second(i, j);
381  }
382  }
383  innerTransformations.push_back(aTrans);
384  }
385  theDimension.push_back(0);
386  theDimension.push_back(1);
388  construct(); // construct (composed) trajectory
389  }
390 #endif
391 
393  }
394 
396  bool GblTrajectory::isValid() const {
397  return constructOK;
398  }
399 
401  unsigned int GblTrajectory::getNumPoints() const {
402  return numAllPoints;
403  }
404 
406 
410 
411  constructOK = false;
412  fitOK = false;
413  unsigned int aLabel = 0;
414  if (numAllPoints < 2) {
415  std::cout << " GblTrajectory construction failed: too few GblPoints "
416  << std::endl;
417  return;
418  }
419  // loop over trajectories
420  numTrajectories = thePoints.size();
421  //std::cout << " numTrajectories: " << numTrajectories << ", " << innerTransformations.size() << std::endl;
422  for (unsigned int iTraj = 0; iTraj < numTrajectories; ++iTraj) {
423  std::vector<GblPoint>::iterator itPoint;
424  for (itPoint = thePoints[iTraj].begin();
425  itPoint < thePoints[iTraj].end(); ++itPoint) {
426  numLocals = std::max(numLocals, itPoint->getNumLocals());
427  numMeasurements += itPoint->hasMeasurement();
428  itPoint->setLabel(++aLabel);
429  }
430  }
431  defineOffsets();
432  calcJacobians();
433  try {
434  prepare();
435  } catch (std::overflow_error &e) {
436  std::cout << " GblTrajectory construction failed: " << e.what()
437  << std::endl;
438  return;
439  }
440 
441  constructOK = true;
442  // number of fit parameters
445  }
446 
448 
453 
454  // loop over trajectories
455  for (unsigned int iTraj = 0; iTraj < numTrajectories; ++iTraj) {
456  // first point is offset
457  thePoints[iTraj].front().setOffset(numOffsets++);
458  // intermediate scatterers are offsets
459  std::vector<GblPoint>::iterator itPoint;
460  for (itPoint = thePoints[iTraj].begin() + 1;
461  itPoint < thePoints[iTraj].end() - 1; ++itPoint) {
462  if (itPoint->hasScatterer()) {
463  itPoint->setOffset(numOffsets++);
464  } else {
465  itPoint->setOffset(-numOffsets);
466  }
467  }
468  // last point is offset
469  thePoints[iTraj].back().setOffset(numOffsets++);
470  }
471  }
472 
475 
476  Matrix5d scatJacobian;
477  // loop over trajectories
478  for (unsigned int iTraj = 0; iTraj < numTrajectories; ++iTraj) {
479  // forward propagation (all)
480  GblPoint* previousPoint = &thePoints[iTraj].front();
481  unsigned int numStep = 0;
482  std::vector<GblPoint>::iterator itPoint;
483  for (itPoint = thePoints[iTraj].begin() + 1;
484  itPoint < thePoints[iTraj].end(); ++itPoint) {
485  if (numStep == 0) {
486  scatJacobian = itPoint->getP2pJacobian();
487  } else {
488  scatJacobian = itPoint->getP2pJacobian() * scatJacobian;
489  }
490  numStep++;
491  itPoint->addPrevJacobian(scatJacobian); // iPoint -> previous scatterer
492  if (itPoint->getOffset() >= 0) {
493  previousPoint->addNextJacobian(scatJacobian); // lastPoint -> next scatterer
494  numStep = 0;
495  previousPoint = &(*itPoint);
496  }
497  }
498  // backward propagation (without scatterers)
499  for (itPoint = thePoints[iTraj].end() - 1;
500  itPoint > thePoints[iTraj].begin(); --itPoint) {
501  if (itPoint->getOffset() >= 0) {
502  scatJacobian = itPoint->getP2pJacobian();
503  continue; // skip offsets
504  }
505  itPoint->addNextJacobian(scatJacobian); // iPoint -> next scatterer
506  scatJacobian = scatJacobian * itPoint->getP2pJacobian();
507  }
508  }
509  }
510 
512 
520  std::pair<std::vector<unsigned int>, MatrixXd> GblTrajectory::getJacobian(
521  int aSignedLabel) const {
522 
523  unsigned int nDim = theDimension.size();
524  unsigned int nCurv = numCurvature;
525  unsigned int nLocals = numLocals;
526  unsigned int nBorder = nCurv + nLocals;
527  unsigned int nParBRL = nBorder + 2 * nDim;
528  unsigned int nParLoc = nLocals + 5;
529  std::vector<unsigned int> anIndex;
530  anIndex.reserve(nParBRL);
531  MatrixXd aJacobian(nParLoc, nParBRL);
532  aJacobian.setZero();
533 
534  unsigned int aLabel = abs(aSignedLabel);
535  unsigned int firstLabel = 1;
536  unsigned int lastLabel = 0;
537  unsigned int aTrajectory = 0;
538  // loop over trajectories
539  for (unsigned int iTraj = 0; iTraj < numTrajectories; ++iTraj) {
540  aTrajectory = iTraj;
541  lastLabel += numPoints[iTraj];
542  if (aLabel <= lastLabel)
543  break;
544  if (iTraj < numTrajectories - 1)
545  firstLabel += numPoints[iTraj];
546  }
547  int nJacobian; // 0: prev, 1: next
548  // check consistency of (index, direction)
549  if (aSignedLabel > 0) {
550  nJacobian = 1;
551  if (aLabel >= lastLabel) {
552  aLabel = lastLabel;
553  nJacobian = 0;
554  }
555  } else {
556  nJacobian = 0;
557  if (aLabel <= firstLabel) {
558  aLabel = firstLabel;
559  nJacobian = 1;
560  }
561  }
562  const GblPoint aPoint = thePoints[aTrajectory][aLabel - firstLabel];
563  std::array<unsigned int, 5> labDer;
564  Matrix5d matDer;
565  getFitToLocalJacobian(labDer, matDer, aPoint, 5, nJacobian);
566 
567  // from local parameters
568  for (unsigned int i = 0; i < nLocals; ++i) {
569  aJacobian(i + 5, i) = 1.0;
570  anIndex.push_back(i + 1);
571  }
572  // from trajectory parameters
573  unsigned int iCol = nLocals;
574  for (unsigned int i = 0; i < 5; ++i) {
575  if (labDer[i] > 0) {
576  anIndex.push_back(labDer[i]);
577  for (unsigned int j = 0; j < 5; ++j) {
578  aJacobian(j, iCol) = matDer(j, i);
579  }
580  ++iCol;
581  }
582  }
583  return std::make_pair(anIndex, aJacobian);
584  }
585 
587 
596  void GblTrajectory::getFitToLocalJacobian(std::array<unsigned int, 5>& anIndex,
597  Matrix5d &aJacobian, const GblPoint &aPoint, unsigned int measDim,
598  unsigned int nJacobian) const {
599 
600  unsigned int nDim = theDimension.size();
601  unsigned int nCurv = numCurvature;
602  unsigned int nLocals = numLocals;
603 
604  int nOffset = aPoint.getOffset();
605 
606  aJacobian.setZero();
607  if (nOffset < 0) // need interpolation
608  {
609  Matrix2d prevW, prevWJ, nextW, nextWJ, matN;
610  Vector2d prevWd, nextWd;
611  aPoint.getDerivatives(0, prevW, prevWJ, prevWd); // W-, W- * J-, W- * d-
612  aPoint.getDerivatives(1, nextW, nextWJ, nextWd); // W+, W+ * J+, W+ * d+
613  const Matrix2d sumWJ(prevWJ + nextWJ);
614  matN = sumWJ.inverse(); // N = (W- * J- + W+ * J+)^-1
615  // derivatives for u_int
616  const Matrix2d prevNW(matN * prevW); // N * W-
617  const Matrix2d nextNW(matN * nextW); // N * W+
618  const Vector2d prevNd(matN * prevWd); // N * W- * d-
619  const Vector2d nextNd(matN * nextWd); // N * W+ * d+
620 
621  unsigned int iOff = nDim * (-nOffset - 1) + nLocals + nCurv + 1; // first offset ('i' in u_i)
622 
623  // local offset
624  if (nCurv > 0) {
625  aJacobian.block<2, 1>(3, 0) = -prevNd - nextNd; // from curvature
626  anIndex[0] = nLocals + 1;
627  }
628  aJacobian.block<2, 2>(3, 1) = prevNW; // from 1st Offset
629  aJacobian.block<2, 2>(3, 3) = nextNW; // from 2nd Offset
630  for (unsigned int i = 0; i < nDim; ++i) {
631  anIndex[1 + theDimension[i]] = iOff + i;
632  anIndex[3 + theDimension[i]] = iOff + nDim + i;
633  }
634 
635  // local slope and curvature
636  if (measDim > 2) {
637  // derivatives for u'_int
638  const Matrix2d prevWPN(nextWJ * prevNW); // W+ * J+ * N * W-
639  const Matrix2d nextWPN(prevWJ * nextNW); // W- * J- * N * W+
640  const Vector2d prevWNd(nextWJ * prevNd); // W+ * J+ * N * W- * d-
641  const Vector2d nextWNd(prevWJ * nextNd); // W- * J- * N * W+ * d+
642  if (nCurv > 0) {
643  aJacobian(0, 0) = 1.0;
644  aJacobian.block<2, 1>(1, 0) = prevWNd - nextWNd; // from curvature
645  }
646  aJacobian.block<2, 2>(1, 1) = -prevWPN; // from 1st Offset
647  aJacobian.block<2, 2>(1, 3) = nextWPN; // from 2nd Offset
648  }
649  } else { // at point
650  // anIndex must be sorted
651  // forward : iOff2 = iOff1 + nDim, index1 = 1, index2 = 3
652  // backward: iOff2 = iOff1 - nDim, index1 = 3, index2 = 1
653  unsigned int iOff1 = nDim * nOffset + nCurv + nLocals + 1; // first offset ('i' in u_i)
654  unsigned int index1 = 3 - 2 * nJacobian; // index of first offset
655  unsigned int iOff2 = iOff1 + nDim * (nJacobian * 2 - 1); // second offset ('i' in u_i)
656  unsigned int index2 = 1 + 2 * nJacobian; // index of second offset
657  // local offset
658  aJacobian(3, index1) = 1.0; // from 1st Offset
659  aJacobian(4, index1 + 1) = 1.0;
660  for (unsigned int i = 0; i < nDim; ++i) {
661  anIndex[index1 + theDimension[i]] = iOff1 + i;
662  }
663 
664  // local slope and curvature
665  if (measDim > 2) {
666  Matrix2d matW, matWJ;
667  Vector2d vecWd;
668  aPoint.getDerivatives(nJacobian, matW, matWJ, vecWd); // W, W * J, W * d
669  double sign = (nJacobian > 0) ? 1. : -1.;
670  if (nCurv > 0) {
671  aJacobian(0, 0) = 1.0;
672  aJacobian.block<2, 1>(1, 0) = -sign * vecWd; // from curvature
673  anIndex[0] = nLocals + 1;
674  }
675  aJacobian.block<2, 2>(1, index1) = -sign * matWJ; // from 1st Offset
676  aJacobian.block<2, 2>(1, index2) = sign * matW; // from 2nd Offset
677  for (unsigned int i = 0; i < nDim; ++i) {
678  anIndex[index2 + theDimension[i]] = iOff2 + i;
679  }
680  }
681  }
682  }
683 
685 
691  void GblTrajectory::getFitToKinkJacobian(std::array<unsigned int, 7>& anIndex,
692  Matrix27d &aJacobian, const GblPoint &aPoint) const {
693 
694  unsigned int nDim = theDimension.size();
695  unsigned int nCurv = numCurvature;
696  unsigned int nLocals = numLocals;
697 
698  int nOffset = aPoint.getOffset();
699 
700  aJacobian.setZero();
701 
702  Matrix2d prevW, prevWJ, nextW, nextWJ;
703  Vector2d prevWd, nextWd;
704  aPoint.getDerivatives(0, prevW, prevWJ, prevWd); // W-, W- * J-, W- * d-
705  aPoint.getDerivatives(1, nextW, nextWJ, nextWd); // W-, W- * J-, W- * d-
706  const Matrix2d sumWJ(prevWJ + nextWJ); // W- * J- + W+ * J+
707  const Vector2d sumWd(prevWd + nextWd); // W+ * d+ + W- * d-
708 
709  unsigned int iOff = (nOffset - 1) * nDim + nCurv + nLocals + 1; // first offset ('i' in u_i)
710 
711  // local offset
712  if (nCurv > 0) {
713  aJacobian.block<2, 1>(0, 0) = -sumWd; // from curvature
714  anIndex[0] = nLocals + 1;
715  }
716  aJacobian.block<2, 2>(0, 1) = prevW; // from 1st Offset
717  aJacobian.block<2, 2>(0, 3) = -sumWJ; // from 2nd Offset
718  aJacobian.block<2, 2>(0, 5) = nextW; // from 1st Offset
719  for (unsigned int i = 0; i < nDim; ++i) {
720  anIndex[1 + theDimension[i]] = iOff + i;
721  anIndex[3 + theDimension[i]] = iOff + nDim + i;
722  anIndex[5 + theDimension[i]] = iOff + nDim * 2 + i;
723  }
724  }
725 
727 
741  unsigned int GblTrajectory::getResults(int aSignedLabel, VectorXd &localPar,
742  MatrixXd &localCov) const {
743  if (not fitOK)
744  return 1;
745  std::pair<std::vector<unsigned int>, MatrixXd> indexAndJacobian =
746  getJacobian(aSignedLabel);
747  unsigned int nParBrl = indexAndJacobian.first.size();
748  VectorXd aVec(nParBrl); // compressed vector
749  for (unsigned int i = 0; i < nParBrl; ++i) {
750  aVec[i] = theVector(indexAndJacobian.first[i] - 1);
751  }
752  MatrixXd aMat = theMatrix.getBlockMatrix(indexAndJacobian.first); // compressed matrix
753  localPar = indexAndJacobian.second * aVec;
754  localCov = indexAndJacobian.second * aMat
755  * indexAndJacobian.second.adjoint();
756  return 0;
757  }
758 
760 
772  unsigned int GblTrajectory::getMeasResults(unsigned int aLabel,
773  unsigned int &numData, VectorXd &aResiduals, VectorXd &aMeasErrors,
774  VectorXd &aResErrors, VectorXd &aDownWeights) {
775  numData = 0;
776  if (not fitOK)
777  return 1;
778 
779  unsigned int firstData = measDataIndex[aLabel - 1]; // first data block with measurement
780  numData = measDataIndex[aLabel] - firstData; // number of data blocks
781  for (unsigned int i = 0; i < numData; ++i) {
782  getResAndErr(firstData + i, (aLabel != skippedMeasLabel), aResiduals(i),
783  aMeasErrors(i), aResErrors(i), aDownWeights(i));
784  }
785  return 0;
786  }
787 
789 
801  unsigned int GblTrajectory::getScatResults(unsigned int aLabel,
802  unsigned int &numData, VectorXd &aResiduals, VectorXd &aMeasErrors,
803  VectorXd &aResErrors, VectorXd &aDownWeights) {
804  numData = 0;
805  if (not fitOK)
806  return 1;
807 
808  unsigned int firstData = scatDataIndex[aLabel - 1]; // first data block with scatterer
809  numData = scatDataIndex[aLabel] - firstData; // number of data blocks
810  for (unsigned int i = 0; i < numData; ++i) {
811  getResAndErr(firstData + i, true, aResiduals(i), aMeasErrors(i),
812  aResErrors(i), aDownWeights(i));
813  }
814  return 0;
815  }
816 
817 #ifdef GBL_EIGEN_SUPPORT_ROOT
818 
833  unsigned int GblTrajectory::getResults(int aSignedLabel, TVectorD &localPar,
834  TMatrixDSym &localCov) const {
835  if (not fitOK)
836  return 1;
837  std::pair<std::vector<unsigned int>, MatrixXd> indexAndJacobian =
838  getJacobian(aSignedLabel);
839  unsigned int nParBrl = indexAndJacobian.first.size();
840  VectorXd aVec(nParBrl); // compressed vector
841  for (unsigned int i = 0; i < nParBrl; ++i) {
842  aVec[i] = theVector(indexAndJacobian.first[i] - 1);
843  }
844  MatrixXd aMat = theMatrix.getBlockMatrix(indexAndJacobian.first); // compressed matrix
845  VectorXd aLocalPar = indexAndJacobian.second * aVec;
846  MatrixXd aLocalCov = indexAndJacobian.second * aMat
847  * indexAndJacobian.second.adjoint();
848  // convert to ROOT
849  unsigned int nParOut = localPar.GetNrows();
850  for (unsigned int i = 0; i < nParOut; ++i) {
851  localPar[i] = aLocalPar(i);
852  for (unsigned int j = 0; j < nParOut; ++j) {
853  localCov(i, j) = aLocalCov(i, j);
854  }
855  }
856  return 0;
857  }
858 
860 
872  unsigned int GblTrajectory::getMeasResults(unsigned int aLabel,
873  unsigned int &numData, TVectorD &aResiduals, TVectorD &aMeasErrors,
874  TVectorD &aResErrors, TVectorD &aDownWeights) {
875  numData = 0;
876  if (not fitOK)
877  return 1;
878 
879  unsigned int firstData = measDataIndex[aLabel - 1]; // first data block with measurement
880  numData = measDataIndex[aLabel] - firstData; // number of data blocks
881  for (unsigned int i = 0; i < numData; ++i) {
882  getResAndErr(firstData + i, (aLabel != skippedMeasLabel), aResiduals[i],
883  aMeasErrors[i], aResErrors[i], aDownWeights[i]);
884  }
885  return 0;
886  }
887 
889 
901  unsigned int GblTrajectory::getScatResults(unsigned int aLabel,
902  unsigned int &numData, TVectorD &aResiduals, TVectorD &aMeasErrors,
903  TVectorD &aResErrors, TVectorD &aDownWeights) {
904  numData = 0;
905  if (not fitOK)
906  return 1;
907 
908  unsigned int firstData = scatDataIndex[aLabel - 1]; // first data block with scatterer
909  numData = scatDataIndex[aLabel] - firstData; // number of data blocks
910  for (unsigned int i = 0; i < numData; ++i) {
911  getResAndErr(firstData + i, true, aResiduals[i], aMeasErrors[i],
912  aResErrors[i], aDownWeights[i]);
913  }
914  return 0;
915  }
916 
917 #endif
918 
920 
924  unsigned int GblTrajectory::getLabels(std::vector<unsigned int> &aLabelList) {
925  if (not constructOK)
926  return 1;
927 
928  unsigned int aLabel = 0;
929  unsigned int nPoint = thePoints[0].size();
930  aLabelList.resize(nPoint);
931  for (unsigned i = 0; i < nPoint; ++i) {
932  aLabelList[i] = ++aLabel;
933  }
934  return 0;
935  }
936 
938 
943  std::vector<std::vector<unsigned int> > &aLabelList) {
944  if (not constructOK)
945  return 1;
946 
947  unsigned int aLabel = 0;
948  aLabelList.resize(numTrajectories);
949  for (unsigned int iTraj = 0; iTraj < numTrajectories; ++iTraj) {
950  unsigned int nPoint = thePoints[iTraj].size();
951  aLabelList[iTraj].resize(nPoint);
952  for (unsigned i = 0; i < nPoint; ++i) {
953  aLabelList[iTraj][i] = ++aLabel;
954  }
955  }
956  return 0;
957  }
958 
960 
970  void GblTrajectory::getResAndErr(unsigned int aData, bool used,
971  double &aResidual, double &aMeasError, double &aResError,
972  double &aDownWeight) {
973 
974  double aMeasVar;
975  unsigned int numLocal;
976  unsigned int* indLocal;
977  double* derLocal;
978  theData[aData].getResidual(aResidual, aMeasVar, aDownWeight, numLocal,
979  indLocal, derLocal);
980  VectorXd aVec(numLocal); // compressed vector of derivatives
981  for (unsigned int j = 0; j < numLocal; ++j) {
982  aVec[j] = derLocal[j];
983  }
984  MatrixXd aMat = theMatrix.getBlockMatrix(numLocal, indLocal); // compressed (covariance) matrix
985  double aFitVar = aVec.transpose() * aMat * aVec; // variance from track fit
986  aMeasError = sqrt(aMeasVar); // error of measurement
987  if (used)
988  aResError = (aFitVar < aMeasVar ? sqrt(aMeasVar - aFitVar) : 0.); // error of (biased) residual
989  else
990  aResError = sqrt(aMeasVar + aFitVar); // error of (unbiased) residual
991  }
992 
995  unsigned int nBorder = numCurvature + numLocals;
997  theMatrix.resize(numParameters, nBorder);
998  double aValue, aWeight;
999  unsigned int* indLocal;
1000  double* derLocal;
1001  unsigned int numLocal;
1002 
1003  std::vector<GblData>::iterator itData;
1004  for (itData = theData.begin(); itData < theData.end(); ++itData) {
1005  // skipped (internal) measurement ?
1006  if (itData->getLabel() == skippedMeasLabel
1007  && itData->getType() == InternalMeasurement)
1008  continue;
1009  itData->getLocalData(aValue, aWeight, numLocal, indLocal, derLocal);
1010  for (unsigned int j = 0; j < numLocal; ++j) {
1011  theVector(indLocal[j] - 1) += derLocal[j] * aWeight * aValue;
1012  }
1013  theMatrix.addBlockMatrix(aWeight, numLocal, indLocal, derLocal);
1014  }
1015  }
1016 
1018 
1022  unsigned int nDim = theDimension.size();
1023  // upper limit
1024  unsigned int maxData = numMeasurements + nDim * (numOffsets - 2)
1025  + externalSeed.rows();
1026  theData.reserve(maxData);
1027  measDataIndex.resize(numAllPoints + 3); // include external seed and measurements
1028  scatDataIndex.resize(numAllPoints + 1);
1029  unsigned int nData = 0;
1030  std::vector<MatrixXd> innerTransDer;
1031  std::vector<std::array<unsigned int, 5> > innerTransLab;
1032  // composed trajectory ?
1033  if (numInnerTrans > 0) {
1034  //std::cout << "composed trajectory" << std::endl;
1035  for (unsigned int iTraj = 0; iTraj < numTrajectories; ++iTraj) {
1036  // innermost point
1037  GblPoint* innerPoint = &thePoints[iTraj].front();
1038  // transformation fit to local track parameters
1039  std::array<unsigned int, 5> firstLabels;
1040  Matrix5d matFitToLocal, matLocalToFit;
1041  getFitToLocalJacobian(firstLabels, matFitToLocal, *innerPoint, 5);
1042  // transformation local track to fit parameters
1043  matLocalToFit = matFitToLocal.inverse();
1044  // transformation external to fit parameters at inner (first) point
1045  innerTransDer.emplace_back(matLocalToFit * innerTransformations[iTraj]);
1046  innerTransLab.push_back(firstLabels);
1047  }
1048  }
1049  Matrix5d matP; // measurements
1050  std::vector<GblPoint>::iterator itPoint;
1051  // limit the scope of proDer:
1052  {
1053  // transform for external parameters
1054  Eigen::Matrix<double, Eigen::Dynamic, 5,
1055  Eigen::ColMajor /* default */, 5, 5> proDer;
1056  // loop over trajectories
1057  for (unsigned int iTraj = 0; iTraj < numTrajectories; ++iTraj) {
1058  for (itPoint = thePoints[iTraj].begin();
1059  itPoint < thePoints[iTraj].end(); ++itPoint) {
1060  Vector5d aMeas, aPrec;
1061  unsigned int nLabel = itPoint->getLabel();
1062  unsigned int measDim = itPoint->hasMeasurement();
1063  if (measDim) {
1064  const MatrixXd localDer = itPoint->getLocalDerivatives();
1066  itPoint->getNumGlobals());
1067  MatrixXd transDer;
1068  itPoint->getMeasurement(matP, aMeas, aPrec);
1069  double minPrecision = itPoint->getMeasPrecMin();
1070  unsigned int iOff = 5 - measDim; // first active component
1071  std::array<unsigned int, 5> labDer;
1072  Matrix5d matDer, matPDer;
1073  unsigned int nJacobian =
1074  (itPoint < thePoints[iTraj].end() - 1) ? 1 : 0; // last point needs backward propagation
1075  getFitToLocalJacobian(labDer, matDer, *itPoint, measDim,
1076  nJacobian);
1077  if (measDim > 2) {
1078  matPDer = matP * matDer;
1079  } else { // 'shortcut' for position measurements
1080  matPDer.block<2, 5>(3, 0) = matP.block<2, 2>(3, 3)
1081  * matDer.block<2, 5>(3, 0);
1082  }
1083 
1084  if (numInnerTrans > 0) {
1085  // transform for external parameters
1086  proDer.resize(measDim, Eigen::NoChange);
1087  // match parameters
1088  unsigned int ifirst = 0;
1089  unsigned int ilabel = 0;
1090  while (ilabel < 5) {
1091  if (labDer[ilabel] > 0) {
1092  while (innerTransLab[iTraj][ifirst]
1093  != labDer[ilabel] and ifirst < 5) {
1094  ++ifirst;
1095  }
1096  if (ifirst >= 5) {
1097  labDer[ilabel] -= 2 * nDim * (iTraj + 1); // adjust label
1098  } else {
1099  // match
1100  labDer[ilabel] = 0; // mark as related to external parameters
1101  for (unsigned int k = iOff; k < 5; ++k) {
1102  proDer(k - iOff, ifirst) = matPDer(k,
1103  ilabel);
1104  }
1105  }
1106  }
1107  ++ilabel;
1108  }
1109  transDer.resize(measDim, numCurvature);
1110  transDer = proDer * innerTransDer[iTraj];
1111  }
1112  for (unsigned int i = iOff; i < 5; ++i) {
1113  if (aPrec(i) > minPrecision) {
1114  GblData aData(nLabel, InternalMeasurement, aMeas(i),
1115  aPrec(i), iTraj,
1116  itPoint - thePoints[iTraj].begin());
1117  aData.addDerivatives(i, labDer, matPDer, iOff, localDer,
1118  numLocals, transDer);
1119  theData.emplace_back(aData);
1120  nData++;
1121  }
1122  }
1123 
1124  }
1125  measDataIndex[nLabel] = nData;
1126  }
1127  }
1128  } // end of scope for proDer
1129 
1130  Matrix2d matT; // measurements
1131  // limit the scope of proDer:
1132  {
1133  // transform for external parameters
1134  Eigen::Matrix<double, Eigen::Dynamic, 5,
1135  Eigen::ColMajor /* default */, 5, 5> proDer;
1136  scatDataIndex[0] = nData;
1137  scatDataIndex[1] = nData;
1138  // loop over trajectories
1139  for (unsigned int iTraj = 0; iTraj < numTrajectories; ++iTraj) {
1140  for (itPoint = thePoints[iTraj].begin() + 1;
1141  itPoint < thePoints[iTraj].end() - 1; ++itPoint) {
1142  Vector2d aMeas, aPrec;
1143  unsigned int nLabel = itPoint->getLabel();
1144  if (itPoint->hasScatterer()) {
1145  itPoint->getScatterer(matT, aMeas, aPrec);
1146  MatrixXd transDer;
1147  std::array<unsigned int, 7> labDer;
1148  Matrix27d matDer, matTDer;
1149  getFitToKinkJacobian(labDer, matDer, *itPoint);
1150  matTDer = matT * matDer;
1151  if (numInnerTrans > 0) {
1152  // transform for external parameters
1153  proDer.resize(nDim, Eigen::NoChange);
1154  // match parameters
1155  unsigned int ifirst = 0;
1156  unsigned int ilabel = 0;
1157  while (ilabel < 7) {
1158  if (labDer[ilabel] > 0) {
1159  while (innerTransLab[iTraj][ifirst]
1160  != labDer[ilabel] and ifirst < 5) {
1161  ++ifirst;
1162  }
1163  if (ifirst >= 5) {
1164  labDer[ilabel] -= 2 * nDim * (iTraj + 1); // adjust label
1165  } else {
1166  // match
1167  labDer[ilabel] = 0; // mark as related to external parameters
1168  for (unsigned int k = 0; k < nDim; ++k) {
1169  proDer(k, ifirst) = matTDer(k, ilabel);
1170  }
1171  }
1172  }
1173  ++ilabel;
1174  }
1175  transDer.resize(nDim, numCurvature);
1176  transDer = proDer * innerTransDer[iTraj];
1177  }
1178  for (unsigned int i = 0; i < nDim; ++i) {
1179  unsigned int iDim = theDimension[i];
1180  if (aPrec(iDim) > 0.) {
1181  GblData aData(nLabel, InternalKink, aMeas(iDim),
1182  aPrec(iDim), iTraj,
1183  itPoint - thePoints[iTraj].begin());
1184  aData.addDerivatives(iDim, labDer, matTDer, numLocals,
1185  transDer);
1186  theData.emplace_back(aData);
1187  nData++;
1188  }
1189  }
1190  }
1191  scatDataIndex[nLabel] = nData;
1192  }
1193  scatDataIndex[thePoints[iTraj].back().getLabel()] = nData;
1194  }
1195  }
1196 
1197  // external seed
1198  if (externalPoint > 0) {
1199  std::pair<std::vector<unsigned int>, MatrixXd> indexAndJacobian =
1201  std::vector<unsigned int> externalSeedIndex = indexAndJacobian.first;
1202  std::vector<double> externalSeedDerivatives(externalSeedIndex.size());
1203  SelfAdjointEigenSolver<MatrixXd> externalSeedEigen(externalSeed);
1204  VectorXd valEigen = externalSeedEigen.eigenvalues();
1205  MatrixXd vecEigen = externalSeedEigen.eigenvectors();
1206  vecEigen = vecEigen.transpose() * indexAndJacobian.second;
1207  for (int i = 0; i < externalSeed.rows(); ++i) {
1208  if (valEigen(i) > 0.) {
1209  for (int j = 0; j < externalSeed.cols(); ++j) {
1210  externalSeedDerivatives[j] = vecEigen(i, j);
1211  }
1212  GblData aData(externalPoint, ExternalSeed, 0., valEigen(i));
1213  aData.addDerivatives(externalSeedIndex,
1214  externalSeedDerivatives);
1215  theData.emplace_back(std::move(aData));
1216  nData++;
1217  }
1218  }
1219  }
1220  measDataIndex[numAllPoints + 1] = nData;
1221  // external measurements
1222  unsigned int nExt = externalMeasurements.rows();
1223  if (nExt > 0) {
1224  std::vector<unsigned int> index(numCurvature);
1225  std::vector<double> derivatives(numCurvature);
1226  for (unsigned int iExt = 0; iExt < nExt; ++iExt) {
1227  for (unsigned int iCol = 0; iCol < numCurvature; ++iCol) {
1228  index[iCol] = iCol + 1;
1229  derivatives[iCol] = externalDerivatives(iExt, iCol);
1230  }
1232  externalPrecisions(iExt));
1233  aData.addDerivatives(index, derivatives);
1234  theData.emplace_back(std::move(aData));
1235  nData++;
1236  }
1237  }
1238  measDataIndex[numAllPoints + 2] = nData;
1239  }
1240 
1243  std::vector<GblData>::iterator itData;
1244  for (itData = theData.begin(); itData < theData.end(); ++itData) {
1245  itData->setPrediction(theVector);
1246  }
1247  }
1248 
1250 
1253  double GblTrajectory::downWeight(unsigned int aMethod) {
1254  double aLoss = 0.;
1255  std::vector<GblData>::iterator itData;
1256  for (itData = theData.begin(); itData < theData.end(); ++itData) {
1257  aLoss += (1. - itData->setDownWeighting(aMethod));
1258  }
1259  return aLoss;
1260  }
1261 
1263 
1275  unsigned int GblTrajectory::fit(double &Chi2, int &Ndf, double &lostWeight,
1276  const std::string& optionList, unsigned int aLabel) {
1277  const double normChi2[4] = { 1.0, 0.8737, 0.9326, 0.8228 };
1278  const std::string methodList = "TtHhCc";
1279 
1280  Chi2 = 0.;
1281  Ndf = -1;
1282  lostWeight = 0.;
1283  if (not constructOK)
1284  return 10;
1285 
1286  unsigned int aMethod = 0;
1287  skippedMeasLabel = aLabel;
1288 
1290  lostWeight = 0.;
1291  unsigned int ierr = 0;
1292  try {
1293 
1295  predict();
1296 
1297  for (unsigned int i = 0; i < optionList.size(); ++i) // down weighting iterations
1298  {
1299  size_t aPosition = methodList.find(optionList[i]);
1300  if (aPosition != std::string::npos) {
1301  aMethod = aPosition / 2 + 1;
1302  lostWeight = downWeight(aMethod);
1305  predict();
1306  }
1307  }
1308  Ndf = -numParameters;
1309  Chi2 = 0.;
1310  for (unsigned int i = 0; i < theData.size(); ++i) {
1311  // skipped (internal) measurement ?
1312  if (theData[i].getLabel() == skippedMeasLabel
1313  && theData[i].getType() == InternalMeasurement)
1314  continue;
1315  Chi2 += theData[i].getChi2();
1316  Ndf++;
1317  }
1318  Chi2 /= normChi2[aMethod];
1319  fitOK = true;
1320 
1321  } catch (int e) {
1322  std::cout << " GblTrajectory::fit exception " << e << std::endl;
1323  ierr = e;
1324  }
1325  return ierr;
1326  }
1327 
1330  double aValue;
1331  double aErr;
1332  unsigned int aTraj;
1333  unsigned int aPoint;
1334  unsigned int aRow;
1335  unsigned int numLocal;
1336  unsigned int* labLocal;
1337  double* derLocal;
1338  std::vector<int> labGlobal;
1339  std::vector<double> derGlobal;
1340 
1341  if (not constructOK)
1342  return;
1343 
1344  // data: measurements, kinks and external seed
1345  labGlobal.reserve(maxNumGlobals);
1346  derGlobal.reserve(maxNumGlobals);
1347  std::vector<GblData>::iterator itData;
1348  for (itData = theData.begin(); itData != theData.end(); ++itData) {
1349  itData->getAllData(aValue, aErr, numLocal, labLocal, derLocal, aTraj,
1350  aPoint, aRow);
1351  if (itData->getType() == InternalMeasurement)
1352  thePoints[aTraj][aPoint].getGlobalLabelsAndDerivatives(aRow,
1353  labGlobal, derGlobal);
1354  else
1355  labGlobal.resize(0);
1356  aMille.addData(aValue, aErr, numLocal, labLocal, derLocal, labGlobal,
1357  derGlobal);
1358  }
1359  aMille.writeRecord();
1360  }
1361 
1363 
1367  if (numInnerTrans) {
1368  std::cout << "Composed GblTrajectory, " << numInnerTrans
1369  << " subtrajectories" << std::endl;
1370  } else {
1371  std::cout << "Simple GblTrajectory" << std::endl;
1372  }
1373  if (theDimension.size() < 2) {
1374  std::cout << " 2D-trajectory" << std::endl;
1375  }
1376  std::cout << " Number of GblPoints : " << numAllPoints
1377  << std::endl;
1378  std::cout << " Number of points with offsets: " << numOffsets << std::endl;
1379  std::cout << " Number of fit parameters : " << numParameters
1380  << std::endl;
1381  std::cout << " Number of measurements : " << numMeasurements
1382  << std::endl;
1383  if (externalMeasurements.rows()) {
1384  std::cout << " Number of ext. measurements : "
1385  << externalMeasurements.rows() << std::endl;
1386  }
1387  if (externalPoint) {
1388  std::cout << " Label of point with ext. seed: " << externalPoint
1389  << std::endl;
1390  }
1391  if (constructOK) {
1392  std::cout << " Constructed OK " << std::endl;
1393  }
1394  if (fitOK) {
1395  std::cout << " Fitted OK " << std::endl;
1396  }
1397  if (level > 0) {
1398  IOFormat CleanFmt(4, 0, ", ", "\n", "[", "]");
1399  if (numInnerTrans) {
1400  std::cout << " Inner transformations" << std::endl;
1401  for (unsigned int i = 0; i < numInnerTrans; ++i) {
1402  std::cout << innerTransformations[i].format(CleanFmt)
1403  << std::endl;
1404  }
1405  }
1406  if (externalMeasurements.rows()) {
1407  std::cout << " External measurements" << std::endl;
1408  std::cout << " Measurements:" << std::endl;
1409  std::cout << externalMeasurements.format(CleanFmt) << std::endl;
1410  std::cout << " Precisions:" << std::endl;
1411  std::cout << externalPrecisions.format(CleanFmt) << std::endl;
1412  std::cout << " Derivatives:" << std::endl;
1413  std::cout << externalDerivatives.format(CleanFmt) << std::endl;
1414  }
1415  if (externalPoint) {
1416  std::cout << " External seed:" << std::endl;
1417  std::cout << externalSeed.format(CleanFmt) << std::endl;
1418  }
1419  if (fitOK) {
1420  std::cout << " Fit results" << std::endl;
1421  std::cout << " Parameters:" << std::endl;
1422  theVector.print();
1423  std::cout << " Covariance matrix (bordered band part):"
1424  << std::endl;
1426  }
1427  }
1428  }
1429 
1431 
1434  void GblTrajectory::printPoints(unsigned int level) {
1435  std::cout << "GblPoints " << std::endl;
1436  for (unsigned int iTraj = 0; iTraj < numTrajectories; ++iTraj) {
1437  std::vector<GblPoint>::iterator itPoint;
1438  for (itPoint = thePoints[iTraj].begin();
1439  itPoint < thePoints[iTraj].end(); ++itPoint) {
1440  itPoint->printPoint(level);
1441  }
1442  }
1443  }
1444 
1447  std::cout << "GblData blocks " << std::endl;
1448  std::vector<GblData>::iterator itData;
1449  for (itData = theData.begin(); itData < theData.end(); ++itData) {
1450  itData->printData();
1451  }
1452  }
1453 
1454 }
void printMatrix() const
Print bordered band matrix.
size
Write out results.
std::vector< unsigned int > scatDataIndex
mapping points to data blocks from scatterers
void addBlockMatrix(double aWeight, const std::vector< unsigned int > *anIndex, const std::vector< double > *aVector)
Add symmetric block matrix.
void construct()
Construct trajectory from list of points.
void addNextJacobian(const Matrix5d &aJac)
Define jacobian to next scatterer (by GBLTrajectory constructor)
Definition: GblPoint.cc:487
Data (block) for independent scalar measurement.
Definition: GblData.h:55
void buildLinearEquationSystem()
Build linear equation system from data (blocks).
void printTrajectory(unsigned int level=0)
Print GblTrajectory.
unsigned int numParameters
Number of fit parameters.
void writeRecord()
Write record to file.
Definition: MilleBinary.cc:113
bool constructOK
Trajectory has been successfully constructed (ready for fit/output)
std::vector< unsigned int > measDataIndex
mapping points to data blocks from measurements
Eigen::MatrixXd getBlockMatrix(const std::vector< unsigned int > anIndex) const
Retrieve symmetric block matrix.
bool fitOK
Trajectory has been successfully fitted (results are valid)
std::vector< unsigned int > theDimension
List of active dimensions (0=u1, 1=u2) in fit.
void addDerivatives(unsigned int iRow, const std::array< unsigned int, 5 > &labDer, const Matrix5d &matDer, unsigned int iOff, const Eigen::MatrixBase< LocalDerivative > &derLocal, unsigned int nLocal, const Eigen::MatrixBase< TrafoDerivative > &derTrans)
Add derivatives from measurement.
Definition: GblData.h:126
unsigned int getLabels(std::vector< unsigned int > &aLabelList)
Get (list of) labels of points on (simple) valid trajectory.
void print() const
Print vector.
Definition: VMatrix.cc:298
void predict()
Calculate predictions for all points.
int getOffset() const
Retrieve offset for point.
Definition: GblPoint.cc:459
unsigned int numOffsets
Number of (points with) offsets on trajectory.
Eigen::VectorXd externalMeasurements
CLHEP::HepMatrix Matrix
Definition: matutil.h:65
bool isValid() const
Retrieve validity of trajectory.
unsigned int getResults(int aSignedLabel, Eigen::VectorXd &localPar, Eigen::MatrixXd &localCov) const
Get fit results at point.
U second(std::pair< T, U > const &p)
unsigned int fit(double &Chi2, int &Ndf, double &lostWeight, const std::string &optionList="", unsigned int aLabel=0)
Perform fit of (valid) trajectory.
unsigned int externalPoint
Label of external point (or 0)
unsigned int numInnerTrans
Number of inner transformations to external parameters.
void calcJacobians()
Calculate Jacobians to previous/next scatterer from point to point ones.
unsigned int numMeasurements
Total number of measurements.
void getFitToKinkJacobian(std::array< unsigned int, 7 > &anIndex, Matrix27d &aJacobian, const GblPoint &aPoint) const
Get jacobian for transformation from (trajectory) fit to kink parameters at point.
T sqrt(T t)
Definition: SSEVec.h:18
void resize(const unsigned int nRows)
Resize vector.
Definition: VMatrix.cc:262
Eigen::MatrixXd externalSeed
Precision (inverse covariance matrix) of external seed.
void defineOffsets()
Define offsets from list of points.
void addData(double aMeas, double aErr, unsigned int numLocal, unsigned int *indLocal, double *derLocal, const std::vector< int > &labGlobal, const std::vector< double > &derGlobal)
Add data block to (end of) record.
Definition: MilleBinary.cc:72
void prepare()
Prepare fit for simple or composed trajectory.
unsigned int getScatResults(unsigned int aLabel, unsigned int &numRes, Eigen::VectorXd &aResiduals, Eigen::VectorXd &aMeasErrors, Eigen::VectorXd &aResErrors, Eigen::VectorXd &aDownWeights)
Get (kink) residuals from fit at point for scatterer.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void printData()
Print GblData blocks for trajectory.
Namespace for the general broken lines package.
double downWeight(unsigned int aMethod)
Down-weight all points.
#define end
Definition: vmac.h:37
std::vector< unsigned int > numPoints
Number of points on (sub)trajectory.
int k[5][pyjets_maxn]
std::pair< std::vector< unsigned int >, Eigen::MatrixXd > getJacobian(int aSignedLabel) const
Get jacobian for transformation from fit to track parameters at point.
void getFitToLocalJacobian(std::array< unsigned int, 5 > &anIndex, Matrix5d &aJacobian, const GblPoint &aPoint, unsigned int measDim, unsigned int nJacobian=1) const
Get (part of) jacobian for transformation from (trajectory) fit to track parameters at point...
void milleOut(MilleBinary &aMille)
Write valid trajectory to Millepede-II binary file.
std::vector< GblData > theData
List of data blocks.
Eigen::Matrix< double, 5, 1 > Vector5d
Definition: GblPoint.h:47
Eigen::Matrix< double, 2, 7 > Matrix27d
Definition: GblData.h:44
unsigned int maxNumGlobals
Max. number of global labels/derivatives per point.
void resize(unsigned int nSize, unsigned int nBorder=1, unsigned int nBand=5)
Resize bordered band matrix.
std::vector< Eigen::MatrixXd > innerTransformations
Transformations at innermost points of.
Eigen::MatrixXd externalDerivatives
void getDerivatives(int aDirection, Eigen::Matrix2d &matW, Eigen::Matrix2d &matWJ, Eigen::Vector2d &vecWd) const
Retrieve derivatives of local track model.
Definition: GblPoint.cc:501
std::vector< std::vector< GblPoint > > thePoints
(list of) List of points on trajectory
BorderedBandMatrix theMatrix
(Bordered band) matrix of linear equation system
VVector theVector
Vector of linear equation system.
unsigned int numTrajectories
Number of trajectories (in composed trajectory)
unsigned int getMeasResults(unsigned int aLabel, unsigned int &numRes, Eigen::VectorXd &aResiduals, Eigen::VectorXd &aMeasErrors, Eigen::VectorXd &aResErrors, Eigen::VectorXd &aDownWeights)
Get residuals from fit at point for measurement.
#define begin
Definition: vmac.h:30
GblTrajectory(const std::vector< GblPoint > &aPointList, bool flagCurv=true, bool flagU1dir=true, bool flagU2dir=true)
Create new (simple) trajectory from list of points.
unsigned int numLocals
Total number of (additional) local parameters.
virtual ~GblTrajectory()
Millepede-II (binary) record.
Definition: MilleBinary.h:68
unsigned int getNumPoints() const
Retrieve number of point from trajectory.
unsigned int numAllPoints
Number of all points on trajectory.
void getResAndErr(unsigned int aData, bool used, double &aResidual, double &aMeadsError, double &aResError, double &aDownWeight)
Get residual and errors from data block.
void printPoints(unsigned int level=0)
Print GblPoints on trajectory.
void solveAndInvertBorderedBand(const VVector &aRightHandSide, VVector &aSolution)
Solve linear equation system, partially calculate inverse.
Definition: Chi2.h:17
unsigned int numCurvature
Number of curvature parameters (0 or 1) or external parameters.
def move(src, dest)
Definition: eostools.py:510
unsigned int skippedMeasLabel
Label of point with measurements skipped in fit (for unbiased residuals) (or 0)
const Matrix5d & getP2pJacobian() const
Retrieve point-to-(previous)point jacobian.
Definition: GblPoint.cc:464
Eigen::VectorXd externalPrecisions
Point on trajectory.
Definition: GblPoint.h:66
Eigen::Matrix< double, 5, 5 > Matrix5d
Definition: GblData.h:43