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 <gsl/gsl_linalg.h> 00029 #include <gsl/gsl_ieee_utils.h> 00030 00031 #include <o2scl/misc.h> 00032 #include <o2scl/array.h> 00033 #include <o2scl/ovector_tlate.h> 00034 #include <o2scl/omatrix_tlate.h> 00035 #include <o2scl/uvector_tlate.h> 00036 #include <o2scl/umatrix_tlate.h> 00037 #include <o2scl/test_mgr.h> 00038 #include <o2scl/linear_solver.h> 00039 00040 #ifndef DOXYGENP 00041 namespace o2scl { 00042 #endif 00043 00044 /// Function class for ode_it_solve 00045 template<class vec_t=o2scl::ovector_base> class ode_it_funct { 00046 00047 public: 00048 00049 ode_it_funct() {} 00050 00051 virtual ~ode_it_funct() {} 00052 00053 /// Using \c x and \c y, return the value of function number \c ieq 00054 virtual double operator()(size_t ieq, double x, vec_t &y) { 00055 return 0.0; 00056 } 00057 00058 }; 00059 00060 /// Function pointer for ode_it_solve 00061 template<class vec_t=o2scl::ovector_base> 00062 class ode_it_funct_fptr : public ode_it_funct<vec_t> { 00063 00064 public: 00065 00066 virtual ~ode_it_funct_fptr() {} 00067 00068 /// Create using a function pointer 00069 ode_it_funct_fptr(double (*fp)(size_t, double, vec_t &)) { 00070 fptr=fp; 00071 } 00072 00073 /// Using \c x and \c y, return the value of function number \c ieq 00074 virtual double operator()(size_t ieq, double x, vec_t &y) { 00075 return fptr(ieq,x,y); 00076 } 00077 00078 protected: 00079 00080 ode_it_funct_fptr() {} 00081 00082 /// The function pointer 00083 double (*fptr)(size_t ieq, double x, vec_t &y); 00084 }; 00085 00086 /// Member function pointer for ode_it_solve 00087 template<class tclass, class vec_t=o2scl::ovector_base> 00088 class ode_it_funct_mfptr : public ode_it_funct<vec_t> { 00089 public: 00090 00091 virtual ~ode_it_funct_mfptr() {} 00092 00093 /// Create using a class instance and member function 00094 ode_it_funct_mfptr(tclass *tp, 00095 double (tclass::*fp)(size_t, double, vec_t &)) { 00096 tptr=tp; 00097 fptr=fp; 00098 } 00099 00100 /// Using \c x and \c y, return the value of function number \c ieq 00101 virtual double operator()(size_t ieq, double x, vec_t &y) { 00102 return (*tptr.*fptr)(ieq,x,y); 00103 } 00104 00105 protected: 00106 00107 ode_it_funct_mfptr() {} 00108 00109 /// The class pointer 00110 tclass *tptr; 00111 00112 /// The member function pointer 00113 double (tclass::*fptr)(size_t ieq, double x, vec_t &y); 00114 00115 }; 00116 00117 /** \brief ODE solver using a generic linear solver to solve 00118 finite-difference equations 00119 00120 \todo Implement as a child of ode_bv_solve ? 00121 \todo Max and average tolerance? 00122 \todo partial correction option? 00123 */ 00124 template <class func_t, class vec_t, class mat_t, class matrix_row_t, 00125 class solver_vec_t, class solver_mat_t> class ode_it_solve { 00126 00127 public: 00128 00129 ode_it_solve() { 00130 h=1.0e-4; 00131 niter=30; 00132 tolf=1.0e-8; 00133 verbose=0; 00134 } 00135 00136 virtual ~ode_it_solve() {} 00137 00138 /// Set level of output (default 0) 00139 int verbose; 00140 00141 /// Stepsize for finite differencing (default \f$ 10^{-4} \f$) 00142 double h; 00143 00144 /// Tolerance (default \f$ 10^{-8} \f$) 00145 double tolf; 00146 00147 /// Maximum number of iterations (default 30) 00148 size_t niter; 00149 00150 /// Set the linear solver 00151 int set_solver(o2scl_linalg::linear_solver<solver_vec_t,solver_mat_t> &ls) { 00152 solver=&ls; 00153 return 0; 00154 } 00155 00156 /// Solve \c derivs with boundary conditions \c left and \c right 00157 int solve(size_t ngrid, size_t neq, size_t nbleft, vec_t &x, mat_t &y, 00158 func_t &derivs, func_t &left, func_t &right, 00159 solver_mat_t &mat, solver_vec_t &rhs, solver_vec_t &dy) { 00160 00161 // Store the functions for simple derivatives 00162 fd=&derivs; 00163 fl=&left; 00164 fr=&right; 00165 00166 // Variable index 00167 size_t ix; 00168 00169 // Number of RHS boundary conditions 00170 size_t nbright=neq-nbleft; 00171 00172 // Number of variables 00173 size_t nvars=ngrid*neq; 00174 00175 bool done=false; 00176 for(size_t it=0;done==false && it<niter;it++) { 00177 00178 ix=0; 00179 00180 mat.set_all(0.0); 00181 00182 // Construct the entries corresponding to the LHS boundary 00183 for(size_t i=0;i<nbleft;i++) { 00184 matrix_row_t yk(y,0); 00185 rhs[ix]=-left(i,x[0],yk); 00186 for(size_t j=0;j<neq;j++) { 00187 int rxa=mat.set(ix,j,fd_left(i,j,x[0],yk)); 00188 } 00189 ix++; 00190 } 00191 00192 // Construct the matrix entries for the internal points 00193 for(size_t k=0;k<ngrid-1;k++) { 00194 size_t kp1=k+1; 00195 double tx=(x[kp1]+x[k])/2.0; 00196 double dx=x[kp1]-x[k]; 00197 matrix_row_t yk(y,k); 00198 matrix_row_t ykp1(y,k+1); 00199 00200 for(size_t i=0;i<neq;i++) { 00201 00202 rhs[ix]=y[k][i]-y[kp1][i]+(x[kp1]-x[k])* 00203 (derivs(i,tx,ykp1)+derivs(i,tx,yk))/2.0; 00204 00205 size_t lhs=k*neq; 00206 for(size_t j=0;j<neq;j++) { 00207 int rxb=mat.set(ix,lhs+j,-fd_derivs(i,j,tx,yk)*dx/2.0); 00208 int rxc=mat.set(ix,lhs+j+neq,-fd_derivs(i,j,tx,ykp1)*dx/2.0); 00209 if (i==j) { 00210 int rxd=mat.set(ix,lhs+j,mat(ix,lhs+j)-1.0); 00211 int rxe=mat.set(ix,lhs+j+neq,mat(ix,lhs+j+neq)+1.0); 00212 } 00213 } 00214 00215 ix++; 00216 00217 } 00218 } 00219 00220 // Construct the entries corresponding to the RHS boundary 00221 for(size_t i=0;i<nbright;i++) { 00222 matrix_row_t ylast(y,ngrid-1); 00223 size_t lhs=neq*(ngrid-1); 00224 00225 rhs[ix]=-right(i,x[ngrid-1],ylast); 00226 00227 for(size_t j=0;j<neq;j++) { 00228 int rxf=mat.set(ix,lhs+j,fd_right(i,j,x[ngrid-1],ylast)); 00229 } 00230 00231 ix++; 00232 00233 } 00234 00235 // Compute correction by calling the linear solver 00236 00237 if (verbose>3) { 00238 for(size_t i=0;i<nvars;i++) { 00239 for(size_t j=0;j<nvars;j++) { 00240 std::cout << mat[i][j] << " "; 00241 } 00242 std::cout << std::endl; 00243 } 00244 for(size_t i=0;i<nvars;i++) { 00245 std::cout << rhs[i] << std::endl; 00246 } 00247 } 00248 00249 int ret=solver->solve(ix,mat,rhs,dy); 00250 00251 if (verbose>3) { 00252 for(size_t i=0;i<nvars;i++) { 00253 std::cout << dy[i] << std::endl; 00254 } 00255 char ch; 00256 std::cin >> ch; 00257 } 00258 00259 // Apply correction and compute its size 00260 00261 double res=0.0; 00262 ix=0; 00263 00264 for(size_t igrid=0;igrid<ngrid;igrid++) { 00265 for(size_t ieq=0;ieq<neq;ieq++) { 00266 y[igrid][ieq]+=dy[ix]; 00267 res+=dy[ix]*dy[ix]; 00268 ix++; 00269 } 00270 } 00271 00272 if (verbose>0) { 00273 // Since we're in the o2scl namespace, we explicitly 00274 // specify std::sqrt() here 00275 std::cout << it << " " << std::sqrt(res) << " " 00276 << tolf << std::endl; 00277 } 00278 00279 // If the correction has become small enough, we're done 00280 if (std::sqrt(res)<=tolf) done=true; 00281 } 00282 00283 if (done==false) { 00284 O2SCL_ERR_RET("Exceeded number of iterations in solve().", 00285 o2scl::gsl_emaxiter); 00286 } 00287 00288 return 0; 00289 } 00290 00291 protected: 00292 00293 /// \name Storage for functions 00294 //@{ 00295 ode_it_funct<vec_t> *fl, *fr, *fd; 00296 //@} 00297 00298 /// Solver 00299 o2scl_linalg::linear_solver<solver_vec_t,solver_mat_t> *solver; 00300 00301 /// Compute the derivatives of the LHS boundary conditions 00302 double fd_left(size_t ieq, size_t ivar, double x, vec_t &y) { 00303 double ret, dydx; 00304 00305 y[ivar]+=h; 00306 ret=(*fl)(ieq,x,y); 00307 00308 y[ivar]-=h; 00309 ret-=(*fl)(ieq,x,y); 00310 00311 ret/=h; 00312 return ret; 00313 } 00314 00315 /// Compute the derivatives of the RHS boundary conditions 00316 double fd_right(size_t ieq, size_t ivar, double x, vec_t &y) { 00317 double ret, dydx; 00318 00319 y[ivar]+=h; 00320 ret=(*fr)(ieq,x,y); 00321 00322 y[ivar]-=h; 00323 ret-=(*fr)(ieq,x,y); 00324 00325 ret/=h; 00326 return ret; 00327 } 00328 00329 /// Compute the finite-differenced part of the differential equations 00330 double fd_derivs(size_t ieq, size_t ivar, double x, vec_t &y) { 00331 double ret, dydx; 00332 00333 y[ivar]+=h; 00334 ret=(*fd)(ieq,x,y); 00335 00336 y[ivar]-=h; 00337 ret-=(*fd)(ieq,x,y); 00338 00339 ret/=h; 00340 00341 return ret; 00342 } 00343 00344 }; 00345 00346 #ifndef DOXYGENP 00347 } 00348 #endif 00349 00350 #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