Object-oriented Scientific Computing Library: Version 0.910
jacobian.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_JACOBIAN_H
00024 #define O2SCL_JACOBIAN_H
00025 
00026 #include <string>
00027 #include <o2scl/omatrix_tlate.h>
00028 #include <o2scl/mm_funct.h>
00029 #include <o2scl/gsl_deriv.h>
00030 #include <o2scl/columnify.h>
00031 
00032 #ifndef DOXYGENP
00033 namespace o2scl {
00034 #endif
00035 
00036   /** \brief Base for a square Jacobian where J is 
00037       computed at x given y=f(x) [abstract base]
00038 
00039       Compute
00040       \f[
00041       J_{ij} = \frac{\partial f_i}{\partial x_j}
00042       \f]
00043 
00044       The \c vec_t objects in operator() could have been written to be
00045       \c const, but they are not \c const so that they can be used as
00046       temporary workspace. They are typically restored to their
00047       original values before operator() exits. The Jacobian is
00048       stored in the order \c J[i][j], i.e. the rows have index
00049       \c i and the columns with index \c j.
00050   */
00051   template<class vec_t=ovector_base, 
00052     class mat_t=omatrix_base> class jac_funct {
00053 
00054     public:  
00055 
00056     jac_funct() {}
00057 
00058     virtual ~jac_funct() {}
00059     
00060     /** \brief The operator()
00061      */
00062     virtual int operator()(size_t nv, vec_t &x, vec_t &y, mat_t &j)=0;
00063     
00064 #ifndef DOXYGENP
00065 
00066     private:
00067 
00068     jac_funct(const jac_funct &);
00069     jac_funct& operator=(const jac_funct&);
00070 
00071 #endif
00072 
00073   };
00074 
00075   /** \brief Function pointer to jacobian
00076    */
00077   template<class vec_t=ovector_base, 
00078     class mat_t=omatrix_base> class jac_funct_fptr : 
00079   public jac_funct<vec_t,mat_t> {
00080     
00081     public:
00082     
00083     /** \brief Specify the function pointer
00084      */
00085     jac_funct_fptr(int (*fp)(size_t nv, vec_t &x, vec_t &y, mat_t &j)) {
00086       fptr=fp;
00087     }
00088     
00089     virtual ~jac_funct_fptr();
00090     
00091     /** \brief The operator()
00092      */
00093     virtual int operator()(size_t nv, vec_t &x, vec_t &y, mat_t &j) {
00094       return fptr(nv,x,y,j);
00095     }
00096 
00097 #ifndef DOXYGENP
00098 
00099     protected:
00100   
00101     jac_funct_fptr() {};
00102     
00103     /// Function pointer
00104     int (*fptr)(size_t nv, vec_t &x, vec_t &y, mat_t &j);
00105     
00106     private:
00107 
00108     jac_funct_fptr(const jac_funct_fptr &);
00109     jac_funct_fptr& operator=(const jac_funct_fptr&);
00110 
00111 #endif
00112 
00113   };
00114 
00115   /** \brief Member function pointer to a Jacobian
00116    */
00117   template <class tclass, class vec_t=ovector_base,
00118     class mat_t=omatrix_base> class jac_funct_mfptr : 
00119   public jac_funct<vec_t,mat_t> {
00120 
00121     public:
00122   
00123     /** \brief Specify the member function pointer
00124      */
00125     jac_funct_mfptr(tclass *tp, int (tclass::*fp)
00126                     (size_t nv, vec_t &x, vec_t &y,
00127                      mat_t &j)) {
00128       tptr=tp;
00129       fptr=fp;
00130     }
00131   
00132     virtual ~jac_funct_mfptr() {};
00133   
00134     /** \brief The operator()
00135      */
00136     virtual int operator()(size_t nv, vec_t &x, vec_t &y, mat_t &j) {
00137       return (*tptr.*fptr)(nv,x,y,j);
00138     }
00139   
00140 #ifndef DOXYGEN_INTERNAL
00141 
00142     protected:
00143   
00144     /// Member function pointer
00145     int (tclass::*fptr)(size_t nv, vec_t &x, vec_t &y, mat_t &j);
00146 
00147     /// Class pointer
00148     tclass *tptr;
00149 
00150 #endif
00151 
00152 #ifndef DOXYGENP
00153   
00154     private:
00155 
00156     jac_funct_mfptr(const jac_funct_mfptr &);
00157     jac_funct_mfptr& operator=(const jac_funct_mfptr&);
00158 
00159 #endif
00160 
00161   };
00162 
00163   /** \brief Const member function pointer to a Jacobian
00164    */
00165   template <class tclass, class vec_t=ovector_base,
00166     class mat_t=omatrix_base> class jac_funct_cmfptr : 
00167   public jac_funct<vec_t,mat_t> {
00168 
00169     public:
00170   
00171     /** \brief Specify the member function pointer
00172      */
00173     jac_funct_cmfptr(tclass *tp, int (tclass::*fp)
00174                      (size_t nv, vec_t &x, vec_t &y,
00175                       mat_t &j) const) {
00176       tptr=tp;
00177       fptr=fp;
00178     }
00179   
00180     virtual ~jac_funct_cmfptr() {};
00181   
00182     /** \brief The operator()
00183      */
00184     virtual int operator()(size_t nv, vec_t &x, vec_t &y, mat_t &j) {
00185       return (*tptr.*fptr)(nv,x,y,j);
00186     }
00187   
00188 #ifndef DOXYGEN_INTERNAL
00189 
00190     protected:
00191   
00192     /// Member function pointer
00193     int (tclass::*fptr)(size_t nv, vec_t &x, vec_t &y, mat_t &j) const;
00194 
00195     /// Class pointer
00196     tclass *tptr;
00197 
00198 #endif
00199 
00200 #ifndef DOXYGENP
00201   
00202     private:
00203 
00204     jac_funct_cmfptr(const jac_funct_cmfptr &);
00205     jac_funct_cmfptr& operator=(const jac_funct_cmfptr&);
00206 
00207 #endif
00208 
00209   };
00210 
00211   /** \brief Base for providing a numerical jacobian [abstract base]
00212 
00213   This is provides a Jacobian which is numerically determined
00214   by differentiating a user-specified function (typically 
00215   of the form of \ref mm_funct). 
00216   */
00217   template<class func_t=mm_funct<>, class vec_t=ovector_base, 
00218     class mat_t=omatrix_base> class jacobian : 
00219   public jac_funct<vec_t,mat_t> {
00220     
00221     public:
00222     
00223     jacobian() {
00224     };
00225     
00226     virtual ~jacobian() {};
00227     
00228     /// Set the function to compute the Jacobian of
00229     virtual int set_function(func_t &f) {
00230       func=&f;
00231       return 0;
00232     }
00233     
00234     /** \brief Evaluate the Jacobian \c j at point \c y(x)
00235      */
00236     virtual int operator()(size_t nv, vec_t &x, vec_t &y, mat_t &j)=0;
00237     
00238 #ifndef DOXYGEN_INTERNAL
00239     
00240     protected:
00241     
00242     /// A pointer to the user-specified function
00243     func_t *func;
00244     
00245     private:
00246     
00247     jacobian(const jacobian &);
00248     jacobian& operator=(const jacobian&);
00249     
00250 #endif
00251     
00252   };
00253   
00254   /** \brief Simple automatic Jacobian
00255 
00256   This class computes a numerical Jacobian by finite differencing.
00257   The stepsize is chosen to be \f$ h_j = \mathrm{epsrel}~x_j \f$ or
00258   \f$ h_j = \mathrm{epsmin} \f$ if \f$ \mathrm{epsrel}~x_j <
00259   \mathrm{epsmin} \f$.
00260       
00261   This is nearly equivalent to the GSL method for computing
00262   Jacobians as in \c multiroots/fdjac.c. To obtain the GSL
00263   behavior, set \ref epsrel to \c GSL_SQRT_DBL_EPSILON and set
00264   \ref epsmin to zero. The \ref gsl_mroot_hybrids class sets \ref
00265   epsrel to \c GSL_SQRT_DBL_EPSILON in its constructor,
00266   but does not set \ref epsmin to zero.
00267 
00268   This class does not separately check the vector and matrix sizes
00269   to ensure they are commensurate. 
00270   */
00271   template<class func_t=mm_funct<>, class vec_t=ovector_base, 
00272     class mat_t=omatrix_base, class alloc_vec_t=ovector, 
00273     class alloc_t=ovector_alloc> class simple_jacobian :
00274   public jacobian<func_t,vec_t,mat_t> {
00275 
00276 #ifndef DOXYGEN_INTERNAL
00277 
00278     protected:
00279 
00280     /// For memory allocation
00281     alloc_t ao;
00282 
00283     /// Function values
00284     alloc_vec_t f;
00285 
00286     /// Function arguments
00287     alloc_vec_t xx;
00288 
00289     /// Size of allocated memory
00290     size_t mem_size;
00291 
00292 #endif
00293 
00294     public:
00295     
00296     simple_jacobian() {
00297       epsrel=1.0e-4;
00298       epsmin=1.0e-15;
00299       err_nonconv=true;
00300       mem_size=0;
00301     }
00302 
00303     virtual ~simple_jacobian() {
00304       ao.free(f);
00305     }
00306   
00307     /** \brief The relative stepsize for finite-differencing 
00308         (default \f$ 10^{-4} \f$ )
00309     */
00310     double epsrel;
00311     
00312     /// The minimum stepsize (default \f$ 10^{-15} \f$)
00313     double epsmin;
00314 
00315     /// If true, call the error handler if the routine does not "converge"
00316     bool err_nonconv;
00317 
00318     /** \brief The operator()
00319      */
00320     virtual int operator()(size_t nv, vec_t &x, vec_t &y, mat_t &jac) {
00321       
00322       size_t i,j;
00323       double h,temp;
00324       bool success=true;
00325 
00326       if (mem_size!=nv) {
00327         if (mem_size>0) ao.free(f);
00328         ao.allocate(f,nv);
00329         ao.allocate(xx,nv);
00330         mem_size=nv;
00331       }
00332       
00333       vector_copy(nv,x,xx);
00334 
00335       for (j=0;j<nv;j++) {
00336         
00337         h=epsrel*fabs(x[j]);
00338         if (fabs(h)<=epsmin) h=epsrel;
00339         
00340         xx[j]=x[j]+h;
00341         (*this->func)(nv,xx,f);
00342         xx[j]=x[j];
00343 
00344         // This is the equivalent of GSL's test of
00345         // gsl_vector_isnull(&col.vector)
00346 
00347         bool nonzero=false;
00348         for (i=0;i<nv;i++) {
00349           temp=(f[i]-y[i])/h;
00350           if (temp!=0.0) nonzero=true;
00351           jac[i][j]=temp;
00352         }
00353         if (nonzero==false) success=false;
00354 
00355       }
00356       
00357       if (success==false) {
00358         O2SCL_CONV2_RET("At least one row of the Jacobian is zero ",
00359                         "in simple_jacobian::operator().",gsl_esing,
00360                         this->err_nonconv);
00361       }
00362       return 0;
00363     }
00364   };
00365 
00366   /** \brief A direct calculation of the jacobian using a \ref deriv object
00367       
00368   Note that it is most often wasteful to use this Jacobian in a
00369   root-finding routine and using more approximate Jacobians is
00370   more efficient. This class is mostly useful for demonstration
00371   purposes.
00372   */
00373   template<class func_t=mm_funct<>, class vec_t=ovector_base, 
00374     class mat_t=omatrix_base> class exact_jacobian : 
00375   public jacobian<func_t,vec_t,mat_t> {
00376     
00377     public:
00378     
00379     exact_jacobian() {
00380       def_deriv.h=1.0e-4;
00381       dptr=&def_deriv;
00382     }
00383     
00384     /** \brief Parameter structure for passing information
00385         
00386     This class is primarily useful for specifying derivatives
00387     for using the jacobian::set_deriv() function.
00388 
00389     \comment
00390     This type needs to be publicly available so that the
00391     user can properly specify a base 1-dimensional derivative
00392     object. 
00393     \endcomment
00394     */
00395     typedef struct {
00396       /// The number of variables
00397       size_t nv;
00398       /// The current x value
00399       size_t xj;
00400       /// The current y value
00401       size_t yi;
00402       /// The x vector
00403       vec_t *x;
00404       /// The y vector
00405       vec_t *y;
00406     } ej_parms;
00407 
00408     /// The default derivative object
00409     gsl_deriv<funct> def_deriv;
00410   
00411     /// Set the derivative object
00412     int set_deriv(deriv<funct> &de) {
00413       dptr=&de;
00414       return 0;
00415     }
00416     
00417     /** \brief The operator()
00418      */
00419     virtual int operator()(size_t nv, vec_t &x, vec_t &y, mat_t &jac) {
00420 
00421       size_t i,j;
00422       double h,temp;
00423 
00424       ej_parms ejp;
00425       ejp.nv=nv;
00426       ejp.x=&x;
00427       ejp.y=&y;
00428     
00429       funct_mfptr_param<exact_jacobian,ej_parms> 
00430       dfnp(this,&exact_jacobian::dfn,ejp);
00431       
00432       for (j=0;j<nv;j++) {
00433         ejp.xj=j;
00434         for (i=0;i<nv;i++) {
00435           ejp.yi=i;
00436           double tmp=(*ejp.x)[j];
00437           jac[i][j]=dptr->calc(tmp,dfnp);
00438           (*ejp.x)[j]=tmp;
00439         }
00440       }
00441 
00442       return 0;
00443     }
00444 
00445 #ifndef DOXYGEN_INTERNAL
00446 
00447     protected:
00448 
00449     /// Pointer to the derivative object
00450     deriv<funct> *dptr;
00451     
00452     /// Function for the derivative object
00453     double dfn(double x, ej_parms &ejp) {
00454       (*ejp.x)[ejp.xj]=x;
00455       (*this->func)(ejp.nv,*ejp.x,*ejp.y);
00456       return (*ejp.y)[ejp.yi];
00457     }
00458 
00459 #endif
00460 
00461   };
00462   
00463 #ifndef DOXYGENP
00464 }
00465 #endif
00466 
00467 #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.