12 const ROOT::Math::IMultiGenFunction *gFunction;
14 const ROOT::Math::IMultiGradFunction *gGradFunction;
19 double minfunction(
const std::vector<double> & x){
22 return (*gFunction)(x.data());
25 TVectorD mingradfunction(TVectorD y){
26 unsigned int size = y.GetNoElements();
27 const double * yy = y.GetMatrixArray();
29 gGradFunction->Gradient(yy,z);
38 RMinimizer::RMinimizer(Option_t *method){
40 if (fMethod.empty() || fMethod==
"Migrad") fMethod=
"BFGS";
44 unsigned int RMinimizer::NCalls()
const {
return gNCalls; }
47 bool RMinimizer::Minimize() {
50 (gFunction)= ObjFunction();
51 (gGradFunction) = GradObjFunction();
56 ROOT::R::TRInterface &r=ROOT::R::TRInterface::Instance();
58 r[
"minfunction"] = ROOT::R::TRFunctionExport(minfunction);
59 r[
"mingradfunction"] = ROOT::R::TRFunctionExport(mingradfunction);
60 r[
"method"] = fMethod.c_str();
61 std::vector<double> stepSizes(StepSizes(), StepSizes()+NDim());
62 std::vector<double> values(X(), X()+NDim());
65 r[
"stepsizes"] = stepSizes;
66 r[
"initialparams"] = values;
69 bool optimxloaded = FALSE;
70 r[
"optimxloaded"] = optimxloaded;
71 r.Execute(
"optimxloaded<-library(optimx,logical.return=TRUE)");
73 int ibool = r.Eval(
"optimxloaded");
74 if (ibool==1) optimxloaded=kTRUE;
80 if (optimxloaded==kTRUE) {
83 cmd = TString::Format(
"result <- optimx( initialparams, minfunction,method='%s',control = list(ndeps=stepsizes,maxit=%d,trace=%d,abstol=%e),hessian=TRUE)",fMethod.c_str(),MaxIterations(),PrintLevel(),Tolerance());
87 cmd = TString::Format(
"result <- optimx( initialparams, minfunction,mingradfunction, method='%s', control = list(ndeps=stepsizes,maxit=%d,trace=%d,abstol=%e),hessian=TRUE)",fMethod.c_str(),MaxIterations(),PrintLevel(),Tolerance());
96 cmd = TString::Format(
"result <- optim( initialparams, minfunction,method='%s',control = list(ndeps=stepsizes,maxit=%d,trace=%d,abstol=%e),hessian=TRUE)",fMethod.c_str(),MaxIterations(),PrintLevel(),Tolerance());
100 cmd = TString::Format(
"result <- optim( initialparams, minfunction,mingradfunction, method='%s', control = list(ndeps=stepsizes,maxit=%d,trace=%d,abstol=%e),hessian=TRUE)",fMethod.c_str(),MaxIterations(),PrintLevel(),Tolerance());
104 std::cout <<
"Calling R with command " << cmd << std::endl;
105 r.Execute(cmd.Data());
110 r.Execute(
"par<-coef(result)");
112 r.Execute(
"hess<-attr(result,\"details\")[,\"nhatend\"]");
114 r.Execute(
"hess<-sapply(hess,function(x) x)");
116 r.Execute(
"hess<-matrix(hess,c(ndim,ndim))");
118 r.Execute(
"cov<-solve(hess)");
120 r.Execute(
"errors<-sqrt(abs(diag(cov)))");
125 r.Execute(
"par<-result$par");
126 r.Execute(
"hess<-result$hessian");
127 r.Execute(
"cov<-solve(hess)");
128 r.Execute(
"errors<-sqrt(abs(diag(cov)))");
133 std::vector<double> vectorPar = r[
"par"];
138 TMatrixD cm = r[
"cov"];
141 std::vector<double> err = r[
"errors"];
144 TMatrixD hm = r[
"hess"];
147 fCovMatrix.ResizeTo(ndim,ndim);
148 fHessMatrix.ResizeTo(ndim,ndim);
155 const double *min=vectorPar.data();
157 SetMinValue((*gFunction)(min));
158 std::cout<<
"Value at minimum ="<<MinValue()<<std::endl;
164 double RMinimizer::CovMatrix(
unsigned int i,
unsigned int j)
const {
165 unsigned int ndim = NDim();
166 if (fCovMatrix==0)
return 0;
167 if (i > ndim || j > ndim)
return 0;
168 return fCovMatrix[i][j];
175 double RMinimizer::HessMatrix(
unsigned int i,
unsigned int j)
const {
176 unsigned int ndim = NDim();
177 if (fHessMatrix==0)
return 0;
178 if (i > ndim || j > ndim)
return 0;
179 return fHessMatrix[i][j];