Object-oriented Scientific Computing Library: Version 0.910
hh_base.h
Go to the documentation of this file.
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 /* linalg/hh.c
00024  * 
00025  * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007 Gerard Jungman, Brian Gough
00026  * 
00027  * This program is free software; you can redistribute it and/or modify
00028  * it under the terms of the GNU General Public License as published by
00029  * the Free Software Foundation; either version 3 of the License, or (at
00030  * your option) any later version.
00031  * 
00032  * This program is distributed in the hope that it will be useful, but
00033  * WITHOUT ANY WARRANTY; without even the implied warranty of
00034  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00035  * General Public License for more details.
00036  * 
00037  * You should have received a copy of the GNU General Public License
00038  * along with this program; if not, write to the Free Software
00039  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 
00040  * 02110-1301, USA.
00041  */
00042 /** \file hh_base.h
00043     \brief File for householder solver
00044 */
00045 #include <o2scl/err_hnd.h>
00046 #include <o2scl/cblas.h>
00047 #include <o2scl/permutation.h>
00048 #include <o2scl/vec_arith.h>
00049 
00050 #ifdef DOXYGENP
00051 namespace o2scl_linalg {
00052 #endif
00053 
00054   /** \brief Solve linear system after Householder decomposition
00055   */
00056   template<class mat_t, class vec_t, class vec2_t>
00057     int HH_solve(size_t n, mat_t &A, const vec_t &b, vec2_t &x) {
00058     for(size_t i=0;i<n;i++) O2SCL_IX(x,i)=O2SCL_IX(b,i);
00059     return HH_svx(n,n,A,x);
00060   }
00061 
00062   /** \brief Solve a linear system after Householder decomposition in place
00063 
00064       \future Handle memory allocation like the tridiagonal 
00065       functions.
00066   */
00067   template<class mat_t, class vec_t>
00068     int HH_svx(size_t N, size_t M, mat_t &A, vec_t &x) {
00069     size_t i, j, k;
00070   
00071     double *d=(double *)malloc(N*sizeof(double));
00072     if (d == 0) {
00073       O2SCL_ERR_RET("Could not allocate memory for workspace for HH_svx().",
00074                     o2scl::gsl_enomem);
00075     }
00076 
00077     /* Perform Householder transformation. */
00078 
00079     for (i = 0; i < N; i++) {
00080       const double aii = O2SCL_IX2(A,i,i);
00081       double alpha;
00082       double f;
00083       double ak;
00084       double max_norm = 0.0;
00085       double r = 0.0;
00086 
00087       for (k = i; k < M; k++) {
00088         double aki = O2SCL_IX2(A,k,i);
00089         r += aki * aki;
00090       }
00091 
00092       if (r == 0.0) {
00093         /* Rank of matrix is less than size1. */
00094         free(d);
00095         O2SCL_ERR_RET("Matrix is rank deficient in HH_svx().",
00096                       o2scl::gsl_esing);
00097       }
00098 
00099       alpha = sqrt (r) * GSL_SIGN (aii);
00100 
00101       ak = 1.0 / (r + alpha * aii);
00102       O2SCL_IX2(A,i,i)=aii+alpha;
00103       
00104       d[i] = -alpha;
00105 
00106       for (k = i + 1; k < N; k++) {
00107         double norm = 0.0;
00108         f = 0.0;
00109         for (j = i; j < M; j++) {
00110           double ajk = O2SCL_IX2(A,j,k);
00111           double aji = O2SCL_IX2(A,j,i);
00112           norm += ajk * ajk;
00113           f += ajk * aji;
00114         }
00115         max_norm = GSL_MAX (max_norm, norm);
00116 
00117         f *= ak;
00118 
00119         for (j = i; j < M; j++) {
00120           double ajk = O2SCL_IX2(A,j,k);
00121           double aji = O2SCL_IX2(A,j,i);
00122           O2SCL_IX2(A,j,k)=ajk-f*aji;
00123         }
00124       }
00125 
00126       if (fabs (alpha) < 2.0 * GSL_DBL_EPSILON * sqrt (max_norm)) {
00127         /* Apparent singularity. */
00128         free(d);
00129         O2SCL_ERR_RET("Apparent singularity in HH_svx().",o2scl::gsl_esing);
00130       }
00131 
00132       /* Perform update of RHS. */
00133 
00134       f = 0.0;
00135       for (j = i; j < M; j++) {
00136         f += O2SCL_IX(x,j)*O2SCL_IX2(A,j,i);
00137       }
00138       f *= ak;
00139       for (j = i; j < M; j++) {
00140         double xj = O2SCL_IX(x,j);
00141         double aji = O2SCL_IX2(A,j,i);
00142         O2SCL_IX(x,j)=xj - f * aji;
00143       }
00144     }
00145 
00146     /* Perform back-substitution. */
00147     
00148     for (i = N; i-- > 0; ) {
00149       double xi = O2SCL_IX(x,i);
00150       double sum = 0.0;
00151       for (k = i + 1; k < N; k++) {
00152         sum += O2SCL_IX2(A,i,k)*O2SCL_IX(x,k);
00153       }
00154     
00155       O2SCL_IX(x,i)=(xi - sum) / d[i];
00156     }
00157 
00158     free(d);
00159 
00160     return o2scl::gsl_success;
00161   }
00162   
00163 #ifdef DOXYGENP
00164 }
00165 #endif
 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.