00001 #pragma ident "$Id: SRI.cpp 2293 2010-02-12 18:14:16Z btolman $"
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00049
00050
00051 #include <string>
00052 #include <vector>
00053 #include <algorithm>
00054 #include <ostream>
00055
00056 #include "SRI.hpp"
00057 #include "Namelist.hpp"
00058 #include "logstream.hpp"
00059
00060 #include "StringUtils.hpp"
00061
00062 using namespace std;
00063
00064 namespace gpstk
00065 {
00066 using namespace StringUtils;
00067
00068
00069
00070 const Matrix<double> SRINullMatrix;
00071
00072
00073
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
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
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
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
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
00138
00139
00140
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
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
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
00180
00181
00182
00183
00184
00185
00186
00187
00188
00189
00190
00191
00192
00193
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217
00218
00219
00220
00221
00222
00223
00224
00225
00226
00227
00228
00229
00230
00231
00232
00233
00234
00235
00236
00237
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
00254 Namelist N0(names);
00255 for(i=1; i<=n; i++) {
00256 for(j=1; j<=m; j++) {
00257 if(NL.labels[n-i] == N0.labels[m-j]) {
00258 N0.swap(m-i,m-j);
00259 break;
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
00269 Sleft = *this;
00270 Sleft.permute(N0);
00271
00272
00273 SRI S1(NL);
00274 S1.R = Matrix<double>(Sleft.R,m-n,m-n,n,n);
00275
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
00292
00293 SRI& SRI::operator+=(const Namelist& NL)
00294 throw(MatrixException,VectorException)
00295 {
00296 try {
00297 Namelist B(names);
00298
00299
00300 B |= NL;
00301
00302 SRI A(B);
00303
00304
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
00322
00323
00324
00325
00326
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;
00334
00335
00336 Namelist add(NL);
00337 add ^= keep;
00338 SRI Sdrop;
00339
00340 split(keep,Sdrop);
00341
00342 *this += add;
00343 this->permute(NL);
00344 }
00345 catch(MatrixException& me) {
00346 GPSTK_RETHROW(me);
00347 }
00348 catch(VectorException& ve) {
00349 GPSTK_RETHROW(ve);
00350 }
00351 }
00352
00353
00354
00355
00356
00357
00358 SRI& SRI::append(const SRI& S)
00359 throw(MatrixException,VectorException)
00360 {
00361 try {
00362
00363 if((names & S.names).size() > 0) {
00364 Exception e("Cannot append duplicate names");
00365 GPSTK_THROW(e);
00366 }
00367
00368
00369 const int I(names.size());
00370 *this += S.names;
00371
00372
00373 if(I+S.names.size() != names.size()) {
00374 Exception e("Append failed");
00375 GPSTK_THROW(e);
00376 }
00377
00378
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
00397
00398 SRI& SRI::operator+=(const SRI& S)
00399 throw(MatrixException,VectorException)
00400 {
00401 try {
00402 Namelist all(names);
00403 all |= S.names;
00404
00405
00406
00407
00408
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
00416
00417 for(j=0; j<m; j++) {
00418
00419
00420 k = all.index(names.labels[j]);
00421 if(k == -1) {
00422 MatrixException me("Algorithm error 1");
00423 GPSTK_THROW(me);
00424 }
00425
00426
00427 for(i=0; i<=j; i++) A(i,k) = R(i,j);
00428
00429 A(j,n) = Z(j);
00430 }
00431
00432
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
00443 Householder<double> HA;
00444 HA(A);
00445
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
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
00482
00483 void SRI::zeroOne(const unsigned int n)
00484 throw()
00485 {
00486 if(n >= R.rows())
00487 return;
00488
00489
00490
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
00498
00499
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
00521
00522
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
00538
00539
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
00555
00556
00557
00558
00559
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
00580 Matrix<double> Ti(T);
00581 if(&invT == &SRINullMatrix)
00582 Ti = inverseSVD(T);
00583 else
00584 Ti = invT;
00585
00586
00587 Matrix<double> B = T * R * Ti;
00588 Vector<double> Q = T * Z;
00589
00590
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
00605
00606
00607
00608
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
00621 Matrix<double> A = R * invT;
00622
00623 Householder<double> HA;
00624 HA(A);
00625 R = HA.A;
00626 }
00627
00628
00629
00630
00631
00632
00633
00634
00635
00636
00637
00638
00639
00640
00641
00642
00643
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;
00657 G(in,0) = 1.0;
00658 G = R * G;
00659 for(i=0; i<n; i++) {
00660 A(ns+i,0) = -G(i,0);
00661 for(j=0; j<n; j++)
00662 if(i<=j) A(ns+i,ns+j) = R(i,j);
00663 A(ns+i,ns+n) = Z(i);
00664 }
00665
00666
00667 Householder<double> HA;
00668 HA(A);
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
00683
00684
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
00696 for(i=0; i<in; i++) Z(i) -= R(i,in)*bias;
00697
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
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
00732
00733
00734
00735
00736
00737
00738
00739
00740
00741
00742
00743
00744
00745
00746
00747
00748
00749
00750
00751
00752
00753
00754
00755
00756
00757
00758
00759
00760
00761
00762 unsigned int i,j,k;
00763
00764 vector<int> indx;
00765 vector<double> value;
00766 for(i=0; i<dropNL.size(); i++) {
00767 int in = names.index(dropNL.getName(i));
00768 if(in > -1) {
00769 indx.push_back(in);
00770 value.push_back(values_in(i));
00771 }
00772
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
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
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;
00794
00795 Ztmp(k) = Z(i);
00796 for(j=i; j<n; j++) Rtmp(k,j) = R(i,j);
00797 k++;
00798 }
00799
00800
00801 Z = Ztmp;
00802
00803
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;
00809
00810 for(i=0; i<=j; i++) R(i,k) = Rtmp(i,j);
00811 k++;
00812 }
00813
00814
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
00830
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));
00873 Vector<double> apZ(apR*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);
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
00903
00904
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--) {
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++)
00922 sum -= R(i,j)*X(j);
00923 X(i) = sum/R(i,i);
00924 }
00925 }
00926
00927
00928
00929
00930
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
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
00967
00968
00969 os << LM;
00970
00971 return os;
00972 }
00973
00974 }
00975
00976