Logo ROOT  
Reference Guide
GSLMultiFitFunctionWrapper.h
Go to the documentation of this file.
1// @(#)root/mathmore:$Id$
2// Authors: L. Moneta Dec 2006
3
4 /**********************************************************************
5 * *
6 * Copyright (c) 2004 ROOT Foundation, CERN/PH-SFT *
7 * *
8 * This library is free software; you can redistribute it and/or *
9 * modify it under the terms of the GNU General Public License *
10 * as published by the Free Software Foundation; either version 2 *
11 * of the License, or (at your option) any later version. *
12 * *
13 * This library is distributed in the hope that it will be useful, *
14 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
15 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU *
16 * General Public License for more details. *
17 * *
18 * You should have received a copy of the GNU General Public License *
19 * along with this library (see file COPYING); if not, write *
20 * to the Free Software Foundation, Inc., 59 Temple Place, Suite *
21 * 330, Boston, MA 02111-1307 USA, or contact the author. *
22 * *
23 **********************************************************************/
24
25// Header file for class GSLMultiMinFunctionWrapper
26//
27// Created by: moneta at Sat Nov 13 14:54:41 2004
28//
29// Last update: Sat Nov 13 14:54:41 2004
30//
31#ifndef ROOT_Math_GSLMultiFitFunctionWrapper
32#define ROOT_Math_GSLMultiFitFunctionWrapper
33
34#include "gsl/gsl_multifit.h"
35
37
38
39#include <cassert>
40
41namespace ROOT {
42namespace Math {
43
44
45
46 typedef double ( * GSLMultiFitFPointer ) ( const gsl_vector *, void *, gsl_vector *);
47 typedef void ( * GSLMultiFitDfPointer ) ( const gsl_vector *, void *, gsl_matrix *);
48 typedef void ( * GSLMultiFitFdfPointer ) ( const gsl_vector *, void *, gsl_vector *, gsl_matrix *);
49
50
51/**
52 wrapper to a multi-dim function withtout derivatives for multi-dimensional
53 minimization algorithm
54
55 @ingroup MultiMin
56*/
57
59
60public:
61
63 {
64 fFunc.f = 0;
65 fFunc.df = 0;
66 fFunc.fdf = 0;
67 fFunc.n = 0;
68 fFunc.p = 0;
69 fFunc.params = 0;
70#if GSL_MAJOR_VERSION > 1
71 fFunc.nevalf = 0;
72 fFunc.nevaldf = 0;
73#endif
74 }
75
76
77 /// Fill gsl function structure from a C++ function iterator and size and number of residuals
78 template<class FuncVector>
79 void SetFunction(const FuncVector & f, unsigned int nres, unsigned int npar ) {
80 const void * p = &f;
81 assert (p != 0);
85 fFunc.n = nres;
86 fFunc.p = npar;
87 fFunc.params = const_cast<void *>(p);
88 }
89
90 gsl_multifit_function_fdf * GetFunc() { return &fFunc; }
91
92
93 private:
94
95 gsl_multifit_function_fdf fFunc;
96
97};
98
99
100
101} // namespace Math
102} // namespace ROOT
103
104#endif /* ROOT_Math_GSLMultiMinFunctionWrapper */
double
Definition: Converters.cxx:921
#define f(i)
Definition: RSha256.hxx:104
typedef void((*Func_t)())
Class for adapting a C++ functor class to C function pointers used by GSL MultiFit Algorithm The temp...
wrapper to a multi-dim function withtout derivatives for multi-dimensional minimization algorithm
void SetFunction(const FuncVector &f, unsigned int nres, unsigned int npar)
Fill gsl function structure from a C++ function iterator and size and number of residuals.
Namespace for new Math classes and functions.
double(* GSLMultiFitFPointer)(const gsl_vector *, void *, gsl_vector *)
void(* GSLMultiFitDfPointer)(const gsl_vector *, void *, gsl_matrix *)
void(* GSLMultiFitFdfPointer)(const gsl_vector *, void *, gsl_vector *, gsl_matrix *)
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...
Definition: StringConv.hxx:21