00001 #pragma ident "$Id: SRIFilter.cpp 3140 2012-06-18 15:03:02Z susancummins $" 00002 00013 //============================================================================ 00014 // 00015 // This file is part of GPSTk, the GPS Toolkit. 00016 // 00017 // The GPSTk is free software; you can redistribute it and/or modify 00018 // it under the terms of the GNU Lesser General Public License as published 00019 // by the Free Software Foundation; either version 2.1 of the License, or 00020 // any later version. 00021 // 00022 // The GPSTk is distributed in the hope that it will be useful, 00023 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00024 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00025 // GNU Lesser General Public License for more details. 00026 // 00027 // You should have received a copy of the GNU Lesser General Public 00028 // License along with GPSTk; if not, write to the Free Software Foundation, 00029 // Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110, USA 00030 // 00031 // Copyright 2004, The University of Texas at Austin 00032 // 00033 //============================================================================ 00034 00035 //============================================================================ 00036 // 00037 //This software developed by Applied Research Laboratories at the University of 00038 //Texas at Austin, under contract to an agency or agencies within the U.S. 00039 //Department of Defense. The U.S. Government retains all rights to use, 00040 //duplicate, distribute, disclose, or release this software. 00041 // 00042 //Pursuant to DoD Directive 523024 00043 // 00044 // DISTRIBUTION STATEMENT A: This software has been approved for public 00045 // release, distribution is unlimited. 00046 // 00047 //============================================================================= 00048 00049 //------------------------------------------------------------------------------------ 00050 #include "SRIFilter.hpp" 00051 #include "RobustStats.hpp" 00052 #include "StringUtils.hpp" 00053 00054 //------------------------------------------------------------------------------------ 00055 // TD 00056 using namespace std; 00057 00058 namespace gpstk 00059 { 00060 00061 using namespace StringUtils; 00062 00063 //------------------------------------------------------------------------------------ 00064 // empty constructor 00065 SRIFilter::SRIFilter(void) 00066 throw() 00067 { 00068 defaults(); 00069 } 00070 00071 //------------------------------------------------------------------------------------ 00072 // constructor given the dimension N. 00073 SRIFilter::SRIFilter(const unsigned int N) 00074 throw() 00075 { 00076 defaults(); 00077 R = Matrix<double>(N,N,0.0); 00078 Z = Vector<double>(N,0.0); 00079 names = Namelist(N); 00080 } 00081 00082 //------------------------------------------------------------------------------------ 00083 // constructor given a Namelist, its dimension determines the SRI dimension. 00084 SRIFilter::SRIFilter(const Namelist& NL) 00085 throw() 00086 { 00087 defaults(); 00088 if(NL.size() <= 0) return; 00089 R = Matrix<double>(NL.size(),NL.size(),0.0); 00090 Z = Vector<double>(NL.size(),0.0); 00091 names = NL; 00092 } 00093 00094 //------------------------------------------------------------------------------------ 00095 // explicit constructor - throw if the dimensions are inconsistent. 00096 SRIFilter::SRIFilter(const Matrix<double>& Rin, 00097 const Vector<double>& Zin, 00098 const Namelist& NLin) 00099 throw(MatrixException) 00100 { 00101 defaults(); 00102 if(Rin.rows() != Rin.cols() || 00103 Rin.rows() != Zin.size() || 00104 Rin.rows() != NLin.size()) { 00105 MatrixException me("Invalid input dimensions: R is " 00106 + asString<int>(Rin.rows()) + "x" 00107 + asString<int>(Rin.cols()) + ", Z has length " 00108 + asString<int>(Zin.size()) + ", and NL has length " 00109 + asString<int>(NLin.size()) 00110 ); 00111 GPSTK_THROW(me); 00112 } 00113 R = Rin; 00114 Z = Zin; 00115 names = NLin; 00116 } 00117 00118 //------------------------------------------------------------------------------------ 00119 // operator= 00120 SRIFilter& SRIFilter::operator=(const SRIFilter& right) 00121 throw() 00122 { 00123 R = right.R; 00124 Z = right.Z; 00125 names = right.names; 00126 //valid = right.valid; 00127 return *this; 00128 } 00129 00130 //------------------------------------------------------------------------------------ 00131 // SRIF (Kalman) measurement update, or least squares update 00132 // Returns (whitened) residuals in D 00133 void SRIFilter::measurementUpdate(const Matrix<double>& H, 00134 Vector<double>& D, 00135 const Matrix<double>& CM) 00136 throw(MatrixException,VectorException) 00137 { 00138 if(H.cols() != R.cols() || H.rows() != D.size() || 00139 (&CM != &SRINullMatrix && (CM.rows() != D.size() || CM.cols() != D.size())) ) { 00140 00141 string msg("\nInvalid input dimensions:\n SRI is "); 00142 msg += asString<int>(R.rows()) + "x" 00143 + asString<int>(R.cols()) + ",\n Partials is " 00144 + asString<int>(H.rows()) + "x" 00145 + asString<int>(H.cols()) + ",\n Data has length " 00146 + asString<int>(D.size()); 00147 if(&CM != &SRINullMatrix) msg += ",\n and Cov is " 00148 + asString<int>(CM.rows()) + "x" 00149 + asString<int>(CM.cols()); 00150 00151 MatrixException me(msg); 00152 GPSTK_THROW(me); 00153 } 00154 try { 00155 Matrix<double> P(H); 00156 Cholesky<double> Ch; 00157 00158 // whiten partials and data 00159 if(&CM != &SRINullMatrix) { 00160 Matrix<double> L; 00161 Ch(CM); 00162 L = inverse(Ch.L); 00163 P = L * P; 00164 D = L * D; 00165 } 00166 00167 // update *this with the whitened information 00168 SrifMU(R, Z, P, D); 00169 00170 // un-whiten residuals 00171 if(&CM != &SRINullMatrix) { 00172 D = Ch.L * D; 00173 } 00174 } 00175 catch(MatrixException& me) { GPSTK_RETHROW(me); } 00176 catch(VectorException& ve) { GPSTK_RETHROW(ve); } 00177 } 00178 00179 //------------------------------------------------------------------------------------ 00180 // SRIF (Kalman) time update see SrifTU for doc. 00181 void SRIFilter::timeUpdate(Matrix<double>& PhiInv, 00182 Matrix<double>& Rw, 00183 Matrix<double>& G, 00184 Vector<double>& Zw, 00185 Matrix<double>& Rwx) 00186 throw(MatrixException) 00187 { 00188 try { SrifTU(R, Z, PhiInv, Rw, G, Zw, Rwx); } 00189 catch(MatrixException& me) { GPSTK_RETHROW(me); } 00190 } 00191 00192 //------------------------------------------------------------------------------------ 00193 // SRIF (Kalman) smoother update see SrifSU for doc. 00194 void SRIFilter::smootherUpdate(Matrix<double>& Phi, 00195 Matrix<double>& Rw, 00196 Matrix<double>& G, 00197 Vector<double>& Zw, 00198 Matrix<double>& Rwx) 00199 throw(MatrixException) 00200 { 00201 try { SrifSU(R, Z, Phi, Rw, G, Zw, Rwx); } 00202 catch(MatrixException& me) { GPSTK_RETHROW(me); } 00203 } 00204 00205 //------------------------------------------------------------------------------------ 00206 void SRIFilter::DMsmootherUpdate(Matrix<double>& P, 00207 Vector<double>& X, 00208 Matrix<double>& Phinv, 00209 Matrix<double>& Rw, 00210 Matrix<double>& G, 00211 Vector<double>& Zw, 00212 Matrix<double>& Rwx) 00213 throw(MatrixException) 00214 { 00215 try { SrifSU_DM(P, X, Phinv, Rw, G, Zw, Rwx); } 00216 catch(MatrixException& me) { GPSTK_RETHROW(me); } 00217 } 00218 00219 //------------------------------------------------------------------------------------ 00220 // reset the computation, i.e. remove all stored information 00221 void SRIFilter::zeroAll(void) 00222 { 00223 SRI::zeroAll(); 00224 } 00225 00226 //------------------------------------------------------------------------------------ 00227 // reset the computation, i.e. remove all stored information, and 00228 // optionally change the dimension. If N is not input, the 00229 // dimension is not changed. 00230 // @param N new SRIFilter dimension (optional). 00231 void SRIFilter::Reset(const int N) 00232 { 00233 if(N > 0 && N != R.rows()) { 00234 R.resize(N,N,0.0); 00235 Z.resize(N,0.0); 00236 } 00237 else 00238 SRI::zeroAll(N); 00239 } 00240 00241 //------------------------------------------------------------------------------------ 00242 // private beyond this 00243 //------------------------------------------------------------------------------------ 00244 00245 //------------------------------------------------------------------------------------ 00246 // Kalman time update. 00247 // This routine uses the Householder transformation to propagate the SRIFilter 00248 // state and covariance through a time step. 00249 // Input: 00250 // R a priori square root information (SRI) matrix (an n by n 00251 // upper triangular matrix) 00252 // Z a priori SRIF state vector, of length n (state is X, Z = R*X). 00253 // PhiInv Inverse of state transition matrix, an n by n matrix. 00254 // PhiInv is destroyed on output. 00255 // Rw a priori square root information matrix for the process 00256 // noise, an ns by ns upper triangular matrix 00257 // G The n by ns matrix associated with process noise. The 00258 // process noise covariance is G*Q*transpose(G) where inverse(Q) 00259 // is transpose(Rw)*Rw. G is destroyed on output. 00260 // Zw a priori 'state' associated with the process noise, 00261 // a vector with ns elements. Usually set to zero by 00262 // the calling routine (for unbiased process noise). 00263 // Rwx An ns by n matrix which is set to zero by this routine 00264 // but is used for output. 00265 // 00266 // Output: 00267 // The updated square root information matrix and SRIF state (R,Z) and 00268 // the matrices which are used in smoothing: Rw, Zw, Rwx. 00269 // Note that PhiInv and G are trashed, and that Rw and Zw are modified. 00270 // 00271 // Return values: 00272 // SrifTU returns void, but throws exceptions if the input matrices 00273 // or vectors have incompatible dimensions or incorrect types. 00274 // 00275 // Method: 00276 // This SRIF time update method treats the process noise and mapping 00277 // information as a separate data equation, and applies a Householder 00278 // transformation to the (appended) equations to solve for an updated 00279 // state. Thus there is another 'state' variable associated with 00280 // whatever state variables have process noise. The matrix G relates 00281 // the process noise variables to the regular state variables, and 00282 // appears in the term GQG(trans) of the covariance. If all n state 00283 // variables have process noise, then ns=n and G is an n by n matrix. 00284 // Since some (or all) of the state variables may not have process 00285 // noise, ns may be zero. [Bierman ftnt pg 122 seems to indicate that 00286 // variables with zero process noise can be handled by ns=n & setting a 00287 // column of G=0. But note that the case of the matrix G=0 is the 00288 // same as ns=0, because the first ns columns would be zero below the 00289 // diagonal in that case anyway, so the HH transformation would be 00290 // null.] 00291 // For startup, all of the a priori information and state arrays may 00292 // be zero. That is, "no information" would imply that R and Z are zero, 00293 // as well as Rw and Zw. A priori information (covariance) and state 00294 // are handled by setting P = inverse(R)*transpose(inverse((R)), Z = R*X. 00295 // There are three ways to handle non-zero process noise covariance. 00296 // (1) If Q is the (known) a priori process noise covariance Q, then 00297 // set Q=Rw(-1)*Rw(-T), and G=1. 00298 // (2) Transform process noise covariance matrix to UDU form, Q=UDU, 00299 // then set G=U and Rw = (D)**-1/2. 00300 // (3) Take the sqrt of process noise covariance matrix Q, then set 00301 // G=this sqrt and Rw = 1. [2 and 3 have been tested.] 00302 // The routine applies a Householder transformation to a large 00303 // matrix formed by concatenation of the input matricies. Two preliminary 00304 // steps are to form Rd = R*PhiInv (stored in PhiInv) and -Rd*G (stored in 00305 // G) by matrix multiplication, and to set Rwx to the zero matrix. 00306 // Then the Householder transformation is applied to the following 00307 // matrix, dimensions are shown in (): 00308 // _ (ns) (n) (1) _ _ _ 00309 // (ns) | Rw 0 Zw | ==> | Rw Rwx Zw | 00310 // (n) | -Rd*G Rd Z | ==> | 0 R Z | . 00311 // - - - - 00312 // The SRI matricies R and Rw remain upper triangular. 00313 // 00314 // For the programmer: after Rwx is set to zero, G is made into 00315 // -Rd*G and PhiInv is made into R*PhiInv, the transformation is applied 00316 // to the matrix: 00317 // _ (ns) (n) (1) _ 00318 // (ns) | Rw Rwx Zw | 00319 // (n) | G PhiInv Z | 00320 // - - 00321 // then the (upper triangular) matrix R is copied out of PhiInv into R. 00322 // ------------------------------------------------------------------- 00323 // The matrix Rwx is related to the sensitivity of the state 00324 // estimate to the unmodeled parameters in Zw. The sensitivity matrix 00325 // is Sen = -inverse(Rw)*Rwx, 00326 // where perturbation in model X = 00327 // Sen * diagonal(a priori sigmas of parameter uncertainties). 00328 // ------------------------------------------------------------------- 00329 // The quantities Rw, Rwx and Zw on output are to be saved and used 00330 // in the sqrt information fixed interval smoother (SRIS), during the 00331 // backward filter process. 00332 // ------------------------------------------------------------------- 00333 // Ref: Bierman, G.J. "Factorization Methods for Discrete Sequential 00334 // Estimation," Academic Press, 1977, pg 121. 00335 // ------------------------------------------------------------------- 00336 template <class T> 00337 void SRIFilter::SrifTU(Matrix<T>& R, 00338 Vector<T>& Z, 00339 Matrix<T>& PhiInv, 00340 Matrix<T>& Rw, 00341 Matrix<T>& G, 00342 Vector<T>& Zw, 00343 Matrix<T>& Rwx) 00344 throw(MatrixException) 00345 { 00346 const T EPS=-T(1.e-200); 00347 unsigned int n=R.rows(),ns=Rw.rows(); 00348 unsigned int i,j,k; 00349 T sum, beta, delta, dum; 00350 00351 if(PhiInv.rows() < n || PhiInv.cols() < n || 00352 G.rows() < n || G.cols() < ns || 00353 R.cols() != n || 00354 Rwx.rows() < ns || Rwx.cols() < n || 00355 Z.size() < n || Zw.size() < ns) { 00356 MatrixException me("Invalid input dimensions:\n R is " 00357 + asString<int>(R.rows()) + "x" 00358 + asString<int>(R.cols()) + ", Z has length " 00359 + asString<int>(Z.size()) + "\n PhiInv is " 00360 + asString<int>(PhiInv.rows()) + "x" 00361 + asString<int>(PhiInv.cols()) + "\n Rw is " 00362 + asString<int>(Rw.rows()) + "x" 00363 + asString<int>(Rw.cols()) + "\n G is " 00364 + asString<int>(G.rows()) + "x" 00365 + asString<int>(G.cols()) + "\n Zw has length " 00366 + asString<int>(Zw.size()) + "\n Rwx is " 00367 + asString<int>(Rwx.rows()) + "x" 00368 + asString<int>(Rwx.cols()) 00369 ); 00370 GPSTK_THROW(me); 00371 } 00372 00373 try { 00374 // initialize 00375 Rwx = T(0); 00376 PhiInv = R * PhiInv; // set PhiInv = Rd = R*PhiInv 00377 G = -PhiInv * G; 00378 // fixed Matrix problem - unary minus should not return an l-value 00379 //G = -(PhiInv * G); // set G = -Rd*G 00380 00381 // temp 00382 //Matrix <T> A; 00383 //A = (Rw || Rwx || Zw) && (G || PhiInv || Z); 00384 //cout << "SrifTU - :\n" << fixed << setw(10) << setprecision(5) << A << endl; 00385 00386 //--------------------------------------------------------------- 00387 for(j=0; j<ns; j++) { // loop over first ns columns 00388 sum = T(0); 00389 for(i=0; i<n; i++) // rows of -Rd*G 00390 sum += G(i,j)*G(i,j); 00391 dum = Rw(j,j); 00392 sum += dum*dum; 00393 sum = (dum > T(0) ? -T(1) : T(1)) * ::sqrt(sum); 00394 delta = dum - sum; 00395 Rw(j,j) = sum; 00396 00397 beta = sum * delta; 00398 if(beta > EPS) continue; 00399 beta = T(1)/beta; 00400 00401 // apply jth Householder transformation 00402 // to submatrix below and right of (j,j) 00403 if(j+1 < ns) { // apply to G 00404 for(k=j+1; k<ns; k++) { // columns to right of diagonal 00405 sum = delta * Rw(j,k); 00406 for(i=0; i<n; i++) // rows of G 00407 sum += G(i,j)*G(i,k); 00408 if(sum == T(0)) continue; 00409 sum *= beta; 00410 Rw(j,k) += sum*delta; 00411 for(i=0; i<n; i++) // rows of G again 00412 G(i,k) += sum * G(i,j); 00413 } 00414 } 00415 00416 // apply jth Householder transformation 00417 // to Rwx and PhiInv 00418 for(k=0; k<n; k++) { // columns of Rwx and PhiInv 00419 sum = delta * Rwx(j,k); 00420 for(i=0; i<n; i++) // rows of PhiInv and G 00421 sum += PhiInv(i,k) * G(i,j); 00422 if(sum == T(0)) continue; 00423 sum *= beta; 00424 Rwx(j,k) += sum*delta; 00425 for(i=0; i<n; i++) // rows of PhiInv and G 00426 PhiInv(i,k) += sum * G(i,j); 00427 } // end loop over columns of Rwx and PhiInv 00428 00429 // apply jth Householder transformation 00430 // to Zw and Z 00431 sum = delta * Zw(j); 00432 for(i=0; i<n; i++) // rows of G and elements of Z 00433 sum += Z(i) * G(i,j); 00434 if(sum == T(0)) continue; 00435 sum *= beta; 00436 Zw(j) += sum * delta; 00437 for(i=0; i<n; i++) // rows of G and elements of Z 00438 Z(i) += sum * G(i,j); 00439 } // end loop over first ns columns 00440 00441 //--------------------------------------------------------------- 00442 for(j=0; j<n; j++) { // loop over columns of Rwx and PhiInv 00443 sum = T(0); 00444 for(i=j+1; i<n; i++) // rows of PhiInv 00445 sum += PhiInv(i,j)*PhiInv(i,j); 00446 dum = PhiInv(j,j); 00447 sum += dum*dum; 00448 sum = (dum > T(0) ? -T(1) : T(1)) * ::sqrt(sum); 00449 delta = dum - sum; 00450 PhiInv(j,j) = sum; 00451 beta = sum*delta; 00452 if(beta > EPS) continue; 00453 beta = T(1)/beta; 00454 00455 // apply jth Householder transformation to columns of PhiInv on row j 00456 for(k=j+1; k<n; k++) { // columns of PhiInv 00457 sum = delta * PhiInv(j,k); 00458 for(i=j+1; i<n; i++) 00459 sum += PhiInv(i,j)*PhiInv(i,k); 00460 if(sum == T(0)) continue; 00461 sum *= beta; 00462 PhiInv(j,k) += sum*delta; 00463 for(i=j+1; i<n; i++) 00464 PhiInv(i,k) += sum * PhiInv(i,j); 00465 } 00466 00467 // apply jth Householder transformation to Z 00468 sum = delta *Z(j); 00469 for(i=j+1; i<n; i++) 00470 sum += Z(i) * PhiInv(i,j); 00471 if(sum == T(0)) continue; 00472 sum *= beta; 00473 Z(j) += sum*delta; 00474 for(i=j+1; i<n; i++) 00475 Z(i) += sum * PhiInv(i,j); 00476 } // end loop over cols of Rwx and PhiInv 00477 00478 // temp 00479 //A = (Rw || Rwx || Zw) && (G || PhiInv || Z); 00480 //cout << "SrifTU + :\n" << fixed << setw(10) << setprecision(5) << A << endl; 00481 00482 // copy transformed R out of PhiInv 00483 for(j=0; j<n; j++) 00484 for(i=0; i<=j; i++) 00485 R(i,j) = PhiInv(i,j); 00486 } 00487 catch(MatrixException& me) { GPSTK_RETHROW(me); } 00488 } // end SrifTU 00489 00490 //------------------------------------------------------------------------------------ 00491 // Kalman smoother update. 00492 // This routine uses the Householder transformation to propagate the SRIF 00493 // state and covariance through a smoother (backward filter) step. 00494 // Input: 00495 // R A priori square root information (SRI) matrix (an N by N 00496 // upper triangular matrix) 00497 // z a priori SRIF state vector, an N vector (state is x, z = R*x). 00498 // Phi State transition matrix, an N by N matrix. Phi is destroyed on output. 00499 // Rw A priori square root information matrix for the process 00500 // noise, an Ns by Ns upper triangular matrix (which has 00501 // Ns(Ns+1)/2 elements). 00502 // G The N by Ns matrix associated with process noise. The 00503 // process noise covariance is GQGtrans where Qinverse 00504 // is Rw(trans)*Rw. G is destroyed on output. 00505 // Zw A priori 'state' associated with the process noise, 00506 // a vector with Ns elements. Zw is destroyed on output. 00507 // Rwx An Ns by N matrix. Rwx is destroyed on output. 00508 // 00509 // The inputs Rw,Zw,Rwx are the output of the SRIF time update, and these and 00510 // Phi and G are associated with the same timestep. 00511 // 00512 // Output: 00513 // The updated square root information matrix and SRIF smoothed state (R,z). 00514 // All other inputs are trashed. 00515 // 00516 // Return values: 00517 // SrifSU returns void, but throws exceptions if the input matrices 00518 // or vectors have incompatible dimensions or incorrect types. 00519 // 00520 // Method: 00521 // The fixed interval square root information smoother (SRIS) is 00522 // composed of two Kalman filters, one identical with the square root 00523 // information filter (SRIF), the other similar but operating on the 00524 // data in reverse order and combining the current (smoothed) state 00525 // with elements output by the SRIF in its forward run and saved. 00526 // Thus a smoother is composed of a forward filter which saves all of 00527 // its output, followed by a backward filter which makes use of that 00528 // saved information. 00529 // This form of the SRIF backward filter algorithm is equivalent to the 00530 // Dyer-McReynolds SRIS algorithm, which uses less computer resources, but 00531 // propagates the state and covariance rather than the SRI (R,z). (As always, 00532 // at any point the state X and covariance P are related to the SRI by 00533 // X = R^-1 * z , P = R^-1 * R^-T.) 00534 // For startup of the backward filter, the state after the final 00535 // measurement update of the SRIF is given another time update, the 00536 // output of which is identified with the a priori values for the 00537 // backward filter. Backward filtering proceeds from there, the N+1st 00538 // point, toward the first point. 00539 // 00540 // In this implementation of the backward filter, the Householder 00541 // transformation is applied to the following matrix 00542 // (dimensions are shown in ()): 00543 // 00544 // _ (Ns) (N) (1) _ _ _ 00545 // (Ns) | Rw+Rwx*G Rwx*Phi Zw | ==> | Rw Rwx Zw | 00546 // (N) | R*G R*Phi z | ==> | 0 R z | . 00547 // - - - - 00548 // The SRI matricies R and Rw remain upper triangular. 00549 // 00550 // For the programmer: First create an NsXNs matrix A, then 00551 // Rw+Rwx*G -> A, Rwx*Phi -> Rwx, R*Phi -> Phi, and R*G -> G, and 00552 // the transformation is applied to the matrix: 00553 // 00554 // _ (Ns) (N) (1) _ 00555 // (Ns) | A Rwx Zw | 00556 // (N) | G Phi z | 00557 // - - 00558 // then the (upper triangular) matrix R is copied out of Phi into R. 00559 // 00560 // Ref: Bierman, G.J. "Factorization Methods for Discrete Sequential 00561 // Estimation," Academic Press, 1977, pg 216. 00562 template <class T> 00563 void SRIFilter::SrifSU(Matrix<T>& R, 00564 Vector<T>& Z, 00565 Matrix<T>& Phi, 00566 Matrix<T>& Rw, 00567 Matrix<T>& G, 00568 Vector<T>& Zw, 00569 Matrix<T>& Rwx) 00570 throw(MatrixException) 00571 { 00572 unsigned int N=R.rows(),Ns=Rw.rows(); 00573 00574 if(Phi.rows() < N || Phi.cols() < N || 00575 G.rows() < N || G.cols() < Ns || 00576 R.cols() != N || 00577 Rwx.rows() < Ns || Rwx.cols() < N || 00578 Z.size() < N || Zw.size() < Ns) { 00579 MatrixException me("Invalid input dimensions:\n R is " 00580 + asString<int>(R.rows()) + "x" 00581 + asString<int>(R.cols()) + ", Z has length " 00582 + asString<int>(Z.size()) + "\n Phi is " 00583 + asString<int>(Phi.rows()) + "x" 00584 + asString<int>(Phi.cols()) + "\n Rw is " 00585 + asString<int>(Rw.rows()) + "x" 00586 + asString<int>(Rw.cols()) + "\n G is " 00587 + asString<int>(G.rows()) + "x" 00588 + asString<int>(G.cols()) + "\n Zw has length " 00589 + asString<int>(Zw.size()) + "\n Rwx is " 00590 + asString<int>(Rwx.rows()) + "x" 00591 + asString<int>(Rwx.cols()) 00592 ); 00593 GPSTK_THROW(me); 00594 } 00595 00596 const T EPS=-T(1.e-200); 00597 int i, j, k; 00598 T sum, beta, delta, diag; 00599 00600 try { 00601 // Rw+Rwx*G -> A 00602 Matrix<T> A; 00603 A = Rw + Rwx*G; 00604 Rwx = Rwx * Phi; 00605 Phi = R * Phi; 00606 G = R * G; 00607 00608 //----------------------------------------- 00609 // HouseHolder Transformation 00610 00611 // Loop over first Ns columns 00612 for(j=0; j<Ns; j++) { // columns of A 00613 sum = T(0); 00614 for(i=j+1; i<Ns; i++) { // rows i below diagonal in A 00615 sum += A(i,j) * A(i,j); 00616 } 00617 for(i=0; i<N; i++) { // all rows i in G 00618 sum += G(i,j) * G(i,j); 00619 } 00620 00621 diag = A(j,j); 00622 sum += diag*diag; 00623 sum = (diag > T(0) ? -T(1) : T(1)) * ::sqrt(sum); 00624 delta = diag - sum; 00625 A(j,j) = sum; 00626 beta = sum*delta; 00627 if(beta > EPS) continue; 00628 beta = T(1)/beta; 00629 00630 // apply jth HH trans to submatrix below and right of j,j 00631 for(k=j+1; k<Ns; k++) { // cols to right of diag 00632 sum = delta * A(j,k); 00633 for(i=j+1; i<Ns; i++) { // rows of A below diagonal 00634 sum += A(i,j)*A(i,k); 00635 } 00636 for(i=0; i<N; i++) { // all rows of G 00637 sum += G(i,j)*G(i,k); 00638 } 00639 if(sum == T(0)) continue; 00640 sum *= beta; 00641 //------------------------------------------ 00642 A(j,k) += sum*delta; 00643 00644 for(i=j+1; i<Ns; i++) { // rows of A > j (same loops again) 00645 A(i,k) += sum * A(i,j); 00646 } 00647 for(i=0; i<N; i++) { // all rows of G (again) 00648 G(i,k) += sum * G(i,j); 00649 } 00650 } 00651 00652 // apply jth HH trans to Rwx and Phi sub-matrices 00653 for(k=0; k<N; k++) { // all columns of Rwx / Phi 00654 sum = delta * Rwx(j,k); 00655 for(i=j+1; i<Ns; i++) { // rows of Rwx below j 00656 sum += A(i,j) * Rwx(i,k); 00657 } 00658 for(i=0; i<N; i++) { // all rows of Phi 00659 sum += G(i,j) * Phi(i,k); 00660 } 00661 if(sum == T(0)) continue; 00662 sum *= beta; 00663 Rwx(j,k) += sum*delta; 00664 for(i=j+1; i<Ns; i++) { // rows of Rwx below j (again) 00665 Rwx(i,k) += sum * A(i,j); 00666 } 00667 for(i=0; i<N; i++) { // all rows of Phi (again) 00668 Phi(i,k) += sum * G(i,j); 00669 } 00670 } 00671 00672 // apply jth HH trans to Zw and Z 00673 sum = delta * Zw(j); 00674 for(i=j+1; i<Ns; i++) { // rows (elements) of Zw below j 00675 sum += A(i,j) * Zw(i); 00676 } 00677 for(i=0; i<N; i++) { // all rows (elements) of Z 00678 sum += Z(i) * G(i,j); 00679 } 00680 if(sum == T(0)) continue; 00681 sum *= beta; 00682 Zw(j) += sum*delta; 00683 for(i=j+1; i<Ns; i++) { // rows of Zw below j (again) 00684 Zw(i) += sum * A(i,j); 00685 } 00686 for(i=0; i<N; i++) { // all rows of Z (again) 00687 Z(i) += sum * G(i,j); 00688 } 00689 } 00690 00691 // Loop over columns past the Ns block: all of Rwx and Phi 00692 for(j=0; j<N; j++) { 00693 sum = T(0); 00694 for(i=j+1; i<N; i++) { // rows of Phi at and below j 00695 sum += Phi(i,j) * Phi(i,j); 00696 } 00697 diag = Phi(j,j); 00698 sum += diag*diag; 00699 sum = (diag > T(0) ? -T(1) : T(1)) * ::sqrt(sum); 00700 delta = diag - sum; 00701 Phi(j,j) = sum; 00702 beta = sum*delta; 00703 if(beta > EPS) continue; 00704 beta = T(1)/beta; 00705 00706 // apply HH trans to Phi sub-block below and right of j,j 00707 for(k=j+1; k<N; k++) { // columns k > j 00708 sum = delta * Phi(j,k); 00709 for(i=j+1; i<N; i++) { // rows below j 00710 sum += Phi(i,j) * Phi(i,k); 00711 } 00712 if(sum == T(0)) continue; 00713 sum *= beta; 00714 Phi(j,k) += sum*delta; 00715 for(i=j+1; i<N; i++) { // rows below j (again) 00716 Phi(i,k) += sum * Phi(i,j); 00717 } 00718 } 00719 // Now apply to the Z column 00720 sum = delta * Z(j); 00721 for(i=j+1; i<N; i++) { // rows of Z below j 00722 sum += Z(i) * Phi(i,j); 00723 } 00724 if(sum == T(0)) continue; 00725 sum *= beta; 00726 Z(j) += sum*delta; 00727 for(i=j+1; i<N; i++) { // rows of Z below j (again) 00728 Z(i) += sum * Phi(i,j); 00729 } 00730 } 00731 //------------------------------ 00732 // Transformation finished 00733 00734 //------------------------------------- 00735 // copy transformed R out of Phi into R 00736 R = T(0); 00737 for(j=0; j<N; j++) { 00738 for(i=0; i<=j; i++) { 00739 R(i,j) = Phi(i,j); 00740 } 00741 } 00742 } 00743 catch(Exception& e) { GPSTK_RETHROW(e); } 00744 } // end SrifSU 00745 00746 //------------------------------------------------------------------------------ 00747 // Covariance/State version of the Kalman smoother update (Dyer-McReynolds). 00748 // This routine implements the Dyer-McReynolds form of the state and covariance 00749 // recursions which constitute the backward filter of the Square Root 00750 // Information Smoother. 00751 // 00752 // Input: (assume N and Ns are greater than zero) 00753 // Vector X(N) A priori state, derived from SRI (R*X=Z) 00754 // Matrix P(N,N) A priori covariance, derived from SRI (P=R^-1*R^-T) 00755 // Matrix Rw(Ns,Ns) Process noise covariance (UT), output of SRIF TU 00756 // Matrix Rwx(Ns,N) PN 'cross term', output of SRIF TU 00757 // Vector Zw(Ns) Process noise state, output of SRIF TU 00758 // Matrix Phinv(N,N) Inverse of state transition, saved at SRIF TU 00759 // Matrix G(N,Ns) Noise coupling matrix, saved at SRIF TU 00760 // Output: 00761 // Updated X and P. The other inputs are trashed. 00762 // 00763 // Method: 00764 // The fixed interval square root information smoother (SRIS) is 00765 // composed of two Kalman filters, one identical with the square root 00766 // information filter (SRIF), the other similar but operating on the 00767 // data in reverse order and combining the current (smoothed) state 00768 // with elements output by the SRIF in its forward run and saved. 00769 // Thus a smoother is composed of a forward filter which saves all of 00770 // its output, followed by a backward filter which makes use of that 00771 // saved information. 00772 // This form of the SRIS algorithm is equivalent to the SRIS backward 00773 // filter Householder transformation algorithm, but uses less computer 00774 // resources. It is not necessary to update both the state and the 00775 // covariance, although doing both at once is less expensive than 00776 // doing them separately. (This routine does both.) 00777 // For startup of the backward filter, the state after the final 00778 // measurement update of the SRIF is given another time update, the 00779 // output of which is identified with the a priori values for the 00780 // backward filter. Backward filtering proceeds from there, the N+1st 00781 // point, toward the first point. 00782 // 00783 // Ref: Bierman, G.J. "Factorization Methods for Discrete Sequential 00784 // Estimation," Academic Press, 1977, pg 216. 00785 template <class T> 00786 void SRIFilter::SrifSU_DM(Matrix<T>& P, 00787 Vector<T>& X, 00788 Matrix<T>& Phinv, 00789 Matrix<T>& Rw, 00790 Matrix<T>& G, 00791 Vector<T>& Zw, 00792 Matrix<T>& Rwx) 00793 throw(MatrixException) 00794 { 00795 unsigned int N=P.rows(),Ns=Rw.rows(); 00796 00797 if(P.cols() != P.rows() || 00798 X.size() != N || 00799 Rwx.cols() != N || 00800 Zw.size() != Ns || 00801 Rwx.rows() != Ns || Rwx.cols() != N || 00802 Phinv.rows() != N || Phinv.cols() != N || 00803 G.rows() != N || G.cols() != Ns ) { 00804 MatrixException me("Invalid input dimensions:\n P is " 00805 + asString<int>(P.rows()) + "x" 00806 + asString<int>(P.cols()) + ", X has length " 00807 + asString<int>(X.size()) + "\n Phinv is " 00808 + asString<int>(Phinv.rows()) + "x" 00809 + asString<int>(Phinv.cols()) + "\n Rw is " 00810 + asString<int>(Rw.rows()) + "x" 00811 + asString<int>(Rw.cols()) + "\n G is " 00812 + asString<int>(G.rows()) + "x" 00813 + asString<int>(G.cols()) + "\n Zw has length " 00814 + asString<int>(Zw.size()) + "\n Rwx is " 00815 + asString<int>(Rwx.rows()) + "x" 00816 + asString<int>(Rwx.cols()) 00817 ); 00818 GPSTK_THROW(me); 00819 } 00820 00821 try { 00822 G = G * inverse(Rw); 00823 Matrix<T> F; 00824 F = ident<T>(N) + G*Rwx; 00825 // update X 00826 Vector<T> C; 00827 C = F*X - G*Zw; 00828 X = Phinv * C; 00829 // update P 00830 P = F*P*transpose(F) + G*transpose(G); 00831 P = Phinv*P*transpose(Phinv); 00832 } 00833 catch(Exception& e) { GPSTK_RETHROW(e); } 00834 } // end SrifSU_DM 00835 00836 // Modification for case with control vector: Xj+1 = Phi*Xj + Gwj + u 00837 template <class T> 00838 void DMsmootherUpdateWithControl(Matrix<double>& P, 00839 Vector<double>& X, 00840 Matrix<double>& Phinv, 00841 Matrix<double>& Rw, 00842 Matrix<double>& G, 00843 Vector<double>& Zw, 00844 Matrix<double>& Rwx, 00845 Vector<double>& U) 00846 throw(MatrixException) 00847 { 00848 unsigned int N=P.rows(),Ns=Rw.rows(); 00849 00850 if(P.cols() != P.rows() || 00851 X.size() != N || 00852 Rwx.cols() != N || 00853 Zw.size() != Ns || 00854 Rwx.rows() != Ns || Rwx.cols() != N || 00855 Phinv.rows() != N || Phinv.cols() != N || 00856 G.rows() != N || G.cols() != Ns || 00857 U.size() != N) { 00858 MatrixException me("Invalid input dimensions:\n P is " 00859 + asString<int>(P.rows()) + "x" 00860 + asString<int>(P.cols()) + ", X has length " 00861 + asString<int>(X.size()) + "\n Phinv is " 00862 + asString<int>(Phinv.rows()) + "x" 00863 + asString<int>(Phinv.cols()) + "\n Rw is " 00864 + asString<int>(Rw.rows()) + "x" 00865 + asString<int>(Rw.cols()) + "\n G is " 00866 + asString<int>(G.rows()) + "x" 00867 + asString<int>(G.cols()) + "\n Zw has length " 00868 + asString<int>(Zw.size()) + "\n Rwx is " 00869 + asString<int>(Rwx.rows()) + "x" 00870 + asString<int>(Rwx.cols()) + "\n U has length " 00871 + asString<int>(U.size()) 00872 ); 00873 GPSTK_THROW(me); 00874 } 00875 00876 try { 00877 G = G * inverse(Rw); 00878 Matrix<T> F; 00879 F = ident<T>(N) + G*Rwx; 00880 // update X 00881 Vector<T> C; 00882 C = F*X - G*Zw - U; 00883 X = Phinv * C; 00884 // update P 00885 P = F*P*transpose(F) + G*transpose(G); 00886 P = Phinv*P*transpose(Phinv); 00887 P += outer(U,U); 00888 } 00889 catch(Exception& e) { GPSTK_RETHROW(e); } 00890 } // end DMsmootherUpdateWithControl 00891 00892 //------------------------------------------------------------------------------------ 00893 } // end namespace gpstk
1.3.9.1