00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 2008, 2009, 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 #ifndef O2SCL_LINEAR_SOLVER_H 00024 #define O2SCL_LINEAR_SOLVER_H 00025 00026 #include <o2scl/ovector_tlate.h> 00027 #include <o2scl/omatrix_tlate.h> 00028 #include <o2scl/permutation.h> 00029 #include <o2scl/lu.h> 00030 #include <o2scl/qr.h> 00031 #include <o2scl/hh.h> 00032 #include <o2scl/array.h> 00033 00034 namespace o2scl_linalg { 00035 00036 /** 00037 \brief A generic solver for the linear system \f$ A x = b \f$ 00038 [abstract base] 00039 00040 A generic solver for dense linear systems. 00041 00042 Those writing production level code should consider calling 00043 LAPACK directly using \o2 objects as described in the \ref 00044 linalg_section section of the User's Guide. 00045 00046 \future The test code uses a Hilbert matrix, which is known 00047 to be ill-conditioned, especially for the larger sizes. This 00048 should probably be changed. 00049 */ 00050 template<class vec_t, class mat_t> class linear_solver { 00051 public: 00052 00053 virtual ~linear_solver() {} 00054 00055 /// Solve square linear system \f$ A x = b \f$ of size \c n 00056 virtual int solve(size_t n, mat_t &a, vec_t &b, vec_t &x)=0; 00057 }; 00058 00059 /// Generic linear solver using LU decomposition 00060 template <class vec_t, class mat_t> class linear_solver_lu : 00061 public linear_solver<vec_t, mat_t> { 00062 00063 public: 00064 00065 /// Solve square linear system \f$ A x = b \f$ of size \c n 00066 virtual int solve(size_t n, mat_t &A, 00067 vec_t &b, vec_t &x) { 00068 int sig; 00069 o2scl::permutation p(n); 00070 LU_decomp(n,A,p,sig); 00071 LU_solve(n,A,p,b,x); 00072 return 0; 00073 }; 00074 00075 virtual ~linear_solver_lu() {} 00076 00077 }; 00078 00079 /// Generic linear solver using QR decomposition 00080 template <class vec_t, class mat_t> class linear_solver_qr : 00081 public linear_solver<vec_t, mat_t> { 00082 00083 public: 00084 00085 /// Solve square linear system \f$ A x = b \f$ of size \c n 00086 virtual int solve(size_t n, mat_t &A, 00087 vec_t &b, vec_t &x) { 00088 double *tau=o2scl::new_array<double>(n); 00089 QR_decomp(n,n,A,tau); 00090 QR_solve(n,A,tau,b,x); 00091 o2scl::delete_array(tau); 00092 return 0; 00093 }; 00094 00095 virtual ~linear_solver_qr() {} 00096 00097 }; 00098 00099 /// Generic Householder linear solver 00100 template <class vec_t, class mat_t> class linear_solver_hh : 00101 public linear_solver<vec_t, mat_t> { 00102 00103 public: 00104 00105 /// Solve square linear system \f$ A x = b \f$ of size \c n 00106 virtual int solve(size_t n, mat_t &A, 00107 vec_t &b, vec_t &x) { 00108 return HH_solve(n,A,b,x); 00109 }; 00110 00111 virtual ~linear_solver_hh() {} 00112 00113 }; 00114 00115 /// GSL solver by LU decomposition 00116 #ifdef DOXYGENP 00117 class gsl_solver_LU : public linear_solver 00118 #else 00119 class gsl_solver_LU : public linear_solver<o2scl::ovector_base, 00120 o2scl::omatrix_base> 00121 #endif 00122 { 00123 00124 public: 00125 00126 /// Solve square linear system \f$ A x = b \f$ of size \c n 00127 virtual int solve(size_t n, o2scl::omatrix_base &A, 00128 o2scl::ovector_base &b, 00129 o2scl::ovector_base &x) { 00130 gsl_permutation *P=gsl_permutation_alloc(n); 00131 int s; 00132 gsl_linalg_LU_decomp((gsl_matrix *)&A,P,&s); 00133 gsl_linalg_LU_solve((gsl_matrix *)&A,P,(gsl_vector *)&b, 00134 (gsl_vector *)&x); 00135 gsl_permutation_free(P); 00136 00137 return 0; 00138 }; 00139 00140 virtual ~gsl_solver_LU() {} 00141 00142 }; 00143 00144 /// GSL solver by QR decomposition 00145 #ifdef DOXYGENP 00146 class gsl_solver_QR : public linear_solver 00147 #else 00148 class gsl_solver_QR : public linear_solver<o2scl::ovector_base, 00149 o2scl::omatrix_base> 00150 #endif 00151 { 00152 00153 public: 00154 00155 /// Solve square linear system \f$ A x = b \f$ of size \c n 00156 virtual int solve(size_t n, o2scl::omatrix_base &A, 00157 o2scl::ovector_base &b, 00158 o2scl::ovector_base &x) { 00159 o2scl::ovector tau(n); 00160 gsl_linalg_QR_decomp((gsl_matrix *)&A,(gsl_vector *)&tau); 00161 gsl_linalg_QR_solve((gsl_matrix *)&A,(gsl_vector *)&tau, 00162 (gsl_vector *)&b,(gsl_vector *)&x); 00163 return 0; 00164 }; 00165 00166 virtual ~gsl_solver_QR() {} 00167 00168 }; 00169 00170 /// GSL Householder solver 00171 #ifdef DOXYGENP 00172 class gsl_solver_HH : public linear_solver 00173 #else 00174 class gsl_solver_HH : public linear_solver<o2scl::ovector_base, 00175 o2scl::omatrix_base> 00176 #endif 00177 { 00178 00179 public: 00180 00181 /// Solve square linear system \f$ A x = b \f$ of size \c n 00182 virtual int solve(size_t n, o2scl::omatrix_base &A, 00183 o2scl::ovector_base &b, 00184 o2scl::ovector_base &x) { 00185 return gsl_linalg_HH_solve((gsl_matrix *)(&A),(gsl_vector *)(&b), 00186 (gsl_vector *)(&x)); 00187 }; 00188 00189 virtual ~gsl_solver_HH() {} 00190 00191 }; 00192 00193 } 00194 00195 #endif
Documentation generated with Doxygen and provided under the GNU Free Documentation License. See License Information for details.
Project hosting provided by
,
O2scl Sourceforge Project Page