![]() |
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_MULTISHOOT_H 00025 #define O2SCL_ODE_BV_MULTISHOOT_H 00026 00027 #include <string> 00028 #include <o2scl/ovector_tlate.h> 00029 #include <o2scl/adapt_step.h> 00030 #include <o2scl/gsl_astep.h> 00031 #include <o2scl/ode_iv_solve.h> 00032 #include <o2scl/gsl_mroot_hybrids.h> 00033 00034 #ifndef DOXYGENP 00035 namespace o2scl { 00036 #endif 00037 00038 /** \brief Multishooting 00039 00040 This class is experimental. 00041 00042 \todo Improve documentation a little and create 00043 testing code 00044 */ 00045 template <class func_t=ode_funct<>, 00046 class vec_t=ovector_base, class alloc_vec_t=ovector, 00047 class alloc_t=ovector_alloc,class vec_int_t=ovector_int_base, 00048 class mat_t=omatrix> class ode_bv_multishoot { 00049 00050 public: 00051 00052 ode_bv_multishoot(){ 00053 oisp=&def_ois; 00054 mrootp=&def_mroot; 00055 } 00056 00057 virtual ~ode_bv_multishoot() { 00058 } 00059 00060 /* Solve the boundary value problem on a mesh */ 00061 virtual int solve(vec_t &mesh, int &n_func, vec_t &y_start, 00062 func_t &left_b, func_t &right_b, func_t &extra_b, 00063 func_t &derivs, vec_t &x_save, mat_t &y_save) { 00064 00065 /* Make copies of the input for later access */ 00066 this->l_mesh=&mesh; 00067 this->l_n_func=&n_func; 00068 this->l_y_start=&y_start; 00069 this->l_left_b=&left_b; 00070 this->l_right_b=&right_b; 00071 this->l_extra_b=&extra_b; 00072 this->l_derivs=&derivs; 00073 this->l_x_save=&x_save; 00074 this->l_y_save=&y_save; 00075 00076 /* vector of variables */ 00077 int nsolve=y_start.size(); 00078 ovector sx(nsolve),sy(nsolve); 00079 sx=y_start; 00080 00081 /* Equation solver */ 00082 mm_funct_mfptr<ode_bv_multishoot<func_t,vec_t, 00083 alloc_vec_t,alloc_t,vec_int_t> > 00084 mfm(this,&ode_bv_multishoot<func_t,vec_t, 00085 alloc_vec_t,alloc_t,vec_int_t>::solve_fun); 00086 00087 /* Run multishooting and save at the last step */ 00088 this->save=false; 00089 int ret=this->mrootp->msolve(nsolve,sx,mfm); 00090 this->save=true; 00091 solve_fun(nsolve,sx,sy); 00092 00093 return ret; 00094 } 00095 00096 /* Set initial value solver */ 00097 int set_iv(ode_iv_solve<func_t,vec_t,alloc_vec_t,alloc_t> &ois) { 00098 oisp=&ois; 00099 return 0; 00100 } 00101 00102 /* Set equation solver */ 00103 int set_mroot(mroot<mm_funct<> > &root) { 00104 mrootp=&root; 00105 return 0; 00106 } 00107 00108 /* Default initial value solver */ 00109 ode_iv_solve<func_t,vec_t,alloc_vec_t,alloc_t> def_ois; 00110 00111 /* Default equation solver */ 00112 gsl_mroot_hybrids<mm_funct<> > def_mroot; 00113 00114 #ifndef DOXYGEN_INTERNAL 00115 00116 protected: 00117 00118 /// The initial value solver 00119 ode_iv_solve<func_t,vec_t,alloc_vec_t,alloc_t> *oisp; 00120 00121 /// The equation solver 00122 gsl_mroot_hybrids<mm_funct<> > *mrootp; 00123 00124 /// Desc 00125 vec_t *l_mesh; 00126 /// Desc 00127 vec_t *l_y_start; 00128 /// Desc 00129 func_t *l_left_b; 00130 /// Desc 00131 func_t *l_right_b; 00132 /// Desc 00133 func_t *l_extra_b; 00134 /// Desc 00135 func_t *l_derivs; 00136 /// Desc 00137 int *l_n_func; 00138 /// Desc 00139 vec_t *l_x_save; 00140 /// Desc 00141 mat_t *l_y_save; 00142 /// Desc 00143 bool save; 00144 00145 /// Function to solve 00146 int solve_fun(size_t nv, const vec_t &sx, vec_t &sy) { 00147 00148 double xa,xb=0.0,h; 00149 ovector y((*this->l_n_func)),y2((*this->l_n_func)); 00150 00151 /* We update y_start in order that derivs knows 00152 all the values of parameters */ 00153 for(size_t i=0;i<(*this->l_y_start).size();i++) { 00154 (*this->l_y_start)[i]=sx[i]; 00155 } 00156 00157 /* A loop on each subinterval */ 00158 for(size_t k=0;k<(*this->l_mesh).size()-1;k++) { 00159 00160 xa=(*this->l_mesh)[k]; 00161 xb=(*this->l_mesh)[k+1]; 00162 h=(xb-xa)/100.0; 00163 00164 /* We load function's value at the left point of the sub-interval */ 00165 if (k==0) { 00166 (*this->l_left_b)(xa,(*this->l_n_func),sx,y); 00167 } else { 00168 for(int i=0;i<(*this->l_n_func);i++) 00169 y[i]=sx[i+(*this->l_n_func)*(1+k)]; 00170 } 00171 00172 if (this->save) { 00173 00174 /* iv_solver if we save */ 00175 int ngrid=((*this->l_x_save).size()-1)/((*this->l_mesh).size()-1)+1; 00176 ovector xxsave(ngrid); 00177 omatrix yysave(ngrid,(*this->l_n_func)); 00178 00179 if (k!=((*this->l_mesh).size()-2)) { 00180 xb=(*this->l_mesh)[k+1]-((*this->l_mesh)[k+1]- 00181 (*this->l_mesh)[k])/ngrid; 00182 } 00183 00184 this->oisp->solve_grid(xa,xb,h,(*this->l_n_func),y,ngrid, 00185 xxsave,yysave,(*this->l_derivs)); 00186 00187 for(int i=0;i<ngrid;i++) { 00188 (*this->l_x_save)[i+k*(ngrid)]=xxsave[i]; 00189 for(int j=0;j<(*this->l_n_func);j++) { 00190 (*this->l_y_save)[i+k*(ngrid)][j]=yysave[i][j]; 00191 } 00192 } 00193 00194 } else { 00195 00196 /* iv_solver if we don't save */ 00197 this->oisp->solve_final_value 00198 (xa,xb,h,(*this->l_n_func),y,y2,(*this->l_derivs)); 00199 00200 } 00201 00202 /* Then we load values at the end of sub-interval */ 00203 if (k==(*this->l_mesh).size()-2) { 00204 (*this->l_right_b)(xb,(*this->l_n_func),sx,y); 00205 } else { 00206 for(int i=0;i<(*this->l_n_func);i++) { 00207 y[i]=sx[i+(*this->l_n_func)*(2+k)]; 00208 } 00209 } 00210 00211 /* Now we take the difference */ 00212 for(int i=0;i<(*this->l_n_func);i++) { 00213 //sy[i+k*(*this->l_n_func)]=(y2[i]-y[i])/y[i]; 00214 sy[i+k*(*this->l_n_func)]= y2[i]-y[i]; 00215 } 00216 00217 } 00218 00219 /* Then load Extra boundary condition */ 00220 (*this->l_extra_b)(xb,(*this->l_n_func),sx,y); 00221 for(int i=0;i<(*this->l_n_func);i++) { 00222 sy[i+(int((*this->l_mesh).size()-1))*(*this->l_n_func)]=y[i]; 00223 } 00224 00225 return 0; 00226 } 00227 00228 #endif 00229 00230 }; 00231 00232 #ifndef DOXYGENP 00233 } 00234 #endif 00235 00236 #endif
Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).