20#ifndef CCGSL_MULTIMIN_FUNCTION_FDF_HPP
21#define CCGSL_MULTIMIN_FUNCTION_FDF_HPP
27#if __cplusplus <= 199711L
115 typedef double (*function_t)(
const vector&,
void*);
119 typedef void (*function_df_t)(
const vector&,
void*,
vector&);
123 typedef void (*function_fdf_t)(
const vector&,
void*,
double&,
vector&);
129 struct FDF :
public base_F {
135 FDF( T& t ) : t( t ){
136 xv.wrap_gsl_vector_without_ownership( 0 );
137 dfv.wrap_gsl_vector_without_ownership( 0 );
146 static double fn( gsl_vector
const* x,
void* params ){
147 FDF<T>* fdft =
reinterpret_cast<FDF<T>*
>( params );
148 fdft->xv.ccgsl_pointer =
const_cast<gsl_vector*
>( x );
149 return fdft->t.f( fdft->xv ); }
156 static void dfn( gsl_vector
const* x,
void* params, gsl_vector*
df ){
157 FDF<T>* fdft =
reinterpret_cast<FDF<T>*
>( params );
158 fdft->xv.ccgsl_pointer =
const_cast<gsl_vector*
>( x );
159 fdft->dfv.ccgsl_pointer =
df;
160 return fdft->t.df( fdft->xv, fdft->dfv ); }
168 static void fdfn( gsl_vector
const* x,
void* params,
double* f, gsl_vector*
df ){
169 FDF<T>* fdft =
reinterpret_cast<FDF<T>*
>( params );
170 fdft->xv.ccgsl_pointer =
const_cast<gsl_vector*
>( x );
171 fdft->dfv.ccgsl_pointer =
df;
172 return fdft->t.fdf( fdft->xv, *f, fdft->dfv ); }
234 function_fdf_constructor<T>( *
this, t );
275#ifdef __GXX_EXPERIMENTAL_CXX0X__
283 std::swap(
df, v.df );
284 std::swap(
fdf, v.fdf );
285 std::swap( params, v.params );
286 std::swap(
count, v.count );
287 v.function =
nullptr;
298 std::swap(
df, v.df );
299 std::swap(
fdf, v.fdf );
300 std::swap( params, v.params );
301 std::swap(
count, v.count );
327 fdf.function =
new function_fdf::FDF<T>( t );
329 fdf.f = &function_fdf::FDF<T>::fn;
330 fdf.df = &function_fdf::FDF<T>::dfn;
331 fdf.fdf = &function_fdf::FDF<T>::fdfn;
332 fdf.params =
fdf.function;
334 fdf.count =
new size_t;
342 fdf.function = v.function;
fdf.count = v.count;
fdf.f = v.f;
fdf.n = v.n;
fdf.df = v.df;
343 fdf.fdf = v.fdf;
fdf.params = v.params;
if(
fdf.count != 0 ) ++*
fdf.count;
350 fdf.f = v.f;
fdf.n = v.n;
fdf.df = v.df;
fdf.fdf = v.fdf;
fdf.params = v.params; }
361#ifdef __GXX_EXPERIMENTAL_CXX0X__
362 return std::move( fn );
Class that extends gsl_function_fdf so that it can be constructed from arbitrary function objects.
Class that extends gsl_multimin_function_fdf so that it can be constructed from arbitrary function ob...
~function_fdf()
The destructor unshares any shared resource.
function_fdf(T &t)
Construct from an object that implements gsl::multimin::function_fdf::Concept.
function_fdf(gsl_multimin_function_fdf &v)
Could construct from a gsl_multimin_function_fdf.
function_fdf(function_fdf const &v)
The copy constructor.
size_t * count
The shared reference count: used for copying this.
function_fdf()
The default constructor is only really useful for assigning to.
function_fdf(function_fdf &&v)
Move constructor.
function_fdf & operator=(function_fdf &&v)
Move operator.
base_F * function
This gives something for params to point to.
friend void function_fdf_constructor(function_fdf &fdf, T &t)
function_fdf & operator=(function_fdf const &v)
The assignment operator.
Class that extends gsl_multimin_function so that it can be constructed from arbitrary function object...
function()
The default constructor is only really useful for assigning to.
This class handles vector objects as shared handles.
int df(double const h, gsl_multifit_nlinear_fdtype const fdtype, vector const &x, vector const &wts, gsl::multifit::nlinear::function_fdf &fdf, vector const &f, matrix &J, vector &work)
C++ version of gsl_multifit_nlinear_df().
gsl_multilarge_nlinear_fdf fdf
Typedef for shorthand.
int gradient(vector const &g, double epsabs)
C++ version of gsl_multimin_test_gradient().
void function_fdf_constructor(function_fdf &fdf, T &t)
function_fdf make_function_fdf(T &t)
Make a gsl::multimin::function_fdf from a function object that implements gsl::multimin::function_fdf...
size_t n(workspace const &w)
C++ version of gsl_rstat_n().
gsl_sf_result result
Typedef for gsl_sf_result.
The gsl package creates an interface to the GNU Scientific Library for C++.
This is an empty abstract base class.
virtual void fdf(const vector &x, double *result, vector &gradient)=0
The function and derivative.
virtual double f(const vector &x)=0
The function.
size_t size() const
The dimension of the function argument.
virtual void df(const vector &x, vector &gradient)=0
The derivative.