Object-oriented Scientific Computing Library: Version 0.910
ode_bv_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_BV_SOLVE_H
00024 #define O2SCL_ODE_BV_SOLVE_H
00025 
00026 #include <string>
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   /** \brief Base class for boundary-value ODE solvers
00038 
00039       This class is experimental.
00040   */
00041   class ode_bv_solve {
00042 
00043     public:
00044 
00045     ode_bv_solve() {
00046       verbose=0;
00047     }
00048 
00049     virtual ~ode_bv_solve() {}
00050     
00051     /** \name Values for index arrays */
00052     //@{
00053     /// Unknown on both the left and right boundaries
00054     static const int unk=0;
00055     /// Known on the right boundary
00056     static const int right=1;
00057     /// Known on the left boundary
00058     static const int left=2;
00059     /// Known on both the left and right boundaries
00060     static const int both=3;
00061     //@}
00062 
00063     /// Set output level
00064     int verbose;
00065     
00066   };
00067   
00068   /** \brief Solve boundary-value ODE problems by shooting from one 
00069       boundary to the other
00070 
00071       This class is experimental.
00072 
00073       Default template arguments
00074       - \c func_t - \ref ode_funct 
00075       - \c vec_t - \ref ovector_base 
00076       - \c alloc_vec_t - \ref ovector 
00077       - \c alloc_t  - \ref ovector_alloc 
00078       - \c vec_int_t - \ref ovector_int_base
00079   */
00080   template<class func_t=ode_funct<>, class vec_t=ovector_base, 
00081     class alloc_vec_t=ovector, class alloc_t=ovector_alloc, 
00082     class vec_int_t=ovector_int_base> class ode_bv_shoot : 
00083   public ode_bv_solve {
00084     
00085     public:
00086     
00087     ode_bv_shoot() {
00088       oisp=&def_ois;
00089       mrootp=&def_mroot;
00090       mem_size=0;
00091     }
00092     
00093     virtual ~ode_bv_shoot() {
00094       if (mem_size>0) {
00095         ao.free(sy);
00096         ao.free(sy2);
00097         ao.free(syerr);
00098         ao.free(sdydx);
00099       }
00100     }
00101 
00102     /// Allocate internal storage
00103     void allocate(size_t n) {
00104       if (n!=mem_size) {
00105         ao.allocate(sy,n);
00106         ao.allocate(sy2,n);
00107         ao.allocate(syerr,n);
00108         ao.allocate(sdydx,n);
00109         mem_size=n;
00110       }
00111       return;
00112     }
00113     
00114     /** \brief Solve the boundary-value problem and store the 
00115         solution
00116 
00117         Given the \c n initial values of the functions in \c ystart,
00118         this function integrates the ODEs specified in \c derivs over
00119         the interval from \c x0 to \c x1 with an initial stepsize of
00120         \c h. The final values of the function are given in \c yend,
00121         the derivatives in \c dydx_end, and the associated errors are
00122         given in \c yerr. The initial values of \c yend and \c yerr
00123         are ignored.
00124 
00125     */
00126     int solve_final_value(double x0, double x1, double h, size_t n,
00127                           vec_t &ystart, vec_t &yend, vec_int_t &index,
00128                           vec_t &yerr, vec_t &dydx_end, func_t &derivs) {
00129       
00130       // Get pointers to inputs for access by solve function
00131       this->l_index=&index;
00132       this->l_ystart=&ystart;
00133       this->l_yend=&yend;
00134       this->l_yerr=&yerr;
00135       this->l_dydx_end=&dydx_end;
00136       this->l_x0=x0;
00137       this->l_x1=x1;
00138       this->l_h=h;
00139       this->l_derivs=&derivs;
00140       this->l_n=n;
00141   
00142       // lhs_unks is the number of unknowns on the LHS
00143       // rhs_conts is the number of constraints on the RHS
00144       int lhs_unks=0, rhs_conts=0;
00145 
00146       // Count the number of variables we'll need to solve for
00147       for (size_t i=0;i<n;i++) {
00148         if (index[i]<2) lhs_unks++;
00149         if (index[i]%2==1) rhs_conts++;
00150       }
00151 
00152       // Make sure that the boundary conditions make sense
00153       if (lhs_unks!=rhs_conts) {
00154         O2SCL_ERR2_RET("Incorrect boundary conditions in ",
00155                        "ode_bv_shoot::solve()",gsl_einval);
00156       } 
00157 
00158       // Make space for the solution
00159       ovector tx(lhs_unks);
00160 
00161       // Copy initial guesses from ystart
00162       int j=0;
00163       for (size_t i=0;i<n;i++) {
00164         if (index[i]<2) {
00165           tx[j]=ystart[i];
00166           j++;
00167         }
00168       }
00169 
00170       allocate(n);
00171   
00172       // The function object
00173       mm_funct_mfptr<ode_bv_shoot<func_t,vec_t,alloc_vec_t,
00174       alloc_t,vec_int_t> > 
00175       mfm(this,&ode_bv_shoot<func_t,vec_t,
00176           alloc_vec_t,alloc_t,vec_int_t>::solve_fun);
00177 
00178       // Solve
00179       int ret=this->mrootp->msolve(lhs_unks,tx,mfm);
00180       if (ret!=0) {
00181         O2SCL_ERR("Solver failed in ode_bv_shoot::solve().",ret);
00182       }
00183 
00184       // Copy the solution back to ystart
00185       j=0;
00186       for (size_t i=0;i<n;i++) {
00187         if (index[i]<2) {
00188           ystart[i]=tx[j];
00189           j++;
00190         }
00191       }
00192 
00193       return 0;
00194     }
00195 
00196     /** \brief Solve the boundary-value problem and store the 
00197         solution
00198      */
00199     template<class mat_t, class mat_row_t> 
00200     int solve_store(double x0, double x1, double h, size_t n,
00201                     vec_t &ystart, vec_t &yend, 
00202                     vec_int_t &index, size_t &n_sol, vec_t &x_sol, 
00203                     mat_t &y_sol, mat_t &yerr_sol, mat_t &dydx_sol, 
00204                     func_t &derivs) {
00205 
00206       mat_row_t yerr(yerr_sol,0);
00207       mat_row_t dydx_end(dydx_sol,0);
00208 
00209       // Solve for the final value
00210       solve_final_value(x0,x1,h,n,ystart,yend,index,yerr,dydx_end,derivs);
00211 
00212       // Copy ystart to y_sol[0]
00213       for(size_t i=0;i<n;i++) {
00214         y_sol[0][i]=ystart[i];
00215       }
00216 
00217       // Evaluate the function one more time to create the table
00218       oisp->solve_store<mat_t,mat_row_t>(this->l_x0,this->l_x1,this->l_h,
00219                                          this->l_n,n_sol,x_sol,
00220                                          y_sol,yerr_sol,dydx_sol,derivs);
00221 
00222       // Copy the stored solution back to ystart and yend
00223       for (size_t i=0;i<n;i++) {
00224         if (index[i]%2==0) {
00225           yend[i]=y_sol[n_sol-1][i];
00226         }
00227         ystart[i]=y_sol[0][i];
00228       }
00229   
00230       return 0;
00231     }
00232     
00233     /** \brief Set initial value solver
00234      */
00235     int set_iv(ode_iv_solve<func_t,vec_t,alloc_vec_t,alloc_t> &ois) {
00236       oisp=&ois;
00237       return 0;
00238     }
00239     
00240     /** \brief Set the equation solver */
00241     int set_mroot(mroot<mm_funct<> > &root) {
00242       mrootp=&root;
00243       return 0;
00244     }
00245 
00246     /// The default initial value solver
00247     ode_iv_solve<func_t,vec_t,alloc_vec_t,alloc_t> def_ois;
00248 
00249     /// The default equation solver
00250     gsl_mroot_hybrids<mm_funct<> > def_mroot;
00251 
00252 #ifndef DOXYGEN_INTERNAL
00253 
00254     protected:
00255     
00256     /// The solver for the initial value problem
00257     ode_iv_solve<func_t,vec_t,alloc_vec_t,alloc_t> *oisp;
00258 
00259     /// The equation solver
00260     mroot<mm_funct<> > *mrootp;
00261 
00262     /// The index defining the boundary conditions
00263     vec_int_t *l_index;
00264 
00265     /// Storage for the starting vector
00266     vec_t *l_ystart;
00267 
00268     /// Storage for the ending vector
00269     vec_t *l_yend;
00270 
00271     /// Storage for the starting vector
00272     vec_t *l_yerr;
00273 
00274     /// Storage for the ending vector
00275     vec_t *l_dydx_end;
00276 
00277     /// Storage for the starting point
00278     double l_x0;
00279 
00280     /// Storage for the ending abcissa
00281     double l_x1;
00282 
00283     /// Storage for the stepsize
00284     double l_h;
00285 
00286     /// The functions to integrate
00287     func_t *l_derivs;
00288 
00289     /// The number of functions
00290     size_t l_n;
00291 
00292     /// \name Temporary storage for \ref solve_fun()
00293     //@{
00294     alloc_vec_t sy, sy2, syerr, sdydx;
00295     //@}
00296     
00297     /// Vector allocation object
00298     alloc_t ao;
00299 
00300     /// Size of recent allocation
00301     size_t mem_size;
00302 
00303     /// The shooting function to be solved by the multidimensional solver
00304     int solve_fun(size_t nv, const vec_t &tx, vec_t &ty) {
00305 
00306       int j;
00307 
00308       // Create leftmost point from combination of 
00309       // starting array and proposed solution
00310       j=0;
00311       for(size_t i=0;i<this->l_n;i++) {
00312         if ((*this->l_index)[i]<2) {
00313           sy[i]=tx[j];
00314           j++;
00315         } else {
00316           sy[i]=(*this->l_ystart)[i];
00317         }
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(this->l_x0,this->l_x1,this->l_h,
00324                               this->l_n,sy,sy2,*this->l_yerr,
00325                               *this->l_dydx_end,*this->l_derivs);
00326   
00327       j=0;
00328       for(size_t i=0;i<this->l_n;i++) {
00329         if ((*this->l_index)[i]%2==1) {
00330           // Construct the equations from the rightmost point
00331           if ((*this->l_yend)[i]==0.0) {
00332             ty[j]=sy2[i];
00333           } else {
00334             ty[j]=(sy2[i]-(*this->l_yend)[i])/(*this->l_yend)[i];
00335           }
00336           j++;
00337         } else {
00338           // Otherwise copy the final values from y2 to *l_yend
00339           (*this->l_yend)[i]=sy2[i];
00340         }
00341       }
00342   
00343       return 0;
00344     }
00345 
00346 #endif
00347 
00348   };
00349 
00350   /** \brief Solve boundary-value ODE problems by shooting from one 
00351       boundary to the other on a grid
00352 
00353       This class is experimental.
00354 
00355       Default template arguments
00356       - \c func_t - \ref ode_funct 
00357       - \c vec_t - \ref ovector_base 
00358       - \c alloc_vec_t - \ref ovector 
00359       - \c alloc_t  - \ref ovector_alloc 
00360       - \c vec_int_t - \ref ovector_int_base
00361   */
00362   template<class mat_t=omatrix, class mat_row_t=omatrix_row, 
00363     class func_t=ode_funct<>, class vec_t=ovector_base, 
00364     class alloc_vec_t=ovector, class alloc_t=ovector_alloc, 
00365     class vec_int_t=ovector_int_base> class ode_bv_shoot_grid : 
00366   public ode_bv_shoot<func_t,vec_t,alloc_vec_t,alloc_t,vec_int_t> {
00367     
00368     public:
00369     
00370     ode_bv_shoot_grid() {
00371     }
00372     
00373     virtual ~ode_bv_shoot_grid() {
00374     }
00375     
00376     /// Desc
00377     int solve_grid(double x0, double x1, double h, size_t n, 
00378                    vec_t &ystart, vec_t &yend, vec_int_t &index,
00379                    size_t nsol, vec_t &xsol, mat_t &ysol, 
00380                    mat_t &err_sol, mat_t &dydx_sol, func_t &derivs) {
00381 
00382       // lhs_unks is the number of unknowns on the LHS
00383       // rhs_conts is the number of constraints on the RHS
00384       int lhs_unks=0, rhs_conts=0;
00385 
00386       // Count the number of variables we'll need to solve for
00387       for (size_t i=0;i<n;i++) {
00388         if (index[i]<2) lhs_unks++;
00389         if (index[i]%2==1) rhs_conts++;
00390       }
00391 
00392       // Make sure that the boundary conditions make sense
00393       if (lhs_unks!=rhs_conts) {
00394         O2SCL_ERR2_RET("Incorrect boundary conditions in ",
00395                        "ode_bv_shoot_grid::solve_grid()",gsl_einval);
00396       } 
00397 
00398       // Get pointers to inputs for access by solve function
00399       this->l_n=n;
00400       this->l_index=&index;
00401       this->l_derivs=&derivs;
00402       this->l_ystart=&ystart;
00403       this->l_yend=&yend;
00404       this->l_h=h;
00405       l_nsol=nsol;
00406       l_xsol=&xsol;
00407       l_ysol=&ysol;
00408       l_dydxsol=&dydx_sol;
00409       l_errsol=&err_sol;
00410 
00411       // Make space for the solution
00412       ovector tx(lhs_unks);
00413 
00414       // Copy initial guesses from ysol
00415       int j=0;
00416       for (size_t i=0;i<n;i++) {
00417         if (index[i]<2) {
00418           tx[j]=ystart[i];
00419           j++;
00420         }
00421       }
00422 
00423       // Allocate internal storage
00424       this->allocate(n);
00425 
00426       // The function object
00427       mm_funct_mfptr<ode_bv_shoot_grid<mat_t,mat_row_t,
00428       func_t,vec_t,alloc_vec_t,alloc_t,vec_int_t> > 
00429       mfm(this,&ode_bv_shoot_grid<mat_t,mat_row_t,func_t,vec_t,
00430           alloc_vec_t,alloc_t,vec_int_t>::solve_grid_fun);
00431   
00432       // Solve
00433       int ret=this->mrootp->msolve(lhs_unks,tx,mfm);
00434       if (ret!=0) {
00435         O2SCL_ERR("Solver failed in ode_bv_shoot_grid::solve_grid().",ret);
00436       }
00437 
00438       // Copy the solution back to ysol
00439       j=0;
00440       for (size_t i=0;i<n;i++) {
00441         if (index[i]<2) {
00442           ysol[0][i]=tx[j];
00443           j++;
00444         }
00445       }
00446 
00447       return 0;
00448     }
00449 
00450 #ifndef DOXYGEN_INTERNAL
00451 
00452     protected:
00453     
00454     /// Desc
00455     size_t l_nsol;
00456 
00457     /// Desc
00458     vec_t *l_xsol;
00459 
00460     /// Desc
00461     mat_t *l_ysol;
00462 
00463     /// Desc
00464     mat_t *l_dydxsol;
00465 
00466     /// Desc
00467     mat_t *l_errsol;
00468 
00469     /// The shooting function to be solved by the multidimensional solver    
00470     int solve_grid_fun(size_t nv, const vec_t &tx, vec_t &ty) {
00471 
00472       int j;
00473 
00474       // Create leftmost point from combination of 
00475       // starting array and proposed solution
00476       j=0;
00477       for(size_t i=0;i<this->l_n;i++) {
00478         if ((*this->l_index)[i]<2) {
00479           (*l_ysol)[0][i]=tx[j];
00480           j++;
00481         }
00482       }
00483 
00484       // Perform the solution
00485       this->oisp->template solve_grid<mat_t,mat_row_t>
00486         (this->l_h,this->l_n,l_nsol,*l_xsol,*l_ysol,*l_errsol,
00487          *l_dydxsol,*this->l_derivs);
00488 
00489       // The last vector
00490       mat_row_t yend2(*l_ysol,l_nsol-1);
00491 
00492       j=0;
00493       for(size_t i=0;i<this->l_n;i++) {
00494         if ((*this->l_index)[i]%2==1) {
00495           // Construct the equations from the rightmost point
00496           if ((*this->l_yend)[i]==0.0) {
00497             ty[j]=yend2[i];
00498           } else {
00499             ty[j]=((*this->l_yend)[i]-yend2[i])/(*this->l_yend)[i];
00500           }
00501           j++;
00502         }
00503       }
00504 
00505       return 0;
00506     }
00507 
00508 #endif
00509 
00510   };
00511 
00512 #ifndef DOXYGENP
00513 }
00514 #endif
00515 
00516 #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.