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_BV_SOLVE_H 00024 #define O2SCL_ODE_BV_SOLVE_H 00025 #include <string> 00026 #include <o2scl/collection.h> 00027 #include <o2scl/ovector_tlate.h> 00028 #include <o2scl/adapt_step.h> 00029 #include <o2scl/gsl_astep.h> 00030 #include <o2scl/ode_iv_solve.h> 00031 #include <o2scl/gsl_mroot_hybrids.h> 00032 00033 #ifndef DOXYGENP 00034 namespace o2scl { 00035 #endif 00036 00037 /** 00038 \brief Solve boundary-value ODE problems with constant 00039 boundary conditions [abstract base] 00040 */ 00041 template<class param_t, class func_t=ode_funct<param_t>, 00042 class vec_t=ovector_base, 00043 class alloc_vec_t=ovector, class alloc_t=ovector_alloc, 00044 class vec_int_t=ovector_int_base> class ode_bv_solve { 00045 00046 public: 00047 00048 ode_bv_solve() { 00049 verbose=0; 00050 } 00051 00052 virtual ~ode_bv_solve() {} 00053 00054 /** 00055 \brief Solve the boundary-value problem 00056 */ 00057 virtual int solve(double x0, double x1, double h, size_t n, 00058 vec_t &ystart, vec_t ¥d, vec_int_t &index, 00059 param_t &pa, func_t &derivs)=0; 00060 00061 /** \name Values for the index array */ 00062 //@{ 00063 /// Unknown on both the left and right boundaries 00064 static const int unk=0; 00065 /// Known on the right boundary 00066 static const int right=1; 00067 /// Known on the left boundary 00068 static const int left=2; 00069 /// Known on both the left and right boundaries 00070 static const int both=3; 00071 //@} 00072 00073 /// Set output level 00074 int verbose; 00075 00076 }; 00077 00078 /** 00079 \brief Solve boundary-value ODE problems by shooting 00080 00081 \future Create a solution table as in ode_iv_solve 00082 */ 00083 template<class param_t, class func_t=ode_funct<param_t>, 00084 class vec_t=ovector_base, class alloc_vec_t=ovector, 00085 class alloc_t=ovector_alloc, class vec_int_t=ovector_int_base> 00086 class ode_bv_shoot : 00087 public ode_bv_solve<param_t,func_t,vec_t,alloc_vec_t,alloc_t,vec_int_t> { 00088 00089 public: 00090 00091 ode_bv_shoot() { 00092 oisp=&def_ois; 00093 mrootp=&def_mroot; 00094 } 00095 00096 virtual ~ode_bv_shoot() {} 00097 00098 /** 00099 \brief Solve the boundary-value problem 00100 */ 00101 virtual int solve(double x0, double x1, double h, size_t n, 00102 vec_t &ystart, vec_t ¥d, vec_int_t &index, 00103 param_t &pa, func_t &derivs) { 00104 00105 // Make copies of the inputs for later access 00106 this->l_index=&index; 00107 this->l_ystart=&ystart; 00108 this->l_yend=¥d; 00109 this->l_x0=x0; 00110 this->l_x1=x1; 00111 this->l_h=h; 00112 this->l_derivs=&derivs; 00113 this->l_n=n; 00114 00115 int nsolve=0, nsolve2=0, j; 00116 00117 // Count the number of variables we'll need to solve for: 00118 for (size_t i=0;i<n;i++) { 00119 if (index[i]<2) nsolve++; 00120 if (index[i]%2==1) nsolve2++; 00121 } 00122 00123 // Make sure that the boundary conditions make sense: 00124 if (nsolve!=nsolve2) { 00125 O2SCL_ERR2_RET("Incorrect boundary conditions in ", 00126 "ode_bv_shoot::solve()",gsl_einval); 00127 } 00128 00129 // Make space for the solution 00130 ovector sx(nsolve), sy(nsolve); 00131 00132 // Copy initial guesses from ystart 00133 j=0; 00134 for (size_t i=0;i<n;i++) { 00135 if (index[i]<2) { 00136 sx[j]=ystart[i]; 00137 j++; 00138 } 00139 } 00140 00141 // Solve 00142 mm_funct_mfptr<ode_bv_shoot<param_t,func_t,vec_t,alloc_vec_t, 00143 alloc_t,vec_int_t>,param_t > 00144 mfm(this,&ode_bv_shoot<param_t,func_t,vec_t, 00145 alloc_vec_t,alloc_t,vec_int_t>::solve_fun); 00146 int ret=this->mrootp->msolve(nsolve,sx,pa,mfm); 00147 solve_fun(nsolve,sx,sy,pa); 00148 00149 // Copy the solution back to ystart and yend 00150 j=0; 00151 for (size_t i=0;i<n;i++) { 00152 if (index[i]<2) { 00153 ystart[i]=sx[j]; 00154 j++; 00155 } 00156 } 00157 00158 if (ret!=0) { 00159 O2SCL_ERR("Solver failed in ode_bv_shoot::solve().",ret); 00160 } 00161 return ret; 00162 } 00163 00164 /** 00165 \brief Set initial value solver 00166 */ 00167 int set_iv(ode_iv_solve<param_t,func_t,vec_t,alloc_vec_t,alloc_t> &ois) { 00168 oisp=&ois; 00169 return 0; 00170 } 00171 00172 /** \brief Set the equation solver */ 00173 int set_mroot(mroot<param_t,mm_funct<param_t> > &root) { 00174 mrootp=&root; 00175 return 0; 00176 } 00177 00178 /// The default initial value solver 00179 ode_iv_solve<param_t,func_t,vec_t,alloc_vec_t,alloc_t> def_ois; 00180 00181 /// The default equation solver 00182 gsl_mroot_hybrids<param_t,mm_funct<param_t> > def_mroot; 00183 00184 #ifndef DOXYGEN_INTERNAL 00185 00186 protected: 00187 00188 /// The solver for the initial value problem 00189 ode_iv_solve<param_t,func_t,vec_t,alloc_vec_t,alloc_t> *oisp; 00190 00191 /// The equation solver 00192 mroot<param_t,mm_funct<param_t> > *mrootp; 00193 00194 /// The index defining the boundary conditions 00195 vec_int_t *l_index; 00196 00197 /// Storage for the starting vector 00198 vec_t *l_ystart; 00199 00200 /// Storage for the ending vector 00201 vec_t *l_yend; 00202 00203 /// Storage for the starting point 00204 double l_x0; 00205 00206 /// Storage for the ending abcissa 00207 double l_x1; 00208 00209 /// Storage for the stepsize 00210 double l_h; 00211 00212 /// The functions to integrate 00213 func_t *l_derivs; 00214 00215 /// The number of functions 00216 size_t l_n; 00217 00218 /// The shooting function to be solved by the multidimensional solver 00219 int solve_fun(size_t nv, const vec_t &sx, vec_t &sy, 00220 param_t &pa) { 00221 int j; 00222 ovector y(this->l_n), y2(this->l_n); 00223 00224 // Create leftmost point from combination of 00225 // starting array and proposed solution 00226 j=0; 00227 for(size_t i=0;i<this->l_n;i++) { 00228 if ((*this->l_index)[i]<2) { 00229 y[i]=sx[j]; 00230 j++; 00231 } else { 00232 y[i]=(*this->l_ystart)[i]; 00233 } 00234 } 00235 00236 // Shoot across 00237 this->oisp->solve_final_value(this->l_x0,this->l_x1,this->l_h, 00238 this->l_n,y,y2,pa,*this->l_derivs); 00239 00240 j=0; 00241 for(size_t i=0;i<this->l_n;i++) { 00242 if ((*this->l_index)[i]%2==1) { 00243 // Construct the equations from the rightmost point 00244 if ((*this->l_yend)[i]==0.0) { 00245 sy[j]=y2[i]; 00246 } else { 00247 sy[j]=(y2[i]-(*this->l_yend)[i])/(*this->l_yend)[i]; 00248 } 00249 j++; 00250 } else { 00251 // Otherwise copy the final values from y2 to *l_yend 00252 (*this->l_yend)[i]=y2[i]; 00253 } 00254 } 00255 00256 return 0; 00257 } 00258 00259 #endif 00260 00261 }; 00262 00263 #ifndef DOXYGENP 00264 } 00265 #endif 00266 00267 #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