00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 2008, 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_MROOT_H 00024 #define O2SCL_MROOT_H 00025 00026 #include <o2scl/ovector_tlate.h> 00027 #include <o2scl/mm_funct.h> 00028 #include <o2scl/jacobian.h> 00029 00030 #ifndef DOXYGENP 00031 namespace o2scl { 00032 #endif 00033 00034 /** 00035 \brief Multidimensional root-finding [abstract base] 00036 00037 <b>The template parameters:</b> 00038 The template parameter \c func_t specifies the functions to 00039 solve and should be a class containing a definition 00040 \code 00041 func_t::operator()(size_t nv, const vec_t &x, vec_t &y, param_t &pa); 00042 \endcode 00043 where \c y is the value of the functions at \c x with parameter \c pa 00044 and \c x and \c y are a array-like classes defining \c operator[] 00045 of size \c nv. If the Jacobian matrix is to be specified by the 00046 user, then the parameter \c jfunc_t specifies the jacobian 00047 and should contain the definition 00048 \code 00049 jfunc_t::operator(size_t nv, vec_t &x, vec_t &y, 00050 mat_t &j, param_t &pa); 00051 \endcode 00052 where \c x is the independent variables, \c y is the array of 00053 function values, and \c j is the Jacobian matrix. This template 00054 parameter can be ignored if only the function msolve() will be used. 00055 00056 \warning Many of the routines assume that the scale of the 00057 functions and their variables is of order unity. The solution 00058 routines may lose precision if this is not the case. 00059 00060 There is an example for the usage of the multidimensional solver 00061 classes given in <tt>examples/ex_mroot.cpp</tt>, see \ref 00062 ex_mroot_sect . 00063 */ 00064 template<class param_t, class func_t, class vec_t=ovector_view, 00065 class jfunc_t=jac_funct<param_t,vec_t,omatrix_view> > class mroot { 00066 00067 public: 00068 00069 mroot() { 00070 ntrial=100; 00071 tolf=1.0e-8; 00072 verbose=0; 00073 tolx=1.0e-12; 00074 last_ntrial=0; 00075 } 00076 00077 virtual ~mroot() {} 00078 00079 /// The maximum value of the functions for success (default 1.0e-8) 00080 double tolf; 00081 00082 /// The minimum allowable stepsize (default 1.0e-12) 00083 double tolx; 00084 00085 /// Output control (default 0) 00086 int verbose; 00087 00088 /// Maximum number of iterations (default 100) 00089 int ntrial; 00090 00091 /// The number of iterations for in the most recent minimization 00092 int last_ntrial; 00093 00094 /// Return the type, \c "mroot". 00095 virtual const char *type() { return "mroot"; } 00096 00097 /// Solve \c func using \c x as an initial guess, returning \c x. 00098 virtual int msolve(size_t n, vec_t &x, param_t &pa, func_t &func)=0; 00099 00100 /** \brief Solve \c func with derivatives \c dfunc using \c x as 00101 an initial guess, returning \c x. 00102 00103 By default, this function just calls msolve() and ignores the 00104 last argument. 00105 */ 00106 virtual int msolve_de(size_t n, vec_t &x, param_t &pa, func_t &func, 00107 jfunc_t &dfunc) { 00108 return msolve(n,x,pa,func); 00109 } 00110 00111 /** 00112 \brief Print out iteration information. 00113 00114 Depending on the value of the variable verbose, this prints out 00115 the iteration information. If verbose=0, then no information is 00116 printed, while if verbose>1, then after each iteration, the 00117 present values of x and y are output to std::cout along with the 00118 iteration number. If verbose>=2 then each iteration waits for a 00119 character. 00120 00121 This is implemented as a template class using a new vector type 00122 because sometimes the internal vector class is distinct from 00123 the user-specified vector class (e.g. in \ref gsl_mroot_hybrids). 00124 */ 00125 template<class vec2_t, class vec3_t> 00126 int print_iter(size_t n, const vec2_t &x, const vec3_t &y, 00127 int iter, double value=0.0, double limit=0.0, 00128 std::string comment="") 00129 { 00130 if (verbose<=0) return o2scl::gsl_success; 00131 00132 char ch; 00133 00134 std::cout << comment << " Iteration: " << iter << std::endl; 00135 for(size_t i=0;i<n;i++) { 00136 std::cout.width(3); 00137 std::cout << i; 00138 if (x[i]>=0.0) { 00139 std::cout << " " << x[i]; 00140 } else { 00141 std::cout << " " << x[i]; 00142 } 00143 if (y[i]>=0.0) { 00144 std::cout << " " << y[i] << std::endl; 00145 } else { 00146 std::cout << " " << y[i] << std::endl; 00147 } 00148 } 00149 std::cout << "Val: " << value << " Lim: " << limit << std::endl; 00150 if (verbose>1) { 00151 std::cout << "Press a key and type enter to continue. "; 00152 std::cin >> ch; 00153 } 00154 00155 return o2scl::gsl_success; 00156 00157 } 00158 00159 }; 00160 00161 #ifndef DOXYGENP 00162 00163 template<> int io_tlate<mroot<void *,mm_funct<void *> > >::input 00164 (cinput *co, in_file_format *ins, mroot<void *, mm_funct<void *> > *ro); 00165 template<> int io_tlate<mroot<void *,mm_funct<void *> > >::output 00166 (coutput *co, out_file_format *outs, mroot<void *, mm_funct<void *> > *ro); 00167 template<> const char *io_tlate<mroot<void *,mm_funct<void *> > >::type(); 00168 00169 #endif 00170 00171 typedef io_tlate<mroot<void *,mm_funct<void *> > > mroot_io_type; 00172 00173 #ifndef DOXYGENP 00174 } 00175 #endif 00176 00177 #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