Object-oriented Scientific Computing Library: Version 0.910
cholesky_base.h
00001 /*
00002   -------------------------------------------------------------------
00003 
00004   Copyright (C) 2006-2012, Andrew W. Steiner
00005 
00006   This file is part of O2scl.
00007 
00008   O2scl is free software; you can redistribute it and/or modify
00009   it under the terms of the GNU General Public License as published by
00010   the Free Software Foundation; either version 3 of the License, or
00011   (at your option) any later version.
00012 
00013   O2scl 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 General Public License for more details.
00017 
00018   You should have received a copy of the GNU General Public License
00019   along with O2scl. If not, see <http://www.gnu.org/licenses/>.
00020 
00021   -------------------------------------------------------------------
00022 */
00023 /* Cholesky Decomposition
00024  *
00025  * Copyright (C) 2000 Thomas Walter
00026  *
00027  * 03 May 2000: Modified for GSL by Brian Gough
00028  * 29 Jul 2005: Additions by Gerard Jungman
00029  * Copyright (C) 2000,2001, 2002, 2003, 2005, 2007 Brian Gough, 
00030  * Gerard Jungman
00031  *
00032  * This is free software; you can redistribute it and/or modify it
00033  * under the terms of the GNU General Public License as published by the
00034  * Free Software Foundation; either version 3, or (at your option) any
00035  * later version.
00036  *
00037  * This source is distributed in the hope that it will be useful, but
00038  * WITHOUT ANY WARRANTY; without even the implied warranty of
00039  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
00040  * General Public License for more details.
00041  */
00042 
00043 /** \brief Compute the in-place Cholesky decomposition of a symmetric
00044     positive-definite square matrix
00045     
00046     On input, the upper triangular part of A is ignored (only the 
00047     lower triangular part and diagonal are used). On output,
00048     the diagonal and lower triangular part contain the matrix L
00049     and the upper triangular part contains L^T. 
00050     
00051     If the matrix is not positive-definite, the error handler 
00052     will be called. 
00053 */
00054 template<class mat_t> void cholesky_decomp(const size_t M, mat_t &A) {
00055   
00056   size_t i,j,k;
00057 
00058   /* [GSL] Do the first 2 rows explicitly. It is simple, and faster.
00059      And one can return if the matrix has only 1 or 2 rows.
00060   */
00061 
00062   double A_00=O2SCL_IX2(A,0,0);
00063   
00064   // AWS: The GSL version stores GSL_NAN in L_00 and then throws
00065   // an error if A_00 <= 0. We throw the error first and then
00066   // the square root should always be safe?
00067 
00068   if (A_00<=0) {
00069     O2SCL_ERR2("Matrix not positive definite (A[0][0]<=0) in ",
00070                "cholesky_decomp().",o2scl::gsl_einval);
00071   }
00072   
00073   double L_00=sqrt(A_00);
00074   O2SCL_IX2(A,0,0)=L_00;
00075   
00076   if (M>1) {
00077     double A_10=O2SCL_IX2(A,1,0);
00078     double A_11=O2SCL_IX2(A,1,1);
00079           
00080     double L_10=A_10/L_00;
00081     double diag=A_11-L_10*L_10;
00082     
00083     if (diag<=0) {
00084       O2SCL_ERR2("Matrix not positive definite (diag<=0 for 2x2) in ",
00085                  "cholesky_decomp().",o2scl::gsl_einval);
00086     }
00087     double L_11=sqrt(diag);
00088 
00089     O2SCL_IX2(A,1,0)=L_10;
00090     O2SCL_IX2(A,1,1)=L_11;
00091   }
00092       
00093   for (k=2;k<M;k++) {
00094     double A_kk=O2SCL_IX2(A,k,k);
00095           
00096     for (i=0;i<k;i++) {
00097       double sum=0.0;
00098 
00099       double A_ki=O2SCL_IX2(A,k,i);
00100       double A_ii=O2SCL_IX2(A,i,i);
00101 
00102       // AWS: Should change to use a form of ddot() here
00103       if (i>0) {
00104         sum=0.0;
00105         for(j=0;j<i;j++) {
00106           sum+=O2SCL_IX2(A,i,j)*O2SCL_IX2(A,k,j);
00107         }
00108       }
00109 
00110       A_ki=(A_ki-sum)/A_ii;
00111       O2SCL_IX2(A,k,i)=A_ki;
00112     } 
00113 
00114     {
00115       double sum=dnrm2_subrow(A,k,0,k);
00116       double diag=A_kk-sum*sum;
00117 
00118       if(diag<=0) {
00119         O2SCL_ERR2("Matrix not positive definite (diag<=0) in ",
00120                    "cholesky_decomp().",o2scl::gsl_einval);
00121       }
00122 
00123       double L_kk=sqrt(diag);
00124       
00125       O2SCL_IX2(A,k,k)=L_kk;
00126     }
00127   }
00128 
00129   /* [GSL] Now copy the transposed lower triangle to the upper
00130      triangle, the diagonal is common.
00131   */
00132       
00133   for (i=1;i<M;i++) {
00134     for (j=0;j<i;j++) {
00135       double A_ij=O2SCL_IX2(A,i,j);
00136       O2SCL_IX2(A,j,i)=A_ij;
00137     }
00138   } 
00139   
00140   return;
00141 }
00142 
00143 /** \brief Solve a symmetric positive-definite linear system after a 
00144     Cholesky decomposition
00145 
00146     Given the Cholesky decomposition of a matrix A in \c LLT, 
00147     this function solves the system <tt>A*x=b</tt>. 
00148 */
00149 template<class mat_t, class vec_t, class vec2_t>
00150   void cholesky_solve(const size_t N, const mat_t &LLT, 
00151                       const vec_t &b, vec2_t &x) {
00152   
00153   // [GSL] Copy x <- b 
00154   vector_copy(N,b,x);
00155   
00156   // [GSL] Solve for c using forward-substitution, L c=b 
00157   o2scl_cblas::dtrsv(o2scl_cblas::o2cblas_RowMajor,
00158                      o2scl_cblas::o2cblas_Lower, 
00159                      o2scl_cblas::o2cblas_NoTrans,
00160                      o2scl_cblas::o2cblas_NonUnit,N,N,LLT,x);
00161   
00162   // [GSL] Perform back-substitution,U x=c 
00163   o2scl_cblas::dtrsv(o2scl_cblas::o2cblas_RowMajor,
00164                      o2scl_cblas::o2cblas_Upper,
00165                      o2scl_cblas::o2cblas_NoTrans,
00166                      o2scl_cblas::o2cblas_NonUnit,N,N,LLT,x);
00167   
00168   return;
00169 }
00170 
00171 /** \brief Solve a linear system in place using a Cholesky decomposition
00172  */
00173 template<class mat_t, class vec_t> 
00174   void cholesky_svx(const size_t N, const mat_t &LLT, vec_t &x) {
00175   
00176   // [GSL] Solve for c using forward-substitution, L c=b 
00177   o2scl_cblas::dtrsv(o2scl_cblas::o2cblas_RowMajor,
00178                      o2scl_cblas::o2cblas_Lower,
00179                      o2scl_cblas::o2cblas_NoTrans,
00180                      o2scl_cblas::o2cblas_NonUnit,N,N,LLT,x);
00181   
00182   // [GSL] Perform back-substitution, U x=c 
00183   o2scl_cblas::dtrsv(o2scl_cblas::o2cblas_RowMajor,
00184                      o2scl_cblas::o2cblas_Upper,
00185                      o2scl_cblas::o2cblas_NoTrans,
00186                      o2scl_cblas::o2cblas_NonUnit,N,N,LLT,x);
00187   
00188   return;
00189 }
00190 
00191 /** \brief Compute the inverse of a symmetric positive definite matrix
00192     given the Cholesky decomposition
00193 
00194     Given a Cholesky decomposition produced by \ref cholesky_decomp(),
00195     this function returns the inverse of that matrix in \c LLT.
00196 */
00197 template<class mat_t> void cholesky_invert(const size_t N, mat_t &LLT) {
00198   
00199   size_t i, j;
00200   double sum;
00201 
00202   // [GSL] invert the lower triangle of LLT
00203   for (i=0;i<N;++i) {
00204     
00205     j=N-i-1;
00206     
00207     O2SCL_IX2(LLT,j,j)=1.0/O2SCL_IX2(LLT,j,j);
00208     double ajj=-O2SCL_IX2(LLT,j,j);
00209     
00210     if (j<N-1) {
00211 
00212       // This section is just the equivalent of dtrmv() for a part of
00213       // the matrix LLT.
00214       {
00215         
00216         size_t ix=N-j-2;
00217         for (size_t ii=N-j-1;ii>0 && ii--;) {
00218           double temp=0.0;
00219           const size_t j_min=0;
00220           const size_t j_max=ii;
00221           size_t jx=j_min;
00222           for (size_t jj=j_min;jj<j_max;jj++) {
00223             temp+=O2SCL_IX2(LLT,jx+j+1,j)*O2SCL_IX2(LLT,ii+j+1,jj+j+1);
00224             jx++;
00225           }
00226           O2SCL_IX2(LLT,ix+j+1,j)=temp+O2SCL_IX2(LLT,ix+j+1,j)*
00227             O2SCL_IX2(LLT,ii+j+1,ii+j+1);
00228           ix--;
00229         }
00230 
00231       }
00232 
00233       o2scl_cblas::dscal_subcol(LLT,j+1,j,N,ajj);
00234 
00235     }
00236   }
00237   
00238   /*
00239     [GSL] The lower triangle of LLT now contains L^{-1}. Now compute
00240     A^{-1}=L^{-t} L^{-1}
00241     
00242     The (ij) element of A^{-1} is column i of L^{-1} dotted into
00243     column j of L^{-1}
00244   */
00245   
00246   for (i=0;i<N;++i) {
00247 
00248     for (j=i+1;j<N;++j) {
00249 
00250       // [GSL] Compute Ainv_{ij}=sum_k Linv_{ki} Linv_{kj}.
00251 
00252       // AWS: Should change to use a form of ddot() here
00253       sum=0.0;
00254       for(size_t k=j;k<N;k++) {
00255         sum+=O2SCL_IX2(LLT,k,i)*O2SCL_IX2(LLT,k,j);
00256       }
00257       
00258       // [GSL] Store in upper triangle
00259       O2SCL_IX2(LLT,i,j)=sum;
00260     }
00261     
00262     // [GSL] now compute the diagonal element
00263     
00264     // AWS: Should change to use a form of ddot() here
00265     sum=0.0;
00266     for(size_t k=i;k<N;k++) {
00267       sum+=O2SCL_IX2(LLT,k,i)*O2SCL_IX2(LLT,k,i);
00268     }
00269 
00270     O2SCL_IX2(LLT,i,i)=sum;
00271   }
00272   
00273   // [GSL] Copy the transposed upper triangle to the lower triangle 
00274 
00275   for (j=1;j<N;j++) {
00276     for (i=0;i<j;i++) {
00277       double A_ij=O2SCL_IX2(LLT,i,j);
00278       O2SCL_IX2(LLT,j,i)=A_ij;
00279     }
00280   } 
00281   
00282   return;
00283 }
00284 
00285 template<class mat_t, class vec_t>
00286   int cholesky_decomp_unit(const size_t N, mat_t &A, vec_t &D) {
00287   
00288   size_t i, j;
00289   
00290   // [GSL] Initial Cholesky decomposition
00291   int stat_chol=cholesky_decomp(N,A);
00292   
00293   // [GSL] Calculate D from diagonal part of initial Cholesky
00294   for(i=0;i<N;++i) {
00295     const double C_ii=O2SCL_IX2(A,i,i);
00296     O2SCL_IX(D,i)=C_ii*C_ii;
00297   }
00298   
00299   // [GSL] Multiply initial Cholesky by 1/sqrt(D) on the right 
00300   for(i=0;i<N;++i) {
00301     for(j=0;j<N;++j) {
00302       O2SCL_IX2(A,i,j)=O2SCL_IX2(A,i,j)/sqrt(O2SCL_IX(D,j));
00303     }
00304   }
00305   
00306   /* [GSL] Because the initial Cholesky contained both L and
00307      transpose(L), the result of the multiplication is not symmetric
00308      anymore; but the lower triangle _is_ correct. Therefore we
00309      reflect it to the upper triangle and declare victory.
00310   */
00311   for(i=0;i<N;++i) {
00312     for(j=i+1;j<N;++j) {
00313       O2SCL_IX2(A,i,j)=O2SCL_IX2(A,j,i);
00314     }
00315   }
00316   
00317   return stat_chol;
00318 }
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines

Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).

Get Object-oriented Scientific Computing
Lib at SourceForge.net. Fast, secure and Free Open Source software
downloads.