8 #include "gsl/gsl_vector.h"
9 #include "gsl/gsl_multimin.h"
10 #include "gsl/gsl_math.h"
11 #include "gsl/gsl_matrix.h"
12 #include "gsl/gsl_linalg.h"
13 #include "gsl/gsl_blas.h"
14 #include "gsl/gsl_errno.h"
21 using namespace CLHEP;
33 (
const FuncMinimum::GenFunc&
func ,
34 FuncMinimum::Arg& arg )
39 const size_t N = func.dimensionality () ;
41 for(
size_t i = 0 ;
i <
N ; ++
i )
43 Genfun::GENFUNCTION
fun = func.partial(
i);
44 m_grad.push_back (fun.clone());
73 double fun_gsl (
const gsl_vector* v ,
78 const FuncMinimum::GenFunc& eq = *(local->
equation());
81 for (
unsigned int i = 0;
i < v->size; ++
i) {
82 arg[
i] = gsl_vector_get (v,
i);
89 void dfun_gsl (
const gsl_vector* v ,
void * params,
97 for (
unsigned int i = 0;
i < v->size; ++
i) {
98 arg[
i] = gsl_vector_get (v,
i);
102 for(
unsigned int i = 0 ;
i < df->size ; ++
i )
104 Genfun::GENFUNCTION f = *(grad[
i]);
105 gsl_vector_set ( df,
i, f(arg) );
113 void fdfun_gsl (
const gsl_vector* v ,
118 *f = fun_gsl( v , params );
119 dfun_gsl ( v , params, df);
135 gsl_vector_view vect = gsl_vector_view_array ( &arg[0] ,
137 FuncMinimumMisc local (func, arg);
139 gsl_multimin_function_fdf
function;
141 function.f = &fun_gsl;
142 function.df = &dfun_gsl;
143 function.fdf = &fdfun_gsl;
144 function.n = vect.vector.size;
145 function.params = (
void*) &local;
149 const gsl_multimin_fdfminimizer_type *T = m_type ;
151 gsl_multimin_fdfminimizer *
s;
153 s = gsl_multimin_fdfminimizer_alloc ( T, vect.vector.size);
155 gsl_multimin_fdfminimizer_set ( s, &
function,
156 &vect.vector, m_step_size, m_tol);
158 for( iter = 0 ; iter < m_max_iter ; ++iter )
160 status = gsl_multimin_fdfminimizer_iterate (s);
165 (
"Error from gsl_multimin_fdfminimizer_iterate '"
169 status = gsl_multimin_test_gradient (s->gradient,
173 if ( status != GSL_CONTINUE ) {
break; }
176 for (
unsigned int i = 0;
i < vect.vector.size; ++
i)
178 gsl_vector_set (&vect.vector,
i, gsl_vector_get (s->x,
i));
181 if (status == GSL_SUCCESS)
184 <<
"We stopped in the method on the " << iter
185 <<
" iteration (we have maximum " << m_max_iter
186 <<
" iterations)" <<
endmsg;
188 msgStream() <<
"The Euclidean norm of gradient = "
189 << gsl_blas_dnrm2 (s->gradient)
190 <<
" by the absolute tolerance = "
191 << m_norm_gradient <<
endmsg;
193 else if (status == GSL_CONTINUE && iter <= m_max_iter )
195 return Error (
"Method finished with '"
201 return Error (
"Method finished with '" +
206 gsl_multimin_fdfminimizer_free (s);
210 return Error (
"Method finished with '"
233 return Error (
"MINIMUM IS NOT FOUND. StatusCode = '"
239 HepSymMatrix cov(arg.dimension(), 0);
240 for (
unsigned int i = 0;
i < arg.dimension(); ++
i)
242 Genfun::GENFUNCTION f = func.partial(
i);
243 for (
unsigned int j =
i; j < arg.dimension(); ++j)
245 Genfun::GENFUNCTION fij = f.partial(j);
246 cov(
i+1, j+1) = 0.5 * fij(arg);
251 covar = cov.inverse(inv);
255 (
"Matrix of Error is not complete successful");
271 return Error (
"Could not initialize base class GaudiTool", sc);
277 m_type = gsl_multimin_fdfminimizer_conjugate_fr ;
279 <<
"Minimization algorithm to be used: "
280 <<
"'gsl_multimin_fdfminimizer_conjugate_fr'"
285 m_type = gsl_multimin_fdfminimizer_conjugate_pr ;
287 <<
"Minimization algorithm to be used: "
288 <<
"'gsl_multimin_fdfminimizer_conjugate_pr'"
293 m_type = gsl_multimin_fdfminimizer_vector_bfgs ;
295 <<
"Minimization algorithm to be used: " <<
296 "'gsl_multimin_fdfminimizer_vector_bfgs'" <<
endmsg;
298 else if (
"steepest_descent" ==
m_algType )
300 m_type = gsl_multimin_fdfminimizer_steepest_descent ;
302 <<
"Minimization algorithm to be used: "
303 <<
"'gsl_multimin_fdfminimizer_steepest_descent'"
320 return Error(
"Could not finalize base class GaudiTool", sc);
unsigned long getCode() const
Get the status code by value.
Genfun::Argument Arg
Argument of function "GenFunc" (.
const Arg & argument() const
const GenFunc * equation() const
const gsl_multimin_fdfminimizer_type * m_type
bool isFailure() const
Test for a status code of FAILURE.
#define DECLARE_COMPONENT(type)
StatusCode minimum(const GenFunc &func, Arg &arg) const override
Find minimum of the function "GenFunc".
Genfun::AbsFunction GenFunc
Function which we minimize (.
This class is used for returning status codes from appropriate routines.
Definition of the basic interface.
The simplest concrete implementation of IFuncMinimum interface.
const Gradient & gradient() const
MsgStream & debug() const
shortcut for the method msgStream(MSG::DEBUG)
CLHEP::HepSymMatrix Covariance
Covariance matrix (matrix of error) (.
double fun(const std::vector< double > &x)
StatusCode finalize() override
FuncMinimum()=delete
default constructor is private
MsgStream & endmsg(MsgStream &s)
MsgStream Modifier: endmsg. Calls the output method of the MsgStream.
StatusCode initialize() override
Overriding initialize.