00001
00002 #ifdef VCL_NEEDS_PRAGMA_INTERFACE
00003 #pragma implementation
00004 #endif
00005
00006
00007
00008
00009
00010
00011
00012
00013 #include "vnl_ldl_cholesky.h"
00014 #include <vcl_cmath.h>
00015 #include <vcl_cassert.h>
00016 #include <vcl_iostream.h>
00017 #include <vnl/algo/vnl_netlib.h>
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029 vnl_ldl_cholesky::vnl_ldl_cholesky(vnl_matrix<double> const & M, Operation mode):
00030 L_(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_ldl_cholesky: WARNING: unsymmetric: " << M << vcl_endl;
00037 }
00038
00039 if (mode != estimate_condition) {
00040
00041 v3p_netlib_dpofa_(L_.data_block(), &n, &n, &num_dims_rank_def_);
00042 if (mode == verbose && num_dims_rank_def_ != 0)
00043 vcl_cerr << "vnl_ldl_cholesky: " << num_dims_rank_def_ << " dimensions of non-posdeffness\n";
00044 } else {
00045 vnl_vector<double> nullvector(n);
00046 v3p_netlib_dpoco_(L_.data_block(), &n, &n, &rcond_, nullvector.data_block(), &num_dims_rank_def_);
00047 if (num_dims_rank_def_ != 0)
00048 vcl_cerr << "vnl_ldl_cholesky: rcond=" << rcond_ << " so " << num_dims_rank_def_ << " dimensions of non-posdeffness\n";
00049 }
00050
00051
00052
00053 d_.set_size(n);
00054
00055
00056 vnl_vector<double> sqrt_d(n);
00057
00058 for (int i=0; i<n; ++i)
00059 {
00060 sqrt_d[i]=L_(i,i);
00061 d_[i]=sqrt_d[i]*sqrt_d[i];
00062 }
00063
00064
00065 for (int i=0; i<n; ++i)
00066 {
00067 double *row = L_[i];
00068 for (int j=0; j<i; ++j) row[j]/=sqrt_d[j];
00069 row[i]=1.0;
00070 for (int j=i+1; j<n; ++j) row[j]=0.0;
00071 }
00072 }
00073
00074
00075 inline double dot(const double* v1, const double* v2, unsigned n)
00076 {
00077 double sum=0.0;
00078 for (unsigned i=0;i<n;++i) sum+= v1[i]*v2[i];
00079 return sum;
00080 }
00081
00082 inline double dot(const double* v1, unsigned s, const double* v2, unsigned n)
00083 {
00084 double sum=0.0;
00085 for (unsigned i=0;i<n;++i,v1+=s) sum+= (*v1)*v2[i];
00086 return sum;
00087 }
00088
00089
00090
00091 void vnl_ldl_cholesky::solve_lx(vnl_vector<double>& x)
00092 {
00093 unsigned n = d_.size();
00094 for (unsigned i=1;i<n;++i)
00095 x[i] -= dot(L_[i],x.data_block(),i);
00096 }
00097
00098
00099
00100
00101 void vnl_ldl_cholesky::inplace_solve(double* x) const
00102 {
00103 unsigned n = d_.size();
00104
00105 for (unsigned i=1;i<n;++i)
00106 x[i] -= dot(L_[i],x,i);
00107
00108
00109 for (unsigned i=0;i<n;++i) x[i]/=d_[i];
00110
00111
00112 const double* L_data = &L_(n-1,n-2);
00113 const double* x_data = &x[n-1];
00114 unsigned c=1;
00115 for (int i=n-2;i>=0;--i,L_data-=(n+1),--x_data,++c)
00116 {
00117 x[i] -= dot(L_data,n,x_data,c);
00118 }
00119 }
00120
00121
00122
00123
00124 double vnl_ldl_cholesky::xt_m_inv_x(const vnl_vector<double>& x) const
00125 {
00126 unsigned n=d_.size();
00127 assert(x.size()==n);
00128 vnl_vector<double> y=x;
00129
00130 double sum = y[0]*y[0]/d_[0];
00131 for (unsigned i=1;i<n;++i)
00132 {
00133 y[i] -= dot(L_[i],y.data_block(),i);
00134 sum += y[i]*y[i]/d_[i];
00135 }
00136 return sum;
00137 }
00138
00139
00140
00141 double vnl_ldl_cholesky::xt_m_x(const vnl_vector<double>& x) const
00142 {
00143 unsigned n=d_.size();
00144 assert(x.size()==n);
00145 double sum=0.0;
00146 const double* xd = x.data_block();
00147 const double* L_col = L_.data_block();
00148 unsigned c=n;
00149 for (unsigned i=0;i<n;++i,++xd,L_col+=(n+1),--c)
00150 {
00151 double xLi = dot(L_col,n,xd,c);
00152 sum+= xLi*xLi*d_[i];
00153 }
00154 return sum;
00155 }
00156
00157
00158
00159
00160
00161 void vnl_ldl_cholesky::solve(vnl_vector<double> const& b,
00162 vnl_vector<double>* xp) const
00163 {
00164 assert(b.size() == d_.size());
00165 *xp = b;
00166 inplace_solve(xp->data_block());
00167 }
00168
00169
00170 vnl_vector<double> vnl_ldl_cholesky::solve(vnl_vector<double> const& b) const
00171 {
00172 assert(b.size() == L_.columns());
00173
00174 vnl_vector<double> ret = b;
00175 solve(b,&ret);
00176 return ret;
00177 }
00178
00179
00180 double vnl_ldl_cholesky::determinant() const
00181 {
00182 unsigned n=d_.size();
00183 double det=1.0;
00184 for (unsigned i=0;i<n;++i) det*=d_[i];
00185 return det;
00186 }
00187
00188
00189
00190
00191
00192
00193
00194 void vnl_ldl_cholesky::rank1_update(const vnl_vector<double>& v)
00195 {
00196 unsigned n = d_.size();
00197 assert(v.size()==n);
00198 double a = 1.0;
00199 vnl_vector<double> w=v;
00200 for (unsigned j=0;j<n;++j)
00201 {
00202 double a2=a+w[j]*w[j]/d_[j];
00203 d_[j]*=a2;
00204 double gamma = w[j]/d_[j];
00205 d_[j]/=a;
00206 a=a2;
00207
00208 for (unsigned p=j+1;p<n;++p)
00209 {
00210 w[p]-=w[j]*L_(p,j);
00211 L_(p,j)+=gamma*w[p];
00212 }
00213 }
00214 }
00215
00216
00217
00218
00219 void vnl_ldl_cholesky::update(const vnl_matrix<double>& W0)
00220 {
00221 unsigned n = d_.size();
00222 assert(W0.rows()==n);
00223 unsigned r = W0.columns();
00224
00225 vnl_matrix<double> W(W0);
00226 vnl_vector<double> a(r,1.0),gamma(r);
00227 for (unsigned j=0;j<n;++j)
00228 {
00229 double* Wj = W[j];
00230 for (unsigned i=0;i<r;++i)
00231 {
00232 double a2=a[i]+Wj[i]*Wj[i]/d_[j];
00233 d_[j]*=a2;
00234 gamma[i]=Wj[i]/d_[j];
00235 d_[j]/=a[i];
00236 a[i]=a2;
00237 }
00238 for (unsigned p=j+1;p<n;++p)
00239 {
00240 double *Wp = W[p];
00241 double *Lp = L_[p];
00242 for (unsigned i=0;i<r;++i)
00243 {
00244 Wp[i]-=Wj[i]*Lp[j];
00245 Lp[j]+=gamma[i]*Wp[i];
00246 }
00247 }
00248 }
00249 }
00250
00251
00252 vnl_matrix<double> vnl_ldl_cholesky::inverse() const
00253 {
00254 if (num_dims_rank_def_) {
00255 vcl_cerr << "vnl_ldl_cholesky: Calling inverse() on rank-deficient matrix\n";
00256 return vnl_matrix<double>();
00257 }
00258
00259 unsigned int n = d_.size();
00260 vnl_matrix<double> R(n,n);
00261 R.set_identity();
00262
00263
00264
00265 for (unsigned int i=0; i<n; ++i)
00266 inplace_solve(R[i]);
00267
00268 return R;
00269 }
00270