Analysis Software
Documentation for sPHENIX simulation software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MeasuredStateOnPlane.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file MeasuredStateOnPlane.cc
1 /* Copyright 2008-2010, Technische Universitaet Muenchen,
2  Authors: Christian Hoeppner & Sebastian Neubert & Johannes Rauch
3 
4  This file is part of GENFIT.
5 
6  GENFIT is free software: you can redistribute it and/or modify
7  it under the terms of the GNU Lesser General Public License as published
8  by the Free Software Foundation, either version 3 of the License, or
9  (at your option) any later version.
10 
11  GENFIT is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU Lesser General Public License for more details.
15 
16  You should have received a copy of the GNU Lesser General Public License
17  along with GENFIT. If not, see <http://www.gnu.org/licenses/>.
18 */
19 
20 #include "MeasuredStateOnPlane.h"
21 #include "AbsTrackRep.h"
22 #include "Exception.h"
23 #include "Tools.h"
24 #include "IO.h"
25 
26 #include <cassert>
27 
28 #include "TDecompChol.h"
29 
30 namespace genfit {
31 
32 void MeasuredStateOnPlane::Print(Option_t*) const {
33  printOut << "genfit::MeasuredStateOnPlane ";
34  printOut << "my address " << this << " my plane's address " << this->sharedPlane_.get() << "; use count: " << sharedPlane_.use_count() << std::endl;
35  printOut << " state vector: "; state_.Print();
36  printOut << " covariance matrix: "; cov_.Print();
37  if (sharedPlane_ != nullptr) {
38  printOut << " defined in plane "; sharedPlane_->Print();
39  TVector3 pos, mom;
40  TMatrixDSym cov(6,6);
41  getRep()->getPosMomCov(*this, pos, mom, cov);
42  printOut << " 3D position: "; pos.Print();
43  printOut << " 3D momentum: "; mom.Print();
44  //printOut << " 6D covariance: "; cov.Print();
45  }
46 }
47 
48 void MeasuredStateOnPlane::blowUpCov(double blowUpFac, bool resetOffDiagonals, double maxVal) {
49 
50  unsigned int dim = cov_.GetNcols();
51 
52  if (resetOffDiagonals) {
53  for (unsigned int i=0; i<dim; ++i) {
54  for (unsigned int j=0; j<dim; ++j) {
55  if (i != j)
56  cov_(i,j) = 0; // reset off-diagonals
57  else
58  cov_(i,j) *= blowUpFac; // blow up diagonals
59  }
60  }
61  }
62  else
63  cov_ *= blowUpFac;
64 
65  // limit
66  if (maxVal > 0.)
67  for (unsigned int i=0; i<dim; ++i) {
68  for (unsigned int j=0; j<dim; ++j) {
69  cov_(i,j) = std::min(cov_(i,j), maxVal);
70  }
71  }
72 
73 }
74 
75 
77  // check if both states are defined in the same plane
78  if (forwardState.getPlane() != backwardState.getPlane()) {
79  Exception e("KalmanFitterInfo::calcAverageState: forwardState and backwardState are not defined in the same plane.", __LINE__,__FILE__);
80  throw e;
81  }
82 
83  // This code is called a lot, so some effort has gone into reducing
84  // the number of temporary objects being constructed.
85 
86 #if 0
87  // For ease of understanding, here's a very explicit implementation
88  // that uses the textbook algorithm:
89  TMatrixDSym fCovInv, bCovInv, smoothed_cov;
90  tools::invertMatrix(forwardState.getCov(), fCovInv);
91  tools::invertMatrix(backwardState.getCov(), bCovInv);
92 
93  tools::invertMatrix(fCovInv + bCovInv, smoothed_cov); // one temporary TMatrixDSym
94 
95  MeasuredStateOnPlane retVal(forwardState);
96  retVal.setState(smoothed_cov*(fCovInv*forwardState.getState() + bCovInv*backwardState.getState())); // four temporary TVectorD's
97  retVal.setCov(smoothed_cov);
98  return retVal;
99 #endif
100 
101  // This is a numerically stable implementation of the averaging
102  // process. We write S1, S2 for the upper diagonal square roots
103  // (Cholesky decompositions) of the covariance matrices, such that
104  // C1 = S1' S1 (transposition indicated by ').
105  //
106  // Then we can write
107  // (C1^-1 + C2^-1)^-1 = (S1inv' S1inv + S2inv' S2inv)^-1
108  // = ( (S1inv', S2inv') . ( S1inv ) )^-1
109  // ( ( S2inv ) )
110  // = ( (R', 0) . Q . Q' . ( R ) )^-1
111  // ( ( 0 ) )
112  // where Q is an orthogonal matrix chosen such that R is upper diagonal.
113  // Since Q'.Q = 1, this reduces to
114  // = ( R'.R )^-1
115  // = R^-1 . (R')^-1.
116  // This gives the covariance matrix of the average and its Cholesky
117  // decomposition.
118  //
119  // In order to get the averaged state (writing x1 and x2 for the
120  // states) we proceed from
121  // C1^-1.x1 + C2^-1.x2 = (S1inv', S2inv') . ( S1inv.x1 )
122  // ( S2inv.x2 )
123  // which by the above can be written as
124  // = (R', 0) . Q . ( S1inv.x1 )
125  // ( S2inv.x2 )
126  // with the same (R, Q) as above.
127  //
128  // The average is then after obvious simplification
129  // average = R^-1 . Q . (S1inv.x1)
130  // (S2inv.x2)
131  //
132  // This is what's implemented below, where we make use of the
133  // tridiagonal shapes of the various matrices when multiplying or
134  // inverting.
135  //
136  // This turns out not only more numerically stable, but because the
137  // matrix operations are simpler, it is also faster than the
138  // straightoforward implementation.
139  //
140  // This is an application of the technique of Golub, G.,
141  // Num. Math. 7, 206 (1965) to the least-squares problem underlying
142  // averaging.
143  TDecompChol d1(forwardState.getCov());
144  bool success = d1.Decompose();
145  TDecompChol d2(backwardState.getCov());
146  success &= d2.Decompose();
147 
148  if (!success) {
149  Exception e("KalmanFitterInfo::calcAverageState: ill-conditioned covariance matrix.", __LINE__,__FILE__);
150  throw e;
151  }
152 
153  int nRows = d1.GetU().GetNrows();
154  assert(nRows == d2.GetU().GetNrows());
155  TMatrixD S1inv, S2inv;
156  tools::transposedInvert(d1.GetU(), S1inv);
157  tools::transposedInvert(d2.GetU(), S2inv);
158 
159  TMatrixD A(2*nRows, nRows);
160  TVectorD b(2 * nRows);
161  double *const bk = b.GetMatrixArray();
162  double *const Ak = A.GetMatrixArray();
163  const double* S1invk = S1inv.GetMatrixArray();
164  const double* S2invk = S2inv.GetMatrixArray();
165  // S1inv and S2inv are lower triangular.
166  for (int i = 0; i < nRows; ++i) {
167  double sum1 = 0;
168  double sum2 = 0;
169  for (int j = 0; j <= i; ++j) {
170  Ak[i*nRows + j] = S1invk[i*nRows + j];
171  Ak[(i + nRows)*nRows + j] = S2invk[i*nRows + j];
172  sum1 += S1invk[i*nRows + j]*forwardState.getState().GetMatrixArray()[j];
173  sum2 += S2invk[i*nRows + j]*backwardState.getState().GetMatrixArray()[j];
174  }
175  bk[i] = sum1;
176  bk[i + nRows] = sum2;
177  }
178 
179  tools::QR(A, b);
180  A.ResizeTo(nRows, nRows);
181 
182  TMatrixD inv;
183  tools::transposedInvert(A, inv);
184  b.ResizeTo(nRows);
185  for (int i = 0; i < nRows; ++i) {
186  double sum = 0;
187  for (int j = i; j < nRows; ++j) {
188  sum += inv.GetMatrixArray()[j*nRows+i] * b[j];
189  }
190  b.GetMatrixArray()[i] = sum;
191  }
192  return MeasuredStateOnPlane(b,
193  TMatrixDSym(TMatrixDSym::kAtA, inv),
194  forwardState.getPlane(),
195  forwardState.getRep(),
196  forwardState.getAuxInfo());
197 }
198 
199 
200 } /* End of namespace genfit */