test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GblPoint.cc
Go to the documentation of this file.
1 /*
2  * GblPoint.cpp
3  *
4  * Created on: Aug 18, 2011
5  * Author: kleinwrt
6  */
7 
9 
11 namespace gbl {
12 
14 
18 GblPoint::GblPoint(const TMatrixD &aJacobian) :
19  theLabel(0), theOffset(0), measDim(0), transFlag(false), measTransformation(), scatFlag(
20  false), localDerivatives(), globalLabels(), globalDerivatives() {
21 
22  for (unsigned int i = 0; i < 5; ++i) {
23  for (unsigned int j = 0; j < 5; ++j) {
24  p2pJacobian(i, j) = aJacobian(i, j);
25  }
26  }
27 }
28 
29 GblPoint::GblPoint(const SMatrix55 &aJacobian) :
30  theLabel(0), theOffset(0), p2pJacobian(aJacobian), measDim(0), transFlag(
31  false), measTransformation(), scatFlag(false), localDerivatives(), globalLabels(), globalDerivatives() {
32 
33 }
34 
36 }
37 
39 
47 void GblPoint::addMeasurement(const TMatrixD &aProjection,
48  const TVectorD &aResiduals, const TVectorD &aPrecision,
49  double minPrecision) {
50  measDim = aResiduals.GetNrows();
51  unsigned int iOff = 5 - measDim;
52  for (unsigned int i = 0; i < measDim; ++i) {
53  measResiduals(iOff + i) = aResiduals[i];
54  measPrecision(iOff + i) = (
55  aPrecision[i] >= minPrecision ? aPrecision[i] : 0.);
56  for (unsigned int j = 0; j < measDim; ++j) {
57  measProjection(iOff + i, iOff + j) = aProjection(i, j);
58  }
59  }
60 }
61 
63 
72 void GblPoint::addMeasurement(const TMatrixD &aProjection,
73  const TVectorD &aResiduals, const TMatrixDSym &aPrecision,
74  double minPrecision) {
75  measDim = aResiduals.GetNrows();
76  TMatrixDSymEigen measEigen(aPrecision);
78  measTransformation = measEigen.GetEigenVectors();
80  transFlag = true;
81  TVectorD transResiduals = measTransformation * aResiduals;
82  TVectorD transPrecision = measEigen.GetEigenValues();
83  TMatrixD transProjection = measTransformation * aProjection;
84  unsigned int iOff = 5 - measDim;
85  for (unsigned int i = 0; i < measDim; ++i) {
86  measResiduals(iOff + i) = transResiduals[i];
87  measPrecision(iOff + i) = (
88  transPrecision[i] >= minPrecision ? transPrecision[i] : 0.);
89  for (unsigned int j = 0; j < measDim; ++j) {
90  measProjection(iOff + i, iOff + j) = transProjection(i, j);
91  }
92  }
93 }
94 
96 
103 void GblPoint::addMeasurement(const TVectorD &aResiduals,
104  const TVectorD &aPrecision, double minPrecision) {
105  measDim = aResiduals.GetNrows();
106  unsigned int iOff = 5 - measDim;
107  for (unsigned int i = 0; i < measDim; ++i) {
108  measResiduals(iOff + i) = aResiduals[i];
109  measPrecision(iOff + i) = (
110  aPrecision[i] >= minPrecision ? aPrecision[i] : 0.);
111  }
112  measProjection = ROOT::Math::SMatrixIdentity();
113 }
114 
116 
124 void GblPoint::addMeasurement(const TVectorD &aResiduals,
125  const TMatrixDSym &aPrecision, double minPrecision) {
126  measDim = aResiduals.GetNrows();
127  TMatrixDSymEigen measEigen(aPrecision);
129  measTransformation = measEigen.GetEigenVectors();
130  measTransformation.T();
131  transFlag = true;
132  TVectorD transResiduals = measTransformation * aResiduals;
133  TVectorD transPrecision = measEigen.GetEigenValues();
134  unsigned int iOff = 5 - measDim;
135  for (unsigned int i = 0; i < measDim; ++i) {
136  measResiduals(iOff + i) = transResiduals[i];
137  measPrecision(iOff + i) = (
138  transPrecision[i] >= minPrecision ? transPrecision[i] : 0.);
139  for (unsigned int j = 0; j < measDim; ++j) {
140  measProjection(iOff + i, iOff + j) = measTransformation(i, j);
141  }
142  }
143 }
144 
146 
150 unsigned int GblPoint::hasMeasurement() const {
151  return measDim;
152 }
153 
155 
160 void GblPoint::getMeasurement(SMatrix55 &aProjection, SVector5 &aResiduals,
161  SVector5 &aPrecision) const {
162  aProjection = measProjection;
163  aResiduals = measResiduals;
164  aPrecision = measPrecision;
165 }
166 
168 
171 void GblPoint::getMeasTransformation(TMatrixD &aTransformation) const {
172  aTransformation.ResizeTo(measDim, measDim);
173  if (transFlag) {
174  aTransformation = measTransformation;
175  } else {
176  aTransformation.UnitMatrix();
177  }
178 }
179 
181 
188 void GblPoint::addScatterer(const TVectorD &aResiduals,
189  const TVectorD &aPrecision) {
190  scatFlag = true;
191  scatResiduals(0) = aResiduals[0];
192  scatResiduals(1) = aResiduals[1];
193  scatPrecision(0) = aPrecision[0];
194  scatPrecision(1) = aPrecision[1];
195  scatTransformation = ROOT::Math::SMatrixIdentity();
196 }
197 
199 
214 void GblPoint::addScatterer(const TVectorD &aResiduals,
215  const TMatrixDSym &aPrecision) {
216  scatFlag = true;
217  TMatrixDSymEigen scatEigen(aPrecision);
218  TMatrixD aTransformation = scatEigen.GetEigenVectors();
219  aTransformation.T();
220  TVectorD transResiduals = aTransformation * aResiduals;
221  TVectorD transPrecision = scatEigen.GetEigenValues();
222  for (unsigned int i = 0; i < 2; ++i) {
223  scatResiduals(i) = transResiduals[i];
224  scatPrecision(i) = transPrecision[i];
225  for (unsigned int j = 0; j < 2; ++j) {
226  scatTransformation(i, j) = aTransformation(i, j);
227  }
228  }
229 }
230 
233  return scatFlag;
234 }
235 
237 
242 void GblPoint::getScatterer(SMatrix22 &aTransformation, SVector2 &aResiduals,
243  SVector2 &aPrecision) const {
244  aTransformation = scatTransformation;
245  aResiduals = scatResiduals;
246  aPrecision = scatPrecision;
247 }
248 
250 
253 void GblPoint::getScatTransformation(TMatrixD &aTransformation) const {
254  aTransformation.ResizeTo(2, 2);
255  if (scatFlag) {
256  for (unsigned int i = 0; i < 2; ++i) {
257  for (unsigned int j = 0; j < 2; ++j) {
258  aTransformation(i, j) = scatTransformation(i, j);
259  }
260  }
261  } else {
262  aTransformation.UnitMatrix();
263  }
264 }
265 
267 
271 void GblPoint::addLocals(const TMatrixD &aDerivatives) {
272  if (measDim) {
273  localDerivatives.ResizeTo(aDerivatives);
274  if (transFlag) {
275  localDerivatives = measTransformation * aDerivatives;
276  } else {
277  localDerivatives = aDerivatives;
278  }
279  }
280 }
281 
283 unsigned int GblPoint::getNumLocals() const {
284  return localDerivatives.GetNcols();
285 }
286 
288 const TMatrixD& GblPoint::getLocalDerivatives() const {
289  return localDerivatives;
290 }
291 
293 
298 void GblPoint::addGlobals(const std::vector<int> &aLabels,
299  const TMatrixD &aDerivatives) {
300  if (measDim) {
301  globalLabels = aLabels;
302  globalDerivatives.ResizeTo(aDerivatives);
303  if (transFlag) {
304  globalDerivatives = measTransformation * aDerivatives;
305  } else {
306  globalDerivatives = aDerivatives;
307  }
308 
309  }
310 }
311 
313 unsigned int GblPoint::getNumGlobals() const {
314  return globalDerivatives.GetNcols();
315 }
316 
318 std::vector<int> GblPoint::getGlobalLabels() const {
319  return globalLabels;
320 }
321 
323 const TMatrixD& GblPoint::getGlobalDerivatives() const {
324  return globalDerivatives;
325 }
326 
328 
331 void GblPoint::setLabel(unsigned int aLabel) {
332  theLabel = aLabel;
333 }
334 
336 unsigned int GblPoint::getLabel() const {
337  return theLabel;
338 }
339 
341 
344 void GblPoint::setOffset(int anOffset) {
345  theOffset = anOffset;
346 }
347 
349 int GblPoint::getOffset() const {
350  return theOffset;
351 }
352 
355  return p2pJacobian;
356 }
357 
359 
363  int ifail = 0;
364 // to optimize: need only two last rows of inverse
365 // prevJacobian = aJac.InverseFast(ifail);
366 // block matrix algebra
367  SMatrix23 CA = aJac.Sub<SMatrix23>(3, 0)
368  * aJac.Sub<SMatrix33>(0, 0).InverseFast(ifail); // C*A^-1
369  SMatrix22 DCAB = aJac.Sub<SMatrix22>(3, 3) - CA * aJac.Sub<SMatrix32>(0, 3); // D - C*A^-1 *B
370  DCAB.InvertFast();
371  prevJacobian.Place_at(DCAB, 3, 3);
372  prevJacobian.Place_at(-DCAB * CA, 3, 0);
373 }
374 
376 
380  nextJacobian = aJac;
381 }
382 
384 
393 void GblPoint::getDerivatives(int aDirection, SMatrix22 &matW, SMatrix22 &matWJ,
394  SVector2 &vecWd) const {
395 
396  SMatrix22 matJ;
397  SVector2 vecd;
398  if (aDirection < 1) {
399  matJ = prevJacobian.Sub<SMatrix22>(3, 3);
400  matW = -prevJacobian.Sub<SMatrix22>(3, 1);
401  vecd = prevJacobian.SubCol<SVector2>(0, 3);
402  } else {
403  matJ = nextJacobian.Sub<SMatrix22>(3, 3);
404  matW = nextJacobian.Sub<SMatrix22>(3, 1);
405  vecd = nextJacobian.SubCol<SVector2>(0, 3);
406  }
407 
408  if (!matW.InvertFast()) {
409  std::cout << " GblPoint::getDerivatives failed to invert matrix: "
410  << matW << std::endl;
411  std::cout
412  << " Possible reason for singular matrix: multiple GblPoints at same arc-length"
413  << std::endl;
414  throw std::overflow_error("Singular matrix inversion exception");
415  }
416  matWJ = matW * matJ;
417  vecWd = matW * vecd;
418 
419 }
420 
422 
425 void GblPoint::printPoint(unsigned int level) const {
426  std::cout << " GblPoint";
427  if (theLabel) {
428  std::cout << ", label " << theLabel;
429  if (theOffset >= 0) {
430  std::cout << ", offset " << theOffset;
431  }
432  }
433  if (measDim) {
434  std::cout << ", " << measDim << " measurements";
435  }
436  if (scatFlag) {
437  std::cout << ", scatterer";
438  }
439  if (transFlag) {
440  std::cout << ", diagonalized";
441  }
442  if (localDerivatives.GetNcols()) {
443  std::cout << ", " << localDerivatives.GetNcols()
444  << " local derivatives";
445  }
446  if (globalDerivatives.GetNcols()) {
447  std::cout << ", " << globalDerivatives.GetNcols()
448  << " global derivatives";
449  }
450  std::cout << std::endl;
451  if (level > 0) {
452  if (measDim) {
453  std::cout << " Measurement" << std::endl;
454  std::cout << " Projection: " << std::endl << measProjection
455  << std::endl;
456  std::cout << " Residuals: " << measResiduals << std::endl;
457  std::cout << " Precision: " << measPrecision << std::endl;
458  }
459  if (scatFlag) {
460  std::cout << " Scatterer" << std::endl;
461  std::cout << " Residuals: " << scatResiduals << std::endl;
462  std::cout << " Precision: " << scatPrecision << std::endl;
463  }
464  if (localDerivatives.GetNcols()) {
465  std::cout << " Local Derivatives:" << std::endl;
466  localDerivatives.Print();
467  }
468  if (globalDerivatives.GetNcols()) {
469  std::cout << " Global Labels:";
470  for (unsigned int i = 0; i < globalLabels.size(); ++i) {
471  std::cout << " " << globalLabels[i];
472  }
473  std::cout << std::endl;
474  std::cout << " Global Derivatives:" << std::endl;
475  globalDerivatives.Print();
476  }
477  std::cout << " Jacobian " << std::endl;
478  std::cout << " Point-to-point " << std::endl << p2pJacobian
479  << std::endl;
480  if (theLabel) {
481  std::cout << " To previous offset " << std::endl << prevJacobian
482  << std::endl;
483  std::cout << " To next offset " << std::endl << nextJacobian
484  << std::endl;
485  }
486  }
487 }
488 
489 }
void getDerivatives(int aDirection, SMatrix22 &matW, SMatrix22 &matWJ, SVector2 &vecWd) const
Retrieve derivatives of local track model.
Definition: GblPoint.cc:393
int i
Definition: DBlmapReader.cc:9
void addMeasurement(const TMatrixD &aProjection, const TVectorD &aResiduals, const TVectorD &aPrecision, double minPrecision=0.)
Add a measurement to a point.
Definition: GblPoint.cc:47
SMatrix55 prevJacobian
Jacobian to previous scatterer (or first measurement)
Definition: GblPoint.h:93
void getScatTransformation(TMatrixD &aTransformation) const
Get scatterer transformation (from diagonalization).
Definition: GblPoint.cc:253
void setLabel(unsigned int aLabel)
Define label of point (by GBLTrajectory constructor)
Definition: GblPoint.cc:331
TMatrixD globalDerivatives
Derivatives of measurement vs additional global (MP-II) parameters.
Definition: GblPoint.h:107
bool transFlag
Transformation exists?
Definition: GblPoint.h:99
bool scatFlag
Scatterer present?
Definition: GblPoint.h:101
void setOffset(int anOffset)
Define offset for point (by GBLTrajectory constructor)
Definition: GblPoint.cc:344
int theOffset
Offset number at point if not negative (else interpolation needed)
Definition: GblPoint.h:91
bool hasScatterer() const
Check for scatterer at a point.
Definition: GblPoint.cc:232
std::vector< int > getGlobalLabels() const
Retrieve global derivatives labels from a point.
Definition: GblPoint.cc:318
unsigned int getLabel() const
Retrieve label of point.
Definition: GblPoint.cc:336
void addLocals(const TMatrixD &aDerivatives)
Add local derivatives to a point.
Definition: GblPoint.cc:271
ROOT::Math::SVector< double, 2 > SVector2
Definition: GblPoint.h:29
unsigned int measDim
Dimension of measurement (1-5), 0 indicates absence of measurement.
Definition: GblPoint.h:95
SVector5 measResiduals
Measurement residuals.
Definition: GblPoint.h:97
int getOffset() const
Retrieve offset for point.
Definition: GblPoint.cc:349
ROOT::Math::SMatrix< double, 5, 5 > SMatrix55
Definition: GblData.h:23
unsigned int hasMeasurement() const
Check for measurement at a point.
Definition: GblPoint.cc:150
const TMatrixD & getGlobalDerivatives() const
Retrieve global derivatives from a point.
Definition: GblPoint.cc:323
void getMeasTransformation(TMatrixD &aTransformation) const
Get measurement transformation (from diagonalization).
Definition: GblPoint.cc:171
ROOT::Math::SVector< double, 5 > SVector5
Definition: GblPoint.h:30
ROOT::Math::SMatrix< double, 3, 2 > SMatrix32
Definition: GblPoint.h:26
SVector2 scatPrecision
Scattering precision (diagonal of inverse covariance matrix)
Definition: GblPoint.h:104
void printPoint(unsigned int level=0) const
Print GblPoint.
Definition: GblPoint.cc:425
ROOT::Math::SMatrix< double, 2 > SMatrix22
Definition: GblPoint.h:22
std::vector< int > globalLabels
Labels of global (MP-II) derivatives.
Definition: GblPoint.h:106
void addScatterer(const TVectorD &aResiduals, const TVectorD &aPrecision)
Add a (thin) scatterer to a point.
Definition: GblPoint.cc:188
int j
Definition: DBlmapReader.cc:9
void addGlobals(const std::vector< int > &aLabels, const TMatrixD &aDerivatives)
Add global derivatives to a point.
Definition: GblPoint.cc:298
unsigned int getNumLocals() const
Retrieve number of local derivatives from a point.
Definition: GblPoint.cc:283
SMatrix55 nextJacobian
Jacobian to next scatterer (or last measurement)
Definition: GblPoint.h:94
SVector5 measPrecision
Measurement precision (diagonal of inverse covariance matrix)
Definition: GblPoint.h:98
GblPoint(const TMatrixD &aJacobian)
Create a point.
Definition: GblPoint.cc:18
TMatrixD measTransformation
Transformation of diagonalization (of meas. precision matrix)
Definition: GblPoint.h:100
void getScatterer(SMatrix22 &aTransformation, SVector2 &aResiduals, SVector2 &aPrecision) const
Retrieve scatterer of a point.
Definition: GblPoint.cc:242
void getMeasurement(SMatrix55 &aProjection, SVector5 &aResiduals, SVector5 &aPrecision) const
Retrieve measurement of a point.
Definition: GblPoint.cc:160
SVector2 scatResiduals
Scattering residuals (initial kinks if iterating)
Definition: GblPoint.h:103
void addNextJacobian(const SMatrix55 &aJac)
Define jacobian to next scatterer (by GBLTrajectory constructor)
Definition: GblPoint.cc:379
TMatrixD localDerivatives
Derivatives of measurement vs additional local (fit) parameters.
Definition: GblPoint.h:105
tuple cout
Definition: gather_cfg.py:145
unsigned int theLabel
Label identifying point.
Definition: GblPoint.h:90
tuple level
Definition: testEve_cfg.py:34
SMatrix55 p2pJacobian
Point-to-point jacobian from previous point.
Definition: GblPoint.h:92
volatile std::atomic< bool > shutdown_flag false
const TMatrixD & getLocalDerivatives() const
Retrieve local derivatives from a point.
Definition: GblPoint.cc:288
SMatrix55 measProjection
Projection from measurement to local system.
Definition: GblPoint.h:96
ROOT::Math::SMatrix< double, 2, 3 > SMatrix23
Definition: GblPoint.h:23
SMatrix22 scatTransformation
Transformation of diagonalization (of scat. precision matrix)
Definition: GblPoint.h:102
const SMatrix55 & getP2pJacobian() const
Retrieve point-to-(previous)point jacobian.
Definition: GblPoint.cc:354
virtual ~GblPoint()
Definition: GblPoint.cc:35
void addPrevJacobian(const SMatrix55 &aJac)
Define jacobian to previous scatterer (by GBLTrajectory constructor)
Definition: GblPoint.cc:362
ROOT::Math::SMatrix< double, 3 > SMatrix33
Definition: GblPoint.h:27
unsigned int getNumGlobals() const
Retrieve number of global derivatives from a point.
Definition: GblPoint.cc:313