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