00001
00002
00003 #include "Math/Minimizer.h"
00004 #include "Math/Factory.h"
00005 #include "Math/Functor.h"
00006
00007 #include "TVirtualFitter.h"
00008
00009 #include "Math/IFunction.h"
00010 #include "Math/Util.h"
00011 #include <cmath>
00012 #include <cassert>
00013
00014 #include <string>
00015 #include <iostream>
00016
00017 #include "TStopwatch.h"
00018 #include "TMatrixD.h"
00019 #include "TVectorD.h"
00020 #include "TRandom3.h"
00021 #include "TMath.h"
00022
00023
00024
00025 int gNCall = 0;
00026 int gNCall2 = 0;
00027 int gNmin = 1000;
00028 int gVerbose = 0;
00029 bool useGradient = false;
00030
00031 bool minos = true;
00032
00033 double gAbsTolerance = 0.005;
00034
00035
00036
00037 typedef void (*FCN)(Int_t &npar, Double_t *gin, Double_t &f, Double_t *u, Int_t flag);
00038
00039
00040
00041
00042
00043 void RosenBrock(Int_t &, Double_t *, Double_t &f, Double_t *par, Int_t )
00044 {
00045 gNCall++;
00046 const Double_t x = par[0];
00047 const Double_t y = par[1];
00048 const Double_t tmp1 = y-x*x;
00049 const Double_t tmp2 = 1-x;
00050 f = 100*tmp1*tmp1+tmp2*tmp2;
00051 }
00052
00053
00054
00055 class RosenBrockFunction : public ROOT::Math::IMultiGenFunction {
00056
00057 public :
00058
00059 virtual ~RosenBrockFunction() {}
00060
00061 unsigned int NDim() const { return 2; }
00062
00063 ROOT::Math::IMultiGenFunction * Clone() const {
00064 return new RosenBrockFunction();
00065 }
00066
00067 const double * TrueMinimum() const {
00068 fTrueMin[0] = 1;
00069 fTrueMin[1] = 1;
00070 return fTrueMin;
00071 }
00072
00073 private:
00074
00075 inline double DoEval (const double * x) const {
00076 #ifdef USE_FREE_FUNC
00077 double f = 0;
00078 int ierr = 0;
00079 int i = 0;
00080 RosenBrock(i,0,f,const_cast<double *>(x),ierr);
00081 return f;
00082 #else
00083 gNCall++;
00084 const Double_t xx = x[0];
00085 const Double_t yy = x[1];
00086 const Double_t tmp1 = yy-xx*xx;
00087 const Double_t tmp2 = 1-xx;
00088 return 100*tmp1*tmp1+tmp2*tmp2;
00089 #endif
00090 }
00091
00092
00093 mutable double fTrueMin[2];
00094 };
00095
00096
00097
00098
00099 class TrigoFletcherFunction : public ROOT::Math::IMultiGradFunction {
00100
00101 public :
00102
00103
00104 TrigoFletcherFunction(unsigned int dim) : fDim(dim) {
00105 double seed = 3;
00106 A.ResizeTo(dim,dim);
00107 B.ResizeTo(dim,dim);
00108 x0.ResizeTo(dim);
00109 sx0.ResizeTo(dim);
00110 cx0.ResizeTo(dim);
00111 sx.ResizeTo(dim);
00112 cx.ResizeTo(dim);
00113 v0.ResizeTo(dim);
00114 v.ResizeTo(dim);
00115 r.ResizeTo(dim);
00116 A.Randomize(-100.,100,seed);
00117 B.Randomize(-100.,100,seed);
00118 for (unsigned int i = 0; i < dim; i++) {
00119 for (unsigned int j = 0; j < dim; j++) {
00120 A(i,j) = int(A(i,j));
00121 B(i,j) = int(B(i,j));
00122 }
00123 }
00124 x0.Randomize(-TMath::Pi(),TMath::Pi(),seed);
00125
00126 for (unsigned int i = 0; i < fDim ; ++i) {
00127 cx0[i] = std::cos(x0[i]);
00128 sx0[i] = std::sin(x0[i]);
00129 }
00130 v0 = A*sx0+B*cx0;
00131 }
00132
00133
00134 unsigned int NDim() const { return fDim; }
00135
00136 ROOT::Math::IMultiGenFunction * Clone() const {
00137 TrigoFletcherFunction * f = new TrigoFletcherFunction(*this);
00138
00139
00140 return f;
00141 }
00142
00143
00144 void StartPoints(double * x, double * s) {
00145 TRandom3 rndm;
00146 const double stepSize = 0.01;
00147 const double deltaAmp = 0.1;
00148 const double pi = TMath::Pi();
00149 for (unsigned int i = 0; i < fDim; ++i) {
00150 double delta = rndm.Uniform(-deltaAmp*pi,deltaAmp*pi);
00151 x[i] = x0(i) + 0.1*delta;
00152 if (x[i] <= - pi) x[i] += 2.*pi;
00153 if (x[i] > pi) x[i] -= 2.*pi;
00154 s[i] = stepSize;
00155 }
00156 }
00157
00158
00159 const double * TrueMinimum() const {
00160 return x0.GetMatrixArray();
00161 }
00162
00163
00164 void Gradient (const double * x, double * g) const {
00165 gNCall2++;
00166
00167 for (unsigned int i = 0; i < fDim ; ++i) {
00168 cx [i] = std::cos(x[i]);
00169 sx [i] = std::sin(x[i]);
00170 }
00171
00172 v = A*sx +B*cx;
00173 r = v0-v;
00174
00175
00176
00177 for (unsigned int i = 0; i < fDim ; ++i) {
00178 g[i] = 0;
00179 for (unsigned int k = 0; k < fDim ; ++k) {
00180 g[i] += 2. * r(k) * ( - A(k,i) * cx(i) + B(k,i) * sx(i) );
00181 }
00182 }
00183
00184 }
00185
00186 #ifdef USE_FDF
00187 void FdF (const double * x, double & f, double * g) const {
00188 gNCall++;
00189
00190 for (unsigned int i = 0; i < fDim ; ++i) {
00191 cx [i] = std::cos(x[i]);
00192 sx [i] = std::sin(x[i]);
00193 }
00194
00195 v = A*sx +B*cx;
00196 r = v0-v;
00197
00198 f = r * r;
00199
00200
00201
00202 for (unsigned int i = 0; i < fDim ; ++i) {
00203 g[i] = 0;
00204 for (unsigned int k = 0; k < fDim ; ++k) {
00205 g[i] += 2. * r(k) * ( - A(k,i) * cx(i) + B(k,i) * sx(i) );
00206 }
00207 }
00208 }
00209 #endif
00210
00211 private:
00212
00213
00214
00215
00216 double DoEval (const double * x) const {
00217 gNCall++;
00218
00219
00220 for (unsigned int i = 0; i < fDim ; ++i) {
00221 cx [i] = std::cos(x[i]);
00222 sx [i] = std::sin(x[i]);
00223 }
00224
00225 v = A*sx +B*cx;
00226 r = v0-v;
00227
00228 return r * r;
00229 }
00230
00231
00232 double DoDerivative (const double * x, unsigned int i ) const {
00233 std::vector<double> g(fDim);
00234 Gradient(x,&g[0]);
00235 return g[i];
00236 }
00237
00238 private:
00239
00240 unsigned int fDim;
00241
00242 TMatrixD A;
00243 TMatrixD B;
00244 TVectorD x0;
00245 mutable TVectorD sx0;
00246 mutable TVectorD cx0;
00247 mutable TVectorD sx;
00248 mutable TVectorD cx;
00249 mutable TVectorD v0;
00250 mutable TVectorD v;
00251 mutable TVectorD r;
00252
00253
00254 };
00255
00256
00257
00258
00259 class ChebyQuadFunction : public ROOT::Math::IMultiGradFunction {
00260
00261 public :
00262
00263 ChebyQuadFunction(unsigned int n) :
00264 fDim(n),
00265 fvec(std::vector<double>(n) ),
00266 fTrueMin(std::vector<double>(n) )
00267 {
00268 }
00269
00270 unsigned int NDim() const { return fDim; }
00271
00272 ROOT::Math::IMultiGenFunction * Clone() const {
00273 return new ChebyQuadFunction(*this);
00274 }
00275
00276 const double * TrueMinimum() const {
00277 return &fTrueMin[0];
00278 }
00279
00280
00281 void StartPoints(double * x, double * s) {
00282 for (unsigned int i = 0; i < fDim; ++i) {
00283 s[i] = 0.01;
00284 x[i] = double(i)/(double(fDim)+1.0);
00285 }
00286 }
00287
00288
00289
00290 void Gradient(const double * x, double * g) const {
00291 gNCall2++;
00292 unsigned int n = fDim;
00293
00294 DoCalculatefi(x);
00295
00296 for (unsigned int j = 0; j < n; ++j) {
00297 g[j] = 0.0;
00298 double t1 = 1.0;
00299 double t2 = 2.0 * x[j] - 1.0;
00300 double t = 2.0 * t2;
00301 double s1 = 0.0;
00302 double s2 = 2.0;
00303 for (unsigned int i = 0; i < n; ++i) {
00304 g[j] += fvec[i] * s2;
00305 double th = 4.0 * t2 + t * s2 - s1;
00306 s1 = s2;
00307 s2 = th;
00308 th = t * t2 - t1;
00309 t1 = t2;
00310 t2 = th;
00311 }
00312 g[j] = 2. * g[j] / double(n);
00313 }
00314
00315
00316 }
00317
00318 private:
00319
00320 double DoEval (const double * x) const {
00321
00322 gNCall++;
00323 DoCalculatefi(x);
00324 double f = 0;
00325 for (unsigned int i = 0; i < fDim; ++i)
00326 f += fvec[i] * fvec[i];
00327
00328 return f;
00329
00330 }
00331
00332 double DoDerivative (const double * x, unsigned int i ) const {
00333 std::vector<double> g(fDim);
00334 Gradient(x,&g[0]);
00335 return g[i];
00336 }
00337
00338 void DoCalculatefi(const double * x) const {
00339
00340 unsigned int n = fDim;
00341 for (unsigned int i = 0; i < n; ++i)
00342 fvec[i] = 0;
00343
00344 for (unsigned int j = 0; j < n; ++j) {
00345 double t1 = 1.0;
00346 double t2 = 2.0 * x[j] - 1.0;
00347 double t = 2.0 * t2;
00348 for (unsigned int i = 0; i < n; ++i) {
00349 fvec[i] += t2;
00350 double th = t * t2 - t1;
00351 t1 = t2;
00352 t2 = th;
00353 }
00354 }
00355
00356
00357 for (unsigned int i = 1; i <= n; ++i) {
00358 int l = i-1;
00359 fvec[l] /= double ( n );
00360 if ( ( i % 2 ) == 0 ) {
00361 fvec[l] += 1.0 / ( double ( i*i ) - 1.0 );
00362 }
00363 }
00364 }
00365
00366 unsigned int fDim;
00367 mutable std::vector<double> fvec;
00368 mutable std::vector<double> fTrueMin;
00369 };
00370
00371
00372
00373 double WoodFunction(const double * par) {
00374 gNCall++;
00375
00376 const double w = par[0];
00377 const double x = par[1];
00378 const double y = par[2];
00379 const double z = par[3];
00380
00381 const double w1 = w-1;
00382 const double x1 = x-1;
00383 const double y1 = y-1;
00384 const double z1 = z-1;
00385 const double tmp1 = x-w*w;
00386 const double tmp2 = z-y*y;
00387
00388 double f = 100*tmp1*tmp1+w1*w1+90*tmp2*tmp2+y1*y1+10.1*(x1*x1+z1*z1)+19.8*x1*z1;
00389 return f;
00390 }
00391
00392
00393 double PowellFunction(const double * par) {
00394 gNCall++;
00395
00396 const double w = par[0];
00397 const double x = par[1];
00398 const double y = par[2];
00399 const double z = par[3];
00400
00401 const double tmp1 = w+10*x;
00402 const double tmp2 = y-z;
00403 const double tmp3 = x-2*y;
00404 const double tmp4 = w-z;
00405
00406 double f = tmp1*tmp1+5*tmp2*tmp2+tmp3*tmp3*tmp3*tmp3+10*tmp4*tmp4*tmp4*tmp4;
00407 return f;
00408 }
00409
00410 double SimpleQuadFunction(const double * par) {
00411 gNCall++;
00412 double x = par[0];
00413 double y = par[1];
00414 double f = x * x + 2 * y * y - x*y;
00415 std::cout << "Quadfunc " << gNCall << "\t" << x << " , " << y << " f = " << f << std::endl;
00416 return f;
00417 }
00418
00419 const double * TrueMinimum(const ROOT::Math::IMultiGenFunction & func) {
00420
00421 const RosenBrockFunction * fRB = dynamic_cast< const RosenBrockFunction *> (&func);
00422 if (fRB != 0) return fRB->TrueMinimum();
00423 const TrigoFletcherFunction * fTF = dynamic_cast< const TrigoFletcherFunction *> (&func);
00424 if (fTF != 0) return fTF->TrueMinimum();
00425
00426
00427 return 0;
00428 }
00429
00430 void printMinimum(const std::vector<double> & x) {
00431 std::cout << "Minimum X values\n";
00432 std::cout << "\t";
00433 int pr = std::cout.precision(12);
00434 unsigned int n = x.size();
00435 for (unsigned int i = 0; i < n; ++i) {
00436 std::cout << x[i];
00437 if ( i != n-1 ) std::cout << " , ";
00438 if ( i > 0 && i % 6 == 0 ) std::cout << "\n\t";
00439 }
00440 std::cout << std::endl;
00441 std::cout.precision(pr);
00442 }
00443
00444 int DoNewMinimization( const ROOT::Math::IMultiGenFunction & func, const double * x0, const double * s0, ROOT::Math::Minimizer * min, double &minval, double &edm, double * minx) {
00445
00446 int iret = 0;
00447
00448 min->SetMaxFunctionCalls(1000000);
00449 min->SetMaxIterations(100000);
00450 min->SetTolerance(gAbsTolerance);
00451 if (func.NDim() >= 10) {
00452 min->SetTolerance(0.01);
00453
00454 }
00455
00456 min->SetPrintLevel(gVerbose);
00457
00458 const ROOT::Math::IMultiGradFunction * gfunc = dynamic_cast<const ROOT::Math::IMultiGradFunction *>(&func);
00459 if (gfunc != 0 && useGradient)
00460 min->SetFunction(*gfunc);
00461 else
00462 min->SetFunction(func);
00463
00464 for (unsigned int i = 0; i < func.NDim(); ++i) {
00465 #ifdef DEBUG
00466 std::cout << "set variable " << i << " to value " << x0[i] << std::endl;
00467 #endif
00468 min->SetVariable(i,"x" + ROOT::Math::Util::ToString(i),x0[i], s0[i]);
00469 }
00470
00471 bool ret = min->Minimize();
00472 minval = min->MinValue();
00473 edm = min->Edm();
00474
00475 if (!ret) {
00476 return -1;
00477 }
00478
00479 const double * xmin = min->X();
00480
00481 bool ok = true;
00482 const double * trueMin = TrueMinimum(func);
00483 if (trueMin != 0) {
00484 for (unsigned int i = 0; i < func.NDim(); ++i)
00485 ok &= (std::fabs(xmin[i]-trueMin[i] ) < gAbsTolerance);
00486 }
00487
00488 if (!ok) iret = -2;
00489 int ncall_migrad = gNCall;
00490
00491
00492 if (minos) {
00493
00494 double el,eu;
00495 for (unsigned int i = 0; i < func.NDim(); ++i) {
00496 ret = min->GetMinosError(i,el,eu);
00497 std::cout << "ncalls " << gNCall << "\t";
00498 if (ret) std::cout << "MINOS error for " << i << " = " << el << " " << eu << std::endl;
00499 else std::cout << "MINOS failed for " << i << std::endl;
00500 }
00501 }
00502
00503 std::cout << "\nMigrad Ncalls:\t " << ncall_migrad << std::endl;
00504 std::cout << "Minos Ncalls :\t " << gNCall - ncall_migrad << std::endl;
00505
00506
00507
00508 std::copy(xmin,xmin+func.NDim(),minx);
00509 min->Clear();
00510
00511 return iret;
00512 }
00513
00514 int DoOldMinimization( FCN func, TVirtualFitter * min, double &minval, double &edm) {
00515
00516 int iret = 0;
00517
00518 assert(min != 0);
00519 min->SetFCN( func );
00520
00521 Double_t arglist[100];
00522 arglist[0] = gVerbose-1;
00523 min->ExecuteCommand("SET PRINT",arglist,1);
00524
00525 min->SetParameter(0,"x",-1.2,0.01,0,0);
00526 min->SetParameter(1,"y", 1.0,0.01,0,0);
00527
00528 arglist[0] = 0;
00529 min->ExecuteCommand("SET NOW",arglist,0);
00530 arglist[0] = 1000;
00531 arglist[1] = gAbsTolerance;
00532
00533 min->ExecuteCommand("MIGRAD",arglist,2);
00534
00535 if (minos) min->ExecuteCommand("MINOS",arglist,0);
00536
00537 Double_t parx,pary;
00538 Double_t we,al,bl;
00539 Char_t parName[32];
00540 min->GetParameter(0,parName,parx,we,al,bl);
00541 min->GetParameter(1,parName,pary,we,al,bl);
00542
00543 bool ok = ( TMath::Abs(parx-1.) < gAbsTolerance &&
00544 TMath::Abs(pary-1.) < gAbsTolerance );
00545
00546
00547 double errdef = 0;
00548 int nvpar, nparx;
00549 min->GetStats(minval,edm,errdef,nvpar,nparx);
00550 if (!ok) iret = -2;
00551
00552 min->Clear();
00553 return iret;
00554
00555 }
00556
00557
00558 int testNewMinimizer( const ROOT::Math::IMultiGenFunction & func, const double * x0, const double * s0, const std::string & minimizer, const std::string & algoType) {
00559
00560 std::cout << "\n************************************************************\n";
00561 std::cout << "\tTest new ROOT::Math::Minimizer\n";
00562 std::cout << "\tMinimizer is " << minimizer << " " << algoType << std::endl;
00563
00564 int iret = 0;
00565 double minval,edm = 0;
00566 std::vector<double> xmin(func.NDim() );
00567
00568 TStopwatch w;
00569 w.Start();
00570
00571 ROOT::Math::Minimizer * min = ROOT::Math::Factory::CreateMinimizer(minimizer, algoType);
00572 if (min == 0) {
00573 std::cout << "Error using minimizer " << minimizer << " " << algoType << std::endl;
00574 return -1;
00575 }
00576
00577
00578 for (int i = 0; i < gNmin; ++i) {
00579 gNCall = 0; gNCall2 = 0;
00580 iret |= DoNewMinimization(func, x0, s0, min,minval,edm,&xmin[0]);
00581 }
00582
00583 w.Stop();
00584 if (iret != 0) std::cout << "\n****** ERROR: Minimization FAILED ! \n";
00585 int pr = std::cout.precision(18);
00586 std::cout << "\nNCalls: \t" << gNCall << " , " << gNCall2
00587 << "\tMinValue: \t" << minval << "\tEdm: \t" << edm; std::cout.precision(pr);
00588 std::cout << "\nTime: \t" << w.RealTime() << " , " << w.CpuTime() << std::endl;
00589 printMinimum(xmin );
00590 std::cout << "\n************************************************************\n";
00591
00592 #ifdef CHECK_WITHMINUIT
00593
00594 if (minimizer == "GSL_BFGS") {
00595 std::cout << "DO Minuit2 from last point\n";
00596 gNCall = 0;
00597 iret |= DoNewMinimization(func, &xmin.front(), s0, "Minuit2","",minval,edm,&xmin[0]);
00598 int pr = std::cout.precision(18);
00599 std::cout << "\nNCalls: \t" << gNCall << "\tMinValue: \t" << minval << "\tEdm: \t" << edm; std::cout.precision(pr);
00600 std::cout << std::endl;
00601 }
00602 #endif
00603
00604 delete min;
00605
00606 return iret;
00607 }
00608
00609
00610 int testOldMinimizer( FCN func, const std::string & fitter, int n=25) {
00611
00612 std::cout << "\n************************************************************\n";
00613 std::cout << "\tTest using TVirtualFitter\n";
00614 std::cout << "\tFitter is " << fitter << std::endl;
00615
00616 int iret = 0;
00617 double minval,edm = 0;
00618
00619 TStopwatch w;
00620 w.Start();
00621
00622 TVirtualFitter::SetDefaultFitter(fitter.c_str());
00623
00624 TVirtualFitter *min = TVirtualFitter::Fitter(0,n);
00625
00626
00627
00628 for (int i = 0; i < gNmin; ++i) {
00629 gNCall = 0;
00630 iret |= DoOldMinimization(func, min,minval,edm);
00631 }
00632
00633 w.Stop();
00634 if (iret != 0) std::cout << "\n****** ERROR: Minimization FAILED ! \n";
00635 int pr = std::cout.precision(18);
00636 std::cout << "\nNCalls: \t" << gNCall << "\tMinValue: \t" << minval << "\tEdm: \t" << edm; std::cout.precision(pr);
00637 std::cout << "\nTime: \t" << w.RealTime() << " , " << w.CpuTime() << std::endl;
00638 std::cout << "\n************************************************************\n";
00639
00640 return iret;
00641 }
00642
00643 int testRosenBrock() {
00644
00645 int iret = 0;
00646
00647
00648 std::cout << "\n************************************************************\n";
00649 std::cout << "\tROSENBROCK function test\n\n";
00650
00651 double s0[2] = {0.01,0.01};
00652
00653
00654 #ifndef DEBUG
00655 gNmin = 2000;
00656 #endif
00657
00658 gNmin = 1;
00659
00660
00661 RosenBrockFunction fRB;
00662 double xRB[2] = { -1.,1.2};
00663
00664 iret |= testNewMinimizer(fRB,xRB,s0,"Minuit","");
00665 iret |= testNewMinimizer(fRB,xRB,s0,"Minuit2","");
00666
00667
00668
00669
00670
00671
00672
00673
00674
00675
00676 return iret;
00677 }
00678
00679
00680 int testChebyQuad() {
00681
00682 int iret = 0;
00683
00684
00685
00686 const int n = 8;
00687 ChebyQuadFunction func(n);
00688
00689 #ifndef DEBUG
00690 gNmin = std::max(1, int(20000/n/n) );
00691 #endif
00692 gNmin = 1;
00693
00694
00695 double s0[n];
00696 double x0[n];
00697 func.StartPoints(x0,s0);
00698
00699 std::cout << "\n************************************************************\n";
00700 std::cout << "\tCHEBYQUAD function test , n = " << n << "\n\n";
00701
00702
00703
00704
00705
00706
00707 double x1[100] = { 0.00712780070646 , 0.0123441993113 , 0.0195428378255 , 0.0283679084192 , 0.0385291131289 , 0.0492202424892 , 0.0591277976178 ,
00708 0.0689433195252 , 0.0791293590525 , 0.088794974369 , 0.0988949579193 , 0.108607151294 , 0.118571075831 ,
00709 0.128605446238 , 0.137918291068 , 0.149177761352 , 0.156665324587 , 0.170851061982 , 0.174688134016 ,
00710 0.192838903364 , 0.193078270803 , 0.209255377225 , 0.217740096779 , 0.225888518345 , 0.241031047421 ,
00711 0.244253844041 , 0.257830449676 , 0.269467652526 , 0.274286498012 , 0.288877029988 , 0.297549406597 ,
00712 0.304950954529 , 0.319230811642 , 0.326387092206 , 0.335229058731 , 0.349178359226 , 0.355905988048 ,
00713 0.365197862755 , 0.379068092603 , 0.385826036925 , 0.394978252826 , 0.408974425717 , 0.415968185065 ,
00714 0.424621041584 , 0.438837361714 , 0.446214149031 , 0.454242324351 , 0.468614308013 , 0.476506553416 ,
00715 0.483916944941 , 0.498229247409 , 0.506794629616 , 0.513736742474 , 0.527712475478 , 0.537073277673 ,
00716 0.543731917673 , 0.557187513963 , 0.567346279639 , 0.57379846397 , 0.586691058785 , 0.597561941009 ,
00717 0.60382873461 , 0.616316037506 , 0.627719652101 , 0.633760038662 , 0.646175283836 , 0.657809344891 ,
00718 0.663569004722 , 0.676314563639 , 0.687674566849 , 0.69332205923 , 0.706839545953 , 0.716907408637 ,
00719 0.723407327715 , 0.738019389561 , 0.744806584048 , 0.754657613362 , 0.769181875619 , 0.772250323489 ,
00720 0.787104833193 , 0.795856360905 , 0.804099304478 , 0.82142178741 , 0.819589601284 , 0.839024540481 ,
00721 0.842457233039 , 0.857393475964 , 0.86408033345 , 0.876329840525 , 0.884867318008 , 0.895744532071 ,
00722 0.905113958163 , 0.915445338697 , 0.925148068352 , 0.935344457785 , 0.945127838313 , 0.955272197168 ,
00723 0.965687518559 , 0.975129521484 , 0.982662007764 };
00724
00725 std::cout << "FUNC " << func(x1) << std::endl;
00726
00727
00728 iret |= testNewMinimizer(func,x0,s0, "Minuit","");
00729 iret |= testNewMinimizer(func,x0,s0, "Minuit2","");
00730
00731
00732
00733
00734
00735 return iret;
00736 }
00737
00738
00739 int testTrigoFletcher() {
00740
00741 int iret = 0;
00742
00743
00744
00745 #ifndef DEBUG
00746 gNmin = 2;
00747 #endif
00748 gNmin = 1;
00749
00750 const int nT = 50;
00751 TrigoFletcherFunction fTrigo(nT);
00752 double sTrigo[nT];
00753 double xTrigo[nT];
00754 fTrigo.StartPoints(xTrigo,sTrigo);
00755
00756 std::cout << "\n************************************************************\n";
00757 std::cout << "\tTRIGONOMETRIC Fletcher function test , n = " << nT << "\n\n";
00758
00759
00760 iret |= testNewMinimizer(fTrigo,xTrigo,sTrigo,"Minuit2","");
00761 iret |= testNewMinimizer(fTrigo,xTrigo,sTrigo,"Minuit","");
00762
00763
00764
00765
00766
00767
00768 return iret;
00769 }
00770
00771
00772 int testWood() {
00773
00774 int iret = 0;
00775
00776
00777
00778
00779
00780
00781 #ifndef DEBUG
00782 gNmin = 1000;
00783 #endif
00784 gNmin = 1;
00785
00786 ROOT::Math::Functor f(&WoodFunction,4);
00787
00788 double x0[4] = { -3, -1, -3, -1 };
00789 double s0[4] = { 0.1, 0.1, 0.1, 0.1};
00790
00791 std::cout << "\n************************************************************\n";
00792 std::cout << "\tWOOD 4 function test \n\n";
00793
00794
00795 iret |= testNewMinimizer(f, x0, s0,"Minuit","");
00796 iret |= testNewMinimizer(f, x0, s0,"Minuit2","");
00797
00798
00799 return iret;
00800 }
00801
00802 int testPowell() {
00803
00804 int iret = 0;
00805
00806
00807
00808
00809 #ifndef DEBUG
00810 gNmin = 1000;
00811 #endif
00812 gNmin = 1;
00813
00814 ROOT::Math::Functor f(&PowellFunction,4);
00815
00816 double x0[4] = { -3, -1, 0, 1 };
00817 double s0[4] = { 0.1, 0.1, 0.1, 0.1};
00818
00819 std::cout << "\n************************************************************\n";
00820 std::cout << "\tPOWELL function test \n\n";
00821
00822
00823 iret |= testNewMinimizer(f, x0, s0,"Minuit","");
00824 iret |= testNewMinimizer(f, x0, s0,"Minuit2","");
00825
00826
00827 return iret;
00828 }
00829
00830
00831 int testQuadFunc() {
00832
00833 int iret = 0;
00834
00835
00836
00837
00838 #ifndef DEBUG
00839 gNmin = 1000;
00840 #endif
00841 gNmin = 1;
00842
00843 ROOT::Math::Functor f(&SimpleQuadFunction,2);
00844
00845 double x0[4] = { -3, -3 };
00846 double s0[4] = { 0.1, 0.1};
00847
00848 std::cout << "\n************************************************************\n";
00849 std::cout << "\tSIMPLE QUAD function test \n\n";
00850
00851
00852 iret |= testNewMinimizer(f, x0, s0,"Minuit","");
00853 iret |= testNewMinimizer(f, x0, s0,"Minuit2","");
00854
00855
00856 return iret;
00857 }
00858
00859
00860 int main() {
00861
00862 int iret = 0;
00863
00864 #ifdef DEBUG
00865 gVerbose = 4;
00866 gNmin = 1;
00867 #endif
00868
00869
00870 iret |= testRosenBrock();
00871 iret |= testChebyQuad();
00872 iret |= testTrigoFletcher();
00873
00874 iret |= testWood();
00875 iret |= testPowell();
00876
00877 iret |= testQuadFunc();
00878
00879
00880
00881 if (iret != 0)
00882 std::cerr << "testMinim :\t FAILED " << std::endl;
00883 else
00884 std::cerr << "testMinim :\t OK " << std::endl;
00885 return iret;
00886
00887 }