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