ccgsl 2.7.2
C++wrappersforGnuScientificLibrary
multimin_function_fdf.hpp
Go to the documentation of this file.
1/*
2 * $Id: multimin_function_fdf.hpp 130 2012-03-27 09:45:00Z jdl3 $
3 * Copyright (C) 2010, 2011, 2012, 2016, 2020, 2021 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_MULTIMIN_FUNCTION_FDF_HPP
21#define CCGSL_MULTIMIN_FUNCTION_FDF_HPP
22
23// For std::pair
24#include<utility>
26// For std::swap
27#if __cplusplus <= 199711L
28#include<algorithm>
29#else
30#include<utility>
31#endif
32
33namespace gsl {
34 namespace multimin {
35 class function_fdf;
36 template<typename T>
71 class function_fdf : public gsl_multimin_function_fdf {
72 public:
80 struct Concept {
86 virtual double f( const vector& x ) = 0;
92 virtual void df( const vector& x, vector& gradient ) = 0;
99 virtual void fdf( const vector& x, double* result, vector& gradient ) = 0;
104 size_t size() const;
105 };
106#ifndef DOXYGEN_SKIP
107 private:
111 struct base_F {
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&);
124 };
128 template<typename T>
129 struct FDF : public base_F {
130 public:
135 FDF( T& t ) : t( t ){
136 xv.wrap_gsl_vector_without_ownership( 0 );
137 dfv.wrap_gsl_vector_without_ownership( 0 );
138 }
139 public:
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 ); }
173 private:
177 gsl::vector xv;
181 gsl::vector dfv;
185 T& t;
186 };
187#endif
191 base_F* function;
195 size_t* count;
196 public:
201 function = 0;
202 f = 0;
203 df = 0;
204 fdf = 0;
205 params = 0;
206 count = 0; // initially nullptr will do
207 }
208 template<typename T>
215 explicit function_fdf( gsl_multimin_function_fdf& v ){
216 // explicitly set function = 0 to indicate nothing to be deleted.
217 function = 0;
218 n = 0;
219 // copy
220 f = v.f;
221 df = v.df;
222 fdf = v.fdf;
223 params = v.params;
224 // just plausibly we could fail to allocate count: no further action needed.
225 count = new size_t;
226 *count = 1; // initially there is just one reference to function
227 }
232 template<typename T>
234 function_fdf_constructor<T>( *this, t );
235 }
236
237 // Copy and move constructors
238 // copy constructor
246 f = v.f;
247 n = v.n;
248 df = v.df;
249 fdf = v.fdf;
250 params = v.params;
251 if( count != 0 ) ++*count; // function_fdf is now shared.
252 }
253 // assignment operator
259 // first, possibly delete anything pointed to by this
260 if( count == 0 or --*count == 0 ){
261 delete function;
262 delete count;
263 }
264 // Then copy
265 function = v.function;
266 f = v.f;
267 n = v.n;
268 df = v.df;
269 fdf = v.fdf;
270 params = v.params;
271 // count = v.count;
272 if( count != 0 ) ++*count; // function_fdf is now shared.
273 return *this;
274 }
275#ifdef __GXX_EXPERIMENTAL_CXX0X__
281 std::swap( f, v.f );
282 std::swap( n, v.n );
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;
288 }
295 std::swap( function, v.function );
296 std::swap( f, v.f );
297 std::swap( n, v.n );
298 std::swap( df, v.df );
299 std::swap( fdf, v.fdf );
300 std::swap( params, v.params );
301 std::swap( count, v.count );
302 return *this;
303 }
304#endif
309 // first, possibly delete anything pointed to by f
310 if( count == 0 or --*count == 0 ){
311 delete function;
312 delete count;
313 }
314 }
315 };
316#ifndef DOXYGEN_SKIP
322 template<typename T>
323 inline void function_fdf_constructor( function_fdf& fdf, T& t ){
324 fdf.function = 0;
325 fdf.count = 0;
326 // try constructing from T
327 fdf.function = new function_fdf::FDF<T>( t );
328 fdf.n = t.size();
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;
333 // just plausibly we could fail to allocate count: no further action needed.
334 fdf.count = new size_t;
335 *fdf.count = 1; // initially there is just one reference to function
336 }
340 template<>
341 inline void function_fdf_constructor( function_fdf& fdf, function_fdf& v ){
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;
344 /* function_fdf is now shared. */ }
348 template<>
349 inline void function_fdf_constructor( function_fdf& fdf, gsl_multimin_function_fdf& v ){
350 fdf.f = v.f; fdf.n = v.n; fdf.df = v.df; fdf.fdf = v.fdf; fdf.params = v.params; }
351#endif
352
358 template<typename T>
360 function_fdf fn( t );
361#ifdef __GXX_EXPERIMENTAL_CXX0X__
362 return std::move( fn );
363#else
364 return fn;
365#endif
366 }
367
368 }
369}
370
371#endif
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.
Definition: vector.hpp:74
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().
Definition: multimin.hpp:402
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().
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
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.