SRI.cpp

Go to the documentation of this file.
00001 #pragma ident "$Id: SRI.cpp 2293 2010-02-12 18:14:16Z btolman $"
00002 
00003 //============================================================================
00004 //
00005 //  This file is part of GPSTk, the GPS Toolkit.
00006 //
00007 //  The GPSTk is free software; you can redistribute it and/or modify
00008 //  it under the terms of the GNU Lesser General Public License as published
00009 //  by the Free Software Foundation; either version 2.1 of the License, or
00010 //  any later version.
00011 //
00012 //  The GPSTk is distributed in the hope that it will be useful,
00013 //  but WITHOUT ANY WARRANTY; without even the implied warranty of
00014 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00015 //  GNU Lesser General Public License for more details.
00016 //
00017 //  You should have received a copy of the GNU Lesser General Public
00018 //  License along with GPSTk; if not, write to the Free Software Foundation,
00019 //  Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00020 //  
00021 //  Copyright 2004, The University of Texas at Austin
00022 //
00023 //============================================================================
00024 
00025 //============================================================================
00026 //
00027 //This software developed by Applied Research Laboratories at the University of
00028 //Texas at Austin, under contract to an agency or agencies within the U.S. 
00029 //Department of Defense. The U.S. Government retains all rights to use,
00030 //duplicate, distribute, disclose, or release this software. 
00031 //
00032 //Pursuant to DoD Directive 523024 
00033 //
00034 // DISTRIBUTION STATEMENT A: This software has been approved for public 
00035 //                           release, distribution is unlimited.
00036 //
00037 //=============================================================================
00038 
00049 // -----------------------------------------------------------------------------------
00050 // system
00051 #include <string>
00052 #include <vector>
00053 #include <algorithm>
00054 #include <ostream>
00055 // geomatics
00056 #include "SRI.hpp"
00057 #include "Namelist.hpp"
00058 #include "logstream.hpp"
00059 // GPSTk
00060 #include "StringUtils.hpp"
00061 
00062 using namespace std;
00063 
00064 namespace gpstk
00065 {
00066 using namespace StringUtils;
00067 
00068    // --------------------------------------------------------------------------------
00069    // used to mark optional input
00070    const Matrix<double> SRINullMatrix;
00071 
00072    //---------------------------------------------------------------------------------
00073    // constructor given the dimension N.
00074    SRI::SRI(const unsigned int N)
00075       throw()
00076    {
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    SRI::SRI(const Namelist& nl)
00085       throw()
00086    {
00087       if(nl.size() <= 0) return;
00088       R = Matrix<double>(nl.size(),nl.size(),0.0);
00089       Z = Vector<double>(nl.size(),0.0);
00090       names = nl;
00091    }
00092 
00093    // --------------------------------------------------------------------------------
00094    // explicit constructor - throw if the dimensions are inconsistent.
00095    SRI::SRI(const Matrix<double>& r,
00096             const Vector<double>& z,
00097             const Namelist& nl)
00098       throw(MatrixException)
00099    {
00100       if(r.rows() != r.cols() || r.rows() != z.size() || r.rows() != nl.size()) {
00101          MatrixException me("Invalid dimensions in explicit SRI constructor:\n R is "
00102                + asString<int>(r.rows()) + "x"
00103                + asString<int>(r.cols()) + ", Z has length "
00104                + asString<int>(z.size()) + " and NL has length "
00105                + asString<int>(nl.size())
00106                );
00107          GPSTK_THROW(me);
00108       }
00109       if(r.rows() <= 0) return;
00110       R = r;
00111       Z = z;
00112       names = nl;
00113    }
00114 
00115    // --------------------------------------------------------------------------------
00116    // copy constructor
00117    SRI::SRI(const SRI& s)
00118       throw()
00119    {
00120       R = s.R;
00121       Z = s.Z;
00122       names = s.names;
00123    }
00124 
00125    // --------------------------------------------------------------------------------
00126    // operator=
00127    SRI& SRI::operator=(const SRI& right)
00128       throw()
00129    {
00130       R = right.R;
00131       Z = right.Z;
00132       names = right.names;
00133       return *this;
00134    }
00135 
00136    // ---------------------------------------------------------------------------
00137    // modify SRIs
00138    // --------------------------------------------------------------------------------
00139    // Permute the SRI elements to match the input Namelist, which may differ with
00140    // the SRI Namelist by AT MOST A PERMUTATION, throw if this is not true.
00141    void SRI::permute(const Namelist& nl)
00142       throw(MatrixException,VectorException)
00143    {
00144       if(identical(names,nl)) return;
00145       if(names != nl) {
00146          MatrixException me("Invalid input: Namelists must be == to permute");
00147          GPSTK_THROW(me);
00148       }
00149 
00150       try {
00151          unsigned int i,j;
00152          // build a permutation matrix
00153          Matrix<double> P(R.rows(),R.rows(),0.0);
00154          for(i=0; i<R.rows(); i++) {
00155             j = nl.index(names.getName(i));
00156             P(j,i) = 1;
00157          }
00158 
00159          Matrix<double> B;
00160          Vector<double> Q;
00161          B = P * R * transpose(P);
00162          Q = P * Z;
00163 
00164          // re-triangularize
00165          R = 0.0;
00166          Z = 0.0;
00167          SrifMU(R,Z,B,Q);
00168          names = nl;
00169       }
00170       catch(MatrixException& me) {
00171          GPSTK_RETHROW(me);
00172       }
00173       catch(VectorException& ve) {
00174          GPSTK_RETHROW(ve);
00175       }
00176    }
00177 
00178    // --------------------------------------------------------------------------------
00179    // Split this SRI (call it S) into two others, S1 and Sleft, where S1 has
00180    // a Namelist identical to the input Namelist (NL); set *this = S1 at the
00181    // end. NL must be a non-empty subset of names, and (names ^ NL) also must
00182    // be non-empty; throw MatrixException if this is not true. The second
00183    // output SRI, Sleft, will have the same names as S, but perhaps permuted.
00184    //
00185    // The routine works by first permuting S so that its Namelist if of the
00186    // form {N2,NL}, where N2 = (names ^ NL); this is possible only if NL is
00187    // a non-trivial subset of names. Then, the rows of S (rows of R and elements
00188    // of Z) naturally separate into the two component SRIs, with zeros in the
00189    // elements of the first SRI which correspond to N2, and those in Sleft
00190    // which correspond to NL.
00191    //
00192    //    Example:    S.name = A B C D E F G and NL = D E F G.
00193    // (Obviously, S may be permuted into such an order whenever this is needed.)
00194    // Note that here the R,Z pair is written in a format reminiscent of the
00195    // set of equations implied by R*X=Z, i.e. 1A+2B+3C+4D+5E+6F+7G=a, etc.
00196    //
00197    //          S (R Z)       =         S1            +         Sleft
00198    // with    names                       NL                  names
00199    //     A B C D E F G           . . . D E F G           A B C D E F G   
00200    //     - - - - - - -  -        - - - - - - -  -        - - - - - - -  -
00201    //     1 2 3 4 5 6 7  a   =    . . . . . . .  .   +    1 2 3 4 5 6 7  a
00202    //       8 9 1 2 3 4  b          . . . . . .  .          8 9 1 2 3 4  b
00203    //         5 6 7 8 9  c            . . . . .  .            5 6 7 8 9  c
00204    //           1 2 3 4  d              1 2 3 4  d              . . . .  d
00205    //             5 6 7  e                5 6 7  e                . . .  e
00206    //               8 9  f                  8 9  f                  . .  f
00207    //                 1  g                    1  g                    .  g
00208    //
00209    // where "." denotes a zero.  The split is simply separating the linear
00210    // equations which make up R*X=Z into two groups; because of the ordering,
00211    // one of the groups of equations (S1) depends only on a particular subset
00212    // of the elements of the state vector, i.e. the elements labelled by the
00213    // Namelist NL.
00214    //
00215    // The equation shown here is an information equation; if the two SRIs S1
00216    // and Sleft were merged again, none of the information would be lost.
00217    // Note that S1 has no dependence on A B C (hence the .'s), and therefore
00218    // its size can be reduced. However S2 still depends on the full names
00219    // Namelist. Sleft is necessarily singular, but S1 is not.
00220    //
00221    // Note that the SRI contains information about both the solution and
00222    // the covariance, i.e. state and noise, and therefore one must be very careful
00223    // in interpreting the results of split and merge (operator+=). [Be especially
00224    // careful about the idea that a merge might be reversible with a split() or
00225    // vice-versa - strictly this is never possible unless the Namelists are
00226    // mutually exclusive - two separate problems.]
00227    //
00228    // For example, suppose two different SRI's, which have some elements in common,
00229    // are merged. The combined SRI will have more information (it can't have less)
00230    // about the common elements, and therefore the solution will be 'better'
00231    // (assuming the underlying model equations for those elements are identical).
00232    // However the noises will also be combined, and the results you get might be
00233    // surprising. Also, note that if you then split the combined SRI again, the
00234    // solution won't change but the noises will be very different; in particular
00235    // the new split part will take all the information with it, so the common states
00236    // will have lower noise than they did in the original SRI.
00237    // See the test program tsri.cpp
00238    //
00239    void SRI::split(const Namelist& NL, SRI& Sleft)
00240       throw(MatrixException,VectorException)
00241    {
00242       try {
00243          Sleft = SRI(0);
00244          unsigned int n,m;
00245          n = NL.size();
00246          m = names.size();
00247          if(n >= m) {
00248             MatrixException me("split: Input Namelist must be a subset of this one");
00249             GPSTK_THROW(me);
00250          }
00251 
00252          unsigned int i,j;
00253             // copy names and permute it so that its end matches NL 
00254          Namelist N0(names);
00255          for(i=1; i<=n; i++) {           // loop (backwards) over names in NL
00256             for(j=1; j<=m; j++) {        // search (backwards) in NO for a match
00257                if(NL.labels[n-i] == N0.labels[m-j]) {  // if found a match
00258                   N0.swap(m-i,m-j);      // then move matching name to end
00259                   break;                 // and go on to next name in NL
00260                }
00261             }
00262             if(j > m) {
00263                MatrixException me("split: Input Namelist is not non-trivial subset");
00264                GPSTK_THROW(me);
00265             }
00266          }
00267 
00268             // copy *this into Sleft, then do the permutation
00269          Sleft = *this;
00270          Sleft.permute(N0);
00271 
00272             // copy parts of Sleft into S1, and then zero out those parts of Sleft
00273          SRI S1(NL);
00274          S1.R = Matrix<double>(Sleft.R,m-n,m-n,n,n);
00275          //S1.Z = Vector<double>(Sleft.Z,m-n,n);
00276          S1.Z.resize(n);
00277          for(i=0; i<n; i++) S1.Z(i) = Sleft.Z(m-n+i);
00278          for(i=m-n; i<m; i++) Sleft.zeroOne(i);
00279 
00280          *this = S1;
00281       }
00282       catch(MatrixException& me) {
00283          GPSTK_RETHROW(me);
00284       }
00285       catch(VectorException& ve) {
00286          GPSTK_RETHROW(ve);
00287       }
00288    }
00289 
00290    // --------------------------------------------------------------------------------
00291    // extend this SRI to include the given Namelist, with no added information;
00292    // names in the input namelist which are not unique are ignored.
00293    SRI& SRI::operator+=(const Namelist& NL)
00294       throw(MatrixException,VectorException)
00295    {
00296       try {
00297          Namelist B(names);
00298             // NB assume that Namelist::operator|=() adds at the _end_
00299             // NB if there are duplicate names, |= will not add them
00300          B |= NL;
00301             // NB assume that this zeros A.R and A.Z
00302          SRI A(B);
00303             // should do this with slices..
00304             // copy into the new SRI
00305          for(unsigned int i=0; i<R.rows(); i++) {
00306             A.Z(i) = Z(i);
00307             for(unsigned int j=0; j<R.cols(); j++) A.R(i,j) = R(i,j);
00308          }
00309          *this = A;
00310          return *this;
00311       }
00312       catch(MatrixException& me) {
00313          GPSTK_RETHROW(me);
00314       }
00315       catch(VectorException& ve) {
00316          GPSTK_RETHROW(ve);
00317       }
00318    }
00319 
00320    // --------------------------------------------------------------------------------
00321    // reshape this SRI to match the input Namelist, by calling other member
00322    // functions, including split(), operator+() and permute()
00323    // Given this SRI and a new Namelist NL, if NL does not match names,
00324    // transform names to match it, using (1) drop elements (this is probably
00325    // optional - you can always keep 'dead' elements), (2) add new elements
00326    // (with zero information), and (3) permute to match NL.
00327    void SRI::reshape(const Namelist& NL)
00328       throw(MatrixException,VectorException)
00329    {
00330       try {
00331          if(names == NL) return;
00332          Namelist keep(names);
00333          keep &= NL;                // keep only those in both names and NL
00334          //Namelist drop(names);    // (drop is unneeded - split would do it)
00335          //drop ^= keep;            // lose those in names but not in keep
00336          Namelist add(NL);
00337          add ^= keep;               // add those in NL but not in keep
00338          SRI Sdrop;                 // make a new SRI to hold the losers
00339          // would like to allow caller access to Sdrop..
00340          split(keep,Sdrop);         // split off only the losers
00341                                     // NB names = drop | keep; drop & keep is empty
00342          *this += add;              // add the new ones
00343          this->permute(NL);         // permute it to match NL
00344       }
00345       catch(MatrixException& me) {
00346          GPSTK_RETHROW(me);
00347       }
00348       catch(VectorException& ve) {
00349          GPSTK_RETHROW(ve);
00350       }
00351    }
00352 
00353    // --------------------------------------------------------------------------------
00354    // append an SRI onto this SRI. Similar to opertor+= but simpler; input SRI is
00355    // simply appended, first using operator+=(Namelist), then filling the new portions
00356    // of R and Z, all without final Householder transformation of result.
00357    // Do not allow a name that is already present to be added: throw.
00358    SRI& SRI::append(const SRI& S)
00359       throw(MatrixException,VectorException)
00360    {
00361       try {
00362          // do not allow duplicates
00363          if((names & S.names).size() > 0) {
00364             Exception e("Cannot append duplicate names");
00365             GPSTK_THROW(e);
00366          }
00367 
00368          // append to names at the end, and to R Z, zero filling
00369          const int I(names.size());
00370          *this += S.names;
00371 
00372          // just in case...to avoid overflow in loop below
00373          if(I+S.names.size() != names.size()) {
00374             Exception e("Append failed");
00375             GPSTK_THROW(e);
00376          }
00377 
00378          // loop over new names, copying data from input into the new SRI
00379          for(int i=0; i<S.names.size(); i++) {
00380             Z(I+i) = S.Z(i);
00381             for(int j=0; j<S.names.size(); j++)
00382                R(I+i,I+j) = S.R(i,j);
00383          }
00384 
00385          return *this;
00386       }
00387       catch(MatrixException& me) {
00388          GPSTK_RETHROW(me);
00389       }
00390       catch(VectorException& ve) {
00391          GPSTK_RETHROW(ve);
00392       }
00393    }
00394 
00395    // --------------------------------------------------------------------------------
00396    // merge this SRI with the given input SRI. ? should this be operator&=() ?
00397    // NB may reorder the names in the resulting Namelist.
00398    SRI& SRI::operator+=(const SRI& S)
00399       throw(MatrixException,VectorException)
00400    {
00401       try {
00402          Namelist all(names);
00403          all |= S.names;      // assumes Namelist::op|= adds unique S.names to _end_
00404 
00405          //all.sort();        // TEMP - for testing with old version
00406 
00407             // stack the (R|Z)'s from both in one matrix;
00408             // all determines the columns, plus last column is for Z
00409          unsigned int i,j,k,n,m,sm;
00410          n = all.labels.size();
00411          m = R.rows();
00412          sm = S.R.rows();
00413          Matrix<double> A(m+sm,n+1,0.0);
00414 
00415             // copy R into A, permuting columns as names differs from all
00416             // loop over columns of R; do Z at the same time using j=row
00417          for(j=0; j<m; j++) {
00418                // find where this column of R goes in A
00419                // (should never throw..)
00420             k = all.index(names.labels[j]);
00421             if(k == -1) {
00422                MatrixException me("Algorithm error 1");
00423                GPSTK_THROW(me);
00424             }
00425 
00426                // copy this col of R into A (R is UT)
00427             for(i=0; i<=j; i++) A(i,k) = R(i,j);
00428                // also the jth element of Z
00429             A(j,n) = Z(j);
00430          }
00431 
00432             // now do the same for S, but put S.R|S.Z below R|Z
00433          for(j=0; j<sm; j++) {
00434             k = all.index(S.names.labels[j]);
00435             if(k == -1) {
00436                MatrixException me("Algorithm error 2");
00437                GPSTK_THROW(me);
00438             }
00439             for(i=0; i<=j; i++) A(m+i,k) = S.R(i,j);
00440             A(m+j,n) = S.Z(j);
00441          }
00442             // now triangularize A and pull out the new R and Z
00443          Householder<double> HA;
00444          HA(A);
00445          // submatrix args are matrix,toprow,topcol,numrows,numcols
00446          R = Matrix<double>(HA.A,0,0,n,n);
00447          Z = Vector<double>(HA.A.colCopy(n));
00448          Z.resize(n);
00449          names = all;
00450 
00451          return *this;
00452       }
00453       catch(MatrixException& me) {
00454          GPSTK_RETHROW(me);
00455       }
00456       catch(VectorException& ve) {
00457          GPSTK_RETHROW(ve);
00458       }
00459    }
00460 
00461    // --------------------------------------------------------------------------------
00462    // merge two SRIs to produce a third. ? should this be operator&() ?
00463    SRI operator+(const SRI& Sleft,
00464                  const SRI& Sright)
00465       throw(MatrixException,VectorException)
00466    {
00467       try {
00468          SRI S(Sleft);
00469          S += Sright;
00470          return S;
00471       }
00472       catch(MatrixException& me) {
00473          GPSTK_RETHROW(me);
00474       }
00475       catch(VectorException& ve) {
00476          GPSTK_RETHROW(ve);
00477       }
00478    }
00479 
00480    // --------------------------------------------------------------------------------
00481    // Zero out the nth row of R and the nth element of Z, removing all
00482    // information about that element.
00483    void SRI::zeroOne(const unsigned int n)
00484       throw()
00485    {
00486       if(n >= R.rows())
00487          return;
00488 
00489       //TD this is not right -- you should permute the element
00490       //to the first row, then zero
00491       for(unsigned int j=n; j<R.cols(); j++) 
00492          R(n,j) = 0.0;
00493       Z(n) = 0.0;
00494    }
00495 
00496    // --------------------------------------------------------------------------------
00497    // Zero out all the first n rows of R and elements of Z, removing all
00498    // information about those elements. Default value of the input is 0,
00499    // meaning zero out the entire SRI.
00500    void SRI::zeroAll(const int n)
00501       throw()
00502    {
00503       if(n <= 0) {
00504          R = 0.0;
00505          Z = 0.0;
00506          return;
00507       }
00508 
00509       if(n >= R.rows())
00510          return;
00511 
00512       for(unsigned int i=0; i<n; i++) {
00513          for(unsigned int j=i; j<R.cols(); j++) 
00514             R(i,j) = 0.0;
00515          Z(i) = 0.0;
00516       }
00517    }
00518 
00519    // --------------------------------------------------------------------------------
00520    // Shift the state vector by a constant vector X0; does not change information
00521    // i.e. let R * X = Z => R * (X-X0) = Z'
00522    // throw on invalid input dimension
00523    void SRI::shift(const Vector<double>& X0)
00524       throw(MatrixException)
00525    {
00526       if(X0.size() != R.cols()) {
00527          MatrixException me("Invalid input dimension: SRI has dimension "
00528                + asString<int>(R.rows()) + " while input has length "
00529                + asString<int>(X0.size())
00530                );
00531          GPSTK_THROW(me);
00532       }
00533       Z = Z - R * X0;
00534    }
00535 
00536    // --------------------------------------------------------------------------------
00537    // Shift the SRI state vector (Z) by a constant vector Z0;
00538    // does not change information. i.e. let Z => Z-Z0
00539    // throw on invalid input dimension
00540    void SRI::shiftZ(const Vector<double>& Z0)
00541       throw(MatrixException)
00542    {
00543       if(Z0.size() != R.cols()) {
00544          MatrixException me("Invalid input dimension: SRI has dimension "
00545                + asString<int>(R.rows()) + " while input has length "
00546                + asString<int>(Z0.size())
00547                );
00548          GPSTK_THROW(me);
00549       }
00550       Z = Z - Z0;
00551    }
00552 
00553    // --------------------------------------------------------------------------------
00554    // Transform this SRI with the transformation matrix T;
00555    // i.e. R -> T * R * inverse(T) and Z -> T * Z. The matrix inverse(T)
00556    // may optionally be supplied as input, otherwise it is computed from
00557    // T. NB names in this SRI are most likely changed; but this routine does
00558    // not change the Namelist. Throw MatrixException if the input has
00559    // the wrong dimension or cannot be inverted.
00560    void SRI::transform(const Matrix<double>& T,
00561                        const Matrix<double>& invT)
00562       throw(MatrixException,VectorException)
00563    {
00564       if(T.rows() != R.rows() ||
00565          T.cols() != R.cols() ||
00566          (&invT != &SRINullMatrix && (invT.rows() != R.rows() ||
00567          invT.cols() != R.cols()))) {
00568             MatrixException me("Invalid input dimension:\n  SRI has dimension "
00569                + asString<int>(R.rows()) + " while T has dimension "
00570                + asString<int>(T.rows()) + "x"
00571                + asString<int>(T.cols()));
00572             if(&invT != &SRINullMatrix) me.addText("\n  and invT has dimension "
00573                            + asString<int>(invT.rows()) + "x"
00574                            + asString<int>(invT.cols()));
00575             GPSTK_THROW(me);
00576       }
00577 
00578       try {
00579             // get the inverse matrix
00580          Matrix<double> Ti(T);
00581          if(&invT == &SRINullMatrix)
00582             Ti = inverseSVD(T);
00583          else
00584             Ti = invT;
00585 
00586             // transform
00587          Matrix<double> B = T * R * Ti;
00588          Vector<double> Q = T * Z;
00589 
00590             // re-triangularize
00591          R = 0.0;
00592          Z = 0.0;
00593          SrifMU(R,Z,B,Q);
00594       }
00595       catch(MatrixException& me) {
00596          GPSTK_RETHROW(me);
00597       }
00598       catch(VectorException& ve) {
00599          GPSTK_RETHROW(ve);
00600       }
00601    }
00602 
00603    // --------------------------------------------------------------------------------
00604    // Transform the state by the transformation matrix T; i.e. X -> T*X,
00605    // without transforming the SRI; this is done by right multiplying R by
00606    // inverse(T), which is the input. Thus R -> R*inverse(T),
00607    // so R*inverse(T)*T*X = Z.  Input is the _inverse_ of the transformation.
00608    // throw MatrixException if input dimensions are wrong.
00609    void SRI::transformState(const Matrix<double>& invT)
00610       throw(MatrixException)
00611    {
00612       if(invT.rows() != R.rows() || invT.cols() != R.rows()) {
00613          MatrixException me("Invalid input dimension: SRI has dimension "
00614             + asString<int>(R.rows()) + " while invT has dimension "
00615             + asString<int>(invT.rows()) + "x"
00616             + asString<int>(invT.cols()));
00617          GPSTK_THROW(me);
00618       }
00619 
00620          // transform
00621       Matrix<double> A = R * invT;
00622          // re-triangularize
00623       Householder<double> HA;
00624       HA(A);
00625       R = HA.A;
00626    }
00627 
00628    // --------------------------------------------------------------------------------
00629    // Decrease the information in this SRI for, or 'Q bump', the element
00630    // with the input index.  This means that the uncertainty and the state
00631    // element given by the index are divided by the input factor q; the
00632    // default input is zero, which means zero out the information (q = infinite).
00633    // A Q bump by factor q is equivalent to 'de-weighting' the element by q.
00634    // No effect if input index is out of range.
00635    //
00636    // Use a specialized form of the time update, with Phi=unity, G(N x 1) = 0
00637    // except 1 for the element (in) getting bumped, and Rw(1 x 1) = 1 / q.
00638    // Note that this bump of the covariance for element k results in
00639    // Cov(k,k) += q (plus, not times!).
00640    // if q is 0, replace q with 1/q, i.e. lose all information, covariance
00641    // goes singular; this is equivalent to (1) permute so that the 'in'
00642    // element is first, (2) zero out the first row of R and the first element
00643    // of Z, (3) permute the first row back to in.
00644    void SRI::Qbump(const unsigned int& in,
00645                    const double& q)
00646       throw(MatrixException,VectorException)
00647    {
00648       try {
00649          if(in >= R.rows()) return;
00650          double factor=0.0;
00651          if(q != 0.0) factor=1.0/q;
00652 
00653          unsigned int ns=1,i,j,n=R.rows();
00654 
00655          Matrix<double> A(n+ns,n+ns+1,0.0), G(n,ns,0.0);
00656          A(0,0) = factor;           // Rw, dimension ns x ns = 1 x 1
00657          G(in,0) = 1.0;
00658          G = R * G;                 // R*Phi*G
00659          for(i=0; i<n; i++) {
00660             A(ns+i,0) = -G(i,0);               //     A =   Rw       0       zw=0
00661             for(j=0; j<n; j++)                 //          -R*Phi*G  R*Phi   Z
00662                if(i<=j) A(ns+i,ns+j) = R(i,j); //
00663             A(ns+i,ns+n) = Z(i);
00664          }
00665 
00666             // triangularize and pull out the new R and Z
00667          Householder<double> HA;                //    A  =  Rw  Rwx  zw
00668          HA(A);                                 //          0    R   z
00669          R = Matrix<double>(HA.A,ns,ns,n,n);
00670          Vector<double> T=HA.A.colCopy(ns+n);
00671          Z = Vector<double>(T,ns,n);
00672       }
00673       catch(MatrixException& me) {
00674          GPSTK_RETHROW(me);
00675       }
00676       catch(VectorException& ve) {
00677          GPSTK_RETHROW(ve);
00678       }
00679    }
00680 
00681    // --------------------------------------------------------------------------------
00682    // Fix the state element with the input index to the input value, and
00683    // collapse the SRI by removing that element.
00684    // No effect if index is out of range.
00685    void SRI::stateFix(const unsigned int& in,
00686                       const double& bias)
00687       throw(MatrixException,VectorException)
00688    {
00689       if(in >= R.rows()) return;
00690 
00691       try {
00692          unsigned int i,j,ii,jj,n=R.rows();
00693          Vector<double> Znew(n-1,0.0);
00694          Matrix<double> Rnew(n-1,n-1,0.0);
00695             // move the X(in) terms to the data vector on the RHS
00696          for(i=0; i<in; i++) Z(i) -= R(i,in)*bias;
00697             // remove row/col in and collapse
00698          for(i=0,ii=0; i<n; i++) {
00699             if(i == in) continue;
00700             Znew(ii) = Z(i);
00701             for(j=i,jj=ii; j<n; j++) {
00702                if(j == in) continue;
00703                Rnew(ii,jj) = R(i,j);
00704                jj++;
00705             }
00706             ii++;
00707          }
00708          R = Rnew;
00709          Z = Znew;
00710          names -= names.labels[in];
00711       }
00712       catch(MatrixException& me) {
00713          GPSTK_RETHROW(me);
00714       }
00715       catch(VectorException& ve) {
00716          GPSTK_RETHROW(ve);
00717       }
00718    }
00719 
00720    // --------------------------------------------------------------------------------
00721    // Vector version of stateFix with several states given in a Namelist.
00722    void SRI::stateFix(const Namelist& dropNL, const Vector<double>& values_in)
00723       throw(MatrixException,VectorException)
00724    {
00725       try {
00726          if(dropNL.size() != values_in.size()) {
00727             VectorException e("Input has inconsistent lengths");
00728             GPSTK_THROW(e);
00729          }
00730 /*
00731          // build a vector of indexes to keep
00732          int i,j;
00733          vector<int> indexes;
00734          for(i=0; i<names.size(); i++) {
00735             j = dropNL.index(names.getName(i)); // index in dropNL, this state
00736             if(j == -1) indexes.push_back(i);// not found in dropNL, so keep
00737          }
00738 
00739          const int n=indexes.size();         // new dimension
00740          if(n == 0) {
00741             Exception e("Cannot drop all states");
00742             GPSTK_THROW(e);
00743          }
00744 
00745          Vector<double> X,newX(n);
00746          Matrix<double> C,newC(n,n);
00747          Namelist newNL;
00748 
00749          double big,small;
00750          getStateAndCovariance(X,C,&small,&big);
00751 
00752          for(i=0; i<n; i++) {
00753             newX(i) = X(indexes[i]);
00754             for(j=0; j<n; j++) newC(i,j) = C(indexes[i],indexes[j]);
00755             newNL += names.getName(indexes[i]);
00756          }
00757 
00758          R = Matrix<double>(inverseUT(upperCholesky(newC)));
00759          Z = Vector<double>(R*newX);
00760          names = newNL;
00761 */
00762          unsigned int i,j,k;
00763             // create a vector of indexes and corresponding values
00764          vector<int> indx;
00765          vector<double> value;
00766          for(i=0; i<dropNL.size(); i++) {
00767             int in = names.index(dropNL.getName(i));   // in must be allowed to be -1
00768             if(in > -1) {
00769                indx.push_back(in);
00770                value.push_back(values_in(i));
00771             }
00772             //else nothing happens
00773          }
00774          const unsigned int m = indx.size();
00775          const unsigned int n = R.rows();
00776          if(m == 0) return;
00777          if(m == n) {
00778             *this = SRI(0);
00779             return;
00780          }
00781             // move the X(in) terms to the data vector on the RHS
00782          for(k=0; k<m; k++)
00783             for(i=0; i<indx[k]; i++)
00784                Z(i) -= R(i,indx[k])*value[k];
00785 
00786             // first remove the rows in indx
00787          bool skip;
00788          Vector<double> Ztmp(n-m,0.0);
00789          Matrix<double> Rtmp(n-m,n,0.0);
00790          for(i=0,k=0; i<n; i++) {
00791             skip = false;
00792             for(j=0; j<m; j++) if(i == indx[j]) { skip=true; break; }
00793             if(skip) continue;      // skip row to be dropped
00794 
00795             Ztmp(k) = Z(i);
00796             for(j=i; j<n; j++) Rtmp(k,j) = R(i,j);
00797             k++;
00798          }
00799 
00800             // Z is now done
00801          Z = Ztmp;
00802 
00803             // now remove columns in indx
00804          R = Matrix<double>(n-m,n-m,0.0);
00805          for(j=0,k=0; j<n; j++) {
00806             skip = false;
00807             for(i=0; i<m; i++) if(j == indx[i]) { skip=true; break; }
00808             if(skip) continue;      // skip col to be dropped
00809 
00810             for(i=0; i<=j; i++) R(i,k) = Rtmp(i,j);
00811             k++;
00812          }
00813 
00814             // remove the names
00815          for(k=0; k<dropNL.size(); k++) {
00816             std::string name(dropNL.getName(k));
00817             names -= name;
00818          }
00819       }
00820       catch(MatrixException& me) {
00821          GPSTK_RETHROW(me);
00822       }
00823       catch(VectorException& ve) {
00824          GPSTK_RETHROW(ve);
00825       }
00826    }
00827 
00828    //---------------------------------------------------------------------------------
00829    // Add a priori or 'constraint' information
00830    // Prefer addAPrioriInformation(inverse(Cov), inverse(Cov)*X);
00831    void SRI::addAPriori(const Matrix<double>& Cov, const Vector<double>& X)
00832       throw(MatrixException)
00833    {
00834       if(Cov.rows() != Cov.cols() || Cov.rows() != R.rows() || X.size() != R.rows()) {
00835          MatrixException me("Invalid input dimensions:\n  SRI has dimension "
00836             + asString<int>(R.rows()) + ",\n  while input is Cov("
00837             + asString<int>(Cov.rows()) + "x"
00838             + asString<int>(Cov.cols()) + ") and X("
00839             + asString<int>(X.size()) + ")."
00840             );
00841          GPSTK_THROW(me);
00842       }
00843 
00844       try {
00845          Matrix<double> InvCov = inverse(Cov);
00846          addAPrioriInformation(InvCov, X);
00847       }
00848       catch(MatrixException& me) {
00849          GPSTK_THROW(me);
00850       }
00851    }
00852 
00853    // --------------------------------------------------------------------------------
00854    void SRI::addAPrioriInformation(const Matrix<double>& InvCov,
00855                                    const Vector<double>& X)
00856       throw(MatrixException)
00857    {
00858       if(InvCov.rows() != InvCov.cols() || InvCov.rows() != R.rows()
00859             || X.size() != R.rows()) {
00860          MatrixException me("Invalid input dimensions:\n  SRI has dimension "
00861             + asString<int>(R.rows()) + ",\n  while input is InvCov("
00862             + asString<int>(InvCov.rows()) + "x"
00863             + asString<int>(InvCov.cols()) + ") and X("
00864             + asString<int>(X.size()) + ")."
00865             );
00866          GPSTK_THROW(me);
00867       }
00868 
00869       try {
00870          Cholesky<double> Ch;
00871          Ch(InvCov);
00872          Matrix<double> apR(transpose(Ch.L));  // R = UT(inv(Cov))
00873          Vector<double> apZ(apR*X);            // Z = R*X
00874          SrifMU(R, Z, apR, apZ);
00875       }
00876       catch(MatrixException& me) {
00877          GPSTK_THROW(me);
00878       }
00879    }
00880 
00881    // --------------------------------------------------------------------------------
00882    void SRI::getConditionNumber(double& small, double& big)
00883       throw(MatrixException)
00884    {
00885       try {
00886          small = big = 0.0;
00887          const int n=R.rows();
00888          if(n == 0) return;
00889          SVD<double> svd;
00890          svd(R);
00891          svd.sort(true);   // now the last s.v. is the smallest
00892          small = svd.S(n-1);
00893          big = svd.S(0);
00894       }
00895       catch(MatrixException& me) {
00896          me.addText("Called by getConditionNumber");
00897          GPSTK_RETHROW(me);
00898       }
00899    }
00900 
00901    // --------------------------------------------------------------------------------
00902    // Compute state without computing covariance. Use the fact that R is upper
00903    // triangular. Throw if and when a zero diagonal element is found; values at larger
00904    // index are still valid. On output *ptr is the largest singular index
00905    void SRI::getState(Vector<double>& X, int *ptr)
00906       throw(MatrixException)
00907    {
00908       const int n=Z.size();
00909       X = Vector<double>(n,0.0);
00910       if(ptr) *ptr = -1;
00911       if(n == 0) return;
00912       int i,j;
00913       for(i=n-1; i>=0; i--) {             // loop over rows, in reverse order
00914          if(R(i,i) == 0.0) {
00915             if(ptr) *ptr = i;
00916             MatrixException me("Singular matrix; zero diagonal element at index "
00917                + asString<int>(i));
00918             GPSTK_THROW(me);
00919          }
00920          double sum=Z(i);
00921          for(j=i+1; j<n; j++)             // sum over columns to right of diagonal
00922             sum -= R(i,j)*X(j);
00923          X(i) = sum/R(i,i);
00924       }
00925    }
00926 
00927    // --------------------------------------------------------------------------------
00928    // get the state X and the covariance matrix C of the state, where
00929    // C = transpose(inverse(R))*inverse(R) and X = inverse(R) * Z.
00930    // Throws MatrixException if R is singular.
00931    void SRI::getStateAndCovariance(Vector<double>& X,
00932                                    Matrix<double>& C,
00933                                    double *ptrSmall,
00934                                    double *ptrBig)
00935       throw(MatrixException,VectorException)
00936    {
00937       try {
00938          Matrix<double> invR;
00939          invR = inverseUT(R,ptrSmall,ptrBig);
00940          C = UTtimesTranspose(invR);
00941          X = invR * Z;
00942       }
00943       catch(MatrixException& me) {
00944          GPSTK_RETHROW(me);
00945       }
00946       catch(VectorException& ve) {
00947          GPSTK_RETHROW(ve);
00948       }
00949    }
00950 
00951    //---------------------------------------------------------------------------------
00952    // output operator
00953    ostream& operator<<(ostream& os, const SRI& S)
00954    {
00955       Namelist NLR(S.names);
00956       Namelist NLC(S.names);
00957       NLC += string("State");
00958       Matrix<double> A;
00959       A = S.R || S.Z;
00960       LabelledMatrix LM(NLR,NLC,A);
00961 
00962       ios_base::fmtflags flags = os.flags();
00963       if(flags & ios_base::scientific) LM.scientific();
00964       LM.setw(os.width());
00965       LM.setprecision(os.precision());
00966       //LM.message("NL");
00967       //LM.linetag("tag");
00968 
00969       os << LM;
00970 
00971       return os;
00972    }
00973 
00974 } // end namespace gpstk
00975 
00976 //------------------------------------------------------------------------------------

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