20 #include "MeasuredStateOnPlane.h"
21 #include "AbsTrackRep.h"
22 #include "Exception.h"
28 #include "TDecompChol.h"
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();
42 printOut <<
" 3D position: "; pos.Print();
43 printOut <<
" 3D momentum: "; mom.Print();
50 unsigned int dim = cov_.GetNcols();
52 if (resetOffDiagonals) {
53 for (
unsigned int i=0; i<dim; ++i) {
54 for (
unsigned int j=0; j<dim; ++j) {
58 cov_(i,j) *= blowUpFac;
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);
78 if (forwardState.getPlane() != backwardState.getPlane()) {
79 Exception e(
"KalmanFitterInfo::calcAverageState: forwardState and backwardState are not defined in the same plane.", __LINE__,__FILE__);
89 TMatrixDSym fCovInv, bCovInv, smoothed_cov;
90 tools::invertMatrix(forwardState.getCov(), fCovInv);
91 tools::invertMatrix(backwardState.getCov(), bCovInv);
93 tools::invertMatrix(fCovInv + bCovInv, smoothed_cov);
96 retVal.setState(smoothed_cov*(fCovInv*forwardState.getState() + bCovInv*backwardState.getState()));
97 retVal.setCov(smoothed_cov);
143 TDecompChol d1(forwardState.getCov());
144 bool success = d1.Decompose();
145 TDecompChol d2(backwardState.getCov());
146 success &= d2.Decompose();
149 Exception e(
"KalmanFitterInfo::calcAverageState: ill-conditioned covariance matrix.", __LINE__,__FILE__);
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);
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();
166 for (
int i = 0; i < nRows; ++i) {
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];
176 bk[i + nRows] = sum2;
180 A.ResizeTo(nRows, nRows);
183 tools::transposedInvert(A, inv);
185 for (
int i = 0; i < nRows; ++i) {
187 for (
int j = i; j < nRows; ++j) {
188 sum += inv.GetMatrixArray()[j*nRows+i] * b[j];
190 b.GetMatrixArray()[i] = sum;
193 TMatrixDSym(TMatrixDSym::kAtA, inv),
194 forwardState.getPlane(),
195 forwardState.getRep(),
196 forwardState.getAuxInfo());
virtual void getPosMomCov(const MeasuredStateOnPlane &state, TVector3 &pos, TVector3 &mom, TMatrixDSym &cov) const =0
Translates MeasuredStateOnPlane into 3D position, momentum and 6x6 covariance.
Exception class for error handling in GENFIT (provides storage for diagnostic information)
#StateOnPlane with additional covariance matrix.
void blowUpCov(double blowUpFac, bool resetOffDiagonals=true, double maxVal=-1.)
Blow up covariance matrix with blowUpFac. Per default, off diagonals are reset to 0 and the maximum v...
Defines for I/O streams used for error and debug printing.
std::ostream printOut
Default stream for output of Print calls.
MeasuredStateOnPlane calcAverageState(const MeasuredStateOnPlane &forwardState, const MeasuredStateOnPlane &backwardState)
Calculate weighted average between two MeasuredStateOnPlanes.