ccgsl 2.7.2
C++wrappersforGnuScientificLibrary
multifit_nlinear_function_fdf.hpp
Go to the documentation of this file.
1/*
2 * $Id
3 * Copyright (C) 2010, 2011, 2012, 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_NLINEAR_FUNCTION_FDF_HPP
21#define CCGSL_MULTIFIT_NLINEAR_FUNCTION_FDF_HPP
22
23// For std::swap
24#if __cplusplus <= 199711L
25#include<algorithm>
26#else
27#include<utility>
28#endif
29#include<gsl/gsl_multifit_nlinear.h>
30#include"vector.hpp"
31#include"matrix.hpp"
32
33namespace gsl {
34 namespace multifit {
35 namespace nlinear {
36 class function_fdf;
37 template<typename T> void function_constructor( function_fdf& f, T& t );
56 class function_fdf : public gsl_multifit_nlinear_fdf {
57 public:
65 struct concept_f {
73 virtual int f( gsl::vector const& x, gsl::vector& f ) = 0;
77 virtual size_t xSize() const = 0;
81 virtual size_t fSize() const = 0;
82 };
90 struct concept_df {
100 virtual int df( gsl::vector const& x, gsl::matrix& J ) = 0;
101 };
109 struct concept_fvv {
119 virtual int fvv( gsl::vector const& x, gsl::vector const& v, gsl::vector& fvv ) = 0;
120 };
121#ifndef DOXYGEN_SKIP
122 struct address_of {
123 typedef int (*f_t)(gsl_vector const*,void*,gsl_vector*);
124 typedef int (*df_t)(gsl_vector const*,void*,gsl_matrix*);
125 typedef int (*fvv_t)(gsl_vector const*,gsl_vector const*,void*,gsl_vector*);
126 // Maybe not needed: f must not be nullptr
127 template<typename U, int (U::*)(gsl::vector const&, gsl::vector&)>
128 struct fSfinae {};
129 template<typename U, int (U::*)(gsl::vector const&, gsl::vector&) const>
130 struct fSfinaeC {};
131 template<typename U> static f_t f(fSfinae<U,&U::f>*){return &F<U>::fn;}
132 template<typename U> static f_t f(fSfinaeC<U,&U::f>*){return &F<U>::fn;}
133 template<typename U> static f_t f(...){return 0;}
134 // address_of<X>::f = nonzero if f implemented
135 template<typename U, int (U::*)(gsl::vector const&, gsl::matrix&)>
136 struct dfSfinae {};
137 template<typename U, int (U::*)(gsl::vector const&, gsl::matrix&) const>
138 struct dfSfinaeC {};
139 template<typename U> static df_t df(dfSfinae<U,&U::df>*){return &F<U>::dfn;}
140 template<typename U> static df_t df(dfSfinaeC<U,&U::df>*){return &F<U>::dfn;}
141 template<typename U> static df_t df(...){return 0;}
142 // address_of<X>::df = nonzero if df implemented
143 template<typename U, int (U::*)(gsl::vector const&, gsl::vector const&, gsl::vector&)>
144 struct fvvSfinae {};
145 template<typename U, int (U::*)(gsl::vector const&, gsl::vector const&, gsl::vector&)
146 const>
147 struct fvvSfinaeC {};
148 template<typename U> static fvv_t fvv(fvvSfinae<U,&U::fvv>*){return &F<U>::fvvn;}
149 template<typename U> static fvv_t fvv(fvvSfinaeC<U,&U::fvv>*){return &F<U>::fvvn;}
150 template<typename U> static fvv_t fvv(...){return 0;}
151 // address_of<X>::fvv = nonzero if fvv implemented
152 };
153 private:
157 struct base_F {
161 typedef int (*function_t)(gsl_vector const*,void*,gsl_vector*);
165 typedef int (*dfunction_t)(gsl_vector const*,void*,gsl_matrix*);
169 typedef int (*fvvfunction_t)(gsl_vector const*,gsl_vector const*,void*,gsl_vector*);
170 };
174 template<typename T>
175 struct F : public base_F {
176 public:
181 F( T& t ) : t( t ){
182 xv.wrap_gsl_vector_without_ownership( 0 );
183 vv.wrap_gsl_vector_without_ownership( 0 );
184 fv.wrap_gsl_vector_without_ownership( 0 );
185 Jv.wrap_gsl_matrix_without_ownership( 0 );
186 }
187 public:
196 static int fn( gsl_vector const* x, void* params, gsl_vector* f ){
197 // share gsl_vector* and gsl::vector contents.
198 F<T>* ft = reinterpret_cast<F<T>*>( params );
199 ft->xv.ccgsl_pointer = const_cast<gsl_vector*>( x );
200 ft->fv.ccgsl_pointer = f;
201 return ft->t.f( ft->xv, ft->fv ); }
210 static int dfn( gsl_vector const* x, void* params, gsl_matrix* J ){
211 // share gsl_vector* and gsl::vector contents.
212 F<T>* ft = reinterpret_cast<F<T>*>( params );
213 ft->xv.ccgsl_pointer = const_cast<gsl_vector*>( x );
214 ft->Jv.ccgsl_pointer = J;
215 return ft->t.df( ft->xv, ft->Jv ); }
225 static int fvvn( gsl_vector const* x, gsl_vector const* v,
226 void* params, gsl_vector* f ){
227 // share gsl_vector* and gsl::vector contents.
228 F<T>* ft = reinterpret_cast<F<T>*>( params );
229 ft->xv.ccgsl_pointer = const_cast<gsl_vector*>( x );
230 ft->vv.ccgsl_pointer = const_cast<gsl_vector*>( v );
231 ft->fv.ccgsl_pointer = f;
232 return ft->t.fvv( ft->xv, ft->vv, ft->fv ); }
233 private:
237 T& t;
241 gsl::vector xv;
245 gsl::vector vv;
249 gsl::vector fv;
253 gsl::matrix Jv;
254 };
258 class Fn : public base_F {
259 public:
265 Fn( int (*const f)(gsl::vector const&,gsl::vector&),
266 int (*const df)(gsl::vector const&,gsl::matrix&),
267 int (*const fvv)(gsl::vector const&,gsl::vector&,gsl::vector&),
268 size_t const xSize, size_t fSize )
269 : f( f ), df( df ), fvv( fvv ), xSize( xSize ), fSize( fSize ){}
273 function_t function(){ return &fn; }
277 dfunction_t dfunction(){ return &dfn; }
281 fvvfunction_t fvvfunction(){ return &fvvn; }
282 private:
286 int (*const f)(gsl::vector const&,gsl::vector&);
290 int (*const df)(gsl::vector const&,gsl::matrix&);
294 int (*const fvv)(gsl::vector const&,gsl::vector&,gsl::vector&);
295 public:
299 size_t const xSize;
303 size_t const fSize;
312 static int fn( gsl_vector const* x, void* params, gsl_vector* f ){
313 // share gsl_vector* and gsl::vector contents.
314 Fn* ft = reinterpret_cast<Fn*>( params );
315 ft->fv.ccgsl_pointer = f;
316 return ft->f( ft->xv, ft->fv ); }
325 static int dfn( gsl_vector const* x, void* params, gsl_matrix* J ){
326 // share gsl_vector* and gsl::vector contents.
327 Fn* ft = reinterpret_cast<Fn*>( params );
328 ft->xv.ccgsl_pointer = const_cast<gsl_vector*>( x );
329 ft->Jv.ccgsl_pointer = J;
330 return ft->df( ft->xv, ft->Jv ); }
340 static int fvvn( gsl_vector const* x, gsl_vector const* v,
341 void* params, gsl_vector* f ){
342 // share gsl_vector* and gsl::vector contents.
343 Fn* ft = reinterpret_cast<Fn*>( params );
344 ft->xv.ccgsl_pointer = const_cast<gsl_vector*>( x );
345 ft->vv.ccgsl_pointer = const_cast<gsl_vector*>( v );
346 ft->fv.ccgsl_pointer = f;
347 return ft->fvv( ft->xv, ft->vv, ft->fv ); }
348 private:
352 gsl::vector xv;
356 gsl::vector vv;
360 gsl::vector fv;
364 gsl::matrix Jv;
365 };
366#endif
370 base_F* func;
374 size_t* count;
375 public:
380 func = 0;
381 f = 0;
382 df = 0;
383 fvv = 0;
384 p = 0;
385 n = 0;
386 params = 0;
387 count = 0; // initially nullptr will do
388 }
389 template<typename T>
397 explicit function_fdf( gsl_multifit_nlinear_fdf& v ){
398 // explicitly set function = 0 to indicate nothing to be deleted.
399 func = 0;
400 // copy
401 f = v.f;
402 df = v.df;
403 fvv = v.fvv;
404 p = v.p;
405 n = v.n;
406 params = v.params;
407 // just plausibly we could fail to allocate count: no further action needed.
408 count = new size_t;
409 *count = 1; // initially there is just one reference to function
410 }
416 template<typename T>
417 function_fdf( T& t ){ function_constructor<T>( *this, t ); }
427 function_fdf( int (* const f)(gsl::vector const&, gsl::vector&),
428 int (* const df)(gsl::vector const&, gsl::matrix&),
429 int (* const fvv)(gsl::vector const&, gsl::vector&, gsl::vector&),
430 size_t const xSize, size_t const fSize ){
431 func = 0;
432 count = 0;
433 func = new function_fdf::Fn( f, df, fvv, xSize, fSize );
434 this->f = (0 == f) ? 0 : Fn::fn;
435 this->df = (0 == df) ? 0 : Fn::dfn;
436 this->fvv = (0 == fvv) ? 0 : Fn::fvvn;
437 this->n = fSize;
438 this->p = xSize;
439 params = func;
440 // just plausibly we could fail to allocate count: no further action needed.
441 count = new size_t;
442 *count = 1; // initially there is just one reference to function
443 }
444
445 // Copy and move constructors
446 // copy constructor
453 function_fdf( function_fdf const& v ) : func( v.func ), count( v.count ){
454 f = v.f;
455 df = v.df;
456 fvv = v.fvv;
457 n = v.n;
458 p = v.p;
459 params = v.params;
460 if( count != 0 ) ++*count; // function is now shared.
461 }
462 // assignment operator
468 // first, possibly delete anything pointed to by this
469 if( count == 0 or --*count == 0 ){
470 delete func;
471 delete count;
472 }
473 // Then copy
474 func = v.func;
475 f = v.f;
476 df = v.df;
477 fvv = v.fvv;
478 n = v.n;
479 p = v.p;
480 params = v.params;
481 // count = v.count;
482 if( count != 0 ) ++*count; // function is now shared.
483 return *this;
484 }
485#ifdef __GXX_EXPERIMENTAL_CXX0X__
490 function_fdf( function_fdf&& v ) : func( v.func ), count( nullptr ){
491 std::swap( f, v.f );
492 std::swap( df, v.df );
493 std::swap( fvv, v.fvv );
494 std::swap( n, v.n );
495 std::swap( p, v.p );
496 std::swap( params, v.params );
497 std::swap( count, v.count );
498 v.func = nullptr;
499 }
506 std::swap( func, v.func );
507 std::swap( f, v.f );
508 std::swap( df, v.df );
509 std::swap( fvv, v.fvv );
510 std::swap( n, v.n );
511 std::swap( p, v.p );
512 std::swap( params, v.params );
513 std::swap( count, v.count );
514 return *this;
515 }
516#endif
521 // first, possibly delete anything pointed to by f
522 if( count == 0 or --*count == 0 ){
523 delete func;
524 delete count;
525 }
526 }
527 };
528#ifndef DOXYGEN_SKIP
536 template<typename T>
537 inline void function_constructor( function_fdf& f, T& t ){
538 f.func = 0;
539 f.count = 0;
540 // try constructing from T
541 f.func = new function_fdf::F<T>( t );
542 f.f = function_fdf::address_of::f<T>(0);
543 f.df = function_fdf::address_of::df<T>(0);
544 f.fvv = function_fdf::address_of::fvv<T>(0);
545 // Checks coplete
546 f.n = t.fSize();
547 f.p = t.xSize();
548 f.params = f.func;
549 // just plausibly we could fail to allocate count: no further action needed.
550 f.count = new size_t;
551 *f.count = 1; // initially there is just one reference to function
552 }
558 template<>
559 inline void function_constructor( function_fdf& f, function_fdf& v ){
560 f.func = v.func; f.count = v.count; f.f = v.f; f.df = v.df; f.fvv = v.fvv;
561 f.n = v.n; f.p = v.p; f.params = v.params; if( f.count != 0 ) ++*f.count;
562 /* function is now shared. */ }
568 template<>
569 inline void function_constructor( function_fdf& f, gsl_multifit_nlinear_fdf& v ){
570 f.f = v.f; f.df = v.df; f.fvv = v.fvv; f.n = v.n; f.p = v.p; f.params = v.params; }
571#endif
577 template<typename T>
579 function_fdf fn( t );
580 return fn;
581 }
582 }
583 }
584}
585
586#endif
Class that extends gsl_function_fdf so that it can be constructed from arbitrary function objects.
size_t * count
The shared reference count: used for copying this.
This class handles matrix objects as shared handles.
Definition: matrix.hpp:72
Class that extends gsl_multifit_nlinear_fdf so that it can be constructed from arbitrary function obj...
function_fdf & operator=(function_fdf const &v)
The assignment operator.
function_fdf & operator=(function_fdf &&v)
Move operator.
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(gsl_multifit_nlinear_fdf &v)
Could construct from a gsl_multifit_nlinear_fdf.
base_F * func
This gives something for params to point to.
friend void function_constructor(function_fdf &, T &)
function_fdf(function_fdf const &v)
The copy constructor.
function_fdf(int(*const f)(gsl::vector const &, gsl::vector &), int(*const df)(gsl::vector const &, gsl::matrix &), int(*const fvv)(gsl::vector const &, gsl::vector &, gsl::vector &), size_t const xSize, size_t const fSize)
Construct from a function with no parameters (but with n function values and arguments).
~function_fdf()
The destructor unshares any shared resource.
function_fdf(function_fdf &&v)
Move constructor.
function_fdf(T &t)
Construct from an object that implements gsl::multifit::nlinear::function_fdf::concept.
This class handles vector objects as shared handles.
Definition: vector.hpp:74
void function_constructor(function_fdf &f, T &t)
function_fdf make_function_fdf(T &t)
Make a gsl::multifit::nlinear::function_fdf from a function object that implements gsl::multifit::nli...
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().
size_t n(workspace const &w)
C++ version of gsl_rstat_n().
Definition: rstat.hpp:299
double F(double phi, double k, mode_t mode)
C++ version of gsl_sf_ellint_F().
Definition: sf_ellint.hpp:138
The gsl package creates an interface to the GNU Scientific Library for C++.
Definition: blas.hpp:34
virtual int df(gsl::vector const &x, gsl::matrix &J)=0
The derivatives (as Jacobian matrix).
virtual int f(gsl::vector const &x, gsl::vector &f)=0
The function.
virtual size_t xSize() const =0
This function should return the number of elements of x in f().
virtual size_t fSize() const =0
This function should return the number of elements of f in f().
virtual int fvv(gsl::vector const &x, gsl::vector const &v, gsl::vector &fvv)=0
The second directional derivatives.