ccgsl 2.7.2
C++wrappersforGnuScientificLibrary
multifit_nlin.hpp
Go to the documentation of this file.
1/*
2 * $Id: multifit_nlin.hpp 309 2013-11-17 16:26:16Z jdl3 $
3 * Copyright (C) 2010, 2011, 2012, 2016, 2020 John D Lamb
4 *
5 * This program is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 2 of the License, or (at
8 * your option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful, but
11 * WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13 * General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program; if not, write to the Free Software
17 * Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
18 */
19
20#ifndef CCGSL_MULTIFIT_NLIN_HPP
21#define CCGSL_MULTIFIT_NLIN_HPP
22#include<gsl/gsl_multifit_nlin.h>
25
26#ifndef DOXYGEN_SKIP
27namespace gsl {
28 namespace multifit {
36 inline int gradient( matrix const& J, vector const& f, vector& g ){
37 return gsl_multifit_gradient( J.get(), f.get(), g.get() ); }
38
46 inline int covar( matrix const& J, double epsrel, matrix& covar ){
47 return gsl_multifit_covar( J.get(), epsrel, covar.get() ); }
48
54 class fsolver {
55 public:
59 typedef gsl_multifit_fsolver_type type;
63 fsolver(){
64 ccgsl_pointer = 0;
65 count = 0; // initially nullptr will do
66 }
67 // Refines random access container
68 // Refines assignable
75 explicit fsolver( type const* T, size_t const n, size_t const p ){
76 ccgsl_pointer = gsl_multifit_fsolver_alloc( T, n, p );
77 // just plausibly we could allocate fsolver but not count
78 try { count = new size_t; } catch( std::bad_alloc& e ){
79 // try to tidy up before rethrowing
80 gsl_multifit_fsolver_free( ccgsl_pointer );
81 throw e;
82 }
83 *count = 1; // initially there is just one reference to ccgsl_pointer
84 }
91 explicit fsolver( gsl_multifit_fsolver* v ){
92 ccgsl_pointer = v;
93 // just plausibly we could fail to allocate count: no further action needed.
94 count = new size_t;
95 *count = 1; // initially there is just one reference to ccgsl_pointer
96 }
97 // copy constructor
102 fsolver( fsolver const& v ){ ccgsl_pointer = v.ccgsl_pointer;
103 count = v.count; if( count != 0 ) ++*count; }
104 // assignment operator
109 fsolver& operator=( fsolver const& v ){
110 // first, possibly delete anything pointed to by this
111 if( count == 0 or --*count == 0 ){
112 if( ccgsl_pointer != 0 ) gsl_multifit_fsolver_free( ccgsl_pointer );
113 delete count;
114 } // Then copy
115 ccgsl_pointer = v.ccgsl_pointer; count = v.count; if( count != 0 ) ++*count; return *this;
116 }
117 // destructor
121 ~fsolver(){
122 if( count == 0 or --*count == 0 ){
123 // could have allocated null pointer
124 if( ccgsl_pointer != 0 ) gsl_multifit_fsolver_free( ccgsl_pointer );
125 delete count;
126 }
127 }
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;
136 }
142 fsolver& operator=( fsolver&& v ){
143 fsolver( std::move( v ) ).swap( *this );
144 return *this;
145 }
146#endif
147 // Refines equality comparable
148 // == operator
155 bool operator==( fsolver const& v ) const { return ccgsl_pointer == v.ccgsl_pointer; }
156 // != operator
163 bool operator!=( fsolver const& v ) const { return not operator==( v ); }
164 // Refines forward container
165 // Refines less than comparable
166 // operator<
175 bool operator<( fsolver const& v ) const { return ccgsl_pointer < v.ccgsl_pointer; }
176 // operator>
185 bool operator>( fsolver const& v ) const { return ccgsl_pointer > v.ccgsl_pointer; }
186 // operator<=
195 bool operator<=( fsolver const& v ) const { return ccgsl_pointer <= v.ccgsl_pointer; }
196 // operator>=
205 bool operator>=( fsolver const& v ) const { return ccgsl_pointer >= v.ccgsl_pointer; }
210 bool empty() const { return ccgsl_pointer == 0; }
211 // swap() --- should work even if sizes don't match
217 void swap( fsolver& v ){
218 std::swap( ccgsl_pointer, v.ccgsl_pointer );
219 std::swap( count, v.count );
220 }
221 private:
225 gsl_multifit_fsolver* ccgsl_pointer;
229 size_t* count;
230 public:
231 // shared reference functions
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__
254 explicit
255#endif
256 operator bool() const { return ccgsl_pointer != 0; }
257
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 );
270 return result;
271 }
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 );
284 return result;
285 }
286
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 );
297 return result;
298 }
299
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 );
314 return result;
315 }
316
322 inline static std::string name( fsolver const& s ){
323 return std::string( gsl_multifit_fsolver_name( s.get() ) ); }
324
332 inline static vector position( fsolver const& s ){
333 vector result( gsl_multifit_fsolver_position( s.get() ) );
334 return result;
335 }
336
337 // Member functions
338
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 );
350 return result;
351 }
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 );
363 return result;
364 }
365
370 int iterate(){
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 );
375 return result;
376 }
377
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 );
391 return result;
392 }
393
398 char const* name(){ return gsl_multifit_fsolver_name( get() ); }
399
406 vector const& position(){
407 p_v.wrap_gsl_vector_without_ownership( gsl_multifit_fsolver_position( get() ) );
408 return p_v;
409 }
410
417 vector const& get_x(){ return x_v; }
424 vector const& get_f(){ return f_v; }
431 vector const& get_dx(){ return dx_v; }
432 private:
436 vector p_v;
440 vector x_v;
444 vector f_v;
448 vector dx_v;
449 public:
450 // Solver types
451 // None yet implemented in GSL
452 };
453
459 class fdfsolver {
460 public:
464 typedef gsl_multifit_fdfsolver_type type;
468 fdfsolver(){
469 ccgsl_pointer = 0;
470 count = 0; // initially nullptr will do
471 }
472 // Refines random access container
473 // Refines assignable
480 explicit fdfsolver( type const* T, size_t const n, size_t const p ){
481 ccgsl_pointer = gsl_multifit_fdfsolver_alloc( T, n, p );
482 // just plausibly we could allocate fdfsolver but not count
483 try { count = new size_t; } catch( std::bad_alloc& e ){
484 // try to tidy up before rethrowing
485 gsl_multifit_fdfsolver_free( ccgsl_pointer );
486 throw e;
487 }
488 *count = 1; // initially there is just one reference to ccgsl_pointer
489 }
496 explicit fdfsolver( gsl_multifit_fdfsolver* v ){
497 ccgsl_pointer = v;
498 // just plausibly we could fail to allocate count: no further action needed.
499 count = new size_t;
500 *count = 1; // initially there is just one reference to ccgsl_pointer
501 }
502 // copy constructor
507 fdfsolver( fdfsolver const& v ){ ccgsl_pointer = v.ccgsl_pointer;
508 count = v.count; if( count != 0 ) ++*count; }
509 // assignment operator
514 fdfsolver& operator=( fdfsolver const& v ){
515 // first, possibly delete anything pointed to by this
516 if( count == 0 or --*count == 0 ){
517 if( ccgsl_pointer != 0 ) gsl_multifit_fdfsolver_free( ccgsl_pointer );
518 delete count;
519 } // Then copy
520 ccgsl_pointer = v.ccgsl_pointer; count = v.count; if( count != 0 ) ++*count; return *this;
521 }
522 // destructor
526 ~fdfsolver(){
527 if( count == 0 or --*count == 0 ){
528 // could have allocated null pointer
529 if( ccgsl_pointer != 0 ) gsl_multifit_fdfsolver_free( ccgsl_pointer );
530 delete count;
531 }
532 }
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;
541 }
547 fdfsolver& operator=( fdfsolver&& v ){
548 fdfsolver( std::move( v ) ).swap( *this );
549 return *this;
550 }
551#endif
552 // Refines equality comparable
553 // == operator
560 bool operator==( fdfsolver const& v ) const { return ccgsl_pointer == v.ccgsl_pointer; }
561 // != operator
568 bool operator!=( fdfsolver const& v ) const { return not operator==( v ); }
569 // Refines forward container
570 // Refines less than comparable
571 // operator<
580 bool operator<( fdfsolver const& v ) const { return ccgsl_pointer < v.ccgsl_pointer; }
581 // operator>
590 bool operator>( fdfsolver const& v ) const { return ccgsl_pointer > v.ccgsl_pointer; }
591 // operator<=
600 bool operator<=( fdfsolver const& v ) const { return ccgsl_pointer <= v.ccgsl_pointer; }
601 // operator>=
610 bool operator>=( fdfsolver const& v ) const { return ccgsl_pointer >= v.ccgsl_pointer; }
615 bool empty() const { return ccgsl_pointer == 0; }
616 // swap() --- should work even if sizes don't match
622 void swap( fdfsolver& v ){
623 std::swap( ccgsl_pointer, v.ccgsl_pointer );
624 std::swap( count, v.count );
625 }
626 private:
630 gsl_multifit_fdfsolver* ccgsl_pointer;
634 size_t* count;
635 public:
636 // shared reference functions
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__
659 explicit
660#endif
661 operator bool() const { return ccgsl_pointer != 0; }
662
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 );
675 return result;
676 }
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 );
689 return result;
690 }
691
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 );
702 return result;
703 }
704
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 );
723 return result;
724 }
725
731 inline static std::string name( fdfsolver const& s ){
732 return std::string( gsl_multifit_fdfsolver_name( s.get() ) ); }
733
741 inline static vector position( fdfsolver const& s ){
742 vector result( gsl_multifit_fdfsolver_position( s.get() ) );
743 return result;
744 }
754 inline static int dif_df( vector const& x, vector const& wts,
755 function_fdf* fdf, vector const& f,
756 matrix& J ){
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,
779 matrix& J ){
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() ); }
791
792 // Member functions
793
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 );
805 return result;
806 }
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 );
818 return result;
819 }
820
825 int iterate(){
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 );
830 return result;
831 }
832
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,
845 &info );
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 );
849 return result;
850 }
851
856 char const* name(){ return gsl_multifit_fdfsolver_name( get() ); }
857
864 vector const& position(){
865 p_v.wrap_gsl_vector_without_ownership( gsl_multifit_fdfsolver_position( get() ) );
866 return p_v;
867 }
868
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; }
897
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,
987 int& info ) const {
988 return gsl_multifit_fdfsolver_test( get(), xtol, gtol, ftol, &info ); }
989
990 private:
994 vector p_v;
998 vector x_v;
1002 vector f_v;
1006 vector dx_v;
1010 matrix J_m;
1011 public:
1012 // Solver types
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; }
1027 };
1028
1029 namespace test {
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 ); }
1040
1047 inline int gradient( vector const& g, double epsabs ){
1048 return gsl_multifit_test_gradient( g.get(), epsabs ); }
1049 }
1050
1056 class fdfridge {
1060 typedef gsl_multifit_fdfsolver_type type;
1061 public:
1065 fdfridge(){
1066 ccgsl_pointer = 0;
1067 count = 0; // initially nullptr will do
1068 }
1069 // Refines random access container
1070 // Refines assignable
1077 explicit fdfridge( type const* T, size_t const n, size_t const p ){
1078 ccgsl_pointer = gsl_multifit_fdfridge_alloc( T, n, p );
1079 // just plausibly we could allocate fdfridge but not count
1080 try { count = new size_t; } catch( std::bad_alloc& e ){
1081 // try to tidy up before rethrowing
1082 gsl_multifit_fdfridge_free( ccgsl_pointer );
1083 throw e;
1084 }
1085 *count = 1; // initially there is just one reference to ccgsl_pointer
1086 }
1093 explicit fdfridge( gsl_multifit_fdfridge* v ){
1094 ccgsl_pointer = v;
1095 // just plausibly we could fail to allocate count: no further action needed.
1096 count = new size_t;
1097 *count = 1; // initially there is just one reference to ccgsl_pointer
1098 }
1099 // copy constructor
1104 fdfridge( fdfridge const& v ){ ccgsl_pointer = v.ccgsl_pointer;
1105 count = v.count; if( count != 0 ) ++*count; }
1106 // assignment operator
1111 fdfridge& operator=( fdfridge const& v ){
1112 // first, possibly delete anything pointed to by this
1113 if( count == 0 or --*count == 0 ){
1114 if( ccgsl_pointer != 0 ) gsl_multifit_fdfridge_free( ccgsl_pointer );
1115 delete count;
1116 } // Then copy
1117 ccgsl_pointer = v.ccgsl_pointer; count = v.count; if( count != 0 ) ++*count; return *this;
1118 }
1119 // destructor
1123 ~fdfridge(){
1124 if( count == 0 or --*count == 0 ){
1125 // could have allocated null pointer
1126 if( ccgsl_pointer != 0 ) gsl_multifit_fdfridge_free( ccgsl_pointer );
1127 delete count;
1128 }
1129 }
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;
1138 }
1144 fdfridge& operator=( fdfridge&& v ){
1145 fdfridge( std::move( v ) ).swap( *this );
1146 return *this;
1147 }
1148#endif
1149 // Refines equality comparable
1150 // == operator
1157 bool operator==( fdfridge const& v ) const { return ccgsl_pointer == v.ccgsl_pointer; }
1158 // != operator
1165 bool operator!=( fdfridge const& v ) const { return not operator==( v ); }
1166 // Refines forward container
1167 // Refines less than comparable
1168 // operator<
1177 bool operator<( fdfridge const& v ) const { return ccgsl_pointer < v.ccgsl_pointer; }
1178 // operator>
1187 bool operator>( fdfridge const& v ) const { return ccgsl_pointer > v.ccgsl_pointer; }
1188 // operator<=
1197 bool operator<=( fdfridge const& v ) const { return ccgsl_pointer <= v.ccgsl_pointer; }
1198 // operator>=
1207 bool operator>=( fdfridge const& v ) const { return ccgsl_pointer >= v.ccgsl_pointer; }
1212 bool empty() const { return ccgsl_pointer == 0; }
1213 // swap() --- should work even if sizes don't match
1219 void swap( fdfridge& v ){
1220 std::swap( ccgsl_pointer, v.ccgsl_pointer );
1221 std::swap( count, v.count );
1222 }
1223 private:
1227 gsl_multifit_fdfridge* ccgsl_pointer;
1231 size_t* count;
1232 public:
1233 // shared reference functions
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__
1256 explicit
1257#endif
1258 operator bool() const { return ccgsl_pointer != 0; }
1259
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 ); }
1423 };
1424
1425 }
1426}
1427#endif
1428#endif
void wrap_gsl_vector_without_ownership(gsl_vector *v)
This function is intended mainly for internal use.
Definition: vector.hpp:272
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.
Definition: multilarge.hpp:40
gsl_multilarge_nlinear_fdf fdf
Typedef for shorthand.
int gradient(vector const &g, double epsabs)
C++ version of gsl_multimin_test_gradient().
Definition: multimin.hpp:402
int delta(vector const &dx, vector const &x, double epsabs, double epsrel)
C++ version of gsl_multiroot_test_delta().
Definition: multiroots.hpp:819
double get(quantile_workspace &w)
C++ version of gsl_rstat_quantile_get().
Definition: rstat.hpp:626
size_t n(workspace const &w)
C++ version of gsl_rstat_n().
Definition: rstat.hpp:299
gsl_sf_result result
Typedef for gsl_sf_result.
Definition: sf_result.hpp:30
The gsl package creates an interface to the GNU Scientific Library for C++.
Definition: blas.hpp:34