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_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 param_t, 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, param_t &pa) { 00053 return 0.0; 00054 } 00055 00056 }; 00057 00058 /// Function pointer for ode_it_solve 00059 template<class param_t, class vec_t=ovector_base> 00060 class ode_it_funct_fptr : public ode_it_funct<param_t,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 &, param_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, param_t &pa) { 00073 return fptr(ieq,x,y,pa); 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, param_t &pa); 00082 }; 00083 00084 /// Member function pointer for ode_it_solve 00085 template<class tclass, class param_t, class vec_t=ovector_base> 00086 class ode_it_funct_mfptr : public ode_it_funct<param_t,vec_t> { 00087 public: 00088 00089 virtual ~ode_it_funct_mfptr() {} 00090 00091 /// Create using a class instance and member function 00092 ode_it_funct_mfptr(tclass *tp, 00093 double (tclass::*fp)(size_t, double, 00094 vec_t &, param_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, param_t &pa) { 00101 return (*tptr.*fptr)(ieq,x,y,pa); 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, param_t &pa); 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 param_t, class func_t=ode_it_funct<param_t>, 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 ode_it_solve() { 00135 h=1.0e-4; 00136 niter=30; 00137 tolf=1.0e-8; 00138 verbose=0; 00139 solver=&def_solver; 00140 alpha=1.0; 00141 } 00142 00143 virtual ~ode_it_solve() {} 00144 00145 /// Set level of output (default 0) 00146 int verbose; 00147 00148 /** \brief Stepsize for finite differencing (default \f$ 10^{-4} \f$) 00149 */ 00150 double h; 00151 00152 /// Tolerance (default \f$ 10^{-8} \f$) 00153 double tolf; 00154 00155 /// Maximum number of iterations (default 30) 00156 size_t niter; 00157 00158 /// Set the linear solver 00159 int set_solver(o2scl_linalg::linear_solver<solver_vec_t,solver_mat_t> &ls) { 00160 solver=&ls; 00161 return 0; 00162 } 00163 00164 /// Size of correction to apply (default 1.0) 00165 double alpha; 00166 00167 /** \brief Solve \c derivs with boundary conditions \c left and 00168 \c right 00169 00170 Given a grid of size \c n_grid and \c n_eq differential equations, 00171 solve them by relaxation. The grid is specified in \c x, which 00172 is a vector of size \c n_grid. The differential equations are 00173 given in \c derivs, the boundary conditions on the left hand 00174 side in \c left, and the boundary conditions on the right hand 00175 side in \c right. The number of boundary conditions on the left 00176 hand side is \c nb_left, and the number of boundary conditions on 00177 the right hand side should be <tt>n_eq-nb_left</tt>. The initial 00178 guess for the solution, a matrix of size <tt>[n_grid][n_eq]</tt> 00179 should be given in \c y. Upon success, \c y will contain an 00180 approximate solution of the differential equations. The matrix 00181 \c mat is workspace of size <tt>[n_grid*n_eq][n_grid*n_eq]</tt>, and 00182 the vectors \c rhs and \c y are workspace of size 00183 <tt>[n_grid*n_eq]</tt>. 00184 00185 */ 00186 int solve(size_t n_grid, size_t n_eq, size_t nb_left, vec_t &x, 00187 mat_t &y, func_t &derivs, func_t &left, func_t &right, 00188 solver_mat_t &mat, solver_vec_t &rhs, solver_vec_t &dy, 00189 param_t &pa) { 00190 00191 // Store the functions for simple derivatives 00192 fd=&derivs; 00193 fl=&left; 00194 fr=&right; 00195 00196 // Variable index 00197 size_t ix; 00198 00199 // Number of RHS boundary conditions 00200 size_t nb_right=n_eq-nb_left; 00201 00202 // Number of variables 00203 size_t nvars=n_grid*n_eq; 00204 00205 bool done=false; 00206 for(size_t it=0;done==false && it<niter;it++) { 00207 00208 ix=0; 00209 00210 mat.set_all(0.0); 00211 00212 // Construct the entries corresponding to the LHS boundary. 00213 // This makes the first nb_left rows of the matrix. 00214 for(size_t i=0;i<nb_left;i++) { 00215 matrix_row_t yk(y,0); 00216 rhs[ix]=-left(i,x[0],yk,pa); 00217 for(size_t j=0;j<n_eq;j++) { 00218 int rxa=mat.set(ix,j,fd_left(i,j,x[0],yk,pa)); 00219 } 00220 ix++; 00221 } 00222 00223 // Construct the matrix entries for the internal points 00224 // This loop adds n_grid-1 sets of n_eq rows 00225 for(size_t k=0;k<n_grid-1;k++) { 00226 size_t kp1=k+1; 00227 double tx=(x[kp1]+x[k])/2.0; 00228 double dx=x[kp1]-x[k]; 00229 matrix_row_t yk(y,k); 00230 matrix_row_t ykp1(y,k+1); 00231 00232 for(size_t i=0;i<n_eq;i++) { 00233 00234 rhs[ix]=y[k][i]-y[kp1][i]+(x[kp1]-x[k])* 00235 (derivs(i,tx,ykp1,pa)+derivs(i,tx,yk,pa))/2.0; 00236 00237 size_t lhs=k*n_eq; 00238 for(size_t j=0;j<n_eq;j++) { 00239 int rxb=mat.set(ix,lhs+j,-fd_derivs(i,j,tx,yk,pa)*dx/2.0); 00240 int rxc=mat.set(ix,lhs+j+n_eq,-fd_derivs(i,j,tx,ykp1,pa)*dx/2.0); 00241 if (i==j) { 00242 int rxd=mat.set(ix,lhs+j,mat(ix,lhs+j)-1.0); 00243 int rxe=mat.set(ix,lhs+j+n_eq,mat(ix,lhs+j+n_eq)+1.0); 00244 } 00245 } 00246 00247 ix++; 00248 00249 } 00250 } 00251 00252 // Construct the entries corresponding to the RHS boundary 00253 // This makes the last nb_right rows of the matrix. 00254 for(size_t i=0;i<nb_right;i++) { 00255 matrix_row_t ylast(y,n_grid-1); 00256 size_t lhs=n_eq*(n_grid-1); 00257 00258 rhs[ix]=-right(i,x[n_grid-1],ylast,pa); 00259 00260 for(size_t j=0;j<n_eq;j++) { 00261 int rxf=mat.set(ix,lhs+j,fd_right(i,j,x[n_grid-1],ylast,pa)); 00262 } 00263 00264 ix++; 00265 00266 } 00267 00268 // Compute correction by calling the linear solver 00269 00270 if (verbose>3) { 00271 std::cout << "Matrix: " << std::endl; 00272 for(size_t i=0;i<nvars;i++) { 00273 for(size_t j=0;j<nvars;j++) { 00274 std::cout << mat[i][j] << " "; 00275 } 00276 std::cout << std::endl; 00277 } 00278 std::cout << "Deviations:" << std::endl; 00279 for(size_t i=0;i<nvars;i++) { 00280 std::cout << rhs[i] << std::endl; 00281 } 00282 } 00283 00284 int ret=solver->solve(ix,mat,rhs,dy); 00285 00286 if (verbose>3) { 00287 std::cout << "Corrections:" << std::endl; 00288 for(size_t i=0;i<nvars;i++) { 00289 std::cout << dy[i] << std::endl; 00290 } 00291 std::cout << "Press a key and press enter to continue: " 00292 << std::endl; 00293 char ch; 00294 std::cin >> ch; 00295 } 00296 00297 // Apply correction and compute its size 00298 00299 double res=0.0; 00300 ix=0; 00301 00302 for(size_t igrid=0;igrid<n_grid;igrid++) { 00303 for(size_t ieq=0;ieq<n_eq;ieq++) { 00304 y[igrid][ieq]+=alpha*dy[ix]; 00305 res+=dy[ix]*dy[ix]; 00306 ix++; 00307 } 00308 } 00309 00310 if (verbose>0) { 00311 // Since we're in the o2scl namespace, we explicitly 00312 // specify std::sqrt() here 00313 std::cout << "ode_it_solve: " << it << " " << std::sqrt(res) << " " 00314 << tolf << std::endl; 00315 if (verbose>1) { 00316 char ch; 00317 std::cout << "Press a key and type enter to continue. "; 00318 std::cin >> ch; 00319 } 00320 } 00321 00322 // If the correction has become small enough, we're done 00323 if (std::sqrt(res)<=tolf) done=true; 00324 } 00325 00326 if (done==false) { 00327 O2SCL_ERR_RET("Exceeded number of iterations in solve().", 00328 o2scl::gsl_emaxiter); 00329 } 00330 00331 return 0; 00332 } 00333 00334 /// Default linear solver 00335 o2scl_linalg::linear_solver_hh<solver_vec_t,solver_mat_t> def_solver; 00336 00337 protected: 00338 00339 /// \name Storage for functions 00340 //@{ 00341 ode_it_funct<param_t,vec_t> *fl, *fr, *fd; 00342 //@} 00343 00344 /// Solver 00345 o2scl_linalg::linear_solver<solver_vec_t,solver_mat_t> *solver; 00346 00347 /** \brief Compute the derivatives of the LHS boundary conditions 00348 00349 This function computes \f$ \partial f_{left,\mathrm{ieq}} / \partial 00350 y_{\mathrm{ivar}} \f$ 00351 */ 00352 virtual double fd_left(size_t ieq, size_t ivar, double x, vec_t &y, 00353 param_t &pa) { 00354 double ret, dydx; 00355 00356 y[ivar]+=h; 00357 ret=(*fl)(ieq,x,y,pa); 00358 00359 y[ivar]-=h; 00360 ret-=(*fl)(ieq,x,y,pa); 00361 00362 ret/=h; 00363 return ret; 00364 } 00365 00366 /** \brief Compute the derivatives of the RHS boundary conditions 00367 00368 This function computes \f$ \partial f_{right,\mathrm{ieq}} / \partial 00369 y_{\mathrm{ivar}} \f$ 00370 */ 00371 virtual double fd_right(size_t ieq, size_t ivar, double x, vec_t &y, 00372 param_t &pa) { 00373 double ret, dydx; 00374 00375 y[ivar]+=h; 00376 ret=(*fr)(ieq,x,y,pa); 00377 00378 y[ivar]-=h; 00379 ret-=(*fr)(ieq,x,y,pa); 00380 00381 ret/=h; 00382 return ret; 00383 } 00384 00385 /** \brief Compute the finite-differenced part of the 00386 differential equations 00387 00388 This function computes \f$ \partial f_{\mathrm{ieq}} / \partial 00389 y_{\mathrm{ivar}} \f$ 00390 */ 00391 virtual double fd_derivs(size_t ieq, size_t ivar, double x, vec_t &y, 00392 param_t &pa) { 00393 double ret, dydx; 00394 00395 y[ivar]+=h; 00396 ret=(*fd)(ieq,x,y,pa); 00397 00398 y[ivar]-=h; 00399 ret-=(*fd)(ieq,x,y,pa); 00400 00401 ret/=h; 00402 00403 return ret; 00404 } 00405 00406 }; 00407 00408 #ifndef DOXYGENP 00409 } 00410 #endif 00411 00412 #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