Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TMinuitMinimizer.cxx
Go to the documentation of this file.
1// @(#)root/minuit:$Id$
2// Author: L. Moneta Wed Oct 25 16:28:55 2006
3
4/**********************************************************************
5 * *
6 * Copyright (c) 2006 LCG ROOT Math Team, CERN/PH-SFT *
7 * *
8 * *
9 **********************************************************************/
10
11// Implementation file for class TMinuitMinimizer
12
13#include "TMinuitMinimizer.h"
14#include "Math/IFunction.h"
16
17#include "TMinuit.h"
18#include "TROOT.h"
19#include "TList.h"
20
21#include "TGraph.h" // needed for scan
22#include "TError.h"
23
24#include "TMatrixDSym.h" // needed for inverting the matrix
25
26#include "ThreadLocalStorage.h"
27
28#include <cassert>
29#include <algorithm>
30#include <functional>
31#include <cmath>
32
33////////////////////////////////////////////////////////////////////////////////
34/// \class TMinuitMinimizer
35/// \see Minuit2 for a newer version of this class
36/// TMinuitMinimizer class implementing the ROOT::Math::Minimizer interface
37/// using TMinuit. This class is normally instantiated using the plug-in manager
38/// (plug-in with name Minuit or TMinuit).
39/// In addition the user can choose the minimizer algorithm:
40/// Migrad (the default one), Simplex, or Minimize (combined Migrad + Simplex)
41////////////////////////////////////////////////////////////////////////////////
42
43// initialize the static instances
44
45// Implement a thread local static member
47 TTHREAD_TLS(ROOT::Math::IMultiGenFunction *) fgFunc = nullptr;
48 return fgFunc;
49}
51bool TMinuitMinimizer::fgUsed = false;
52bool TMinuitMinimizer::fgUseStaticMinuit = true; // default case use static Minuit instance
53
55
56
58 fUsed(false),
59 fMinosRun(false),
60 fDim(ndim),
61 fType(type),
62 fMinuit(nullptr)
63{
64 // Constructor for TMinuitMinimier class via an enumeration specifying the minimization
65 // algorithm type. Supported types are : kMigrad, kSimplex, kCombined (a combined
66 // Migrad + Simplex minimization) and kMigradImproved (a Migrad mininimization folloed by an
67 // improved search for global minima). The default type is Migrad (kMigrad).
68
69 // initialize if npar is given
70 if (fDim > 0) InitTMinuit(fDim);
71}
72
73TMinuitMinimizer::TMinuitMinimizer(const char * type, unsigned int ndim ) :
74 fUsed(false),
75 fMinosRun(false),
76 fDim(ndim),
77 fMinuit(nullptr)
78{
79 // constructor from a char * for the algorithm type, used by the plug-in manager
80 // The names supported (case unsensitive) are:
81 // Migrad (default), Simplex, Minimize (for the combined Migrad+ Simplex) and Migrad_imp
82
83 // select type from the string
84 std::string algoname(type);
85 std::transform(algoname.begin(), algoname.end(), algoname.begin(), (int(*)(int)) tolower );
86
88 if (algoname == "simplex") algoType = ROOT::Minuit::kSimplex;
89 if (algoname == "minimize" ) algoType = ROOT::Minuit::kCombined;
90 if (algoname == "migradimproved" ) algoType = ROOT::Minuit::kMigradImproved;
91 if (algoname == "scan" ) algoType = ROOT::Minuit::kScan;
92 if (algoname == "seek" ) algoType = ROOT::Minuit::kSeek;
93
94 fType = algoType;
95
96 // initialize if npar is given
97 if (fDim > 0) InitTMinuit(fDim);
98
99}
100
102{
103 // Destructor implementation.
104 if (fMinuit && !fgUseStaticMinuit) {
105 delete fMinuit;
106 fgMinuit = nullptr;
107 }
108}
109
111 // static method to control usage of global TMinuit instance
112 bool prev = fgUseStaticMinuit;
114 return prev;
115}
116
118
119 // when called a second time check dimension - create only if needed
120 // initialize the minuit instance - recreating a new one if needed
121 if (fMinuit ==nullptr || dim > fMinuit->fMaxpar) {
122
123 // case not using the global instance - recreate it all the time
124 if (fgUseStaticMinuit) {
125
126 // re-use gMinuit as static instance of TMinuit
127 // which can be accessed by the user after minimization
128 // check if fgMinuit is different than gMinuit
129 // case 1: fgMinuit not zero but fgMinuit has been deleted (not in gROOT): set to zero
130 // case 2: fgMinuit not zero and exists in global list : set fgMinuit to gMinuit
131 // case 3: fgMinuit zero - and gMinuit not zero: create a new instance locally to avoid conflict
132 if (fgMinuit != gMinuit) {
133 // if object exists in gROOT remove it to avoid a memory leak
134 if (fgMinuit ) {
135 if (gROOT->GetListOfSpecials()->FindObject(fgMinuit) == nullptr) {
136 // case 1: object does not exists in gROOT - means it has been deleted
137 fgMinuit = nullptr;
138 }
139 else {
140 // case 2: object exists - but gMinuit points to something else
141 // restore gMinuit to the one used before by TMinuitMinimizer
143 }
144 }
145 else {
146 // case 3: avoid reusing existing one - maintain fgMinuit to zero
147 // otherwise we will get a double delete if user deletes externally gMinuit
148 // in this case we will loose gMinuit instance
149// fgMinuit = gMinuit;
150// fgUsed = true; // need to reset in case other gMinuit instance is later used
151 }
152 }
153
154 // check if need to create a new TMinuit instance
155 if (fgMinuit == nullptr) {
156 fgUsed = false;
157 fgMinuit = new TMinuit(dim);
158 }
159 else if (fgMinuit->GetNumPars() != int(dim) ) {
160 delete fgMinuit;
161 fgUsed = false;
162 fgMinuit = new TMinuit(dim);
163 }
164
166 }
167
168 else {
169 // re- create all the time a new instance of TMinuit (fgUseStaticMinuit is false)
170 if (fMinuit) delete fMinuit;
171 fMinuit = new TMinuit(dim);
173 fgUsed = false;
174 }
175
176 } // endif fMinuit ==0 || dim > fMaxpar
177
178 fDim = dim;
179
181
182 // set print level in TMinuit
183 double arglist[1];
184 int ierr= 0;
185 // TMinuit level is shift by 1 -1 means 0;
186 arglist[0] = PrintLevel() - 1;
187 fMinuit->mnexcm("SET PRINT",arglist,1,ierr);
189}
190
191
193 // Set the objective function to be minimized, by passing a function object implement the
194 // basic multi-dim Function interface. In this case the derivatives will be
195 // calculated by Minuit
196 // Here a TMinuit instance is created since only at this point we know the number of parameters
197
198 const bool hasGrad = func.HasGradient();
199
200 fDim = func.NDim();
201
202 // create TMinuit if needed
204
205 // assign to the static pointer (NO Thread safety here)
206 GetGlobalFuncPtr() = const_cast<ROOT::Math::IMultiGenFunction *>(&func);
208
209 double arglist[1];
210 int ierr = 0;
211
212 if(hasGrad) {
213 // set gradient
214 // by default do not check gradient calculation
215 // it cannot be done here, check can be done only after having defined the parameters
216 arglist[0] = 1;
217 fMinuit->mnexcm("SET GRAD",arglist,1,ierr);
218 } else {
219 // switch off gradient calculations
220 fMinuit->mnexcm("SET NOGrad",arglist,0,ierr);
221 }
222}
223
224void TMinuitMinimizer::Fcn( int &, double * , double & f, double * x , int /* iflag */) {
225 // implementation of FCN static function used internally by TMinuit.
226 // Adapt IMultiGenFunction interface to TMinuit FCN static function
227 f = GetGlobalFuncPtr()->operator()(x);
228}
229
230void TMinuitMinimizer::FcnGrad( int &, double * g, double & f, double * x , int iflag ) {
231 // implementation of FCN static function used internally by TMinuit.
232 // Adapt IMultiGradFunction interface to TMinuit FCN static function in the case of user
233 // provided gradient.
235
236 assert(gFunc != nullptr);
237 f = gFunc->operator()(x);
238
239 // calculates also derivatives
240 if (iflag == 2) gFunc->Gradient(x,g);
241}
242
243bool TMinuitMinimizer::SetVariable(unsigned int ivar, const std::string & name, double val, double step) {
244 // set a free variable.
245 if (!CheckMinuitInstance()) return false;
246
248
249 // clear after minimization when setting params
250 if (fUsed) DoClear();
251
252 // check if parameter was defined and in case it was fixed, release it
254
255 int iret = fMinuit->DefineParameter(ivar , name.c_str(), val, step, 0., 0. );
256 return (iret == 0);
257}
258
259bool TMinuitMinimizer::SetLimitedVariable(unsigned int ivar, const std::string & name, double val, double step, double lower, double upper) {
260 // set a limited variable.
261 if (!CheckMinuitInstance()) return false;
262
264
265 // clear after minimization when setting params
266 if (fUsed) DoClear();
267
268 // check if parameter was defined and in case it was fixed, release it
270
271 int iret = fMinuit->DefineParameter(ivar, name.c_str(), val, step, lower, upper );
272 return (iret == 0);
273}
274
275bool TMinuitMinimizer::SetLowerLimitedVariable(unsigned int ivar , const std::string & name , double val , double step , double lower ) {
276 // set a lower limited variable
277 // since is not supported in TMinuit , just use a artificial large value
278 Warning("TMinuitMinimizer::SetLowerLimitedVariable","not implemented - use as upper limit 1.E7 instead of +inf");
279 return SetLimitedVariable(ivar, name, val , step, lower, lower+ 1.E7); // use 1.E7 which will make TMinuit happy
280}
281
282bool TMinuitMinimizer::SetUpperLimitedVariable(unsigned int ivar , const std::string & name , double val , double step , double upper ) {
283 // set a upper limited variable
284 // since is not supported in TMinuit , just use a artificial large negative value
285 Warning("TMinuitMinimizer::SetUpperLimitedVariable","not implemented - use as lower limit -1.E7 instead of -inf");
286 return SetLimitedVariable(ivar, name, val , step, upper -1.E7, upper);
287}
288
289
291 // check instance of fMinuit
292 if (fMinuit == nullptr) {
293 Error("TMinuitMinimizer::CheckMinuitInstance","Invalid TMinuit pointer. Need to call first SetFunction");
294 return false;
295 }
296 return true;
297}
298
299bool TMinuitMinimizer::CheckVarIndex(unsigned int ivar) const {
300 // check index of Variable (assume fMinuit exists)
301 if ((int) ivar >= fMinuit->fNu ) {
302 Error("TMinuitMinimizer::CheckVarIndex","Invalid parameter index");
303 return false;
304 }
305 return true;
306}
307
308
309bool TMinuitMinimizer::SetFixedVariable(unsigned int ivar, const std::string & name, double val) {
310 // set a fixed variable.
311 if (!CheckMinuitInstance()) return false;
312
313 // clear after minimization when setting params
315
316 // clear after minimization when setting params
317 if (fUsed) DoClear();
318
319 // put an arbitrary step (0.1*abs(value) otherwise TMinuit consider the parameter as constant
320 // constant parameters are treated differently (they are ignored inside TMinuit and not considered in the
321 // total list of parameters)
322 double step = ( val != 0) ? 0.1 * std::abs(val) : 0.1;
323 int iret = fMinuit->DefineParameter(ivar, name.c_str(), val, step, 0., 0. );
324 if (iret == 0) iret = fMinuit->FixParameter(ivar);
325 return (iret == 0);
326}
327
328bool TMinuitMinimizer::SetVariableValue(unsigned int ivar, double val) {
329 // set the value of an existing variable
330 // parameter must exist or return false
331 if (!CheckMinuitInstance()) return false;
332
333 double arglist[2];
334 int ierr = 0;
335
336 arglist[0] = ivar+1; // TMinuit starts from 1
337 arglist[1] = val;
338 fMinuit->mnexcm("SET PAR",arglist,2,ierr);
339 return (ierr==0);
340}
341
342bool TMinuitMinimizer::SetVariableStepSize(unsigned int ivar, double step) {
343 // set the step-size of an existing variable
344 // parameter must exist or return false
345 if (!CheckMinuitInstance()) return false;
346 // need to re-implement re-calling mnparm
347 // get first current parameter values and limits
348 double curval,err, lowlim, uplim;
349 int iuint; // internal index
351 fMinuit->mnpout(ivar, name, curval, err, lowlim, uplim,iuint);
352 if (iuint == -1) return false;
353 int iret = fMinuit->DefineParameter(ivar, name, curval, step, lowlim, uplim );
354 return (iret == 0);
355
356}
357
358bool TMinuitMinimizer::SetVariableLowerLimit(unsigned int ivar, double lower ) {
359 // set the limits of an existing variable
360 // parameter must exist or return false
361 Warning("TMinuitMinimizer::SetVariableLowerLimit","not implemented - use as upper limit 1.E30 instead of +inf");
362 return SetVariableLimits(ivar, lower, 1.E30);
363}
364bool TMinuitMinimizer::SetVariableUpperLimit(unsigned int ivar, double upper ) {
365 // set the limits of an existing variable
366 // parameter must exist or return false
367 Warning("TMinuitMinimizer::SetVariableUpperLimit","not implemented - - use as lower limit -1.E30 instead of +inf");
368 return SetVariableLimits(ivar, -1.E30, upper);
369}
370
371bool TMinuitMinimizer::SetVariableLimits(unsigned int ivar, double lower, double upper) {
372 // set the limits of an existing variable
373 // parameter must exist or return false
374
375 if (!CheckMinuitInstance()) return false;
376 // need to re-implement re-calling mnparm
377 // get first current parameter values and limits
378 double curval,err, lowlim, uplim;
379 int iuint; // internal index
381 fMinuit->mnpout(ivar, name, curval, err, lowlim, uplim,iuint);
382 if (iuint == -1) return false;
383 int iret = fMinuit->DefineParameter(ivar, name, curval, err, lower, upper );
384 return (iret == 0);
385
386}
387
388bool TMinuitMinimizer::FixVariable(unsigned int ivar) {
389 // Fix an existing variable
390 if (!CheckMinuitInstance()) return false;
391 if (!CheckVarIndex(ivar)) return false;
392 int iret = fMinuit->FixParameter(ivar);
393 return (iret == 0);
394}
395
396bool TMinuitMinimizer::ReleaseVariable(unsigned int ivar) {
397 // Fix an existing variable
398 if (!CheckMinuitInstance()) return false;
399 if (!CheckVarIndex(ivar)) return false;
400 int iret = fMinuit->Release(ivar);
401 return (iret == 0);
402}
403
404bool TMinuitMinimizer::IsFixedVariable(unsigned int ivar) const {
405 // query if variable is fixed
406 if (!CheckMinuitInstance()) return false;
407 if (!CheckVarIndex(ivar)) return false;
408 return (fMinuit->fNiofex[ivar] == 0 );
409}
410
412 // retrieve variable settings (all set info on the variable)
413 if (!CheckMinuitInstance()) return false;
414 if (!CheckVarIndex(ivar)) return false;
415 double curval,err, lowlim, uplim;
416 int iuint; // internal index
418 fMinuit->mnpout(ivar, name, curval, err, lowlim, uplim,iuint);
419 if (iuint == -1) return false;
420 var.Set(name.Data(), curval, err, lowlim, uplim);
421 if (IsFixedVariable(ivar)) var.Fix();
422 return true;
423}
424
425
426
427std::string TMinuitMinimizer::VariableName(unsigned int ivar) const {
428 // return the variable name
429 if (!CheckMinuitInstance()) return std::string();
430 if (!CheckVarIndex(ivar)) return std::string();
431 return fMinuit->fCpnam[ivar].Data();
432}
433
434int TMinuitMinimizer::VariableIndex(const std::string & ) const {
435 // return variable index
436 Error("TMinuitMinimizer::VariableIndex"," find index of a variable from its name is not implemented in TMinuit");
437 return -1;
438}
439
441 // perform the minimization using the algorithm chosen previously by the user
442 // By default Migrad is used.
443 // Return true if the found minimum is valid and update internal cached values of
444 // minimum values, errors and covariance matrix.
445 // Status of minimizer is set to:
446 // migradResult + 10*minosResult + 100*hesseResult + 1000*improveResult
447
448
449 if (fMinuit == nullptr) {
450 Error("TMinuitMinimizer::Minimize","invalid TMinuit pointer. Need to call first SetFunction and SetVariable");
451 return false;
452 }
453
454
455 // total number of parameter defined in Minuit is fNu
456 if (fMinuit->fNu < static_cast<int>(fDim) ) {
457 Error("TMinuitMinimizer::Minimize","The total number of defined parameters is different than the function dimension, npar = %d, dim = %d",fMinuit->fNu, fDim);
458 return false;
459 }
460
461 int printlevel = PrintLevel();
462
463 // total number of free parameter is 0
464 if (fMinuit->fNpar <= 0) {
465 // retrieve parameters values from TMinuit
467 fMinuit->fAmin = (*GetGlobalFuncPtr())(&fParams.front());
468 if (printlevel > 0) Info("TMinuitMinimizer::Minimize","There are no free parameter - just compute the function value");
469 return true;
470 }
471
472
473 double arglist[10];
474 int ierr = 0;
475
476
477 // set error and print level
478 arglist[0] = ErrorDef();
479 fMinuit->mnexcm("SET Err",arglist,1,ierr);
480
481 arglist[0] = printlevel - 1;
482 fMinuit->mnexcm("SET PRINT",arglist,1,ierr);
483
484 // suppress warning in case Printlevel() == 0
485 if (printlevel == 0) fMinuit->mnexcm("SET NOW",arglist,0,ierr);
486
487 // set precision if needed
488 if (Precision() > 0) {
489 arglist[0] = Precision();
490 fMinuit->mnexcm("SET EPS",arglist,1,ierr);
491 }
492
493 // set strategy
494 int strategy = Strategy();
495 if (strategy >=0 && strategy <=2 ) {
496 arglist[0] = strategy;
497 fMinuit->mnexcm("SET STR",arglist,1,ierr);
498 }
499
500 arglist[0] = MaxFunctionCalls();
501 arglist[1] = Tolerance();
502
503 int nargs = 2;
504
505 switch (fType){
507 // case of Migrad
508 fMinuit->mnexcm("MIGRAD",arglist,nargs,ierr);
509 break;
511 // case of combined (Migrad+ simplex)
512 fMinuit->mnexcm("MINIMIZE",arglist,nargs,ierr);
513 break;
515 // case of Simlex
516 fMinuit->mnexcm("SIMPLEX",arglist,nargs,ierr);
517 break;
519 // case of Scan (scan all parameters with default values)
520 nargs = 0;
521 fMinuit->mnexcm("SCAN",arglist,nargs,ierr);
522 break;
524 // case of Seek (random find minimum in a hypercube around current parameter values
525 // use Tolerance as measures for standard deviation (if < 1) used default value in Minuit ( supposed to be 3)
526 nargs = 1;
527 if (arglist[1] >= 1.) nargs = 2;
528 fMinuit->mnexcm("SEEK",arglist,nargs,ierr);
529 break;
530 default:
531 // default: use Migrad
532 fMinuit->mnexcm("MIGRAD",arglist,nargs,ierr);
533
534 }
535
536 fgUsed = true;
537 fUsed = true;
538
539 fStatus = ierr;
540 int minErrStatus = ierr;
541
542 if (printlevel>2) Info("TMinuitMinimizer::Minimize","Finished to run MIGRAD - status %d",ierr);
543
544 // run improved if needed
545 if (ierr == 0 && fType == ROOT::Minuit::kMigradImproved) {
546 fMinuit->mnexcm("IMPROVE",arglist,1,ierr);
547 fStatus += 1000*ierr;
548 if (printlevel>2) Info("TMinuitMinimizer::Minimize","Finished to run IMPROVE - status %d",ierr);
549 }
550
551
552 // check if Hesse needs to be run
553 // Migrad runs inside it automatically for strategy >=1. Do also
554 // in case improve or other minimizers are used
555 // run Hesse also in case cov matrix has not been computed or has been made pos-def
556 // or a valid error analysis is requested (when IsValidError() == true)
557 if (minErrStatus == 0 && (IsValidError() || ( strategy >=1 && CovMatrixStatus() < 3) ) ) {
558 fMinuit->mnexcm("HESSE",arglist,1,ierr);
559 fStatus += 100*ierr;
560 // here ierr is zero when Hesse fails CovMatrixStatus = 0 or 1 (2 is when was made posdef)
561 if (ierr == 0 && CovMatrixStatus() < 2){
562 fStatus += 100*(CovMatrixStatus()+1);
563 }
564 // should check here cov matrix status??? (if <3 flag error ?)
565 if (printlevel>2) Info("TMinuitMinimizer::Minimize","Finished to run HESSE - status %d",ierr);
566 }
567
568 // retrieve parameters and errors from TMinuit
570
571 if (minErrStatus == 0) {
572
573 // store global min results (only if minimization is OK)
574 // ignore cases when Hesse or IMprove return error different than zero
576
577 // need to re-run Minos again if requested
578 fMinosRun = false;
579
580 return true;
581
582 }
583 return false;
584
585}
586
588 // retrieve from TMinuit minimum parameter values
589 // and errors
590
591 assert(fMinuit != nullptr);
592
593 // get parameter values
594 if (fParams.size() != fDim) fParams.resize( fDim);
595 if (fErrors.size() != fDim) fErrors.resize( fDim);
596 for (unsigned int i = 0; i < fDim; ++i) {
597 fMinuit->GetParameter( i, fParams[i], fErrors[i]);
598 }
599}
600
602 // get covariance error matrix from TMinuit
603 // when some parameters are fixed filled the corresponding rows and column with zero's
604
605 assert(fMinuit != nullptr);
606
607 unsigned int nfree = NFree();
608
609 unsigned int ndim2 = fDim*fDim;
610 if (fCovar.size() != ndim2 ) fCovar.resize(fDim*fDim);
611 if (nfree >= fDim) { // no fixed parameters
612 fMinuit->mnemat(&fCovar.front(), fDim);
613 }
614 else {
615 // case of fixed params need to take care
616 std::vector<double> tmpMat(nfree*nfree);
617 fMinuit->mnemat(&tmpMat.front(), nfree);
618
619 unsigned int l = 0;
620 for (unsigned int i = 0; i < fDim; ++i) {
621
622 if ( fMinuit->fNiofex[i] > 0 ) { // not fixed ?
623 unsigned int m = 0;
624 for (unsigned int j = 0; j <= i; ++j) {
625 if ( fMinuit->fNiofex[j] > 0 ) { //not fixed
626 fCovar[i*fDim + j] = tmpMat[l*nfree + m];
627 fCovar[j*fDim + i] = fCovar[i*fDim + j];
628 m++;
629 }
630 }
631 l++;
632 }
633 }
634
635 }
636}
637
638unsigned int TMinuitMinimizer::NCalls() const {
639 // return total number of function calls
640 if (fMinuit == nullptr) return 0;
641 return fMinuit->fNfcn;
642}
643
645 // return minimum function value
646
647 // use part of code from mnstat
648 if (!fMinuit) return 0;
649 double minval = fMinuit->fAmin;
650 if (minval == fMinuit->fUndefi) return 0;
651 return minval;
652}
653
654double TMinuitMinimizer::Edm() const {
655 // return expected distance from the minimum
656
657 // use part of code from mnstat
658 if (!fMinuit) return -1;
659 if (fMinuit->fAmin == fMinuit->fUndefi || fMinuit->fEDM == fMinuit->fBigedm) return fMinuit->fUp;
660 return fMinuit->fEDM;
661}
662
663unsigned int TMinuitMinimizer::NFree() const {
664 // return number of free parameters
665 if (!fMinuit) return 0;
666 if (fMinuit->fNpar < 0) return 0;
667 return fMinuit->fNpar;
668}
669
670bool TMinuitMinimizer::GetCovMatrix(double * cov) const {
671 // get covariance matrix
672 int covStatus = CovMatrixStatus();
673 if ( fCovar.size() != fDim*fDim || covStatus < 2) {
674 Error("TMinuitMinimizer::GetHessianMatrix","Hessian matrix has not been computed - status %d",covStatus);
675 return false;
676 }
677 std::copy(fCovar.begin(), fCovar.end(), cov);
678 TMatrixDSym cmat(fDim,cov);
679 return true;
680}
681
682bool TMinuitMinimizer::GetHessianMatrix(double * hes) const {
683 // get Hessian - inverse of covariance matrix
684 // just invert it
685 // but need to get the compact form to avoid the zero for the fixed parameters
686 int covStatus = CovMatrixStatus();
687 if ( fCovar.size() != fDim*fDim || covStatus < 2) {
688 Error("TMinuitMinimizer::GetHessianMatrix","Hessian matrix has not been computed - status %d",covStatus);
689 return false;
690 }
691 // case of fixed params need to take care
692 unsigned int nfree = NFree();
693 TMatrixDSym mat(nfree);
694 fMinuit->mnemat(mat.GetMatrixArray(), nfree);
695 // invert the matrix
696 // probably need to check if failed. In that case inverse is equal to original
697 mat.Invert();
698
699 unsigned int l = 0;
700 for (unsigned int i = 0; i < fDim; ++i) {
701 if ( fMinuit->fNiofex[i] > 0 ) { // not fixed ?
702 unsigned int m = 0;
703 for (unsigned int j = 0; j <= i; ++j) {
704 if ( fMinuit->fNiofex[j] > 0 ) { //not fixed
705 hes[i*fDim + j] = mat(l,m);
706 hes[j*fDim + i] = hes[i*fDim + j];
707 m++;
708 }
709 }
710 l++;
711 }
712 }
713 return true;
714}
715// if ( fCovar.size() != fDim*fDim ) return false;
716// TMatrixDSym mat(fDim, &fCovar.front() );
717// std::copy(mat.GetMatrixArray(), mat.GetMatrixArray()+ mat.GetNoElements(), hes);
718// return true;
719// }
720
722 // return status of covariance matrix
723 // status: 0= not calculated at all
724 // 1= approximation only, not accurate
725 // 2= full matrix, but forced positive-definite
726 // 3= full accurate covariance matrix
727
728 // use part of code from mnstat
729 if (!fMinuit) return 0;
730 if (fMinuit->fAmin == fMinuit->fUndefi) return 0;
731 return fMinuit->fISW[1];
732}
733
734double TMinuitMinimizer::GlobalCC(unsigned int i) const {
735 // global correlation coefficient for parameter i
736 if (!fMinuit) return 0;
737 if (!fMinuit->fGlobcc) return 0;
738 if (int(i) >= fMinuit->fNu) return 0;
739 // get internal number in Minuit
740 int iin = fMinuit->fNiofex[i];
741 // index in TMinuit starts from 1
742 if (iin < 1) return 0;
743 return fMinuit->fGlobcc[iin-1];
744}
745
746bool TMinuitMinimizer::GetMinosError(unsigned int i, double & errLow, double & errUp, int ) {
747 // Perform Minos analysis for the given parameter i
748
749 if (fMinuit == nullptr) {
750 Error("TMinuitMinimizer::GetMinosError","invalid TMinuit pointer. Need to call first SetFunction and SetVariable");
751 return false;
752 }
753
754 // check if parameter is fixed
755 if (fMinuit->fNiofex[i] == 0 ) {
756 if (PrintLevel() > 0) Info("TMinuitMinimizer::GetMinosError","Parameter %s is fixed. There are no Minos error to calculate. Ignored.",VariableName(i).c_str());
757 errLow = 0; errUp = 0;
758 return true;
759 }
760
761 double arglist[2];
762 int ierr = 0;
763
764
765 // set error, print level, precision and strategy if they have changed
766 if (fMinuit->fUp != ErrorDef() ) {
767 arglist[0] = ErrorDef();
768 fMinuit->mnexcm("SET Err",arglist,1,ierr);
769 }
770
771 if (fMinuit->fISW[4] != (PrintLevel()-1) ) {
772 arglist[0] = PrintLevel()-1;
773 fMinuit->mnexcm("SET PRINT",arglist,1,ierr);
774 // suppress warning in case Printlevel() == 0
775 if (PrintLevel() == 0) fMinuit->mnexcm("SET NOW",arglist,0,ierr);
776 }
777 if (fMinuit->fIstrat != Strategy() ) {
778 arglist[0] = Strategy();
779 fMinuit->mnexcm("SET STR",arglist,1,ierr);
780 }
781
782 if (Precision() > 0 && fMinuit->fEpsma2 != Precision() ) {
783 arglist[0] = Precision();
784 fMinuit->mnexcm("SET EPS",arglist,1,ierr);
785 }
786
787
788 // syntax of MINOS command is: MINOS [maxcalls] [parno]
789 // if parno = 0 all parameters are done
790 arglist[0] = MaxFunctionCalls();
791 arglist[1] = i+1; // par number starts from 1 in TMInuit
792
793 int nargs = 2;
794 fMinuit->mnexcm("MINOS",arglist,nargs,ierr);
795 bool isValid = (ierr == 0);
796 // check also the status from fCstatu
797 if (isValid && fMinuit->fCstatu != "SUCCESSFUL") {
798 if (fMinuit->fCstatu == "FAILURE" ) {
799 // in this case MINOS failed on all parameters, so it is not valid !
800 ierr = 5;
801 isValid = false;
802 }
803 if (fMinuit->fCstatu == "PROBLEMS") ierr = 6;
804 ierr = 7; // this should be the case UNCHANGED
805 }
806
807 fStatus += 10*ierr;
808 fMinosStatus = ierr;
809
810 fMinosRun = true;
811
812 // retrieve parameters in case a new minimum has been found
813 if (fMinuit->fCstatu == "SUCCESSFUL")
815
816 double errParab = 0;
817 double gcor = 0;
818 // what returns if parameter fixed or constant or at limit ?
819 fMinuit->mnerrs(i,errUp,errLow, errParab, gcor);
820
821 // do not flag errors case of PROBLEMS or UNCHANGED (
822 return isValid;
823
824}
825
827 // reset TMinuit
828
829 fMinuit->mncler();
830
831 //reset the internal Minuit random generator to its initial state
832 double val = 3;
833 int inseed = 12345;
834 fMinuit->mnrn15(val,inseed);
835
836 fUsed = false;
837 fgUsed = false;
838
839}
840
842 // check if a parameter is defined and in case it was fixed released
843 // TMinuit is not able to release free parameters by redefining them
844 // so we need to force the release
845 if (fMinuit == nullptr) return;
846 if (fMinuit->GetNumFixedPars() == 0) return;
847 // check if parameter has already been defined
848 if (int(ivar) >= fMinuit->GetNumPars() ) return;
849
850 // check if parameter is fixed
851 for (int i = 0; i < fMinuit->fNpfix; ++i) {
852 if (fMinuit->fIpfix[i] == ivar+1 ) {
853 // parameter is fixed
854 fMinuit->Release(ivar);
855 return;
856 }
857 }
858
859}
860
861
863 // print-out results using classic Minuit format (mnprin)
864 if (fMinuit == nullptr) return;
865
866 // print minimizer result
867 if (PrintLevel() > 2)
869 else
871}
872
874 // suppress Minuit2 warnings
875 double arglist = 0;
876 int ierr = 0;
877 if (nowarn)
878 fMinuit->mnexcm("SET NOW",&arglist,0,ierr);
879 else
880 fMinuit->mnexcm("SET WAR",&arglist,0,ierr);
881}
882
883
884bool TMinuitMinimizer::Contour(unsigned int ipar, unsigned int jpar, unsigned int &npoints, double * x, double * y) {
885 // contour plot for parameter i and j
886 // need a valid FunctionMinimum otherwise exits
887 if (fMinuit == nullptr) {
888 Error("TMinuitMinimizer::Contour"," invalid TMinuit instance");
889 return false;
890 }
891
892 // set error and print level
893 double arglist[1];
894 int ierr = 0;
895 arglist[0] = ErrorDef();
896 fMinuit->mnexcm("SET Err",arglist,1,ierr);
897
898 arglist[0] = PrintLevel()-1;
899 fMinuit->mnexcm("SET PRINT",arglist,1,ierr);
900
901 // suppress warning in case Printlevel() == 0
902 if (PrintLevel() == 0) fMinuit->mnexcm("SET NOW",arglist,0,ierr);
903
904 // set precision if needed
905 if (Precision() > 0) {
906 arglist[0] = Precision();
907 fMinuit->mnexcm("SET EPS",arglist,1,ierr);
908 }
909
910
911 if (npoints < 4) {
912 Error("TMinuitMinimizer::Contour","Cannot make contour with so few points");
913 return false;
914 }
915 int npfound = 0;
916 // parameter numbers in mncont start from zero
917 fMinuit->mncont( ipar,jpar,npoints, x, y,npfound);
918 if (npfound<4) {
919 // mncont did go wrong
920 Error("TMinuitMinimizer::Contour","Cannot find more than 4 points");
921 return false;
922 }
923 if (npfound!=(int)npoints) {
924 // mncont did go wrong
925 Warning("TMinuitMinimizer::Contour","Returning only %d points ",npfound);
926 npoints = npfound;
927 }
928 return true;
929
930}
931
932bool TMinuitMinimizer::Scan(unsigned int ipar, unsigned int & nstep, double * x, double * y, double xmin, double xmax) {
933 // scan a parameter (variable) around the minimum value
934 // the parameters must have been set before
935 // if xmin=0 && xmax == 0 by default scan around 2 sigma of the error
936 // if the errors are also zero then scan from min and max of parameter range
937 // (if parameters are limited Minuit scan from min and max instead of 2 sigma by default)
938 // (force in that case to use errors)
939
940 // scan is not implemented for TMinuit, the way to return the array is only via the graph
941 if (fMinuit == nullptr) {
942 Error("TMinuitMinimizer::Scan"," invalid TMinuit instance");
943 return false;
944 }
945
946 // case of default xmin and xmax
947 if (xmin >= xmax && (int) ipar < fMinuit->GetNumPars() ) {
948 double val = 0; double err = 0;
950 double xlow = 0; double xup = 0 ;
951 int iuint = 0;
952 // in mnpout index starts from ze
953 fMinuit->mnpout( ipar, name, val, err, xlow, xup, iuint);
954 // redefine 2 sigma for all parameters by default (TMinuit does 1 sigma and range if limited)
955 if (iuint > 0 && err > 0) {
956 xmin = val - 2.*err;
957 xmax = val + 2 * err;
958 }
959 }
960
961 double arglist[4];
962 int ierr = 0;
963
964 arglist[0] = PrintLevel()-1;
965 fMinuit->mnexcm("SET PRINT",arglist,1,ierr);
966 // suppress warning in case Printlevel() == 0
967 if (PrintLevel() == 0) fMinuit->mnexcm("SET NOW",arglist,0,ierr);
968
969 // set precision if needed
970 if (Precision() > 0) {
971 arglist[0] = Precision();
972 fMinuit->mnexcm("SET EPS",arglist,1,ierr);
973 }
974
975 if (nstep == 0) return false;
976 arglist[0] = ipar+1; // TMinuit starts from 1
977 arglist[1] = nstep; //+2; // TMinuit deletes two points
978 int nargs = 2;
979 if (xmax > xmin ) {
980 arglist[2] = xmin;
981 arglist[3] = xmax;
982 nargs = 4;
983 }
984 fMinuit->mnexcm("SCAN",arglist,nargs,ierr);
985 if (ierr) {
986 Error("TMinuitMinimizer::Scan"," Error executing command SCAN");
987 return false;
988 }
989 // get TGraph object
990 TGraph * gr = dynamic_cast<TGraph *>(fMinuit->GetPlot() );
991 if (!gr) {
992 Error("TMinuitMinimizer::Scan"," Error in returned graph object");
993 return false;
994 }
995 nstep = std::min(gr->GetN(), (int) nstep);
996
997
998 std::copy(gr->GetX(), gr->GetX()+nstep, x);
999 std::copy(gr->GetY(), gr->GetY()+nstep, y);
1000 nstep = gr->GetN();
1001 return true;
1002}
1003
1005 // perform calculation of Hessian
1006
1007 if (fMinuit == nullptr) {
1008 Error("TMinuitMinimizer::Hesse","invalid TMinuit pointer. Need to call first SetFunction and SetVariable");
1009 return false;
1010 }
1011
1012
1013 double arglist[10];
1014 int ierr = 0;
1015
1016 // set error and print level
1017 arglist[0] = ErrorDef();
1018 fMinuit->mnexcm("SET ERR",arglist,1,ierr);
1019
1020 int printlevel = PrintLevel();
1021 arglist[0] = printlevel - 1;
1022 fMinuit->mnexcm("SET PRINT",arglist,1,ierr);
1023
1024 // suppress warning in case Printlevel() == 0
1025 if (printlevel == 0) fMinuit->mnexcm("SET NOW",arglist,0,ierr);
1026
1027 // set precision if needed
1028 if (Precision() > 0) {
1029 arglist[0] = Precision();
1030 fMinuit->mnexcm("SET EPS",arglist,1,ierr);
1031 }
1032
1033 arglist[0] = MaxFunctionCalls();
1034
1035 fMinuit->mnexcm("HESSE",arglist,1,ierr);
1036 fStatus += 100*ierr;
1037
1038 if (ierr != 0) return false;
1039
1040 // retrieve results (parameter and error matrix)
1041 // only if result is OK
1042
1045
1046 return true;
1047}
1048
1050 // set debug mode
1051
1052 if (fMinuit == nullptr) {
1053 Error("TMinuitMinimizer::SetDebug","invalid TMinuit pointer. Need to call first SetFunction and SetVariable");
1054 return false;
1055 }
1056 int ierr = 0;
1057 double arglist[1];
1058 arglist[0] = 1;
1059 if (on)
1060 fMinuit->mnexcm("SET DEBUG",arglist,1,ierr);
1061 else
1062 fMinuit->mnexcm("SET NODEBUG",arglist,1,ierr);
1063
1064 return (ierr == 0);
1065}
1066// } // end namespace Fit
1067
1068// } // end namespace ROOT
#define f(i)
Definition RSha256.hxx:104
#define g(i)
Definition RSha256.hxx:105
#define ClassImp(name)
Definition Rtypes.h:382
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
Definition TError.h:125
void Info(const char *location, const char *msgfmt,...)
Use this function for informational messages.
Definition TError.cxx:218
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:185
void Warning(const char *location, const char *msgfmt,...)
Use this function in warning situations.
Definition TError.cxx:229
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void on
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
char name[80]
Definition TGX11.cxx:110
float xmin
float xmax
static ROOT::Math::IMultiGenFunction *& GetGlobalFuncPtr()
R__EXTERN TMinuit * gMinuit
Definition TMinuit.h:271
#define gROOT
Definition TROOT.h:406
Class, describing value, limits and step size of the parameters Provides functionality also to set/re...
void Set(const std::string &name, double value, double step)
set value and name (unlimited parameter)
void Fix()
fix the parameter
Documentation for the abstract class IBaseFunctionMultiDim.
Definition IFunction.h:61
virtual unsigned int NDim() const =0
Retrieve the dimension of the function.
Interface (abstract class) for multi-dimensional functions providing a gradient calculation.
Definition IFunction.h:168
virtual void Gradient(const T *x, T *grad) const
Evaluate all the vector of function derivatives (gradient) at a point x.
Definition IFunction.h:177
double Tolerance() const
absolute tolerance
Definition Minimizer.h:300
unsigned int MaxFunctionCalls() const
max number of function calls
Definition Minimizer.h:294
double Precision() const
precision of minimizer in the evaluation of the objective function ( a value <=0 corresponds to the l...
Definition Minimizer.h:304
int fStatus
status of minimizer
Definition Minimizer.h:371
int Strategy() const
strategy
Definition Minimizer.h:307
double ErrorDef() const
return the statistical scale used for calculate the error is typically 1 for Chi2 and 0....
Definition Minimizer.h:317
bool IsValidError() const
return true if Minimizer has performed a detailed error validation (e.g. run Hesse for Minuit)
Definition Minimizer.h:320
int PrintLevel() const
minimizer configuration parameters
Definition Minimizer.h:291
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
Double_t * GetY() const
Definition TGraph.h:140
Int_t GetN() const
Definition TGraph.h:132
Double_t * GetX() const
Definition TGraph.h:139
TMatrixTSym< Element > & Invert(Double_t *det=nullptr)
Invert the matrix and calculate its determinant Notice that the LU decomposition is used instead of B...
const Element * GetMatrixArray() const override
TMinuitMinimizer class: ROOT::Math::Minimizer implementation based on TMinuit.
bool FixVariable(unsigned int) override
fix an existing variable
TMinuitMinimizer(ROOT::Minuit::EMinimizerType type=ROOT::Minuit::kMigrad, unsigned int ndim=0)
Default constructor.
void RetrieveErrorMatrix()
retrieve error matrix from TMinuit
bool SetFixedVariable(unsigned int, const std::string &, double) override
set fixed variable (override if minimizer supports them )
static TMinuit * fgMinuit
bool SetVariableLimits(unsigned int ivar, double lower, double upper) override
set the limits of an existing variable
static void Fcn(int &, double *, double &f, double *, int)
implementation of FCN for Minuit
bool SetVariable(unsigned int ivar, const std::string &name, double val, double step) override
set free variable
bool ReleaseVariable(unsigned int) override
release an existing variable
static void FcnGrad(int &, double *g, double &f, double *, int)
implementation of FCN for Minuit when user provided gradient is used
ROOT::Minuit::EMinimizerType fType
bool SetVariableLowerLimit(unsigned int, double) override
set the lower-limit of an existing variable
int VariableIndex(const std::string &name) const override
get index of variable given a variable given a name return always -1 .
double MinValue() const override
return minimum function value
bool SetVariableStepSize(unsigned int, double) override
set the step size of an existing variable
void RetrieveParams()
retrieve minimum parameters and errors from TMinuit
~TMinuitMinimizer() override
Destructor (no operations)
std::vector< double > fErrors
void SetFunction(const ROOT::Math::IMultiGenFunction &func) override
set the function to minimize
bool CheckMinuitInstance() const
check TMinuit instance
static bool fgUseStaticMinuit
std::vector< double > fCovar
bool Hesse() override
perform a full calculation of the Hessian matrix for error calculation
bool CheckVarIndex(unsigned int ivar) const
check parameter
bool Minimize() override
method to perform the minimization
bool SetVariableValue(unsigned int, double) override
set the value of an existing variable
int CovMatrixStatus() const override
return status of covariance matrix
void SuppressMinuitWarnings(bool nowarn=true)
suppress the minuit warnings (if called with false will enable them) By default they are suppressed o...
bool SetVariableUpperLimit(unsigned int, double) override
set the upper-limit of an existing variable
double Edm() const override
return expected distance reached from the minimum
static bool UseStaticMinuit(bool on=true)
static function to switch on/off usage of static global TMinuit instance (gMinuit) By default it is u...
void PrintResults() override
return reference to the objective function virtual const ROOT::Math::IGenFunction & Function() const ...
unsigned int NCalls() const override
number of function calls to reach the minimum
void DoReleaseFixParameter(int ivar)
release a parameter that is fixed when it is redefined
bool SetLimitedVariable(unsigned int ivar, const std::string &name, double val, double step, double, double) override
set upper/lower limited variable (override if minimizer supports them )
bool GetMinosError(unsigned int i, double &errLow, double &errUp, int=0) override
minos error for variable i, return false if Minos failed
bool GetVariableSettings(unsigned int, ROOT::Fit::ParameterSettings &) const override
get variable settings in a variable object (like ROOT::Fit::ParamsSettings)
double GlobalCC(unsigned int) const override
global correlation coefficient for variable i
bool IsFixedVariable(unsigned int) const override
query if an existing variable is fixed (i.e.
std::string VariableName(unsigned int ivar) const override
return reference to the objective function virtual const ROOT::Math::IGenFunction & Function() const;
bool SetUpperLimitedVariable(unsigned int ivar, const std::string &name, double val, double step, double upper) override
set upper limit variable (override if minimizer supports them )
bool SetLowerLimitedVariable(unsigned int ivar, const std::string &name, double val, double step, double lower) override
set lower limit variable (override if minimizer supports them )
unsigned int NFree() const override
number of free variables (real dimension of the problem) this is <= Function().NDim() which is the to...
std::vector< double > fParams
void InitTMinuit(int ndim)
initialize the TMinuit instance
bool GetHessianMatrix(double *h) const override
Fill the passed array with the Hessian matrix elements The Hessian matrix is the matrix of the second...
bool Scan(unsigned int i, unsigned int &nstep, double *x, double *y, double xmin=0, double xmax=0) override
scan a parameter i around the minimum.
bool SetDebug(bool on=true)
set debug mode. Return true if setting was successful
bool GetCovMatrix(double *cov) const override
Fill the passed array with the covariance matrix elements if the variable is fixed or const the value...
bool Contour(unsigned int i, unsigned int j, unsigned int &npoints, double *xi, double *xj) override
find the contour points (xi,xj) of the function for parameter i and j around the minimum The contour ...
virtual Int_t GetParameter(Int_t parNo, Double_t &currentValue, Double_t &currentError) const
return parameter value and error
Definition TMinuit.cxx:841
virtual Int_t FixParameter(Int_t parNo)
fix a parameter
Definition TMinuit.cxx:827
virtual Int_t GetNumPars() const
returns the total number of parameters that have been defined as fixed or free.
Definition TMinuit.cxx:872
virtual Int_t GetNumFixedPars() const
returns the number of currently fixed parameters
Definition TMinuit.cxx:855
virtual Int_t Release(Int_t parNo)
release a parameter
Definition TMinuit.cxx:894
Double_t fUndefi
Definition TMinuit.h:60
Int_t fNu
Definition TMinuit.h:130
Double_t * fGlobcc
Definition TMinuit.h:74
virtual void mncler()
Resets the parameter list to UNDEFINED.
Definition TMinuit.cxx:1103
Double_t fUp
Definition TMinuit.h:50
TString * fCpnam
Character to be plotted at the X,Y contour positions.
Definition TMinuit.h:165
Int_t fNfcn
Definition TMinuit.h:145
Double_t fBigedm
Definition TMinuit.h:61
TString fCstatu
Definition TMinuit.h:167
virtual TObject * GetPlot() const
Definition TMinuit.h:200
Int_t fNpfix
Definition TMinuit.h:37
Int_t fISW[7]
Definition TMinuit.h:141
Int_t * fIpfix
Definition TMinuit.h:129
virtual void SetFCN(void(*fcn)(Int_t &, Double_t *, Double_t &f, Double_t *, Int_t))
To set the address of the minimization function.
Definition TMinuit.cxx:920
virtual void mncont(Int_t ke1, Int_t ke2, Int_t nptu, Double_t *xptu, Double_t *yptu, Int_t &ierrf)
Find points along a contour where FCN is minimum.
Definition TMinuit.cxx:1395
Double_t fEpsma2
Definition TMinuit.h:57
virtual void mnpout(Int_t iuext, TString &chnam, Double_t &val, Double_t &err, Double_t &xlolim, Double_t &xuplim, Int_t &iuint) const
Provides the user with information concerning the current status.
Definition TMinuit.cxx:6247
Int_t fIstrat
Definition TMinuit.h:150
virtual void mnemat(Double_t *emat, Int_t ndim)
Calculates the external error matrix from the internal matrix.
Definition TMinuit.cxx:2501
virtual void mnrn15(Double_t &val, Int_t &inseed)
This is a super-portable random number generator.
Definition TMinuit.cxx:6619
virtual void mnerrs(Int_t number, Double_t &eplus, Double_t &eminus, Double_t &eparab, Double_t &gcc)
Utility routine to get MINOS errors.
Definition TMinuit.cxx:2578
Int_t fNpar
Definition TMinuit.h:41
virtual void mnexcm(const char *comand, Double_t *plist, Int_t llist, Int_t &ierflg)
Interprets a command and takes appropriate action.
Definition TMinuit.cxx:2664
Int_t fMaxpar
Definition TMinuit.h:39
Int_t * fNiofex
Definition TMinuit.h:127
virtual Int_t DefineParameter(Int_t parNo, const char *name, Double_t initVal, Double_t initErr, Double_t lowerLimit, Double_t upperLimit)
Define a parameter.
Definition TMinuit.cxx:695
virtual void mnprin(Int_t inkode, Double_t fval)
Prints the values of the parameters at the time of the call.
Definition TMinuit.cxx:6304
Double_t fAmin
Definition TMinuit.h:49
Double_t fEDM
Definition TMinuit.h:51
Basic string class.
Definition TString.h:139
const char * Data() const
Definition TString.h:376
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
TGraphErrors * gr
Definition legend1.C:25
TMarker m
Definition textangle.C:8
TLine l
Definition textangle.C:4