SRIleastSquares.cpp

Go to the documentation of this file.
00001 #pragma ident "$Id: SRIleastSquares.cpp 1894 2009-05-11 20:00:55Z afarris $"
00002 
00003 
00004 //============================================================================
00005 //
00006 //  This file is part of GPSTk, the GPS Toolkit.
00007 //
00008 //  The GPSTk is free software; you can redistribute it and/or modify
00009 //  it under the terms of the GNU Lesser General Public License as published
00010 //  by the Free Software Foundation; either version 2.1 of the License, or
00011 //  any later version.
00012 //
00013 //  The GPSTk is distributed in the hope that it will be useful,
00014 //  but WITHOUT ANY WARRANTY; without even the implied warranty of
00015 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00016 //  GNU Lesser General Public License for more details.
00017 //
00018 //  You should have received a copy of the GNU Lesser General Public
00019 //  License along with GPSTk; if not, write to the Free Software Foundation,
00020 //  Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00021 //
00022 //  Copyright 2004, The University of Texas at Austin
00023 //
00024 //============================================================================
00025 
00026 //============================================================================
00027 //
00028 //This software developed by Applied Research Laboratories at the University of
00029 //Texas at Austin, under contract to an agency or agencies within the U.S. 
00030 //Department of Defense. The U.S. Government retains all rights to use,
00031 //duplicate, distribute, disclose, or release this software. 
00032 //
00033 //Pursuant to DoD Directive 523024 
00034 //
00035 // DISTRIBUTION STATEMENT A: This software has been approved for public 
00036 //                           release, distribution is unlimited.
00037 //
00038 //=============================================================================
00039 
00047 //------------------------------------------------------------------------------------
00048 // GPSTk includes
00049 #include "SRIleastSquares.hpp"
00050 #include "RobustStats.hpp"
00051 #include "StringUtils.hpp"
00052 
00053 //------------------------------------------------------------------------------------
00054 using namespace std;
00055 
00056 namespace gpstk {
00057 using namespace StringUtils;
00058 
00059 //------------------------------------------------------------------------------------
00060 // empty constructor
00061 SRIleastSquares::SRIleastSquares(void) throw()
00062 { defaults(); }
00063 
00064 //------------------------------------------------------------------------------------
00065 // constructor given the dimension N.
00066 SRIleastSquares::SRIleastSquares(const unsigned int N)
00067    throw()
00068 {
00069    defaults();
00070    R = Matrix<double>(N,N,0.0);
00071    Z = Vector<double>(N,0.0);
00072    names = Namelist(N);
00073 }
00074 
00075 //------------------------------------------------------------------------------------
00076 // constructor given a Namelist, its dimension determines the SRI dimension.
00077 SRIleastSquares::SRIleastSquares(const Namelist& NL)
00078    throw()
00079 {
00080    defaults();
00081    if(NL.size() <= 0) return;
00082    R = Matrix<double>(NL.size(),NL.size(),0.0);
00083    Z = Vector<double>(NL.size(),0.0);
00084    names = NL;
00085 }
00086 
00087 //------------------------------------------------------------------------------------
00088 // explicit constructor - throw if the dimensions are inconsistent.
00089 SRIleastSquares::SRIleastSquares(const Matrix<double>& Rin,
00090                      const Vector<double>& Zin,
00091                      const Namelist& NLin)
00092    throw(MatrixException)
00093 {
00094    defaults();
00095    if(Rin.rows() != Rin.cols() ||
00096       Rin.rows() != Zin.size() ||
00097       Rin.rows() != NLin.size()) {
00098       MatrixException me("Invalid input dimensions: R is "
00099          + asString<int>(Rin.rows()) + "x"
00100          + asString<int>(Rin.cols()) + ", Z has length "
00101          + asString<int>(Zin.size()) + ", and NL has length "
00102          + asString<int>(NLin.size())
00103          );
00104       GPSTK_THROW(me);
00105    }
00106    R = Rin;
00107    Z = Zin;
00108    names = NLin;
00109 }
00110 
00111 //------------------------------------------------------------------------------------
00112 // operator=
00113 SRIleastSquares& SRIleastSquares::operator=(const SRIleastSquares& right)
00114    throw()
00115 {
00116    R = right.R;
00117    Z = right.Z;
00118    names = right.names;
00119    iterationsLimit = right.iterationsLimit;
00120    convergenceLimit = right.convergenceLimit;
00121    divergenceLimit = right.divergenceLimit;
00122    doWeight = right.doWeight;
00123    doRobust = right.doRobust;
00124    doLinearize = right.doLinearize;
00125    doSequential = right.doSequential;
00126    doVerbose = right.doVerbose;
00127    valid = right.valid;
00128    number_iterations = right.number_iterations;
00129    number_batches = right.number_batches;
00130    rms_convergence = right.rms_convergence;
00131    condition_number = right.condition_number;
00132    Xsave = right.Xsave;
00133    return *this;
00134 }
00135 
00136 //------------------------------------------------------------------------------------
00137 // SRI least squares update (not the Kalman measurement update).
00138 // Given data and measurement covariance, compute a solution and
00139 // covariance using the appropriate least squares algorithm.
00140 // @param D   Data vector, length M
00141 //               Input:  raw data
00142 //               Output: post-fit residuals
00143 // @param X   Solution vector, length N
00144 //               Input:  nominal solution X0 (zero when doLinearized is false)
00145 //               Output: final solution
00146 // @param Cov Covariance matrix, dimension (N,N)
00147 //               Input:  (If doWeight is true) inverse measurement covariance
00148 //                       or weight matrix(M,M)
00149 //               Output: Solution covariance matrix (N,N)
00150 // @param LSF Pointer to a function which is used to define the equation to be solved.
00151 // LSF arguments are:
00152 //            X  Nominal solution (input)
00153 //            f  Values of the equation f(X) (length M) (output)
00154 //            P  Partials matrix df/dX evaluated at X (dimension M,N) (output)
00155 //        When doLinearize is false, LSF should ignore X and return the (constant)
00156 //        partials matrix in P and zero in f.
00157 // @return 0 ok
00158 //               -1 Problem is underdetermined (M<N) // TD -- naturalized sol?
00159 //               -2 Problem is singular
00160 //               -3 Algorithm failed to converge
00161 //               -4 Algorithm diverged
00162 //
00163 // Reference for robust least squares: Mason, Gunst and Hess,
00164 // "Statistical Design and Analysis of Experiments," Wiley, New York, 1989, pg 593.
00165 //
00166 // Notes on the algorithm:
00167 // Least squares, including linearized (iterative) and sequential processing.
00168 // This class will solve the equation f(X) = D, a vector equation in which the
00169 // solution vector X is of length N, and the data vector D is of length M.
00170 // The function f(X) may be linear, in which case it is of the form
00171 // P*X=D where P is a constant matrix,
00172 // or non-linear, in which case it will be linearized by expanding about a given
00173 // nominal solution X0:
00174 //          df |
00175 //          -- |     * dX = D - f(X0),
00176 //          dX |X=X0
00177 // where dX is defined as (X-X0), the new solution is X, and the partials matrix is
00178 // P=(df/dX)|X=X0. Dimensions are P(M,N)*dX(N) = D(M) - f(X0)(M).
00179 // Linearized problems are iterated until the solution converges (stops changing). 
00180 // 
00181 // The solution may be weighted by a measurement covariance matrix MCov,
00182 // or weight matrix W (in which case MCov = inverse(W)). MCov must be non-singular.
00183 // 
00184 // Options are to make the algorithm linearized (via the boolean input variable
00185 // doLinearize) and/or sequential (doSequential).
00186 // 
00187 //    - linearized. When doLinearize is true, the algorithm solves the linearized
00188 //    version of the measurement equation (see above), rather than the simple
00189 //    linear version P*X=D. Also when doLinearize is true, the code will iterate
00190 //    (repeat until convergence) the linearized algorithm; if you don't want to
00191 //    iterate, set the limit on the number of iterations to zero.
00192 //    NB In this case, a solution must be found for each nominal solution
00193 //    (i.e. the information matrix must be non-singular); otherwise there can be
00194 //    no iteration.
00195 // 
00196 //    - sequential. When doSequential is true, the class will save the accumulated
00197 //    information from all the calls to this routine since the last reset()
00198 //    within the class. This means the resulting solution is determined by ALL the
00199 //    data fed to the class since the last reset(). In this case the data is fed
00200 //    to the algorithm in 'batches', which may be of any size.
00201 // 
00202 //    NB When doLinearize is true, the information stored in the class has a
00203 //    different interpretation than it does in the linear case.
00204 //    Calling Solve(X,Cov) will NOT give the solution vector X, but rather the
00205 //    latest update (X-X0) = (X-Xsave).
00206 // 
00207 //    NB In the linear case, the result you get from sequentially processing
00208 //    a large dataset in many small batches is identical to what you would get
00209 //    by processing all the data in one big batch. This is NOT true in the
00210 //    linearized case, because the information at each batch is dependent on the
00211 //    nominal state. See the next comment.
00212 // 
00213 //    NB Sequential, linearized LS really makes sense only when the state is
00214 //    changing. It is difficult to get a good solution in this case with small
00215 //    batches, because the stored information is dependent on the (final) state
00216 //    solution at each batch. Start with a good nominal state, or with a large
00217 //    batch of data that will produce one.
00218 // 
00219 // The general Least Squares algorithm is:
00220 //  0. set i=0.
00221 //  1. If non-sequential, or if this is the first call, set R=z=0
00222 //     (However doing this prevents you from adding apriori/constraint information)
00223 //  2. Let X = X0 (X0 = initial nominal solution - input). if linear, X0==0.
00224 //  3. Save SRIsave=SRI and X0save=X0                       (SRI is the pair R,z)
00225 //  4. start iteration i here.
00226 //  5. increment the number of iterations i
00227 //  6. Compute partials matrix P and f(X0) by calling LSF(X0,f,P).
00228 //        if linear, LSF returns the constant P and f(X0)=0.
00229 //  7. Set R = SRIsave.R + P(T)*inverse(MCov)*P                 (T means transpose)
00230 //  8. Set z = SRIsave.z + P(T)*inverse(MCov)*(D-f(X0))
00231 //  9. [The measurement equation is now P*DX=d-F(X0)
00232 //        where DX=(X-X0save); in the linear case it is PX = d and DX = X ]
00233 // 10. Solve z = Rx to get
00234 //          Cov = inverse(R)
00235 //       and DX = inverse(R)*z OR
00236 // 11. Set X = X0save + DX
00237 //     [or in the linear case X = DX
00238 // 12. Compute RMS change in X: rms = ||X-X0||/N    (not X-X0save)
00239 // 13. if linear goto quit [else linearized]
00240 // 14. If rms > divergence limit, goto quit(failure).
00241 // 15. If i > 1 and rms < convergence limit, goto quit(success)
00242 // 16. If i (number of iterations) >= iteration limit, goto quit(failure)
00243 // 17. Set X0 = X
00244 // 18. Return to step 5.
00245 // 19. quit: if(sequential and failed) set SRI=SRIsave.
00246 // 
00247 // From the code:
00248 //  1a. Save SRI (i.e. R, Z) in Rapriori, Zapriori
00249 //  2a. If non-sequential, or if this is the first call, set R=z=0 -- DON'T
00250 //  3a. If sequential and not the first call, X = Xsave
00251 //  4a. if linear, X0=0; else X0 is input. Let NominalX = X0
00252 //  5a. set number_iterations = 0
00253 //  6a. start iteration
00254 //  7a. increment number_iterations
00255 //  8a. get partials and f from LSfunc using NominalX
00256 //  9a. if robust, compute weight matrix
00257 // 10a. if number_iterations > 1, restore (R,Z) = (Rapriori,Zapriori)
00258 // 11a. MU : R,Z,Partials,D-f(NominalX),MeasCov(if weighted)
00259 // 12a. Invert to get Xsol  [ Xsol = X-NominalX or, if linear = X]
00260 // 13a. if linearized, add NominalX to Xsol; Xsol now == X = new estimate
00261 // 14a. if linear and not robust, quit here
00262 // 15a. if linearized, compute rms_convergence = RMS(Xsol - NominalX)
00263 // 16a. if robust, recompute weights and define rms_convergence = RMS(old-new wts)
00264 // 17a. failed? if so, and sequential, restore (R,Z) = (Rapriori,Zapriori); quit
00265 // 18a. success? quit
00266 // 19a. if linearized NominalX = Xsol;  if robust NominalX = X
00267 // 20a. iterate - return to 6a.
00268 // 21a. set X = Xsol for return value
00269 // 22a. save X for next time : Xsave = X
00270 //
00271 int SRIleastSquares::dataUpdate(Vector<double>& D,
00272                                 Vector<double>& X,
00273                                 Matrix<double>& Cov,
00274                                 void (LSF)(Vector<double>& X,
00275                                            Vector<double>& f,
00276                                            Matrix<double>& P)) throw(MatrixException)
00277 {
00278    const int M = D.size();
00279    const int N = R.rows();
00280    if(doVerbose) cout << "\nSRIleastSquares::leastSquaresUpdate : M,N are "
00281       << M << "," << N << endl;
00282 
00283    // errors
00284    if(N == 0) {
00285       MatrixException me("Called with zero-sized SRIleastSquares");
00286       GPSTK_THROW(me);
00287    }
00288    if(doLinearize && M < N) {
00289       MatrixException me(
00290             string("When linearizing, problem must not be underdetermined:\n")
00291             + string("   data dimension is ") + asString(M)
00292             + string(" while state dimension is ") + asString(N));
00293       GPSTK_THROW(me);
00294    }
00295    if(doSequential && R.rows() != X.size()) {
00296       MatrixException me("Sequential problem has inconsistent dimensions:\n  SRI is "
00297          + asString<int>(R.rows()) + "x"
00298          + asString<int>(R.cols()) + " while X has length "
00299          + asString<int>(X.size()));
00300       GPSTK_THROW(me);
00301    }
00302    if(doWeight && doRobust) {
00303       MatrixException me("Cannot have doWeight and doRobust both true.");
00304       GPSTK_THROW(me);
00305    }
00306    // TD disallow Robust and Linearized ? why?
00307    // TD disallow Robust and Sequential ? why?
00308 
00309 try {
00310    int i,j,iret;
00311    double big,small;
00312    Vector<double> f(M),Xsol(N),NominalX,Res(M),Wts(M,1.0),OldWts(M,1.0);
00313    Matrix<double> Partials(M,N),MeasCov(M,M);
00314    const Matrix<double> Rapriori(R);
00315    const Vector<double> Zapriori(Z);
00316 
00317    // save measurement covariance matrix
00318    if(doWeight) MeasCov=Cov;
00319 
00320    // NO ... this prevents you from giving it apriori information...
00321    // if the first time, clear the stored information
00322    //if(!doSequential || number_batches==0)
00323    //   zeroAll();
00324 
00325    // if sequential and not the first call, NominalX must be the last solution
00326    if(doSequential && number_batches != 0) X = Xsave;
00327 
00328    // nominal solution
00329    if(!doLinearize) {
00330       if(X.size() != N) X=Vector<double>(N);
00331       X = 0.0;
00332    }
00333    NominalX = X;
00334 
00335    valid = false;
00336    condition_number = 0.0;
00337    rms_convergence = 0.0;
00338    number_iterations = 0;
00339    iret = 0;
00340 
00341    // iteration loop
00342    do {
00343       number_iterations++;
00344 
00345       // call LSF to get f(NominalX) and Partials(NominalX)
00346       LSF(NominalX,f,Partials);
00347 
00348       // Res will be both pre- and post-fit data residuals
00349       Res = D-f;
00350       if(doVerbose) {
00351          cout << "\nSRIleastSquares::leastSquaresUpdate :";
00352          if(doLinearize || doRobust)
00353             cout << " Iteration " << number_iterations;
00354          cout << endl;
00355          LabelledVector LNX(names,NominalX);
00356          LNX.message(" Nominal X:");
00357          cout << LNX << endl;
00358          cout << " Pre-fit data residuals:  "
00359             << fixed << setprecision(6) << Res << endl;
00360       }
00361 
00362       // build measurement covariance matrix for robust LS
00363       if(doRobust) {
00364          MeasCov = 0.0;
00365          for(i=0; i<M; i++) MeasCov(i,i) = 1.0 / (Wts(i)*Wts(i));
00366       }
00367 
00368       // restore apriori information
00369       if(number_iterations > 1) {
00370          R = Rapriori;
00371          Z = Zapriori;
00372       }
00373 
00374       // update information with simple MU
00375       if(doVerbose) {
00376          cout << " Meas Cov:";
00377          for(i=0; i<M; i++) cout << " " << MeasCov(i,i);
00378          cout << endl;
00379          cout << " Partials:\n" << Partials << endl;
00380       }
00381       //if(doRobust || doWeight)
00382       //   measurementUpdate(Partials,Res,MeasCov);
00383       //else
00384       //   measurementUpdate(Partials,Res);
00385       {
00386          Matrix<double> P(Partials);
00387          Cholesky<double> Ch;
00388          if(doRobust || doWeight) {
00389             Ch(MeasCov);
00390             Matrix<double> L = inverse(Ch.L);
00391             P = L * P;
00392             Res = L * Res;
00393          }
00394 
00395          // update with whitened information
00396          SrifMU(R, Z, P, Res);
00397 
00398          // un-whiten the residuals
00399          if(doRobust || doWeight)
00400             Res = Ch.L * Res;
00401       }
00402 
00403       if(doVerbose) {
00404          cout << " Updated information matrix\n" << LabelledMatrix(names,R) << endl;
00405          cout << " Updated information vector\n" << LabelledVector(names,Z) << endl;
00406       }
00407 
00408       // invert
00409       try { getStateAndCovariance(Xsol,Cov,&small,&big); }
00410       catch(SingularMatrixException& sme) {
00411          iret = -2;
00412          break;
00413       }
00414       condition_number = big/small;
00415       if(doVerbose) {
00416          cout << " Condition number: " << scientific << condition_number
00417             << fixed << endl;
00418          cout << " Post-fit data residuals:  "
00419             << fixed << setprecision(6) << Res << endl;
00420       }
00421 
00422       // update X: when linearized, solution = dX
00423       if(doLinearize) {
00424          Xsol += NominalX;
00425       }
00426       if(doVerbose) {
00427          LabelledVector LXsol(names,Xsol);
00428          LXsol.message(" Updated X:");
00429          cout << LXsol << endl;
00430       }
00431 
00432       // linear non-robust is done..
00433       if(!doLinearize && !doRobust) break;
00434 
00435       // test for convergence of linearization
00436       if(doLinearize) {
00437          rms_convergence = RMS(Xsol - NominalX);
00438          if(doVerbose) {
00439             cout << " RMS convergence : "
00440                << scientific << rms_convergence << fixed << endl;
00441          }
00442       }
00443 
00444       // test for convergence of robust weighting, and compute new weights
00445       if(doRobust) {
00446          // must de-weight post-fit residuals
00447          LSF(Xsol,f,Partials);
00448          Res = D-f;
00449 
00450          // compute a new set of weights
00451          double mad,median;
00452          //for(mad=0.0,i=0; i<M; i++)
00453          //   mad += Wts(i)*Res(i)*Res(i);
00454          //mad = sqrt(mad)/sqrt(Robust::TuningA*(M-1));
00455          mad = Robust::MedianAbsoluteDeviation(&(Res[0]),Res.size(),median);
00456 
00457          OldWts = Wts;
00458          for(i=0; i<M; i++) {
00459             if(Res(i) < -RobustTuningT*mad)
00460                Wts(i) = -RobustTuningT*mad/Res(i);
00461             else if(Res(i) > RobustTuningT*mad)
00462                Wts(i) = RobustTuningT*mad/Res(i);
00463             else
00464                Wts(i) = 1.0;
00465          }
00466 
00467          // test for convergence
00468          rms_convergence = RMS(OldWts - Wts);
00469          if(doVerbose) cout << " Convergence: "
00470             << scientific << setprecision(3) << rms_convergence << endl;
00471       }
00472 
00473       // failures
00474       if(rms_convergence > divergenceLimit) iret=-4;
00475       if(number_iterations >= iterationsLimit) iret=-3;
00476       if(iret) {
00477          if(doSequential) {
00478             R = Rapriori;
00479             Z = Zapriori;
00480          }
00481          break;
00482       }
00483 
00484       // success
00485       if(number_iterations > 1 && rms_convergence < convergenceLimit) break;
00486 
00487       // prepare for another iteration
00488       if(doLinearize)
00489          NominalX = Xsol;
00490       if(doRobust)
00491          NominalX = X;
00492 
00493    } while(1); // end iteration loop
00494 
00495    number_batches++;
00496    if(doVerbose) cout << "Return from SRIleastSquares::leastSquaresUpdate\n\n";
00497 
00498    if(iret) return iret;
00499    valid = true;
00500 
00501    // output the solution
00502    Xsave = X = Xsol;
00503 
00504    // put residuals of fit into data vector, or weights if Robust
00505    if(doRobust) D = OldWts;
00506    else         D = Res;
00507 
00508    return iret;
00509 }
00510 catch(Exception& e) { GPSTK_RETHROW(e); }
00511 }
00512 
00513 //------------------------------------------------------------------------------------
00514 // output operator
00515 ostream& operator<<(ostream& os, const SRIleastSquares& srif)
00516 {
00517    Namelist NL(srif.names);
00518    NL += string("State");
00519    Matrix<double> A;
00520    A = srif.R || srif.Z;
00521    LabelledMatrix LM(NL,A);
00522    LM.setw(os.width());
00523    LM.setprecision(os.precision());
00524    os << LM;
00525    return os;
00526 }
00527 
00528 //------------------------------------------------------------------------------------
00529 // reset the computation, i.e. remove all stored information
00530 void SRIleastSquares::zeroAll(void)
00531 {
00532    SRI::zeroAll();
00533    Xsave = 0.0;
00534    number_batches = 0;
00535 }
00536 
00537 //------------------------------------------------------------------------------------
00538 // reset the computation, i.e. remove all stored information, and
00539 // optionally change the dimension. If N is not input, the
00540 // dimension is not changed.
00541 // @param N new SRIleastSquares dimension (optional).
00542 void SRIleastSquares::Reset(const int N) throw(Exception)
00543 {
00544    try {
00545       if(N > 0 && N != R.rows()) {
00546          R.resize(N,N,0.0);
00547          Z.resize(N,0.0);
00548       }
00549       else
00550          SRI::zeroAll(N);
00551       if(N > 0) Xsave.resize(N);
00552       Xsave = 0.0;
00553       number_batches = 0;
00554    }
00555    catch(Exception& e) { GPSTK_RETHROW(e); }
00556 }
00557 
00558 //------------------------------------------------------------------------------------
00559 } // end namespace gpstk

Generated on Wed Feb 8 03:31:03 2012 for GPS ToolKit Software Library by  doxygen 1.3.9.1