20#ifndef CCGSL_MULTIFIT_NLIN_HPP
21#define CCGSL_MULTIFIT_NLIN_HPP
22#include<gsl/gsl_multifit_nlin.h>
36 inline int gradient( matrix
const& J, vector
const& f, vector& g ){
37 return gsl_multifit_gradient( J.get(), f.get(), g.get() ); }
46 inline int covar( matrix
const& J,
double epsrel, matrix&
covar ){
47 return gsl_multifit_covar( J.get(), epsrel,
covar.get() ); }
59 typedef gsl_multifit_fsolver_type
type;
75 explicit fsolver(
type const* T,
size_t const n,
size_t const p ){
76 ccgsl_pointer = gsl_multifit_fsolver_alloc( T,
n, p );
78 try { count =
new size_t; }
catch( std::bad_alloc& e ){
80 gsl_multifit_fsolver_free( ccgsl_pointer );
91 explicit fsolver( gsl_multifit_fsolver* v ){
102 fsolver( fsolver
const& v ){ ccgsl_pointer = v.ccgsl_pointer;
103 count = v.count;
if( count != 0 ) ++*count; }
109 fsolver& operator=( fsolver
const& v ){
111 if( count == 0 or --*count == 0 ){
112 if( ccgsl_pointer != 0 ) gsl_multifit_fsolver_free( ccgsl_pointer );
115 ccgsl_pointer = v.ccgsl_pointer; count = v.count;
if( count != 0 ) ++*count;
return *
this;
122 if( count == 0 or --*count == 0 ){
124 if( ccgsl_pointer != 0 ) gsl_multifit_fsolver_free( ccgsl_pointer );
128#ifdef __GXX_EXPERIMENTAL_CXX0X__
133 fsolver( fsolver&& v ) : ccgsl_pointer( v.ccgsl_pointer ), count( nullptr ){
134 std::swap( count, v.count );
135 v.ccgsl_pointer =
nullptr;
142 fsolver& operator=( fsolver&& v ){
143 fsolver( std::move( v ) ).swap( *
this );
155 bool operator==( fsolver
const& v )
const {
return ccgsl_pointer == v.ccgsl_pointer; }
163 bool operator!=( fsolver
const& v )
const {
return not operator==( v ); }
175 bool operator<( fsolver
const& v )
const {
return ccgsl_pointer < v.ccgsl_pointer; }
185 bool operator>( fsolver
const& v )
const {
return ccgsl_pointer > v.ccgsl_pointer; }
195 bool operator<=( fsolver
const& v )
const {
return ccgsl_pointer <= v.ccgsl_pointer; }
205 bool operator>=( fsolver
const& v )
const {
return ccgsl_pointer >= v.ccgsl_pointer; }
210 bool empty()
const {
return ccgsl_pointer == 0; }
217 void swap( fsolver& v ){
218 std::swap( ccgsl_pointer, v.ccgsl_pointer );
219 std::swap( count, v.count );
225 gsl_multifit_fsolver* ccgsl_pointer;
236 gsl_multifit_fsolver*
get()
const {
return ccgsl_pointer; }
242 bool unique()
const {
return count != 0 and *count == 1; }
247 size_t use_count()
const {
return count == 0 ? 0 : *count; }
253#ifdef __GXX_EXPERIMENTAL_CXX0X__
256 operator bool()
const {
return ccgsl_pointer != 0; }
265 inline static int set( fsolver& s, function* f, vector
const& x ){
266 int const result = gsl_multifit_fsolver_set( s.get(), f, x.get() );
267 s.x_v.wrap_gsl_vector_without_ownership( s.get()->x );
268 s.f_v.wrap_gsl_vector_without_ownership( s.get()->f );
269 s.dx_v.wrap_gsl_vector_without_ownership( s.get()->dx );
279 inline static int set( fsolver& s, function& f, vector
const& x ){
280 int const result = gsl_multifit_fsolver_set( s.get(), &f, x.get() );
281 s.x_v.wrap_gsl_vector_without_ownership( s.get()->x );
282 s.f_v.wrap_gsl_vector_without_ownership( s.get()->f );
283 s.dx_v.wrap_gsl_vector_without_ownership( s.get()->dx );
292 inline static int iterate( fsolver& s ){
293 int const result = gsl_multifit_fsolver_iterate( s.get() );
294 s.x_v.wrap_gsl_vector_without_ownership( s.get()->x );
295 s.f_v.wrap_gsl_vector_without_ownership( s.get()->f );
296 s.dx_v.wrap_gsl_vector_without_ownership( s.get()->dx );
308 inline static int driver( fsolver& s,
size_t const maxiter,
309 double const epsabs,
double const epsrel ){
310 int const result = gsl_multifit_fsolver_driver( s.get(), maxiter, epsabs, epsrel );
311 s.x_v.wrap_gsl_vector_without_ownership( s.get()->x );
312 s.f_v.wrap_gsl_vector_without_ownership( s.get()->x );
313 s.dx_v.wrap_gsl_vector_without_ownership( s.get()->dx );
322 inline static std::string
name( fsolver
const& s ){
323 return std::string( gsl_multifit_fsolver_name( s.get() ) ); }
332 inline static vector
position( fsolver
const& s ){
333 vector
result( gsl_multifit_fsolver_position( s.get() ) );
345 int set( function* f, vector
const& x ){
346 int const result = gsl_multifit_fsolver_set(
get(), f, x.get() );
347 x_v.wrap_gsl_vector_without_ownership(
get()->x );
348 f_v.wrap_gsl_vector_without_ownership(
get()->f );
349 dx_v.wrap_gsl_vector_without_ownership(
get()->dx );
358 int set( function& f, vector
const& x ){
359 int const result = gsl_multifit_fsolver_set(
get(), &f, x.get() );
360 x_v.wrap_gsl_vector_without_ownership(
get()->x );
361 f_v.wrap_gsl_vector_without_ownership(
get()->f );
362 dx_v.wrap_gsl_vector_without_ownership(
get()->dx );
371 int const&
result = gsl_multifit_fsolver_iterate(
get() );
372 x_v.wrap_gsl_vector_without_ownership(
get()->x );
373 f_v.wrap_gsl_vector_without_ownership(
get()->f );
374 dx_v.wrap_gsl_vector_without_ownership(
get()->dx );
385 int driver(
const size_t maxiter,
386 const double epsabs,
const double epsrel ){
387 int const&
result = gsl_multifit_fsolver_driver(
get(), maxiter, epsabs, epsrel );
388 x_v.wrap_gsl_vector_without_ownership(
get()->x );
389 f_v.wrap_gsl_vector_without_ownership(
get()->f );
390 dx_v.wrap_gsl_vector_without_ownership(
get()->dx );
398 char const*
name(){
return gsl_multifit_fsolver_name(
get() ); }
417 vector
const& get_x(){
return x_v; }
424 vector
const& get_f(){
return f_v; }
431 vector
const& get_dx(){
return dx_v; }
464 typedef gsl_multifit_fdfsolver_type
type;
480 explicit fdfsolver(
type const* T,
size_t const n,
size_t const p ){
481 ccgsl_pointer = gsl_multifit_fdfsolver_alloc( T,
n, p );
483 try { count =
new size_t; }
catch( std::bad_alloc& e ){
485 gsl_multifit_fdfsolver_free( ccgsl_pointer );
496 explicit fdfsolver( gsl_multifit_fdfsolver* v ){
507 fdfsolver( fdfsolver
const& v ){ ccgsl_pointer = v.ccgsl_pointer;
508 count = v.count;
if( count != 0 ) ++*count; }
514 fdfsolver& operator=( fdfsolver
const& v ){
516 if( count == 0 or --*count == 0 ){
517 if( ccgsl_pointer != 0 ) gsl_multifit_fdfsolver_free( ccgsl_pointer );
520 ccgsl_pointer = v.ccgsl_pointer; count = v.count;
if( count != 0 ) ++*count;
return *
this;
527 if( count == 0 or --*count == 0 ){
529 if( ccgsl_pointer != 0 ) gsl_multifit_fdfsolver_free( ccgsl_pointer );
533#ifdef __GXX_EXPERIMENTAL_CXX0X__
538 fdfsolver( fdfsolver&& v ) : ccgsl_pointer( v.ccgsl_pointer ), count( nullptr ){
539 std::swap( count, v.count );
540 v.ccgsl_pointer =
nullptr;
547 fdfsolver& operator=( fdfsolver&& v ){
548 fdfsolver( std::move( v ) ).swap( *
this );
560 bool operator==( fdfsolver
const& v )
const {
return ccgsl_pointer == v.ccgsl_pointer; }
568 bool operator!=( fdfsolver
const& v )
const {
return not operator==( v ); }
580 bool operator<( fdfsolver
const& v )
const {
return ccgsl_pointer < v.ccgsl_pointer; }
590 bool operator>( fdfsolver
const& v )
const {
return ccgsl_pointer > v.ccgsl_pointer; }
600 bool operator<=( fdfsolver
const& v )
const {
return ccgsl_pointer <= v.ccgsl_pointer; }
610 bool operator>=( fdfsolver
const& v )
const {
return ccgsl_pointer >= v.ccgsl_pointer; }
615 bool empty()
const {
return ccgsl_pointer == 0; }
622 void swap( fdfsolver& v ){
623 std::swap( ccgsl_pointer, v.ccgsl_pointer );
624 std::swap( count, v.count );
630 gsl_multifit_fdfsolver* ccgsl_pointer;
641 gsl_multifit_fdfsolver*
get()
const {
return ccgsl_pointer; }
647 bool unique()
const {
return count != 0 and *count == 1; }
652 size_t use_count()
const {
return count == 0 ? 0 : *count; }
658#ifdef __GXX_EXPERIMENTAL_CXX0X__
661 operator bool()
const {
return ccgsl_pointer != 0; }
670 inline static int set( fdfsolver& s, function_fdf* f, vector
const& x ){
671 int const result = gsl_multifit_fdfsolver_set( s.get(), f, x.get() );
672 s.x_v.wrap_gsl_vector_without_ownership( s.get()->x );
673 s.f_v.wrap_gsl_vector_without_ownership( s.get()->f );
674 s.dx_v.wrap_gsl_vector_without_ownership( s.get()->dx );
684 inline static int set( fdfsolver& s, function_fdf& f, vector
const& x ){
685 int const result = gsl_multifit_fdfsolver_set( s.get(), &f, x.get() );
686 s.x_v.wrap_gsl_vector_without_ownership( s.get()->x );
687 s.f_v.wrap_gsl_vector_without_ownership( s.get()->f );
688 s.dx_v.wrap_gsl_vector_without_ownership( s.get()->dx );
697 inline static int iterate( fdfsolver& s ){
698 int const result = gsl_multifit_fdfsolver_iterate( s.get() );
699 s.x_v.wrap_gsl_vector_without_ownership( s.get()->x );
700 s.f_v.wrap_gsl_vector_without_ownership( s.get()->f );
701 s.dx_v.wrap_gsl_vector_without_ownership( s.get()->dx );
715 inline static int driver( fdfsolver& s,
size_t const maxiter,
716 double const xtol,
double const gtol,
717 double const ftol,
int& info ){
718 int const result = gsl_multifit_fdfsolver_driver( s.get(), maxiter,
719 xtol, gtol, ftol, &info );
720 s.x_v.wrap_gsl_vector_without_ownership( s.get()->x );
721 s.f_v.wrap_gsl_vector_without_ownership( s.get()->f );
722 s.dx_v.wrap_gsl_vector_without_ownership( s.get()->dx );
731 inline static std::string
name( fdfsolver
const& s ){
732 return std::string( gsl_multifit_fdfsolver_name( s.get() ) ); }
741 inline static vector
position( fdfsolver
const& s ){
742 vector
result( gsl_multifit_fdfsolver_position( s.get() ) );
754 inline static int dif_df( vector
const& x, vector
const& wts,
755 function_fdf*
fdf, vector
const& f,
757 return gsl_multifit_fdfsolver_dif_df( x.get(), wts.get(),
fdf, f.get(), J.get() ); }
766 inline static int dif_fdf( vector
const& x, function_fdf *
fdf, vector& f, matrix& J ){
767 return gsl_multifit_fdfsolver_dif_fdf( x.get(),
fdf, f.get(), J.get() ); }
777 inline static int dif_df( vector
const& x, vector
const& wts,
778 function_fdf&
fdf, vector
const& f,
780 return gsl_multifit_fdfsolver_dif_df( x.get(), wts.get(), &
fdf, f.get(), J.get() ); }
789 inline static int dif_fdf( vector
const& x, function_fdf &
fdf, vector& f, matrix& J ){
790 return gsl_multifit_fdfsolver_dif_fdf( x.get(), &
fdf, f.get(), J.get() ); }
800 int set( function_fdf* f, vector
const& x ){
801 int const result = gsl_multifit_fdfsolver_set(
get(), f, x.get() );
802 x_v.wrap_gsl_vector_without_ownership(
get()->x );
803 f_v.wrap_gsl_vector_without_ownership(
get()->f );
804 dx_v.wrap_gsl_vector_without_ownership(
get()->dx );
813 int set( function_fdf& f, vector
const& x ){
814 int const result = gsl_multifit_fdfsolver_set(
get(), &f, x.get() );
815 x_v.wrap_gsl_vector_without_ownership(
get()->x );
816 f_v.wrap_gsl_vector_without_ownership(
get()->f );
817 dx_v.wrap_gsl_vector_without_ownership(
get()->dx );
826 int const result = gsl_multifit_fdfsolver_iterate(
get() );
827 x_v.wrap_gsl_vector_without_ownership(
get()->x );
828 f_v.wrap_gsl_vector_without_ownership(
get()->f );
829 dx_v.wrap_gsl_vector_without_ownership(
get()->dx );
842 int driver(
size_t const maxiter,
double const xtol,
double const gtol,
843 double const ftol,
int& info ){
844 int const result = gsl_multifit_fdfsolver_driver(
get(), maxiter, xtol, gtol, xtol,
846 x_v.wrap_gsl_vector_without_ownership(
get()->x );
847 f_v.wrap_gsl_vector_without_ownership(
get()->f );
848 dx_v.wrap_gsl_vector_without_ownership(
get()->dx );
856 char const*
name(){
return gsl_multifit_fdfsolver_name(
get() ); }
875 vector
const& get_x(){
return x_v; }
882 vector
const& get_f(){
return f_v; }
889 vector
const& get_dx(){
return dx_v; }
896 matrix
const& get_J(){
return J_m; }
905 inline int wset( function_fdf* f, vector
const& x, vector
const& wts ){
906 return gsl_multifit_fdfsolver_wset(
get(), f, x.get(), wts.get() ); }
914 inline int wset( function_fdf& f, vector
const& x, vector
const& wts ){
915 return gsl_multifit_fdfsolver_wset(
get(), &f, x.get(), wts.get() ); }
921 inline int fdfsolver_jac( matrix& J ){
922 return gsl_multifit_fdfsolver_jac(
get(), J.get() ); }
927 inline vector fdfsolver_residual()
const {
928 return vector( gsl_multifit_fdfsolver_residual(
get() ) ); }
933 inline size_t fdfsolver_niter()
const {
return gsl_multifit_fdfsolver_niter(
get() ); }
942 inline static int eval_wf( multifit::function_fdf&
fdf, vector
const& x,
943 vector
const& wts, vector& y ){
944 return gsl_multifit_eval_wf( &
fdf, x.get(), wts.get(), y.get() ); }
953 inline static int eval_wf( gsl_multifit_function_fdf*
fdf, vector
const& x,
954 vector
const& wts, vector& y ){
955 return gsl_multifit_eval_wf(
fdf, x.get(), wts.get(), y.get() ); }
964 inline static int eval_wdf( multifit::function_fdf&
fdf, vector
const& x,
965 vector
const& wts, matrix& dy ){
966 return gsl_multifit_eval_wdf( &
fdf, x.get(), wts.get(), dy.get() ); }
975 inline static int eval_wdf( gsl_multifit_function_fdf*
fdf, vector
const& x,
976 vector
const& wts, matrix& dy ){
977 return gsl_multifit_eval_wdf(
fdf, x.get(), wts.get(), dy.get() ); }
986 inline int fdfsolver_test(
double const xtol,
double const gtol,
double const ftol,
988 return gsl_multifit_fdfsolver_test(
get(), xtol, gtol, ftol, &info ); }
1017 inline static type const* lmder(){
return gsl_multifit_fdfsolver_lmder; }
1022 inline static type const* lmsder(){
return gsl_multifit_fdfsolver_lmsder; }
1026 inline static type const* lmniel(){
return gsl_multifit_fdfsolver_lmniel; }
1038 inline int delta( vector
const& dx, vector
const& x,
double epsabs,
double epsrel ){
1039 return gsl_multifit_test_delta( dx.get(), x.get(), epsabs, epsrel ); }
1047 inline int gradient( vector
const& g,
double epsabs ){
1048 return gsl_multifit_test_gradient( g.get(), epsabs ); }
1060 typedef gsl_multifit_fdfsolver_type
type;
1077 explicit fdfridge(
type const* T,
size_t const n,
size_t const p ){
1078 ccgsl_pointer = gsl_multifit_fdfridge_alloc( T,
n, p );
1080 try { count =
new size_t; }
catch( std::bad_alloc& e ){
1082 gsl_multifit_fdfridge_free( ccgsl_pointer );
1093 explicit fdfridge( gsl_multifit_fdfridge* v ){
1104 fdfridge( fdfridge
const& v ){ ccgsl_pointer = v.ccgsl_pointer;
1105 count = v.count;
if( count != 0 ) ++*count; }
1111 fdfridge& operator=( fdfridge
const& v ){
1113 if( count == 0 or --*count == 0 ){
1114 if( ccgsl_pointer != 0 ) gsl_multifit_fdfridge_free( ccgsl_pointer );
1117 ccgsl_pointer = v.ccgsl_pointer; count = v.count;
if( count != 0 ) ++*count;
return *
this;
1124 if( count == 0 or --*count == 0 ){
1126 if( ccgsl_pointer != 0 ) gsl_multifit_fdfridge_free( ccgsl_pointer );
1130#ifdef __GXX_EXPERIMENTAL_CXX0X__
1135 fdfridge( fdfridge&& v ) : ccgsl_pointer( v.ccgsl_pointer ), count( nullptr ){
1136 std::swap( count, v.count );
1137 v.ccgsl_pointer =
nullptr;
1144 fdfridge& operator=( fdfridge&& v ){
1145 fdfridge( std::move( v ) ).swap( *
this );
1157 bool operator==( fdfridge
const& v )
const {
return ccgsl_pointer == v.ccgsl_pointer; }
1165 bool operator!=( fdfridge
const& v )
const {
return not operator==( v ); }
1177 bool operator<( fdfridge
const& v )
const {
return ccgsl_pointer < v.ccgsl_pointer; }
1187 bool operator>( fdfridge
const& v )
const {
return ccgsl_pointer > v.ccgsl_pointer; }
1197 bool operator<=( fdfridge
const& v )
const {
return ccgsl_pointer <= v.ccgsl_pointer; }
1207 bool operator>=( fdfridge
const& v )
const {
return ccgsl_pointer >= v.ccgsl_pointer; }
1212 bool empty()
const {
return ccgsl_pointer == 0; }
1219 void swap( fdfridge& v ){
1220 std::swap( ccgsl_pointer, v.ccgsl_pointer );
1221 std::swap( count, v.count );
1227 gsl_multifit_fdfridge* ccgsl_pointer;
1238 gsl_multifit_fdfridge*
get()
const {
return ccgsl_pointer; }
1244 bool unique()
const {
return count != 0 and *count == 1; }
1249 size_t use_count()
const {
return count == 0 ? 0 : *count; }
1255#ifdef __GXX_EXPERIMENTAL_CXX0X__
1258 operator bool()
const {
return ccgsl_pointer != 0; }
1264 inline char const* fdfridge_name()
const {
1265 return gsl_multifit_fdfridge_name(
get() ); }
1270 inline vector fdfridge_position()
const {
1271 return vector( gsl_multifit_fdfridge_position(
get() ) ); }
1276 inline vector fdfridge_residual()
const {
1277 return vector( gsl_multifit_fdfridge_residual(
get() ) ); }
1282 inline size_t fdfridge_niter()
const {
1283 return gsl_multifit_fdfridge_niter(
get() ); }
1291 inline int fdfridge_set( multifit::function_fdf* f, vector
const& x,
double const lambda ){
1292 return gsl_multifit_fdfridge_set(
get(), f, x.get(), lambda ); }
1301 inline int fdfridge_wset( multifit::function_fdf* f, vector
const& x,
double const lambda,
1302 vector
const& wts ){
1303 return gsl_multifit_fdfridge_wset(
get(), f, x.get(), lambda, wts.get() ); }
1311 inline int fdfridge_set2( multifit::function_fdf* f, vector
const& x,
1312 vector
const& lambda ){
1313 return gsl_multifit_fdfridge_set2(
get(), f, x.get(), lambda.get() ); }
1322 inline int fdfridge_wset2( multifit::function_fdf* f, vector
const& x,
1323 vector
const& lambda, vector
const& wts ){
1324 return gsl_multifit_fdfridge_wset2(
get(), f, x.get(), lambda.get(), wts.get() ); }
1332 inline int fdfridge_set3( multifit::function_fdf* f, vector
const& x, matrix
const& L ){
1333 return gsl_multifit_fdfridge_set3(
get(), f, x.get(), L.get() ); }
1342 inline int fdfridge_wset3( multifit::function_fdf* f, vector
const& x, matrix
const& L,
1343 vector
const& wts ){
1344 return gsl_multifit_fdfridge_wset3(
get(), f, x.get(), L.get(), wts.get() ); }
1352 inline int fdfridge_set( multifit::function_fdf& f, vector
const& x,
double const lambda ){
1353 return gsl_multifit_fdfridge_set(
get(), &f, x.get(), lambda ); }
1362 inline int fdfridge_wset( multifit::function_fdf& f, vector
const& x,
double const lambda,
1363 vector
const& wts ){
1364 return gsl_multifit_fdfridge_wset(
get(), &f, x.get(), lambda, wts.get() ); }
1372 inline int fdfridge_set2( multifit::function_fdf& f, vector
const& x,
1373 vector
const& lambda ){
1374 return gsl_multifit_fdfridge_set2(
get(), &f, x.get(), lambda.get() ); }
1383 inline int fdfridge_wset2( multifit::function_fdf& f, vector
const& x,
1384 vector
const& lambda, vector
const& wts ){
1385 return gsl_multifit_fdfridge_wset2(
get(), &f, x.get(), lambda.get(), wts.get() ); }
1393 inline int fdfridge_set3( multifit::function_fdf& f, vector
const& x, matrix
const& L ){
1394 return gsl_multifit_fdfridge_set3(
get(), &f, x.get(), L.get() ); }
1403 inline int fdfridge_wset3( multifit::function_fdf& f, vector
const& x, matrix
const& L,
1404 vector
const& wts ){
1405 return gsl_multifit_fdfridge_wset3(
get(), &f, x.get(), L.get(), wts.get() ); }
1410 inline int fdfridge_iterate(){
return gsl_multifit_fdfridge_iterate(
get() ); }
1420 inline int fdfridge_driver(
size_t const maxiter,
double const xtol,
double const gtol,
1421 double const ftol,
int& info ){
1422 return gsl_multifit_fdfridge_driver(
get(), maxiter, xtol, gtol, ftol, &info ); }
void wrap_gsl_vector_without_ownership(gsl_vector *v)
This function is intended mainly for internal use.
int test(double const xtol, double const gtol, double const ftol, int *info, workspace const &w)
C++ version of gsl_multifit_nlinear_test().
char const * name(workspace const &w)
C++ version of gsl_multifit_nlinear_name().
int iterate(workspace &w)
C++ version of gsl_multifit_nlinear_iterate().
vector position(workspace const &w)
C++ version of gsl_multifit_nlinear_position().
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().
gsl_multilarge_linear_type type
Typedef for shorthand.
gsl_multilarge_nlinear_fdf fdf
Typedef for shorthand.
int gradient(vector const &g, double epsabs)
C++ version of gsl_multimin_test_gradient().
int delta(vector const &dx, vector const &x, double epsabs, double epsrel)
C++ version of gsl_multiroot_test_delta().
double get(quantile_workspace &w)
C++ version of gsl_rstat_quantile_get().
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++.