core/vnl/algo/vnl_cholesky.cxx

Go to the documentation of this file.
00001 // This is core/vnl/algo/vnl_cholesky.cxx
00002 #ifdef VCL_NEEDS_PRAGMA_INTERFACE
00003 #pragma implementation
00004 #endif
00005 //:
00006 // \file
00007 // vnl_cholesky
00008 // \author Andrew W. Fitzgibbon, Oxford RRG
00009 // Created: 08 Dec 96
00010 //
00011 //-----------------------------------------------------------------------------
00012 
00013 #include "vnl_cholesky.h"
00014 #include <vcl_cmath.h> // pow()
00015 #include <vcl_cassert.h>
00016 #include <vcl_iostream.h>
00017 #include <vnl/algo/vnl_netlib.h> // dpofa_(), dposl_(), dpoco_(), dpodi_()
00018 
00019 //: Cholesky decomposition.
00020 // Make cholesky decomposition of M optionally computing
00021 // the reciprocal condition number.  If mode is estimate_condition, the
00022 // condition number and an approximate nullspace are estimated, at a cost
00023 // of a factor of (1 + 18/n).  Here's a table of 1 + 18/n:
00024 // \verbatim
00025 // n:              3      5     10     50    100    500   1000
00026 // slowdown:     7.0    4.6    2.8    1.4   1.18   1.04   1.02
00027 // \endverbatim
00028 
00029 vnl_cholesky::vnl_cholesky(vnl_matrix<double> const & M, Operation mode):
00030   A_(M)
00031 {
00032   long n = M.columns();
00033   assert(n == (int)(M.rows()));
00034   num_dims_rank_def_ = -1;
00035   if (vcl_fabs(M(0,n-1) - M(n-1,0)) > 1e-8) {
00036     vcl_cerr << "vnl_cholesky: WARNING: unsymmetric: " << M << vcl_endl;
00037   }
00038 
00039   if (mode != estimate_condition) {
00040     // Quick factorization
00041     v3p_netlib_dpofa_(A_.data_block(), &n, &n, &num_dims_rank_def_);
00042     if (mode == verbose && num_dims_rank_def_ != 0)
00043       vcl_cerr << "vnl_cholesky: " << num_dims_rank_def_ << " dimensions of non-posdeffness\n";
00044   } else {
00045     vnl_vector<double> nullvector(n);
00046     v3p_netlib_dpoco_(A_.data_block(), &n, &n, &rcond_, nullvector.data_block(), &num_dims_rank_def_);
00047     if (num_dims_rank_def_ != 0)
00048       vcl_cerr << "vnl_cholesky: rcond=" << rcond_ << " so " << num_dims_rank_def_ << " dimensions of non-posdeffness\n";
00049   }
00050 }
00051 
00052 //: Solve least squares problem M x = b.
00053 //  The right-hand-side vcl_vector x may be b,
00054 //  which will give a fractional increase in speed.
00055 void vnl_cholesky::solve(vnl_vector<double> const& b, vnl_vector<double>* x) const
00056 {
00057   assert(b.size() == A_.columns());
00058 
00059   *x = b;
00060   long n = A_.columns();
00061   v3p_netlib_dposl_(A_.data_block(), &n, &n, x->data_block());
00062 }
00063 
00064 //: Solve least squares problem M x = b.
00065 vnl_vector<double> vnl_cholesky::solve(vnl_vector<double> const& b) const
00066 {
00067   assert(b.size() == A_.columns());
00068 
00069   long n = A_.columns();
00070   vnl_vector<double> ret = b;
00071   v3p_netlib_dposl_(A_.data_block(), &n, &n, ret.data_block());
00072   return ret;
00073 }
00074 
00075 //: Compute determinant.
00076 double vnl_cholesky::determinant() const
00077 {
00078   long n = A_.columns();
00079   vnl_matrix<double> I = A_;
00080   double det[2];
00081   long job = 10;
00082   v3p_netlib_dpodi_(I.data_block(), &n, &n, det, &job);
00083   return det[0] * vcl_pow(10.0, det[1]);
00084 }
00085 
00086 // : Compute inverse.  Not efficient.
00087 vnl_matrix<double> vnl_cholesky::inverse() const
00088 {
00089   if (num_dims_rank_def_) {
00090     vcl_cerr << "vnl_cholesky: Calling inverse() on rank-deficient matrix\n";
00091     return vnl_matrix<double>();
00092   }
00093 
00094   long n = A_.columns();
00095   vnl_matrix<double> I = A_;
00096   long job = 01;
00097   v3p_netlib_dpodi_(I.data_block(), &n, &n, 0, &job);
00098 
00099   // Copy lower triangle into upper
00100   for (int i = 0; i < n; ++i)
00101     for (int j = i+1; j < n; ++j)
00102       I(i,j) = I(j,i);
00103 
00104   return I;
00105 }
00106 
00107 //: Return lower-triangular factor.
00108 vnl_matrix<double> vnl_cholesky::lower_triangle() const
00109 {
00110   unsigned n = A_.columns();
00111   vnl_matrix<double> L(n,n);
00112   // Zap upper triangle and transpose
00113   for (unsigned i = 0; i < n; ++i) {
00114     L(i,i) = A_(i,i);
00115     for (unsigned j = i+1; j < n; ++j) {
00116       L(j,i) = A_(j,i);
00117       L(i,j) = 0;
00118     }
00119   }
00120   return L;
00121 }
00122 
00123 
00124 //: Return upper-triangular factor.
00125 vnl_matrix<double> vnl_cholesky::upper_triangle() const
00126 {
00127   unsigned n = A_.columns();
00128   vnl_matrix<double> U(n,n);
00129   // Zap lower triangle and transpose
00130   for (unsigned i = 0; i < n; ++i) {
00131     U(i,i) = A_(i,i);
00132     for (unsigned j = i+1; j < n; ++j) {
00133       U(i,j) = A_(j,i);
00134       U(j,i) = 0;
00135     }
00136   }
00137   return U;
00138 }
00139 

Generated on Sat Nov 22 05:06:19 2008 for core/vnl by  doxygen 1.5.1