![]() |
Object-oriented Scientific Computing Library: Version 0.910
|
00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2008, 2009, 2010, 2011, Julien Garaud and 00005 Andrew W. Steiner 00006 00007 This file is part of O2scl. 00008 00009 O2scl is free software; you can redistribute it and/or modify 00010 it under the terms of the GNU General Public License as published by 00011 the Free Software Foundation; either version 3 of the License, or 00012 (at your option) any later version. 00013 00014 O2scl is distributed in the hope that it will be useful, 00015 but WITHOUT ANY WARRANTY; without even the implied warranty of 00016 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00017 GNU General Public License for more details. 00018 00019 You should have received a copy of the GNU General Public License 00020 along with O2scl. If not, see <http://www.gnu.org/licenses/>. 00021 00022 ------------------------------------------------------------------- 00023 */ 00024 #ifndef O2SCL_ODE_BV_MSHOOT_H 00025 #define O2SCL_ODE_BV_MSHOOT_H 00026 00027 #include <o2scl/ode_bv_solve.h> 00028 00029 #ifndef DOXYGENP 00030 namespace o2scl { 00031 #endif 00032 00033 /** \brief Solve boundary-value ODE problems by multishooting 00034 with a generic nonlinear solver 00035 00036 This class is experimental. 00037 00038 Default template arguments 00039 - \c func_t - \ref ode_funct 00040 - \c mat_t - \ref omatrix_base 00041 - \c vec_t - \ref ovector_base 00042 - \c alloc_vec_t - \ref ovector 00043 - \c alloc_t - \ref ovector_alloc 00044 - \c vec_int_t - \ref ovector_int_base 00045 00046 \future Make a class which performs an iterative linear 00047 solver which uses sparse matrices like ode_it_solve? 00048 */ 00049 template<class func_t=ode_funct<>, class mat_t=omatrix_base, 00050 class vec_t=ovector_base, class alloc_vec_t=ovector, 00051 class alloc_t=ovector_alloc, class vec_int_t=ovector_int_base> 00052 class ode_bv_mshoot : public ode_bv_solve { 00053 00054 public: 00055 00056 ode_bv_mshoot() { 00057 oisp=&def_ois; 00058 mrootp=&def_mroot; 00059 mem_size=0; 00060 } 00061 00062 virtual ~ode_bv_mshoot() { 00063 if (mem_size>0) { 00064 ao.free(sy); 00065 ao.free(sy2); 00066 ao.free(syerr); 00067 ao.free(sdydx); 00068 ao.free(sdydx2); 00069 } 00070 } 00071 00072 /** \brief Solve the boundary-value problem and store the 00073 solution 00074 */ 00075 int solve_final_value(double h, size_t n, size_t n_bound, 00076 vec_t &x_bound, mat_t &y_bound, 00077 vec_int_t &index, func_t &derivs) { 00078 00079 // Make copies of the inputs for later access 00080 this->l_index=&index; 00081 this->l_nbound=n_bound; 00082 this->l_xbound=&x_bound; 00083 this->l_ybound=&y_bound; 00084 this->l_h=h; 00085 this->l_derivs=&derivs; 00086 this->l_n=n; 00087 00088 int lhs_unks=0, rhs_conts=0; 00089 this->l_lhs_unks=lhs_unks; 00090 00091 // Count the number of variables we'll need to solve for 00092 for (size_t i=0;i<n;i++) { 00093 if (index[i]<2) lhs_unks++; 00094 if (index[i]%2==1) rhs_conts++; 00095 } 00096 00097 // Make sure that the boundary conditions make sense 00098 if (lhs_unks!=rhs_conts) { 00099 O2SCL_ERR2_RET("Incorrect boundary conditions in ", 00100 "ode_bv_mshoot::solve()",gsl_einval); 00101 } 00102 00103 // The number of variables to solve for 00104 int n_solve=lhs_unks+n*(n_bound-2); 00105 00106 // Make space for the solution 00107 ovector tx(n_solve); 00108 00109 // Copy initial guess from y_bound 00110 size_t j=0; 00111 for (size_t i=0;i<n;i++) { 00112 if (index[i]<2) { 00113 tx[j]=y_bound[0][i]; 00114 j++; 00115 } 00116 } 00117 for(size_t k=1;k<n_bound-1;k++) { 00118 for(size_t i=0;i<n;i++) { 00119 tx[j]=y_bound[k][i]; 00120 j++; 00121 } 00122 } 00123 00124 // Allocate memory as needed 00125 if (n!=mem_size) { 00126 ao.allocate(sy,n); 00127 ao.allocate(sy2,n); 00128 ao.allocate(syerr,n); 00129 ao.allocate(sdydx,n); 00130 ao.allocate(sdydx2,n); 00131 mem_size=n; 00132 } 00133 00134 // The function object 00135 mm_funct_mfptr<ode_bv_mshoot<func_t,mat_t,vec_t,alloc_vec_t, 00136 alloc_t,vec_int_t> > 00137 mfm(this,&ode_bv_mshoot<func_t,mat_t,vec_t, 00138 alloc_vec_t,alloc_t,vec_int_t>::solve_fun); 00139 00140 // Solve 00141 int ret=this->mrootp->msolve(n_solve,tx,mfm); 00142 if (ret!=0) { 00143 O2SCL_ERR("Solver failed in ode_bv_mshoot::solve().",ret); 00144 } 00145 00146 return ret; 00147 } 00148 00149 /** \brief Solve the boundary-value problem and store the 00150 solution 00151 */ 00152 template<class mat_row_t> 00153 int solve_store(double h, size_t n, size_t n_bound, vec_t &x_bound, 00154 mat_t &y_bound, vec_int_t &index, size_t &n_sol, 00155 vec_t &x_sol, mat_t &y_sol, mat_t &dydx_sol, 00156 mat_t &yerr_sol, func_t &derivs) { 00157 00158 if (n_bound<2) { 00159 O2SCL_ERR2_RET("Not enough boundaries (must be at least two) in ", 00160 "ode_bv_mshoot::solve_store().",gsl_einval); 00161 } 00162 if (n_sol<n_bound) { 00163 O2SCL_ERR2_RET("Not enough room to store boundaries in ", 00164 "ode_bv_mshoot::solve_store().",gsl_einval); 00165 } 00166 00167 // Solve to fix x_bound and y_bound 00168 solve_final_value(h,n,n_bound,x_bound,y_bound,index,derivs); 00169 00170 // Find the indices which correspond to the boundaries 00171 ovector_int inxs(n_bound); 00172 for(size_t k=0;k<n_bound;k++) { 00173 inxs[k]=((size_t)(((double)k)/((double)n_bound)* 00174 (((double)(n_sol))-1.0+1.0e-12))); 00175 std::cout << k << " " << inxs[k] << " " << n_sol << std::endl; 00176 } 00177 // Double check that each interval has some space 00178 for(size_t k=1;k<n_bound-1;k++) { 00179 if (inxs[k]==inxs[k-1] || inxs[k]==inxs[k+1]) { 00180 O2SCL_ERR2_RET("Not enough room to store boundaries in ", 00181 "ode_bv_mshoot::solve_store().",gsl_einval); 00182 } 00183 } 00184 00185 // Now create the table 00186 for(size_t k=0;k<n_bound-1;k++) { 00187 size_t n_sol_tmp=inxs[k+1]; 00188 mat_row_t ystart(y_bound,k); 00189 00190 std::cout << "Old boundaries: " << inxs << std::endl; 00191 00192 oisp->template solve_store<mat_t,mat_row_t> 00193 (x_bound[k],x_bound[k+1],h,n,ystart, 00194 n_sol_tmp,x_sol,y_sol,dydx_sol,yerr_sol,derivs,inxs[k]); 00195 00196 std::cout << "New boundaries: " << n_sol_tmp << std::endl; 00197 00198 // If it didn't use all the space, shift the indexes 00199 // accordingly 00200 if (((int)n_sol_tmp)<inxs[k+1]) { 00201 for(size_t k2=k+1;k2<n_bound;k2++) { 00202 inxs[k2]=((size_t)(((double)k2-k-1)/((double)(n_bound-k-1))* 00203 (((double)(n_sol-n_sol_tmp)))))+n_sol_tmp-1; 00204 } 00205 } 00206 00207 std::cout << "New boundaries: " << inxs << std::endl; 00208 } 00209 00210 n_sol=inxs[n_bound-1]; 00211 00212 return 0; 00213 } 00214 00215 /** \brief Set initial value solver 00216 */ 00217 int set_iv(ode_iv_solve<func_t,vec_t,alloc_vec_t,alloc_t> &ois) { 00218 oisp=&ois; 00219 return 0; 00220 } 00221 00222 /** \brief Set the equation solver */ 00223 int set_mroot(mroot<mm_funct<> > &root) { 00224 mrootp=&root; 00225 return 0; 00226 } 00227 00228 /// The default initial value solver 00229 ode_iv_solve<func_t,vec_t,alloc_vec_t,alloc_t> def_ois; 00230 00231 /// The default equation solver 00232 gsl_mroot_hybrids<mm_funct<> > def_mroot; 00233 00234 #ifndef DOXYGEN_INTERNAL 00235 00236 protected: 00237 00238 /// The solver for the initial value problem 00239 ode_iv_solve<func_t,vec_t,alloc_vec_t,alloc_t> *oisp; 00240 00241 /// The equation solver 00242 mroot<mm_funct<> > *mrootp; 00243 00244 /// The index defining the boundary conditions 00245 vec_int_t *l_index; 00246 00247 /// Storage for the starting vector 00248 vec_t *l_xbound; 00249 00250 /// Storage for the ending vector 00251 mat_t *l_ybound; 00252 00253 /// Storage for the stepsize 00254 double l_h; 00255 00256 /// The functions to integrate 00257 func_t *l_derivs; 00258 00259 /// The number of functions 00260 size_t l_n; 00261 00262 /// The number of boundaries 00263 size_t l_nbound; 00264 00265 /// The number of unknowns on the LHS 00266 size_t l_lhs_unks; 00267 00268 /// \name Temporary storage for \ref solve_fun() 00269 //@{ 00270 alloc_vec_t sy, sy2, syerr, sdydx, sdydx2; 00271 //@} 00272 00273 /// Vector allocation object 00274 alloc_t ao; 00275 00276 /// Size of recent allocation 00277 size_t mem_size; 00278 00279 /// The shooting function to be solved by the multidimensional solver 00280 int solve_fun(size_t nv, const vec_t &tx, vec_t &ty) { 00281 00282 // Counters to index through function parameters tx and ty 00283 size_t j_x=0, j_y=0; 00284 00285 // Shorthand for the boundary specifications 00286 vec_t &xb=*this->l_xbound; 00287 mat_t &yb=*this->l_ybound; 00288 00289 // Set up the boundaries from the values in tx 00290 for(size_t k=0;k<l_nbound-1;k++) { 00291 if (k==0) { 00292 for(size_t i=0;i<this->l_n;i++) { 00293 if ((*this->l_index)[i]<2) { 00294 yb[k][i]=tx[j_x]; 00295 j_x++; 00296 } 00297 } 00298 } else { 00299 for(size_t i=0;i<this->l_n;i++) { 00300 yb[k][i]=tx[j_x]; 00301 j_x++; 00302 } 00303 } 00304 } 00305 00306 // Integrate between all of the boundaries 00307 for(size_t k=0;k<l_nbound-1;k++) { 00308 00309 double x0=xb[k]; 00310 double x1=xb[k+1]; 00311 00312 // Setup the start vector sy. This extra copying from l_ybound 00313 // to sy might be avoided by using a mat_row_t object to send 00314 // l_ybound directly to the solver? 00315 00316 for(size_t i=0;i<this->l_n;i++) { 00317 sy[i]=yb[k][i]; 00318 } 00319 00320 // Shoot across. We specify memory for the derivatives and 00321 // errors to prevent ode_iv_solve from having to continuously 00322 // allocate and reallocate it. 00323 oisp->solve_final_value(x0,x1,this->l_h,this->l_n,sy, 00324 sy2,*this->l_derivs); 00325 00326 if (k!=l_nbound-2) { 00327 // Construct equations for the internal RHS boundaries 00328 for(size_t i=0;i<this->l_n;i++) { 00329 ty[j_y]=sy2[i]-yb[k+1][i]; 00330 j_y++; 00331 } 00332 } else { 00333 // Construct the equations for the rightmost boundary 00334 for(size_t i=0;i<this->l_n;i++) { 00335 if ((*this->l_index)[i]%2==1) { 00336 double yright=yb[k+1][i]; 00337 if (yright==0.0) { 00338 ty[j_y]=sy2[i]; 00339 } else { 00340 ty[j_y]=(sy2[i]-yright)/yright; 00341 } 00342 j_y++; 00343 } 00344 } 00345 } 00346 00347 // End of loop to integrate between all boundaries 00348 } 00349 00350 return 0; 00351 } 00352 00353 #endif 00354 00355 }; 00356 00357 #ifndef DOXYGENP 00358 } 00359 #endif 00360 00361 #endif
Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).