BubbleProfiler
0.3.0
by Peter Athron, Csaba Balazs, Michael Bardsley, Andrew Fowlie, Dylan Harries & Graham White
|
#include <gsl_root_finder.hpp>
Public Types | |
enum | Algorithm { Algorithm::GSL_hybrids, Algorithm::GSL_hybrid, Algorithm::GSL_dnewton } |
using | Function_argument = Eigen::Matrix< double, N, 1 > |
using | Function_value = Eigen::Matrix< double, N, 1 > |
Public Member Functions | |
virtual | ~GSL_root_finder ()=default |
virtual const Function_argument & | get_solution () const override |
void | set_algorithm (Algorithm a) |
void | set_relative_error (double e) |
void | set_absolute_error (double e) |
void | set_use_stepsize_convergence_check (bool flag) |
void | set_max_iterations (std::size_t i) |
![]() | |
virtual | ~Root_finder ()=default |
Root_finder_status | find_root (F &&function, const Eigen::Matrix< double, N, 1 > &guess) |
Protected Types | |
using | Function = std::function< Function_value(const Function_argument &)> |
Protected Member Functions | |
virtual Root_finder_status | solve (const Function &, const Function_argument &) override |
![]() | |
virtual Root_finder_status | solve (const std::function< Eigen::Matrix< double, N, 1 >(const Eigen::Matrix< double, N, 1 > &)> &, const Eigen::Matrix< double, N, 1 > &)=0 |
Private Member Functions | |
const gsl_multiroot_fsolver_type * | get_gsl_solver_type () const |
Static Private Member Functions | |
static int | gsl_function (const gsl_vector *, void *, gsl_vector *) |
Private Attributes | |
std::size_t | n_dims |
Algorithm | solver_algorithm {Algorithm::GSL_hybrids} |
double | rel_err {1.e-6} |
double | abs_err {1.e-6} |
bool | use_stepsize_convergence {false} |
std::size_t | max_iterations {100} |
Function | function {nullptr} |
Function_argument | solution {} |
Definition at line 33 of file gsl_root_finder.hpp.
|
protected |
Definition at line 57 of file gsl_root_finder.hpp.
using BubbleProfiler::GSL_root_finder< N >::Function_argument = Eigen::Matrix<double,N,1> |
Definition at line 35 of file gsl_root_finder.hpp.
using BubbleProfiler::GSL_root_finder< N >::Function_value = Eigen::Matrix<double,N,1> |
Definition at line 36 of file gsl_root_finder.hpp.
|
strong |
Enumerator | |
---|---|
GSL_hybrids | |
GSL_hybrid | |
GSL_dnewton |
Definition at line 38 of file gsl_root_finder.hpp.
|
virtualdefault |
|
private |
Definition at line 187 of file gsl_root_finder.hpp.
References BubbleProfiler::GSL_root_finder< N >::GSL_dnewton, BubbleProfiler::GSL_root_finder< N >::GSL_hybrid, BubbleProfiler::GSL_root_finder< N >::GSL_hybrids, and BubbleProfiler::GSL_root_finder< N >::solver_algorithm.
Referenced by BubbleProfiler::GSL_root_finder< N >::solve().
|
inlineoverridevirtual |
Implements BubbleProfiler::Root_finder< Eigen::Matrix< double, N, 1 > >.
Definition at line 44 of file gsl_root_finder.hpp.
References BubbleProfiler::GSL_root_finder< N >::solution.
|
staticprivate |
Definition at line 76 of file gsl_root_finder.hpp.
Referenced by BubbleProfiler::GSL_root_finder< N >::solve().
|
inline |
Definition at line 50 of file gsl_root_finder.hpp.
References BubbleProfiler::GSL_root_finder< N >::abs_err.
|
inline |
Definition at line 48 of file gsl_root_finder.hpp.
References BubbleProfiler::GSL_root_finder< N >::solver_algorithm.
|
inline |
Definition at line 54 of file gsl_root_finder.hpp.
References BubbleProfiler::GSL_root_finder< N >::max_iterations.
|
inline |
Definition at line 49 of file gsl_root_finder.hpp.
References BubbleProfiler::GSL_root_finder< N >::rel_err.
|
inline |
Definition at line 51 of file gsl_root_finder.hpp.
References BubbleProfiler::GSL_root_finder< N >::use_stepsize_convergence.
|
overrideprotectedvirtual |
Definition at line 112 of file gsl_root_finder.hpp.
References BubbleProfiler::GSL_root_finder< N >::abs_err, BubbleProfiler::FAIL, BubbleProfiler::GSL_root_finder< N >::get_gsl_solver_type(), BubbleProfiler::GSL_root_finder< N >::gsl_function(), BubbleProfiler::make_raii_guard(), BubbleProfiler::GSL_root_finder< N >::max_iterations, BubbleProfiler::GSL_root_finder< N >::n_dims, BubbleProfiler::GSL_root_finder< N >::rel_err, BubbleProfiler::GSL_root_finder< N >::solution, BubbleProfiler::SUCCESS, and BubbleProfiler::GSL_root_finder< N >::use_stepsize_convergence.
|
private |
Definition at line 65 of file gsl_root_finder.hpp.
Referenced by BubbleProfiler::GSL_root_finder< N >::set_absolute_error(), and BubbleProfiler::GSL_root_finder< N >::solve().
|
private |
Definition at line 68 of file gsl_root_finder.hpp.
|
private |
Definition at line 67 of file gsl_root_finder.hpp.
Referenced by BubbleProfiler::GSL_root_finder< N >::set_max_iterations(), and BubbleProfiler::GSL_root_finder< N >::solve().
|
private |
Definition at line 62 of file gsl_root_finder.hpp.
Referenced by BubbleProfiler::GSL_root_finder< N >::solve().
|
private |
Definition at line 64 of file gsl_root_finder.hpp.
Referenced by BubbleProfiler::GSL_root_finder< N >::set_relative_error(), and BubbleProfiler::GSL_root_finder< N >::solve().
|
private |
Definition at line 69 of file gsl_root_finder.hpp.
Referenced by BubbleProfiler::GSL_root_finder< N >::get_solution(), and BubbleProfiler::GSL_root_finder< N >::solve().
|
private |
Definition at line 63 of file gsl_root_finder.hpp.
Referenced by BubbleProfiler::GSL_root_finder< N >::get_gsl_solver_type(), and BubbleProfiler::GSL_root_finder< N >::set_algorithm().
|
private |
Definition at line 66 of file gsl_root_finder.hpp.
Referenced by BubbleProfiler::GSL_root_finder< N >::set_use_stepsize_convergence_check(), and BubbleProfiler::GSL_root_finder< N >::solve().