Object-oriented Scientific Computing Library: Version 0.910
ode_it_solve.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_ODE_IT_SOLVE_H
00024 #define O2SCL_ODE_IT_SOLVE_H
00025 
00026 #include <iostream>
00027 
00028 #include <o2scl/misc.h>
00029 #include <o2scl/array.h>
00030 #include <o2scl/ovector_tlate.h>
00031 #include <o2scl/omatrix_tlate.h>
00032 #include <o2scl/uvector_tlate.h>
00033 #include <o2scl/umatrix_tlate.h>
00034 #include <o2scl/test_mgr.h>
00035 #include <o2scl/linear_solver.h>
00036 
00037 #ifndef DOXYGENP
00038 namespace o2scl {
00039 #endif
00040 
00041   /// Function class for ode_it_solve
00042   template<class vec_t=ovector_base> class ode_it_funct {
00043   
00044   public:
00045   
00046   ode_it_funct() {
00047   }
00048   
00049   virtual ~ode_it_funct() {}
00050   
00051   /// Using \c x and \c y, return the value of function number \c ieq 
00052   virtual double operator()(size_t ieq, double x, vec_t &y) {
00053     return 0.0;
00054   }
00055 
00056   };
00057 
00058   /// Function pointer for ode_it_solve
00059   template<class vec_t=ovector_base> 
00060     class ode_it_funct_fptr : public ode_it_funct<vec_t> {
00061 
00062   public:
00063 
00064   virtual ~ode_it_funct_fptr() {}
00065 
00066   /// Create using a function pointer
00067   ode_it_funct_fptr(double (*fp)(size_t, double, vec_t &)) {
00068     fptr=fp;
00069   }
00070 
00071   /// Using \c x and \c y, return the value of function number \c ieq 
00072   virtual double operator()(size_t ieq, double x, vec_t &y) {
00073     return fptr(ieq,x,y);
00074   }
00075 
00076   protected:
00077 
00078   ode_it_funct_fptr() {}
00079 
00080   /// The function pointer
00081   double (*fptr)(size_t ieq, double x, vec_t &y);
00082 
00083   };
00084 
00085   /// Member function pointer for ode_it_solve
00086   template<class tclass, class vec_t=ovector_base>
00087     class ode_it_funct_mfptr : public ode_it_funct<vec_t> {
00088   public:
00089 
00090   virtual ~ode_it_funct_mfptr() {}
00091 
00092   /// Create using a class instance and member function
00093   ode_it_funct_mfptr(tclass *tp, 
00094                      double (tclass::*fp)(size_t, double, vec_t &)) {
00095     tptr=tp;
00096     fptr=fp;
00097   }
00098 
00099   /// Using \c x and \c y, return the value of function number \c ieq 
00100   virtual double operator()(size_t ieq, double x, vec_t &y) {
00101     return (*tptr.*fptr)(ieq,x,y);
00102   }
00103 
00104   protected:
00105 
00106   ode_it_funct_mfptr() {}
00107 
00108   /// The class pointer
00109   tclass *tptr;
00110 
00111   /// The member function pointer
00112   double (tclass::*fptr)(size_t ieq, double x, vec_t &y);
00113 
00114   };
00115 
00116   /** \brief ODE solver using a generic linear solver to solve 
00117       finite-difference equations
00118 
00119       \note This class does not work with all matrix and vector
00120       types because using operator[][] would require it to be
00121       inefficient for use with sparse matrices.
00122 
00123       \future Create a GSL-like set() and iterate() interface
00124       \future Implement as a child of ode_bv_solve ?
00125       \future Max and average tolerance?
00126   */
00127   template <class func_t=ode_it_funct<>, 
00128     class vec_t=ovector_base, class mat_t=omatrix_base, 
00129     class matrix_row_t=omatrix_row, class solver_vec_t=ovector_base, 
00130     class solver_mat_t=omatrix_base> class ode_it_solve {
00131     
00132   public:
00133 
00134   bool make_mats;
00135   
00136   ode_it_solve() {
00137     h=1.0e-4;
00138     niter=30;
00139     tol_rel=1.0e-8;
00140     verbose=0;
00141     solver=&def_solver;
00142     alpha=1.0;
00143     make_mats=false;
00144   }
00145 
00146   virtual ~ode_it_solve() {}
00147 
00148   /// Set level of output (default 0)
00149   int verbose;
00150 
00151   /** \brief Stepsize for finite differencing (default \f$ 10^{-4} \f$)
00152    */
00153   double h;
00154 
00155   /// Tolerance (default \f$ 10^{-8} \f$)
00156   double tol_rel;
00157   
00158   /// Maximum number of iterations (default 30)
00159   size_t niter;
00160   
00161   /// Set the linear solver
00162   int set_solver(o2scl_linalg::linear_solver<solver_vec_t,solver_mat_t> &ls) {
00163     solver=&ls;
00164     return 0;
00165   }
00166     
00167   /// Size of correction to apply (default 1.0)
00168   double alpha;
00169 
00170   /** \brief Solve \c derivs with boundary conditions \c left and 
00171       \c right
00172 
00173       Given a grid of size \c n_grid and \c n_eq differential equations,
00174       solve them by relaxation. The grid is specified in \c x, which
00175       is a vector of size \c n_grid. The differential equations are
00176       given in \c derivs, the boundary conditions on the left hand
00177       side in \c left, and the boundary conditions on the right hand
00178       side in \c right. The number of boundary conditions on the left
00179       hand side is \c nb_left, and the number of boundary conditions on
00180       the right hand side should be <tt>n_eq-nb_left</tt>. The initial
00181       guess for the solution, a matrix of size <tt>[n_grid][n_eq]</tt>
00182       should be given in \c y. Upon success, \c y will contain an
00183       approximate solution of the differential equations. The matrix
00184       \c mat is workspace of size <tt>[n_grid*n_eq][n_grid*n_eq]</tt>, and
00185       the vectors \c rhs and \c y are workspace of size
00186       <tt>[n_grid*n_eq]</tt>.
00187 
00188   */
00189   int solve(size_t n_grid, size_t n_eq, size_t nb_left, vec_t &x, 
00190             mat_t &y, func_t &derivs, func_t &left, func_t &right,
00191             solver_mat_t &mat, solver_vec_t &rhs, solver_vec_t &dy) {
00192 
00193     // Store the functions for simple derivatives
00194     fd=&derivs;
00195     fl=&left;
00196     fr=&right;
00197     
00198     // Variable index
00199     size_t ix;
00200 
00201     // Number of RHS boundary conditions
00202     size_t nb_right=n_eq-nb_left;
00203 
00204     // Number of variables
00205     size_t nvars=n_grid*n_eq;
00206     
00207     bool done=false;
00208     for(size_t it=0;done==false && it<niter;it++) {
00209       
00210       ix=0;
00211       
00212       mat.set_all(0.0);
00213 
00214       // Construct the entries corresponding to the LHS boundary. 
00215       // This makes the first nb_left rows of the matrix.
00216       for(size_t i=0;i<nb_left;i++) {
00217         matrix_row_t yk(y,0);
00218         rhs[ix]=-left(i,x[0],yk);
00219         for(size_t j=0;j<n_eq;j++) {
00220           int rxa=mat.set(ix,j,fd_left(i,j,x[0],yk));
00221         }
00222         ix++;
00223       }
00224 
00225       // Construct the matrix entries for the internal points
00226       // This loop adds n_grid-1 sets of n_eq rows
00227       for(size_t k=0;k<n_grid-1;k++) {
00228         size_t kp1=k+1;
00229         double tx=(x[kp1]+x[k])/2.0;
00230         double dx=x[kp1]-x[k];
00231         matrix_row_t yk(y,k);
00232         matrix_row_t ykp1(y,k+1);
00233         
00234         for(size_t i=0;i<n_eq;i++) {
00235           
00236           rhs[ix]=y[k][i]-y[kp1][i]+(x[kp1]-x[k])*
00237             (derivs(i,tx,ykp1)+derivs(i,tx,yk))/2.0;
00238           
00239           size_t lhs=k*n_eq;
00240           for(size_t j=0;j<n_eq;j++) {
00241             int rxb=mat.set(ix,lhs+j,-fd_derivs(i,j,tx,yk)*dx/2.0);
00242             int rxc=mat.set(ix,lhs+j+n_eq,-fd_derivs(i,j,tx,ykp1)*dx/2.0);
00243             if (i==j) {
00244               int rxd=mat.set(ix,lhs+j,mat(ix,lhs+j)-1.0);
00245               int rxe=mat.set(ix,lhs+j+n_eq,mat(ix,lhs+j+n_eq)+1.0);
00246             }
00247           }
00248 
00249           ix++;
00250 
00251         }
00252       }
00253       
00254       // Construct the entries corresponding to the RHS boundary
00255       // This makes the last nb_right rows of the matrix.
00256       for(size_t i=0;i<nb_right;i++) {
00257         matrix_row_t ylast(y,n_grid-1);
00258         size_t lhs=n_eq*(n_grid-1);
00259         
00260         rhs[ix]=-right(i,x[n_grid-1],ylast);
00261         
00262         for(size_t j=0;j<n_eq;j++) {
00263           int rxf=mat.set(ix,lhs+j,fd_right(i,j,x[n_grid-1],ylast));
00264         }
00265           
00266         ix++;
00267 
00268       }
00269       
00270       // Compute correction by calling the linear solver
00271 
00272       if (verbose>3) {
00273         std::cout << "Matrix: " << std::endl;
00274         for(size_t i=0;i<nvars;i++) {
00275           for(size_t j=0;j<nvars;j++) {
00276             std::cout << mat[i][j] << " ";
00277           }
00278           std::cout << std::endl;
00279         }
00280         std::cout << "Deviations:" << std::endl;
00281         for(size_t i=0;i<nvars;i++) {
00282           std::cout << rhs[i] << std::endl;
00283         }
00284       }
00285 
00286       if (make_mats) return 0;
00287         
00288       int ret=solver->solve(ix,mat,rhs,dy);
00289 
00290       if (verbose>3) {
00291         std::cout << "Corrections:" << std::endl;
00292         for(size_t i=0;i<nvars;i++) {
00293           std::cout << dy[i] << std::endl;
00294         }
00295         std::cout << "Press a key and press enter to continue: " 
00296                   << std::endl;
00297         char ch;
00298         std::cin >> ch;
00299       }
00300       
00301       // Apply correction and compute its size
00302 
00303       double res=0.0;
00304       ix=0;
00305 
00306       for(size_t igrid=0;igrid<n_grid;igrid++) {
00307         for(size_t ieq=0;ieq<n_eq;ieq++) {
00308           y[igrid][ieq]+=alpha*dy[ix];
00309           res+=dy[ix]*dy[ix];
00310           ix++;
00311         }
00312       }
00313 
00314       if (verbose>0) {
00315         // Since we're in the o2scl namespace, we explicitly
00316         // specify std::sqrt() here
00317         std::cout << "ode_it_solve: " << it << " " << std::sqrt(res) << " " 
00318                   << tol_rel << std::endl;
00319         if (verbose>1) {
00320           char ch;
00321           std::cout << "Press a key and type enter to continue. ";
00322           std::cin >> ch;
00323         }
00324       }
00325       
00326       // If the correction has become small enough, we're done
00327       if (std::sqrt(res)<=tol_rel) done=true;
00328     }
00329 
00330     if (done==false) {
00331       O2SCL_ERR_RET("Exceeded number of iterations in solve().",
00332                     o2scl::gsl_emaxiter);
00333     }
00334 
00335     return 0;
00336   }
00337   
00338   /// Default linear solver
00339   o2scl_linalg::linear_solver_hh<solver_vec_t,solver_mat_t> def_solver;
00340   
00341   protected:
00342   
00343   /// \name Storage for functions
00344   //@{
00345   ode_it_funct<vec_t> *fl, *fr, *fd;
00346   //@}
00347 
00348   /// Solver
00349   o2scl_linalg::linear_solver<solver_vec_t,solver_mat_t> *solver;
00350 
00351   /** \brief Compute the derivatives of the LHS boundary conditions
00352 
00353       This function computes \f$ \partial f_{left,\mathrm{ieq}} / \partial
00354       y_{\mathrm{ivar}} \f$
00355   */
00356   virtual double fd_left(size_t ieq, size_t ivar, double x, vec_t &y) {
00357 
00358     double ret, dydx;
00359     
00360     y[ivar]+=h;
00361     ret=(*fl)(ieq,x,y);
00362     
00363     y[ivar]-=h;
00364     ret-=(*fl)(ieq,x,y);
00365     
00366     ret/=h;
00367     return ret;
00368   }
00369   
00370   /** \brief Compute the derivatives of the RHS boundary conditions
00371         
00372       This function computes \f$ \partial f_{right,\mathrm{ieq}} / \partial
00373       y_{\mathrm{ivar}} \f$
00374   */
00375   virtual double fd_right(size_t ieq, size_t ivar, double x, vec_t &y) {
00376 
00377     double ret, dydx;
00378     
00379     y[ivar]+=h;
00380     ret=(*fr)(ieq,x,y);
00381     
00382     y[ivar]-=h;
00383     ret-=(*fr)(ieq,x,y);
00384     
00385     ret/=h;
00386     return ret;
00387   }
00388   
00389   /** \brief Compute the finite-differenced part of the 
00390       differential equations
00391 
00392       This function computes \f$ \partial f_{\mathrm{ieq}} / \partial
00393       y_{\mathrm{ivar}} \f$
00394   */
00395   virtual double fd_derivs(size_t ieq, size_t ivar, double x, vec_t &y) {
00396 
00397     double ret, dydx;
00398     
00399     y[ivar]+=h;
00400     ret=(*fd)(ieq,x,y);
00401     
00402     y[ivar]-=h;
00403     ret-=(*fd)(ieq,x,y);
00404     
00405     ret/=h;
00406     
00407     return ret;
00408   }
00409   
00410   };
00411 
00412 #ifndef DOXYGENP
00413 }
00414 #endif
00415 
00416 #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.