![]() |
Object-oriented Scientific Computing Library: Version 0.910
|
Solution of one equation in one variable is accomplished by children of the class root.
For one-dimensional solving, use cern_root or gsl_root_brent if you have the root bracketed, or gsl_root_stef if you have the derivative available. If you have neither a bracket or a derivative, you can use cern_mroot_root.
The root base class provides the structure for three different solving methods:
x
x1
and x2
. The values of the function at x1
and x2
should have different signs.x
and the derivative of the function df
.If not all of these three functions are overloaded, then the source code in the root base class is designed to try to automatically provide the solution using the remaining functions. Most of the one-dimensional solving routines, in their original form, are written in the second or third form above. For example, gsl_root_brent is originally a bracketing routine of the form root::solve_bkt(), but calls to either root::solve() or root::solve_de() will attempt to automatically bracket the function given the initial guess that is provided. Of course, it is frequently most efficient to use the solver in the way it was intended.
Solution of more than one equation is accomplished by descendants of the class mroot . The higher-level interface is provided by the function mroot::msolve() .
For multi-dimensional solving, you can use either cern_mroot or gsl_mroot_hybrids. While cern_mroot cannot utilize user-supplied derivatives, gsl_mroot_hybrids can use user-supplied derivative information (as in the GSL hybridsj method) using the function gsl_mroot_hybrids::msolve_de() .
This demonstrates several ways of using the multi-dimensional solvers to solve the equations
/* Example: ex_mroot.cpp ------------------------------------------------------------------- Several ways to use an O2scl solver to solve a simple function */ #include <cmath> #include <o2scl/test_mgr.h> #include <o2scl/mm_funct.h> #include <o2scl/gsl_mroot_hybrids.h> #include <o2scl/cern_mroot.h> using namespace std; using namespace o2scl; int gfn(size_t nv, const ovector_base &x, ovector_base &y) { y[0]=sin(x[1]-0.2); y[1]=sin(x[0]-0.25); return 0; } typedef double arr_t[2]; typedef double mat_t[2][2]; class cl { public: // Store the number of function and derivative evaluations int nf, nd; int mfn(size_t nv, const ovector_base &x, ovector_base &y) { y[0]=sin(x[1]-0.2); y[1]=sin(x[0]-0.25); nf++; return 0; } int operator()(size_t nv, const ovector_base &x, ovector_base &y) { y[0]=sin(x[1]-0.2); y[1]=sin(x[0]-0.25); nf++; return 0; } int mfnd(size_t nv, ovector_base &x, ovector_base &y, omatrix_base &j) { j[0][0]=0.0; j[0][1]=cos(x[1]-0.2); j[1][0]=cos(x[0]-0.25); j[1][1]=0.0; nd++; return 0; } int mfna(size_t nv, const arr_t &x, arr_t &y) { y[0]=sin(x[1]-0.2); y[1]=sin(x[0]-0.25); return 0; } int mfnad(size_t nv, arr_t &x, arr_t &y, mat_t &j) { j[0][0]=0.0; j[0][1]=cos(x[1]-0.2); j[1][0]=cos(x[0]-0.25); j[1][1]=0.0; return 0; } }; int main(void) { cl acl; ovector x(2); double xa[2]; int i; size_t tmp; int r1, r2, r3; bool done; test_mgr t; t.set_output_level(1); /* Using a member function with \ref ovector objects */ mm_funct_mfptr<cl> f1(&acl,&cl::mfn); gsl_mroot_hybrids<> cr1; x[0]=0.5; x[1]=0.5; acl.nf=0; int ret1=cr1.msolve(2,x,f1); cout << "GSL solver (numerical Jacobian): " << endl; cout << "Return value: " << ret1 << endl; cout << "Number of iterations: " << cr1.last_ntrial << endl; cout << "Number of function evaluations: " << acl.nf << endl; cout << endl; t.test_rel(x[0],0.25,1.0e-6,"1a"); t.test_rel(x[1],0.2,1.0e-6,"1b"); /* Using the CERNLIB solver */ cern_mroot<> cr2; x[0]=0.5; x[1]=0.5; acl.nf=0; int ret2=cr2.msolve(2,x,f1); cout << "CERNLIB solver (numerical Jacobian): " << endl; cout << "Return value: " << ret2 << endl; cout << "INFO parameter: " << cr2.get_info() << endl; cout << "Number of function evaluations: " << acl.nf << endl; cout << endl; t.test_rel(x[0],0.25,1.0e-6,"2a"); t.test_rel(x[1],0.2,1.0e-6,"2b"); /* Using a member function with \ref ovector objects, but using the GSL-like interface with set() and iterate(). */ gsl_mroot_hybrids<> cr3; x[0]=0.5; x[1]=0.5; cr3.allocate(2); cr3.set(2,x,f1); done=false; do { r3=cr3.iterate(); double resid=fabs(cr3.f[0])+fabs(cr3.f[1]); if (resid<cr3.tol_rel || r3>0) done=true; } while (done==false); t.test_rel(cr3.x[0],0.25,1.0e-6,"3a"); t.test_rel(cr3.x[1],0.2,1.0e-6,"3b"); cr3.free(); /* Now instead of using the automatic Jacobian, using a user-specified Jacobian. */ jac_funct_mfptr<cl,ovector_base,omatrix_base> j4(&acl,&cl::mfnd); x[0]=0.5; x[1]=0.5; acl.nf=0; acl.nd=0; int ret4=cr1.msolve_de(2,x,f1,j4); cout << "GSL solver (analytic Jacobian): " << endl; cout << "Return value: " << ret4 << endl; cout << "Number of iterations: " << cr1.last_ntrial << endl; cout << "Number of function evaluations: " << acl.nf << endl; cout << "Number of Jacobian evaluations: " << acl.nd << endl; cout << endl; t.test_rel(x[0],0.25,1.0e-6,"4a"); t.test_rel(x[1],0.2,1.0e-6,"4b"); /* Using a user-specified Jacobian and the GSL-like interface */ gsl_mroot_hybrids<> cr5; x[0]=0.5; x[1]=0.5; cr5.allocate(2); cr5.set_de(2,x,f1,j4); done=false; do { r3=cr5.iterate(); double resid=fabs(cr5.f[0])+fabs(cr5.f[1]); if (resid<cr5.tol_rel || r3>0) done=true; } while (done==false); t.test_rel(cr5.x[0],0.25,1.0e-6,"5a"); t.test_rel(cr5.x[1],0.2,1.0e-6,"5b"); cr5.free(); /* Using C-style arrays instead of ovector objects */ mm_funct_mfptr<cl,arr_t> f6(&acl,&cl::mfna); gsl_mroot_hybrids<mm_funct_mfptr<cl,arr_t>,arr_t, arr_t,array_alloc<arr_t> > cr6; xa[0]=0.5; xa[1]=0.5; cr6.msolve(2,xa,f6); t.test_rel(xa[0],0.25,1.0e-6,"6a"); t.test_rel(xa[1],0.2,1.0e-6,"6b"); /* Using the CERNLIB solver with C-style arrays instead of ovector objects */ cern_mroot<mm_funct_mfptr<cl,arr_t>,arr_t, arr_t,array_alloc<arr_t> > cr7; xa[0]=0.5; xa[1]=0.5; cr7.msolve(2,xa,f6); t.test_rel(xa[0],0.25,1.0e-6,"7a"); t.test_rel(xa[1],0.2,1.0e-6,"7b"); /* Using C-style arrays with a user-specified Jacobian */ jac_funct_mfptr<cl,arr_t,mat_t> j8(&acl,&cl::mfnad); gsl_mroot_hybrids<mm_funct_mfptr<cl,arr_t>,arr_t, arr_t,array_alloc<arr_t>,mat_t,mat_t, array_2d_alloc<mat_t>,jac_funct<arr_t,mat_t> > cr8; xa[0]=0.5; xa[1]=0.5; cr8.msolve_de(2,xa,f6,j8); t.test_rel(xa[0],0.25,1.0e-6,"8a"); t.test_rel(xa[1],0.2,1.0e-6,"8b"); /* Using a class with an operator(). Note that there can be only one operator() function in each class. */ gsl_mroot_hybrids<cl> cr9; x[0]=0.5; x[1]=0.5; cr9.msolve(2,x,acl); t.test_rel(x[0],0.25,1.0e-6,"9a"); t.test_rel(x[1],0.2,1.0e-6,"9b"); /* Using a function pointer to a global function. */ typedef int (*gfnt)(size_t, const ovector_base &, ovector_base &); gsl_mroot_hybrids<gfnt> cr10; gfnt f10=&gfn; x[0]=0.5; x[1]=0.5; cr10.msolve(2,x,f10); t.test_rel(x[0],0.25,1.0e-6,"10a"); t.test_rel(x[1],0.2,1.0e-6,"10b"); t.report(); return 0; } // End of example
Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).