Object-oriented Scientific Computing Library: Version 0.910
gsl_bsimp.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 /* ode-initval/bsimp.c
00024  *
00025  * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2004 Gerard Jungman
00026  *
00027  * This program is free software; you can redistribute it and/or modify
00028  * it under the terms of the GNU General Public License as published by
00029  * the Free Software Foundation; either version 3 of the License, or (at
00030  * your option) any later version.
00031  *
00032  * This program is distributed in the hope that it will be useful, but
00033  * WITHOUT ANY WARRANTY; without even the implied warranty of
00034  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00035  * General Public License for more details.
00036  *
00037  * You should have received a copy of the GNU General Public License
00038  * along with this program; if not, write to the Free Software
00039  * Foundation, Inc., 51 Franklin Street, Fifth Floor, 
00040  * Boston, MA 02110-1301, USA.
00041  */
00042 
00043 #ifndef O2SCL_GSL_BSIMP_H
00044 #define O2SCL_GSL_BSIMP_H
00045 
00046 #include <gsl/gsl_math.h>
00047 
00048 #include <o2scl/err_hnd.h>
00049 #include <o2scl/ode_step.h>
00050 #include <o2scl/ode_jac_funct.h>
00051 #include <o2scl/vec_arith.h>
00052 #include <o2scl/lu.h>
00053 #include <o2scl/permutation.h>
00054 
00055 #ifndef DOXYGENP
00056 namespace o2scl {
00057 #endif
00058 
00059   /** \brief Bulirsch-Stoer implicit ODE stepper (GSL)
00060       
00061       Bader-Deuflhard implicit extrapolative stepper (\ref Bader83). 
00062 
00063       \note The variable <tt>h_next</tt> was defined in the original
00064       GSL version has been removed here, as it was unused by the
00065       stepper routine. 
00066 
00067       \note At the moment, this class retains the GSL approach to
00068       handling non-integer return values in the user-specified
00069       derivative function. If the user-specified function returns \c
00070       gsl_efailed, then the stepper attempts to decrease the stepsize
00071       to continue. If the user-specified function returns a non-zero
00072       value other than \c gsl_efailed, or if the Jacobian evaluation
00073       returns any non-zero value, then the stepper aborts and returns
00074       the error value without calling the error handler. This behavior
00075       may change in the future.
00076 
00077       \future Create an example with a stiff diff eq. which requires
00078       this kind of stepper
00079 
00080       \future More detailed documentation about the algorithm
00081 
00082       \future Figure out if this should be a child of ode_step or
00083       adapt_step. The function step_local() is actually its own ODE
00084       stepper and could be reimplemented as an object of type ode_step.
00085 
00086       \future I don't like setting yerr to GSL_POSINF, there should
00087       be a better way to force an adaptive stepper which is calling
00088       this stepper to readjust the stepsize.
00089 
00090       \future The functions deuf_kchoice() and compute_weights() can
00091       be moved out of this header file
00092 
00093       \future Rework internal arrays as uvectors?
00094 
00095       \future Rework linear solver to be amenable to using
00096       a sparse matrix solver
00097   */
00098   template<class func_t, class jac_func_t, 
00099     class vec_t=ovector_base, class alloc_vec_t=ovector, 
00100     class alloc_t=ovector_alloc, class mat_t=omatrix_base,
00101     class alloc_mat_t=omatrix, class mat_alloc_t=omatrix_alloc> 
00102     class gsl_bsimp {
00103 
00104 #ifndef DOXYGEN_INTERAL
00105 
00106   protected:
00107   
00108   /// Size of allocated vectors
00109   size_t dim;
00110 
00111   /// Memory allocator for objects of type \c alloc_vec_t
00112   alloc_t ao;
00113 
00114   /// Memory allocator for objects of type \c alloc_mat_t
00115   mat_alloc_t mo;
00116 
00117   /// Function specifying derivatives
00118   func_t *funcp;
00119 
00120   /// Jacobian
00121   jac_func_t *jfuncp;
00122 
00123   /// Number of entries in the Bader-Deuflhard extrapolation sequence
00124   static const int sequence_count=8;
00125 
00126   /// Largest index of the Bader-Deuflhard extrapolation sequence
00127   static const int sequence_max=7;
00128 
00129   /// Workspace for extrapolation 
00130   gsl_matrix *d;                
00131 
00132   /// Workspace for linear system matrix 
00133   gsl_matrix *a_mat;            
00134 
00135   /** \brief Workspace for extrapolation 
00136         
00137       (This state variable was named 'x' in GSL.)
00138   */
00139   double ex_wk[sequence_max];       
00140 
00141   /// \name State info
00142   //@{
00143   size_t k_current;
00144   size_t k_choice;
00145   double eps;
00146   //@}
00147 
00148   /// \name Workspace for extrapolation step 
00149   //@{
00150   double *yp;
00151   double *y_save;
00152   double *yerr_save;
00153   double *y_extrap_save;
00154   alloc_vec_t y_extrap_sequence;
00155   double *extrap_work;
00156   alloc_vec_t dfdt;
00157   alloc_vec_t y_temp;
00158   alloc_vec_t delta_temp;
00159   double *weight;
00160   alloc_mat_t dfdy;
00161   //@}
00162     
00163   /// \name Workspace for the basic stepper
00164   //@{
00165   alloc_vec_t rhs_temp;
00166   double *delta;
00167   //@}
00168 
00169   /// Order of last step 
00170   size_t order;
00171 
00172   /// Compute weights
00173   int compute_weights(const double y[], double w[], size_t n) {
00174     size_t i;
00175     // w[i] is 1 if y[i] is zero and the absolute value of y[i]
00176     // otherwise
00177     for (i = 0; i < n; i++) {
00178       double u = fabs(y[i]);
00179       w[i] = (u > 0.0) ? u : 1.0;
00180     }
00181     return 0;
00182   }
00183   
00184   /** \brief Calculate a choice for the "order" of the method, using the
00185       Deuflhard criteria. 
00186 
00187       Used in the allocate() function.
00188   */
00189   size_t deuf_kchoice(double eps2, size_t dimension) {
00190 
00191     const double safety_f = 0.25;
00192     const double small_eps = safety_f * eps2;
00193       
00194     double a_work[sequence_count];
00195     double alpha[sequence_max][sequence_max];
00196     
00197     /* Bader-Deuflhard extrapolation sequence */
00198     static const int bd_sequence[sequence_count] =
00199     { 2, 6, 10, 14, 22, 34, 50, 70 };
00200       
00201     int i, k;
00202     
00203     a_work[0] = bd_sequence[0] + 1.0;
00204     
00205     for (k = 0; k < sequence_max; k++) {
00206       a_work[k + 1] = a_work[k] + bd_sequence[k + 1];
00207     }
00208 
00209     for (i = 0; i < sequence_max; i++) {
00210       alpha[i][i] = 1.0;
00211       for (k = 0; k < i; k++) {
00212         const double tmp1 = a_work[k + 1] - a_work[i + 1];
00213         const double tmp2 = (a_work[i + 1] - a_work[0] + 1.0) * 
00214           (2 * k + 1);
00215         alpha[k][i] = pow (small_eps, tmp1 / tmp2);
00216       }
00217     }
00218 
00219     a_work[0] += dimension;
00220 
00221     for (k = 0; k < sequence_max; k++) {
00222       a_work[k + 1] = a_work[k] + bd_sequence[k + 1];
00223     }
00224 
00225     for (k = 0; k < sequence_max - 1; k++) {
00226       if (a_work[k + 2] > a_work[k + 1] * alpha[k][k + 1]) {
00227         break;
00228       }
00229     }
00230 
00231     return k;
00232   }
00233 
00234   /** \brief Polynomial extrapolation
00235 
00236       Compute the step of index <tt>i_step</tt> using polynomial
00237       extrapolation to evaulate functions by fitting a polynomial
00238       to estimates <tt>(x_i,y_i)</tt> and output the result to
00239       <tt>y_0</tt> and <tt>y_0_err</tt>. 
00240 
00241       The index <tt>i_step</tt> begins with zero.
00242         
00243   */
00244   int poly_extrap(gsl_matrix *dloc, const double x[],
00245                   const unsigned int i_step, const double x_i,
00246                   const vec_t &y_i, vec_t &y_0, vec_t &y_0_err, 
00247                   double work[]) {
00248     size_t j, k;
00249       
00250     o2scl::vector_copy(dim,y_i,y_0_err);
00251     o2scl::vector_copy(dim,y_i,y_0);
00252 
00253     if (i_step == 0) {
00254 
00255       for (j = 0; j < dim; j++) {
00256         gsl_matrix_set (dloc, 0, j, y_i[j]);
00257       }
00258 
00259     } else {
00260         
00261       o2scl::vector_copy(dim,y_i,work);
00262         
00263       for (k = 0; k < i_step; k++) {
00264         double deltaloc = 1.0 / (x[i_step - k - 1] - x_i);
00265         const double f1 = deltaloc * x_i;
00266         const double f2 = deltaloc * x[i_step - k - 1];
00267           
00268         for (j = 0; j < dim; j++) {
00269           const double q_kj = gsl_matrix_get (dloc, k, j);
00270           gsl_matrix_set (dloc, k, j, y_0_err[j]);
00271           deltaloc = work[j] - q_kj;
00272           y_0_err[j] = f1 * deltaloc;
00273           work[j] = f2 * deltaloc;
00274           y_0[j] += y_0_err[j];
00275         }
00276       }
00277         
00278       for (j = 0; j < dim; j++) {
00279         gsl_matrix_set(dloc, i_step, j, y_0_err[j]);
00280       }
00281     }
00282     return 0;
00283   }
00284 
00285   /** \brief Basic implicit Bulirsch-Stoer step
00286         
00287       Divide the step <tt>h_total</tt> into <tt>n_step</tt> smaller
00288       steps and do the Bader-Deuflhard semi-implicit iteration. This
00289       function starts at <tt>t0</tt> with function values
00290       <tt>y</tt>, derivatives <tt>yp_loc</tt>, and information from
00291       the Jacobian to compute the final value <tt>y_out</tt>.
00292   */
00293   int step_local(const double t0, const double h_total, 
00294                  const unsigned int n_step, const double y[], 
00295                  const double yp_loc[], const vec_t &dfdt_loc,
00296                  const mat_t &dfdy_loc, vec_t &y_out) {
00297       
00298     double *const w=weight;
00299       
00300     const double h = h_total / n_step;
00301     double t = t0 + h;
00302 
00303     double sum;
00304 
00305     /* This is the factor sigma referred to in equation 3.4 of the
00306        paper.  A relative change in y exceeding sigma indicates a
00307        runaway behavior. According to the authors suitable values
00308        for sigma are >>1.  I have chosen a value of 100*dim. BJG
00309     */
00310     const double max_sum = 100.0 * dim;
00311 
00312     int signum, status;
00313     size_t i, j;
00314     size_t n_inter;
00315 
00316     /* Calculate the matrix for the linear system. */
00317     for (i = 0; i < dim; i++) {
00318       for (j = 0; j < dim; j++) {
00319         gsl_matrix_set(a_mat,i,j,-h*dfdy_loc[i][j]);
00320       }
00321       gsl_matrix_set(a_mat,i,i,gsl_matrix_get(a_mat,i,i)+1.0);
00322     }
00323 
00324     /* LU decomposition for the linear system. */
00325       
00326     o2scl::permutation p_vec(dim);
00327     omatrix *om=(omatrix *)a_mat;
00328     o2scl_linalg::LU_decomp(dim,*om,p_vec,signum);
00329 
00330     /* Compute weighting factors */
00331     compute_weights(y,w,dim);
00332 
00333     /* Initial step. */
00334 
00335     for (i = 0; i < dim; i++) {
00336       y_temp[i]=h*(yp_loc[i]+h*dfdt_loc[i]);
00337     }
00338 
00339     o2scl_linalg::LU_solve(dim,*om,p_vec,y_temp,delta_temp);
00340       
00341     sum = 0.0;
00342     for (i = 0; i < dim; i++) {
00343       const double di = delta_temp[i];
00344       delta[i] = di;
00345       y_temp[i] = y[i] + di;
00346       sum += fabs(di) / w[i];
00347     }
00348     if (sum > max_sum) {
00349       return gsl_efailed;
00350     }
00351 
00352     /* Intermediate steps. */
00353 
00354     status=(*funcp)(t,dim,y_temp,y_out);
00355     if (status) {
00356       return status;
00357     }
00358       
00359     for (n_inter = 1; n_inter < n_step; n_inter++) {
00360 
00361       for (i = 0; i < dim; i++) {
00362         rhs_temp[i] = h*y_out[i]-delta[i];
00363       }
00364         
00365       o2scl_linalg::LU_solve(dim,*om,p_vec,rhs_temp,delta_temp);
00366 
00367       sum = 0.0;
00368       for (i = 0; i < dim; i++) {
00369         delta[i] += 2.0 * delta_temp[i];
00370         y_temp[i] += delta[i];
00371         sum += fabs(delta[i]) / w[i];
00372       }
00373       if (sum > max_sum) {
00374         return gsl_efailed ;
00375       }
00376         
00377       t += h;
00378         
00379       status=(*funcp)(t,dim,y_temp,y_out);
00380       if (status) {
00381         return status;
00382       }
00383     }
00384 
00385 
00386     /* Final step. */
00387 
00388     for (i = 0; i < dim; i++) {
00389       rhs_temp[i]=h*y_out[i]-delta[i];
00390     }
00391       
00392     o2scl_linalg::LU_solve(dim,*om,p_vec,rhs_temp,delta_temp);
00393       
00394     sum = 0.0;
00395     for (i = 0; i < dim; i++) {
00396       y_out[i] = y_temp[i]+delta_temp[i];
00397       sum += fabs(delta_temp[i])/w[i];
00398     }
00399 
00400     if (sum > max_sum) {
00401       return gsl_efailed;
00402     }
00403       
00404     return gsl_success;
00405   }
00406 
00407   /// Allocate memory for a system of size \c n
00408   int allocate(size_t n)  {
00409 
00410     if (dim>0) free();
00411 
00412     dim=n;
00413       
00414     d=gsl_matrix_alloc(sequence_max,n);
00415     a_mat=gsl_matrix_alloc(n,n);
00416       
00417     yp=(double *)malloc(n*sizeof(double));
00418 
00419     // AWS, 12/2/08 - This was added to ensure memory reallocation 
00420     // resets the stepper just like reset() does
00421     for(size_t i=0;i<n;i++) yp[i]=0.0;
00422 
00423     y_save=(double *)malloc(n*sizeof(double));
00424     yerr_save=(double *)malloc(n*sizeof(double));
00425     y_extrap_save=(double *)malloc(n*sizeof(double));
00426     extrap_work=(double *)malloc(n*sizeof(double));
00427     weight=(double *)malloc(n*sizeof(double));
00428 
00429     mo.allocate(dfdy,n,n);
00430     ao.allocate(dfdt,n);
00431     ao.allocate(y_extrap_sequence,n);
00432     ao.allocate(y_temp,n);
00433     ao.allocate(rhs_temp,n);
00434     ao.allocate(delta_temp,n);
00435 
00436     delta=(double *)malloc(n*sizeof(double));
00437 
00438     // This choice of epsilon is not necessarily optimal, it has
00439     // a "FIXME" comment in the original GSL code
00440     size_t k_choice_loc=deuf_kchoice(GSL_SQRT_DBL_EPSILON,n);
00441     k_choice=k_choice_loc;
00442     k_current=k_choice_loc;
00443     order=2*k_choice_loc;
00444       
00445     return 0;
00446   }
00447 
00448   /// Free allocated memory
00449   void free() {
00450     if (dim>0) {
00451       std::free(delta);
00452         
00453       mo.free(dfdy,dim);
00454       ao.free(rhs_temp);
00455       ao.free(dfdt);
00456       ao.free(y_temp);
00457       ao.free(y_extrap_sequence);
00458       ao.free(delta_temp);
00459         
00460       std::free(weight);
00461       std::free(extrap_work);
00462       std::free(y_extrap_save);
00463       std::free(y_save);
00464       std::free(yerr_save);
00465       std::free(yp);
00466         
00467       gsl_matrix_free(a_mat);
00468       gsl_matrix_free(d);
00469       dim=0;
00470     }
00471   }
00472 
00473 #endif
00474     
00475   public:
00476     
00477   gsl_bsimp() {
00478     dim=0;
00479   }
00480       
00481   virtual ~gsl_bsimp() {
00482     if (dim>0) free();
00483   }
00484 
00485   /// Reset stepper
00486   int reset() {
00487     for(size_t i=0;i<dim;i++) yp[i]=0.0;
00488     return gsl_success;
00489   }
00490     
00491   /** \brief Perform an integration step
00492         
00493       Given initial value of the n-dimensional function in \c y and
00494       the derivative in \c dydx (which must generally be computed
00495       beforehand) at the point \c x, take a step of size \c h giving
00496       the result in \c yout, the uncertainty in \c yerr, and the new
00497       derivative in \c dydx_out (at \c x+h) using function \c derivs
00498       to calculate derivatives.  Implementations which do not
00499       calculate \c yerr and/or \c dydx_out do not reference these
00500       variables so that a blank \c vec_t can be given. All of the
00501       implementations allow \c yout=y and \c dydx_out=dydx if
00502       necessary.
00503   */
00504   virtual int step(double x, double h, size_t n, vec_t &y, vec_t &dydx, 
00505                    vec_t &yout, vec_t &yerr, vec_t &dydx_out, 
00506                    func_t &derivs, jac_func_t &jac) {
00507 
00508     int ret;
00509 
00510     funcp=&derivs;
00511     jfuncp=&jac;
00512 
00513     if (n!=dim) allocate(n);
00514 
00515     /* Bader-Deuflhard extrapolation sequence */
00516     static const int bd_sequence[sequence_count] =
00517     { 2, 6, 10, 14, 22, 34, 50, 70 };
00518       
00519     double t_local=x;
00520 
00521     size_t i, k;
00522 
00523     if (h + t_local == t_local) {
00524       // This section is labeled with a "FIXME" comment in the
00525       // original GSL code. I'm not sure why, but an error is
00526       // sensible here.
00527       O2SCL_ERR_RET("Stepsize underflow in gsl_bsimp::step().",
00528                     gsl_eundrflw);
00529     }
00530 
00531     /* Save inputs */
00532     o2scl::vector_copy(dim,y,y_extrap_save);
00533     o2scl::vector_copy(dim,y,y_save);
00534     o2scl::vector_copy(dim,yerr,yerr_save);
00535 
00536     // Copy derivative
00537     o2scl::vector_copy(dim,dydx,yp);
00538       
00539     // Evaluate the Jacobian for the system. */
00540     ret=jac(t_local,dim,y,dfdy,dfdt);
00541     if (ret != gsl_success) {
00542       return ret;
00543     }
00544       
00545     /* Make a series of refined extrapolations, up to the specified
00546        maximum order, which was calculated based on the Deuflhard
00547        criterion in the deuf_kchoice() function (which is called by
00548        allocate() ).
00549     */
00550     for (k = 0; k <= k_current; k++) {
00551 
00552       const unsigned int N = bd_sequence[k];
00553       const double r = (h / N);
00554       const double x_k = r * r;
00555 
00556       // Each step computes a value of y_extrap_sequence,
00557       // using the number of sub-steps dictated by
00558       // the BD sequence
00559       int status=step_local(t_local,h,N,y_extrap_save,yp,
00560                             dfdt,dfdy,y_extrap_sequence);
00561 
00562       if (status == gsl_efailed) {
00563         /* If the local step fails, set the error to infinity in
00564            order to force a reduction in the step size 
00565         */
00566         for (i = 0; i < dim; i++) {
00567           yerr[i] = GSL_POSINF;
00568         }
00569         break;
00570       } else if (status != gsl_success) {
00571         return status;
00572       }
00573         
00574       // Use the information in y_extrap_sequence to compute 
00575       // the new value of y and yerr .
00576       ex_wk[k] = x_k;
00577       poly_extrap(d,ex_wk,k,x_k,y_extrap_sequence,y,yerr,extrap_work);
00578     }
00579 
00580     /* Evaluate dydt_out[]. */
00581 
00582     ret=derivs(t_local+h,dim,y,dydx_out);
00583       
00584     // If we failed, copy the old values back to y and yerr
00585     if (ret != gsl_success) {
00586       o2scl::vector_copy(dim,y_save,y);
00587       o2scl::vector_copy(dim,yerr_save,yerr);
00588       return ret;
00589     }
00590 
00591     return gsl_success;
00592   }
00593     
00594   };
00595   
00596 #ifndef DOXYGENP
00597 }
00598 #endif
00599 
00600 #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.