Skip to content

Commit 66d6d37

Browse files
committed
added scipy.C macro in fit tutorial
1 parent 403d32f commit 66d6d37

File tree

1 file changed

+83
-0
lines changed

1 file changed

+83
-0
lines changed

tutorials/fit/scipy.C

Lines changed: 83 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,83 @@
1+
#include "Math/ScipyMinimizer.h"
2+
#include "Math/MultiNumGradFunction.h"
3+
#include <Fit/ParameterSettings.h>
4+
#include "Math/Functor.h"
5+
#include <string>
6+
#include "Math/MinimizerOptions.h"
7+
#include "TStopwatch.h"
8+
9+
10+
double RosenBrock(const double *xx )
11+
{
12+
const Double_t x = xx[0];
13+
const Double_t y = xx[1];
14+
const Double_t tmp1 = y-x*x;
15+
const Double_t tmp2 = 1-x;
16+
return 100*tmp1*tmp1+tmp2*tmp2;
17+
}
18+
19+
double RosenBrockGrad(const double *x, unsigned int ipar)
20+
{
21+
if (ipar == 0)
22+
return -2 * (1 - x[0]) + 200 * (x[1] - x[0] * x[0]) * (-2 * x[0]);
23+
else
24+
return 200 * (x[1] - x[0] * x[0]);
25+
}
26+
27+
bool RosenBrockHessian(const std::vector<double> &xx, double *hess)
28+
{
29+
const double x = xx[0];
30+
const double y = xx[1];
31+
32+
hess[0] = 1200*x*x - 400*y + 2;
33+
hess[1] = -400*x;
34+
hess[2] = -400*x;
35+
hess[3] = 200;
36+
37+
return true;
38+
}
39+
40+
// methods that requires hessian to work "dogleg", "trust-ncg","trust-exact","trust-krylov"
41+
using namespace std;
42+
int scipy()
43+
{
44+
45+
std::string methods[]={"Nelder-Mead","L-BFGS-B","Powell","CG","BFGS","TNC","COBYLA","SLSQP","trust-constr","Newton-CG", "dogleg", "trust-ncg","trust-exact","trust-krylov"};
46+
TStopwatch t;
47+
for(const std::string &text : methods)
48+
{
49+
ROOT::Math::Experimental::ScipyMinimizer minimizer(text.c_str());
50+
minimizer.SetMaxFunctionCalls(1000000);
51+
minimizer.SetMaxIterations(100000);
52+
minimizer.SetTolerance(1e-3);
53+
minimizer.SetExtraOption("gtol",1e-3);
54+
ROOT::Math::GradFunctor f(&RosenBrock,&RosenBrockGrad,2);
55+
double step[2] = {0.01,0.01};
56+
double variable[2] = { -1.2,1.0};
57+
58+
minimizer.SetFunction(f);
59+
minimizer.SetHessianFunction(RosenBrockHessian);
60+
61+
// variables to be minimized!
62+
minimizer.SetVariable(0,"x",variable[0], step[0]);
63+
minimizer.SetVariable(1,"y",variable[1], step[1]);
64+
minimizer.SetVariableLimits(0, -2.0, 2.0);
65+
minimizer.SetVariableLimits(1, -2.0, 2.0);
66+
67+
t.Reset();
68+
t.Start();
69+
minimizer.Minimize();
70+
t.Stop();
71+
const double *xs = minimizer.X();
72+
cout << "Minimum: f(" << xs[0] << "," << xs[1] << "): "
73+
<< RosenBrock(xs) << endl;
74+
cout << "Cpu Time (sec) = " << t.CpuTime() <<endl<< "Real Time (sec) = " << t.RealTime() << endl;
75+
cout << endl << "===============" << endl;
76+
}
77+
return 0;
78+
}
79+
80+
int main()
81+
{
82+
return scipy();
83+
}

0 commit comments

Comments
 (0)