20#ifndef CCGSL_MULTIFIT_NLINEAR_HPP
21#define CCGSL_MULTIFIT_NLINEAR_HPP
24#include<gsl/gsl_multifit_nlinear.h>
60 explicit workspace( gsl_multifit_nlinear_type
const* T,
61 gsl_multifit_nlinear_parameters
const* params,
62 size_t const n,
size_t const p ){
65 try {
count =
new size_t; }
catch( std::bad_alloc& e ){
79 explicit workspace( gsl_multifit_nlinear_workspace* v ){
118#ifdef __GXX_EXPERIMENTAL_CXX0X__
125 std::swap(
count, v.count );
126 v.ccgsl_pointer =
nullptr;
134 workspace( std::move( v ) ).swap( *
this );
149 ( gsl_multifit_nlinear_workspace* w ){
264#ifdef __GXX_EXPERIMENTAL_CXX0X__
293 class CallbackImplement {
298 typedef void (*function_t)(
size_t const,
void*,gsl_multifit_nlinear_workspace
const*);
312 virtual function_t getFunction() = 0;
317 template<
typename T>
class F :
public base_F {
329 static void function(
size_t const i,
void* params,
330 gsl_multifit_nlinear_workspace
const* w ){
331 T* t =
reinterpret_cast<T*
>( params );
333 W.wrap_gsl_multifit_nlinear_workspace_without_ownership
334 (
const_cast<gsl_multifit_nlinear_workspace*
>( w ) );
335 return t->function( i, W );
341 function_t getFunction(){
return &F<T>::function; }
347 template<
typename T> CallbackImplement( T& t ){
359 CallbackImplement( CallbackImplement
const& v ) : f( v.f ), count( v.count ){
360 if( count != 0 ) ++*count;
367 CallbackImplement& operator=( CallbackImplement
const& v ){
369 if( count == 0 or --*count == 0 ){
376 if( count != 0 ) ++*count;
379#ifdef __GXX_EXPERIMENTAL_CXX0X__
384 CallbackImplement( CallbackImplement&& v ) : f( v.f ), count( nullptr ){
385 std::swap( count, v.count );
393 CallbackImplement& operator=( CallbackImplement&& v ){
395 std::swap( count, v.count );
402 ~CallbackImplement(){
404 if( count == 0 or --*count == 0 ){
413 function_t getFunction(){
return f->getFunction(); }
429 inline gsl_multifit_nlinear_parameters
442 return gsl_multifit_nlinear_init( x.
get(), &
fdf, w.
get() ); }
454 return gsl_multifit_nlinear_winit( x.
get(), wts.
get(), &
fdf, w.
get() ); }
488 void* callback_params =
reinterpret_cast<void*
>( &callback );
489 CallbackImplement callbackImplement( callback );
490 void (*
function)(
size_t const,
void*,gsl_multifit_nlinear_workspace
const*)
491 = callbackImplement.getFunction();
492 return gsl_multifit_nlinear_driver( maxiter, xtol, gtol, ftol,
function,
493 callback_params, &info, w.
get() ); }
498 gsl_multifit_nlinear_workspace
const* w);
517 void* callback_params,
519 return gsl_multifit_nlinear_driver( maxiter, xtol, gtol, ftol,
521 callback_params, info, w.
get() ); }
530 result.wrap_gsl_matrix_without_ownership( gsl_multifit_nlinear_jac( w.
get() ) );
547 result.wrap_gsl_vector_without_ownership( gsl_multifit_nlinear_position( w.
get() ) );
557 result.wrap_gsl_vector_without_ownership( gsl_multifit_nlinear_residual( w.
get() ) );
574 return gsl_multifit_nlinear_rcond( &
rcond, w.
get() ); }
583 return gsl_multifit_nlinear_rcond(
rcond, w.
get() ); }
591 return gsl_multifit_nlinear_trs_name( w.
get() ); }
603 return gsl_multifit_nlinear_eval_f( &
fdf, x.
get(), swts.
get(), y.
get() ); }
618 gsl_multifit_nlinear_fdtype
const fdtype,
621 return gsl_multifit_nlinear_eval_df( x.
get(), f.
get(), swts.
get(), h,
641 return gsl_multifit_nlinear_eval_fvv( h, x.
get(), v.
get(), f.
get(),
643 yvv.
get(), work.
get() ); }
655 return gsl_multifit_nlinear_covar( J.
get(), epsrel,
covar.get() ); }
666 test(
double const xtol,
double const gtol,
double const ftol,
int* info,
668 return gsl_multifit_nlinear_test( xtol, gtol, ftol, info, w.
get() ); }
679 test(
double const xtol,
double const gtol,
double const ftol,
int& info,
681 return gsl_multifit_nlinear_test( xtol, gtol, ftol, &info, w.
get() ); }
695 df(
double const h, gsl_multifit_nlinear_fdtype
const fdtype,
719 return gsl_multifit_nlinear_fdfvv( h, x.
get(), v.
get(), f.
get(), J.
get(),
This class handles matrix objects as shared handles.
gsl_matrix * get()
Get the gsl_matrix.
Class that extends gsl_multifit_function so that it can be constructed from arbitrary function object...
You can create callbacks as a subclass of this class.
virtual void function(size_t const iter, workspace const &w)=0
This is the callback function.
Class that extends gsl_multifit_nlinear_fdf so that it can be constructed from arbitrary function obj...
Workspace for nonlinear linear least squares fitting.
workspace(workspace &&v)
Move constructor.
~workspace()
The destructor only deletes the pointers if count reaches zero.
bool operator>(workspace const &v) const
A container needs to define an ordering for sorting.
gsl_multifit_nlinear_workspace * ccgsl_pointer
The shared pointer.
workspace(gsl_multifit_nlinear_workspace *v)
Could construct from a gsl_multifit_workspace.
size_t use_count() const
Find how many workspace objects share this pointer.
bool operator!=(workspace const &v) const
Two workspace are different if their elements are not identical.
bool operator<=(workspace const &v) const
A container needs to define an ordering for sorting.
workspace()
The default constructor is only really useful for assigning to.
void wrap_gsl_multifit_nlinear_workspace_without_ownership(gsl_multifit_nlinear_workspace *w)
This function is intended mainly for internal use.
workspace(gsl_multifit_nlinear_type const *T, gsl_multifit_nlinear_parameters const *params, size_t const n, size_t const p)
The default constructor creates a new workspace with n elements.
workspace & operator=(workspace const &v)
The assignment operator.
workspace(workspace const &v)
The copy constructor.
bool unique() const
Find if this is the only object sharing the gsl_multifit_nlinear_workspace.
bool empty() const
Find if the workspace is empty.
size_t * count
The shared reference count.
bool operator>=(workspace const &v) const
A container needs to define an ordering for sorting.
gsl_multifit_nlinear_workspace * get() const
Get the gsl_multifit_nlinear_workspace.
workspace & operator=(workspace &&v)
Move operator.
bool operator==(workspace const &v) const
Two workspace are identically equal if their elements are identical.
bool operator<(workspace const &v) const
A container needs to define an ordering for sorting.
void swap(workspace &v)
Swap two workspace objects.
This class handles vector objects as shared handles.
gsl_vector * get()
Get the gsl_vector.
int test(double const xtol, double const gtol, double const ftol, int *info, workspace const &w)
C++ version of gsl_multifit_nlinear_test().
int eval_f(gsl::multifit::nlinear::function_fdf &fdf, vector const &x, vector const &swts, vector &y)
C++ version of gsl_multifit_nlinear_eval_f().
char const * name(workspace const &w)
C++ version of gsl_multifit_nlinear_name().
int eval_fvv(double const h, vector const &x, vector const &v, vector const &f, matrix const &J, vector const &swts, gsl::multifit::nlinear::function_fdf &fdf, vector &yvv, vector &work)
C++ version of gsl_multifit_nlinear_eval_fvv().
gsl_multifit_nlinear_parameters default_parameters()
C++ version of gsl_multifit_nlinear_default_parameters().
int iterate(workspace &w)
C++ version of gsl_multifit_nlinear_iterate().
char const * trs_name(workspace const &w)
C++ version of gsl_multifit_nlinear_trs_name().
double avratio(workspace const &w)
C++ version of gsl_multifit_nlinear_avratio().
vector position(workspace const &w)
C++ version of gsl_multifit_nlinear_position().
int init(vector const &x, gsl::multifit::nlinear::function_fdf &fdf, workspace &w)
C++ version of gsl_multifit_nlinear_init().
int covar(matrix const &J, double const epsrel, matrix &covar)
C++ version of gsl_multifit_nlinear_covar().
int driver(size_t const maxiter, double const xtol, double const gtol, double const ftol, CallbackBase &callback, int &info, workspace &w)
C++ version of gsl_multifit_nlinear_driver().
int winit(vector const &x, vector const &wts, gsl::multifit::nlinear::function_fdf &fdf, workspace &w)
C++ version of gsl_multifit_nlinear_winit().
void(* driver_callback)(const size_t iter, void *params, gsl_multifit_nlinear_workspace const *w)
Function pointer for use by driver.
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().
int rcond(double &rcond, workspace const &w)
C++ version of gsl_multifit_nlinear_rcond().
int fdfvv(double const h, vector const &x, vector const &v, vector const &f, matrix const &J, vector const &swts, gsl::multifit::nlinear::function_fdf &fdf, vector &fvv, vector &work)
C++ version of gsl_multifit_nlinear_fdfvv().
matrix jac(workspace const &w)
C++ version of gsl_multifit_nlinear_jac().
vector residual(workspace const &w)
C++ version of gsl_multifit_nlinear_residual().
int eval_df(vector const &x, vector const &f, vector const &swts, double const h, gsl_multifit_nlinear_fdtype const fdtype, gsl::multifit::nlinear::function_fdf &fdf, matrix &df, vector &work)
C++ version of gsl_multifit_nlinear_eval_df().
size_t niter(workspace const &w)
C++ version of gsl_multifit_nlinear_niter().
gsl_multilarge_nlinear_fdf fdf
Typedef for shorthand.
gsl_multilarge_nlinear_fdtype fdtype
Typedef for shorthand.
size_t n(workspace const &w)
C++ version of gsl_rstat_n().
double F(double phi, double k, mode_t mode)
C++ version of gsl_sf_ellint_F().
gsl_sf_result result
Typedef for gsl_sf_result.
The gsl package creates an interface to the GNU Scientific Library for C++.