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//
35// TMinuitMinimizer class implementing the ROOT::Math::Minimizer interface using
36// TMinuit.
37// This class is normally instantiates using the plug-in manager
38// (plug-in with name Minuit or TMinuit)
39// In addition the user can choose the minimizer algorithm: Migrad (the default one), Simplex, or Minimize (combined Migrad + Simplex)
40//
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(0)
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(0)
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 = 0;
107 }
108}
109
111 Minimizer()
112{
113 // Implementation of copy constructor (it is private).
114}
115
117{
118 // Implementation of assignment operator (private)
119 if (this == &rhs) return *this; // time saving self-test
120 return *this;
121}
122
124 // static method to control usage of global TMinuit instance
125 bool prev = fgUseStaticMinuit;
127 return prev;
128}
129
131
132 // when called a second time check dimension - create only if needed
133 // initialize the minuit instance - recreating a new one if needed
134 if (fMinuit ==0 || dim > fMinuit->fMaxpar) {
135
136 // case not using the global instance - recreate it all the time
137 if (fgUseStaticMinuit) {
138
139 // re-use gMinuit as static instance of TMinuit
140 // which can be accessed by the user after minimization
141 // check if fgMinuit is different than gMinuit
142 // case 1: fgMinuit not zero but fgMinuit has been deleted (not in gROOT): set to zero
143 // case 2: fgMinuit not zero and exists in global list : set fgMinuit to gMinuit
144 // case 3: fgMinuit zero - and gMinuit not zero: create a new instance locally to avoid conflict
145 if (fgMinuit != gMinuit) {
146 // if object exists in gROOT remove it to avoid a memory leak
147 if (fgMinuit ) {
148 if (gROOT->GetListOfSpecials()->FindObject(fgMinuit) == 0) {
149 // case 1: object does not exists in gROOT - means it has been deleted
150 fgMinuit = 0;
151 }
152 else {
153 // case 2: object exists - but gMinuit points to something else
154 // restore gMinuit to the one used before by TMinuitMinimizer
156 }
157 }
158 else {
159 // case 3: avoid reusing existing one - mantain fgMinuit to zero
160 // otherwise we will get a double delete if user deletes externally gMinuit
161 // in this case we will loose gMinuit instance
162// fgMinuit = gMinuit;
163// fgUsed = true; // need to reset in case other gMinuit instance is later used
164 }
165 }
166
167 // check if need to create a new TMinuit instance
168 if (fgMinuit == 0) {
169 fgUsed = false;
170 fgMinuit = new TMinuit(dim);
171 }
172 else if (fgMinuit->GetNumPars() != int(dim) ) {
173 delete fgMinuit;
174 fgUsed = false;
175 fgMinuit = new TMinuit(dim);
176 }
177
179 }
180
181 else {
182 // re- create all the time a new instance of TMinuit (fgUseStaticMinuit is false)
183 if (fMinuit) delete fMinuit;
184 fMinuit = new TMinuit(dim);
186 fgUsed = false;
187 }
188
189 } // endif fMinuit ==0 || dim > fMaxpar
190
191 fDim = dim;
192
194
195 // set print level in TMinuit
196 double arglist[1];
197 int ierr= 0;
198 // TMinuit level is shift by 1 -1 means 0;
199 arglist[0] = PrintLevel() - 1;
200 fMinuit->mnexcm("SET PRINT",arglist,1,ierr);
202}
203
204
206 // Set the objective function to be minimized, by passing a function object implement the
207 // basic multi-dim Function interface. In this case the derivatives will be
208 // calculated by Minuit
209 // Here a TMinuit instance is created since only at this point we know the number of parameters
210
211
212 fDim = func.NDim();
213
214 // create TMinuit if needed
216
217 // assign to the static pointer (NO Thread safety here)
218 GetGlobalFuncPtr() = const_cast<ROOT::Math::IMultiGenFunction *>(&func);
220
221 // switch off gradient calculations
222 double arglist[1];
223 int ierr = 0;
224 fMinuit->mnexcm("SET NOGrad",arglist,0,ierr);
225}
226
228 // Set the objective function to be minimized, by passing a function object implement the
229 // multi-dim gradient Function interface. In this case the function derivatives are provided
230 // by the user via this interface and there not calculated by Minuit.
231
232 fDim = func.NDim();
233
234 // create TMinuit if needed
236
237 // assign to the static pointer (NO Thread safety here)
238 GetGlobalFuncPtr() = const_cast<ROOT::Math::IMultiGradFunction *>(&func);
240
241 // set gradient
242 // by default do not check gradient calculation
243 // it cannot be done here, check can be done only after having defined the parameters
244 double arglist[1];
245 int ierr = 0;
246 arglist[0] = 1;
247 fMinuit->mnexcm("SET GRAD",arglist,1,ierr);
248}
249
250void TMinuitMinimizer::Fcn( int &, double * , double & f, double * x , int /* iflag */) {
251 // implementation of FCN static function used internally by TMinuit.
252 // Adapt IMultiGenFunction interface to TMinuit FCN static function
253 f = GetGlobalFuncPtr()->operator()(x);
254}
255
256void TMinuitMinimizer::FcnGrad( int &, double * g, double & f, double * x , int iflag ) {
257 // implementation of FCN static function used internally by TMinuit.
258 // Adapt IMultiGradFunction interface to TMinuit FCN static function in the case of user
259 // provided gradient.
261
262 assert(gFunc != 0);
263 f = gFunc->operator()(x);
264
265 // calculates also derivatives
266 if (iflag == 2) gFunc->Gradient(x,g);
267}
268
269bool TMinuitMinimizer::SetVariable(unsigned int ivar, const std::string & name, double val, double step) {
270 // set a free variable.
271 if (!CheckMinuitInstance()) return false;
272
274
275 // clear after minimization when setting params
276 if (fUsed) DoClear();
277
278 // check if parameter was defined and in case it was fixed, release it
280
281 int iret = fMinuit->DefineParameter(ivar , name.c_str(), val, step, 0., 0. );
282 return (iret == 0);
283}
284
285bool TMinuitMinimizer::SetLimitedVariable(unsigned int ivar, const std::string & name, double val, double step, double lower, double upper) {
286 // set a limited variable.
287 if (!CheckMinuitInstance()) return false;
288
290
291 // clear after minimization when setting params
292 if (fUsed) DoClear();
293
294 // check if parameter was defined and in case it was fixed, release it
296
297 int iret = fMinuit->DefineParameter(ivar, name.c_str(), val, step, lower, upper );
298 return (iret == 0);
299}
300
301bool TMinuitMinimizer::SetLowerLimitedVariable(unsigned int ivar , const std::string & name , double val , double step , double lower ) {
302 // set a lower limited variable
303 // since is not supported in TMinuit , just use a artificial large value
304 Warning("TMinuitMinimizer::SetLowerLimitedVariable","not implemented - use as upper limit 1.E7 instead of +inf");
305 return SetLimitedVariable(ivar, name, val , step, lower, lower+ 1.E7); // use 1.E7 which will make TMinuit happy
306}
307
308bool TMinuitMinimizer::SetUpperLimitedVariable(unsigned int ivar , const std::string & name , double val , double step , double upper ) {
309 // set a upper limited variable
310 // since is not supported in TMinuit , just use a artificial large negative value
311 Warning("TMinuitMinimizer::SetUpperLimitedVariable","not implemented - use as lower limit -1.E7 instead of -inf");
312 return SetLimitedVariable(ivar, name, val , step, upper -1.E7, upper);
313}
314
315
317 // check instance of fMinuit
318 if (fMinuit == 0) {
319 Error("TMinuitMinimizer::CheckMinuitInstance","Invalid TMinuit pointer. Need to call first SetFunction");
320 return false;
321 }
322 return true;
323}
324
325bool TMinuitMinimizer::CheckVarIndex(unsigned int ivar) const {
326 // check index of Variable (assume fMinuit exists)
327 if ((int) ivar >= fMinuit->fNu ) {
328 Error("TMinuitMinimizer::CheckVarIndex","Invalid parameter index");
329 return false;
330 }
331 return true;
332}
333
334
335bool TMinuitMinimizer::SetFixedVariable(unsigned int ivar, const std::string & name, double val) {
336 // set a fixed variable.
337 if (!CheckMinuitInstance()) return false;
338
339 // clear after minimization when setting params
341
342 // clear after minimization when setting params
343 if (fUsed) DoClear();
344
345 // put an arbitrary step (0.1*abs(value) otherwise TMinuit consider the parameter as constant
346 // constant parameters are treated differently (they are ignored inside TMinuit and not considered in the
347 // total list of parameters)
348 double step = ( val != 0) ? 0.1 * std::abs(val) : 0.1;
349 int iret = fMinuit->DefineParameter(ivar, name.c_str(), val, step, 0., 0. );
350 if (iret == 0) iret = fMinuit->FixParameter(ivar);
351 return (iret == 0);
352}
353
354bool TMinuitMinimizer::SetVariableValue(unsigned int ivar, double val) {
355 // set the value of an existing variable
356 // parameter must exist or return false
357 if (!CheckMinuitInstance()) return false;
358
359 double arglist[2];
360 int ierr = 0;
361
362 arglist[0] = ivar+1; // TMinuit starts from 1
363 arglist[1] = val;
364 fMinuit->mnexcm("SET PAR",arglist,2,ierr);
365 return (ierr==0);
366}
367
368bool TMinuitMinimizer::SetVariableStepSize(unsigned int ivar, double step) {
369 // set the step-size of an existing variable
370 // parameter must exist or return false
371 if (!CheckMinuitInstance()) return false;
372 // need to re-implement re-calling mnparm
373 // get first current parameter values and limits
374 double curval,err, lowlim, uplim;
375 int iuint; // internal index
377 fMinuit->mnpout(ivar, name, curval, err, lowlim, uplim,iuint);
378 if (iuint == -1) return false;
379 int iret = fMinuit->DefineParameter(ivar, name, curval, step, lowlim, uplim );
380 return (iret == 0);
381
382}
383
384bool TMinuitMinimizer::SetVariableLowerLimit(unsigned int ivar, double lower ) {
385 // set the limits of an existing variable
386 // parameter must exist or return false
387 Warning("TMinuitMinimizer::SetVariableLowerLimit","not implemented - use as upper limit 1.E30 instead of +inf");
388 return SetVariableLimits(ivar, lower, 1.E30);
389}
390bool TMinuitMinimizer::SetVariableUpperLimit(unsigned int ivar, double upper ) {
391 // set the limits of an existing variable
392 // parameter must exist or return false
393 Warning("TMinuitMinimizer::SetVariableUpperLimit","not implemented - - use as lower limit -1.E30 instead of +inf");
394 return SetVariableLimits(ivar, -1.E30, upper);
395}
396
397bool TMinuitMinimizer::SetVariableLimits(unsigned int ivar, double lower, double upper) {
398 // set the limits of an existing variable
399 // parameter must exist or return false
400
401 if (!CheckMinuitInstance()) return false;
402 // need to re-implement re-calling mnparm
403 // get first current parameter values and limits
404 double curval,err, lowlim, uplim;
405 int iuint; // internal index
407 fMinuit->mnpout(ivar, name, curval, err, lowlim, uplim,iuint);
408 if (iuint == -1) return false;
409 int iret = fMinuit->DefineParameter(ivar, name, curval, err, lower, upper );
410 return (iret == 0);
411
412}
413
414bool TMinuitMinimizer::FixVariable(unsigned int ivar) {
415 // Fix an existing variable
416 if (!CheckMinuitInstance()) return false;
417 if (!CheckVarIndex(ivar)) return false;
418 int iret = fMinuit->FixParameter(ivar);
419 return (iret == 0);
420}
421
422bool TMinuitMinimizer::ReleaseVariable(unsigned int ivar) {
423 // Fix an existing variable
424 if (!CheckMinuitInstance()) return false;
425 if (!CheckVarIndex(ivar)) return false;
426 int iret = fMinuit->Release(ivar);
427 return (iret == 0);
428}
429
430bool TMinuitMinimizer::IsFixedVariable(unsigned int ivar) const {
431 // query if variable is fixed
432 if (!CheckMinuitInstance()) return false;
433 if (!CheckVarIndex(ivar)) return false;
434 return (fMinuit->fNiofex[ivar] == 0 );
435}
436
438 // retrieve variable settings (all set info on the variable)
439 if (!CheckMinuitInstance()) return false;
440 if (!CheckVarIndex(ivar)) return false;
441 double curval,err, lowlim, uplim;
442 int iuint; // internal index
444 fMinuit->mnpout(ivar, name, curval, err, lowlim, uplim,iuint);
445 if (iuint == -1) return false;
446 var.Set(name.Data(), curval, err, lowlim, uplim);
447 if (IsFixedVariable(ivar)) var.Fix();
448 return true;
449}
450
451
452
453std::string TMinuitMinimizer::VariableName(unsigned int ivar) const {
454 // return the variable name
455 if (!CheckMinuitInstance()) return std::string();
456 if (!CheckVarIndex(ivar)) return std::string();
457 return fMinuit->fCpnam[ivar].Data();
458}
459
460int TMinuitMinimizer::VariableIndex(const std::string & ) const {
461 // return variable index
462 Error("TMinuitMinimizer::VariableIndex"," find index of a variable from its name is not implemented in TMinuit");
463 return -1;
464}
465
467 // perform the minimization using the algorithm chosen previously by the user
468 // By default Migrad is used.
469 // Return true if the found minimum is valid and update internal chached values of
470 // minimum values, errors and covariance matrix.
471 // Status of minimizer is set to:
472 // migradResult + 10*minosResult + 100*hesseResult + 1000*improveResult
473
474
475 if (fMinuit == 0) {
476 Error("TMinuitMinimizer::Minimize","invalid TMinuit pointer. Need to call first SetFunction and SetVariable");
477 return false;
478 }
479
480
481 // total number of parameter defined in Minuit is fNu
482 if (fMinuit->fNu < static_cast<int>(fDim) ) {
483 Error("TMinuitMinimizer::Minimize","The total number of defined parameters is different than the function dimension, npar = %d, dim = %d",fMinuit->fNu, fDim);
484 return false;
485 }
486
487 int printlevel = PrintLevel();
488
489 // total number of free parameter is 0
490 if (fMinuit->fNpar <= 0) {
491 // retrieve parameters values from TMinuit
493 fMinuit->fAmin = (*GetGlobalFuncPtr())(&fParams.front());
494 if (printlevel > 0) Info("TMinuitMinimizer::Minimize","There are no free parameter - just compute the function value");
495 return true;
496 }
497
498
499 double arglist[10];
500 int ierr = 0;
501
502
503 // set error and print level
504 arglist[0] = ErrorDef();
505 fMinuit->mnexcm("SET Err",arglist,1,ierr);
506
507 arglist[0] = printlevel - 1;
508 fMinuit->mnexcm("SET PRINT",arglist,1,ierr);
509
510 // suppress warning in case Printlevel() == 0
511 if (printlevel == 0) fMinuit->mnexcm("SET NOW",arglist,0,ierr);
512
513 // set precision if needed
514 if (Precision() > 0) {
515 arglist[0] = Precision();
516 fMinuit->mnexcm("SET EPS",arglist,1,ierr);
517 }
518
519 // set strategy
520 int strategy = Strategy();
521 if (strategy >=0 && strategy <=2 ) {
522 arglist[0] = strategy;
523 fMinuit->mnexcm("SET STR",arglist,1,ierr);
524 }
525
526 arglist[0] = MaxFunctionCalls();
527 arglist[1] = Tolerance();
528
529 int nargs = 2;
530
531 switch (fType){
533 // case of Migrad
534 fMinuit->mnexcm("MIGRAD",arglist,nargs,ierr);
535 break;
537 // case of combined (Migrad+ simplex)
538 fMinuit->mnexcm("MINIMIZE",arglist,nargs,ierr);
539 break;
541 // case of Simlex
542 fMinuit->mnexcm("SIMPLEX",arglist,nargs,ierr);
543 break;
545 // case of Scan (scan all parameters with default values)
546 nargs = 0;
547 fMinuit->mnexcm("SCAN",arglist,nargs,ierr);
548 break;
550 // case of Seek (random find minimum in a hypercube around current parameter values
551 // use Tolerance as measures for standard deviation (if < 1) used default value in Minuit ( supposed to be 3)
552 nargs = 1;
553 if (arglist[1] >= 1.) nargs = 2;
554 fMinuit->mnexcm("SEEK",arglist,nargs,ierr);
555 break;
556 default:
557 // default: use Migrad
558 fMinuit->mnexcm("MIGRAD",arglist,nargs,ierr);
559
560 }
561
562 fgUsed = true;
563 fUsed = true;
564
565 fStatus = ierr;
566 int minErrStatus = ierr;
567
568 if (printlevel>2) Info("TMinuitMinimizer::Minimize","Finished to run MIGRAD - status %d",ierr);
569
570 // run improved if needed
571 if (ierr == 0 && fType == ROOT::Minuit::kMigradImproved) {
572 fMinuit->mnexcm("IMPROVE",arglist,1,ierr);
573 fStatus += 1000*ierr;
574 if (printlevel>2) Info("TMinuitMinimizer::Minimize","Finished to run IMPROVE - status %d",ierr);
575 }
576
577
578 // check if Hesse needs to be run
579 // Migrad runs inside it automatically for strategy >=1. Do also
580 // in case improve or other minimizers are used
581 // run Hesse also in case cov matrix has not been computed or has been made pos-def
582 // or a valid error analysis is requested (when IsValidError() == true)
583 if (minErrStatus == 0 && (IsValidError() || ( strategy >=1 && CovMatrixStatus() < 3) ) ) {
584 fMinuit->mnexcm("HESSE",arglist,1,ierr);
585 fStatus += 100*ierr;
586 // here ierr is zero when Hesse fails CovMatrixStatus = 0 or 1 (2 is when was made posdef)
587 if (ierr == 0 && CovMatrixStatus() < 2){
588 fStatus += 100*(CovMatrixStatus()+1);
589 }
590 // should check here cov matrix status??? (if <3 flag error ?)
591 if (printlevel>2) Info("TMinuitMinimizer::Minimize","Finished to run HESSE - status %d",ierr);
592 }
593
594 // retrieve parameters and errors from TMinuit
596
597 if (minErrStatus == 0) {
598
599 // store global min results (only if minimization is OK)
600 // ignore cases when Hesse or IMprove return error different than zero
602
603 // need to re-run Minos again if requested
604 fMinosRun = false;
605
606 return true;
607
608 }
609 return false;
610
611}
612
614 // retrieve from TMinuit minimum parameter values
615 // and errors
616
617 assert(fMinuit != 0);
618
619 // get parameter values
620 if (fParams.size() != fDim) fParams.resize( fDim);
621 if (fErrors.size() != fDim) fErrors.resize( fDim);
622 for (unsigned int i = 0; i < fDim; ++i) {
623 fMinuit->GetParameter( i, fParams[i], fErrors[i]);
624 }
625}
626
628 // get covariance error matrix from TMinuit
629 // when some parameters are fixed filled the corresponding rows and column with zero's
630
631 assert(fMinuit != 0);
632
633 unsigned int nfree = NFree();
634
635 unsigned int ndim2 = fDim*fDim;
636 if (fCovar.size() != ndim2 ) fCovar.resize(fDim*fDim);
637 if (nfree >= fDim) { // no fixed parameters
638 fMinuit->mnemat(&fCovar.front(), fDim);
639 }
640 else {
641 // case of fixed params need to take care
642 std::vector<double> tmpMat(nfree*nfree);
643 fMinuit->mnemat(&tmpMat.front(), nfree);
644
645 unsigned int l = 0;
646 for (unsigned int i = 0; i < fDim; ++i) {
647
648 if ( fMinuit->fNiofex[i] > 0 ) { // not fixed ?
649 unsigned int m = 0;
650 for (unsigned int j = 0; j <= i; ++j) {
651 if ( fMinuit->fNiofex[j] > 0 ) { //not fixed
652 fCovar[i*fDim + j] = tmpMat[l*nfree + m];
653 fCovar[j*fDim + i] = fCovar[i*fDim + j];
654 m++;
655 }
656 }
657 l++;
658 }
659 }
660
661 }
662}
663
664unsigned int TMinuitMinimizer::NCalls() const {
665 // return total number of function calls
666 if (fMinuit == 0) return 0;
667 return fMinuit->fNfcn;
668}
669
671 // return minimum function value
672
673 // use part of code from mnstat
674 if (!fMinuit) return 0;
675 double minval = fMinuit->fAmin;
676 if (minval == fMinuit->fUndefi) return 0;
677 return minval;
678}
679
680double TMinuitMinimizer::Edm() const {
681 // return expected distance from the minimum
682
683 // use part of code from mnstat
684 if (!fMinuit) return -1;
685 if (fMinuit->fAmin == fMinuit->fUndefi || fMinuit->fEDM == fMinuit->fBigedm) return fMinuit->fUp;
686 return fMinuit->fEDM;
687}
688
689unsigned int TMinuitMinimizer::NFree() const {
690 // return number of free parameters
691 if (!fMinuit) return 0;
692 if (fMinuit->fNpar < 0) return 0;
693 return fMinuit->fNpar;
694}
695
696bool TMinuitMinimizer::GetCovMatrix(double * cov) const {
697 // get covariance matrix
698 int covStatus = CovMatrixStatus();
699 if ( fCovar.size() != fDim*fDim || covStatus < 2) {
700 Error("TMinuitMinimizer::GetHessianMatrix","Hessian matrix has not been computed - status %d",covStatus);
701 return false;
702 }
703 std::copy(fCovar.begin(), fCovar.end(), cov);
704 TMatrixDSym cmat(fDim,cov);
705 return true;
706}
707
708bool TMinuitMinimizer::GetHessianMatrix(double * hes) const {
709 // get Hessian - inverse of covariance matrix
710 // just invert it
711 // but need to get the compact form to avoid the zero for the fixed parameters
712 int covStatus = CovMatrixStatus();
713 if ( fCovar.size() != fDim*fDim || covStatus < 2) {
714 Error("TMinuitMinimizer::GetHessianMatrix","Hessian matrix has not been computed - status %d",covStatus);
715 return false;
716 }
717 // case of fixed params need to take care
718 unsigned int nfree = NFree();
719 TMatrixDSym mat(nfree);
720 fMinuit->mnemat(mat.GetMatrixArray(), nfree);
721 // invert the matrix
722 // probably need to check if failed. In that case inverse is equal to original
723 mat.Invert();
724
725 unsigned int l = 0;
726 for (unsigned int i = 0; i < fDim; ++i) {
727 if ( fMinuit->fNiofex[i] > 0 ) { // not fixed ?
728 unsigned int m = 0;
729 for (unsigned int j = 0; j <= i; ++j) {
730 if ( fMinuit->fNiofex[j] > 0 ) { //not fixed
731 hes[i*fDim + j] = mat(l,m);
732 hes[j*fDim + i] = hes[i*fDim + j];
733 m++;
734 }
735 }
736 l++;
737 }
738 }
739 return true;
740}
741// if ( fCovar.size() != fDim*fDim ) return false;
742// TMatrixDSym mat(fDim, &fCovar.front() );
743// std::copy(mat.GetMatrixArray(), mat.GetMatrixArray()+ mat.GetNoElements(), hes);
744// return true;
745// }
746
748 // return status of covariance matrix
749 // status: 0= not calculated at all
750 // 1= approximation only, not accurate
751 // 2= full matrix, but forced positive-definite
752 // 3= full accurate covariance matrix
753
754 // use part of code from mnstat
755 if (!fMinuit) return 0;
756 if (fMinuit->fAmin == fMinuit->fUndefi) return 0;
757 return fMinuit->fISW[1];
758}
759
760double TMinuitMinimizer::GlobalCC(unsigned int i) const {
761 // global correlation coefficient for parameter i
762 if (!fMinuit) return 0;
763 if (!fMinuit->fGlobcc) return 0;
764 if (int(i) >= fMinuit->fNu) return 0;
765 // get internal number in Minuit
766 int iin = fMinuit->fNiofex[i];
767 // index in TMinuit starts from 1
768 if (iin < 1) return 0;
769 return fMinuit->fGlobcc[iin-1];
770}
771
772bool TMinuitMinimizer::GetMinosError(unsigned int i, double & errLow, double & errUp, int ) {
773 // Perform Minos analysis for the given parameter i
774
775 if (fMinuit == 0) {
776 Error("TMinuitMinimizer::GetMinosError","invalid TMinuit pointer. Need to call first SetFunction and SetVariable");
777 return false;
778 }
779
780 // check if parameter is fixed
781 if (fMinuit->fNiofex[i] == 0 ) {
782 if (PrintLevel() > 0) Info("TMinuitMinimizer::GetMinosError","Parameter %s is fixed. There are no Minos error to calculate. Ignored.",VariableName(i).c_str());
783 errLow = 0; errUp = 0;
784 return true;
785 }
786
787 double arglist[2];
788 int ierr = 0;
789
790
791 // set error, print level, precision and strategy if they have changed
792 if (fMinuit->fUp != ErrorDef() ) {
793 arglist[0] = ErrorDef();
794 fMinuit->mnexcm("SET Err",arglist,1,ierr);
795 }
796
797 if (fMinuit->fISW[4] != (PrintLevel()-1) ) {
798 arglist[0] = PrintLevel()-1;
799 fMinuit->mnexcm("SET PRINT",arglist,1,ierr);
800 // suppress warning in case Printlevel() == 0
801 if (PrintLevel() == 0) fMinuit->mnexcm("SET NOW",arglist,0,ierr);
802 }
803 if (fMinuit->fIstrat != Strategy() ) {
804 arglist[0] = Strategy();
805 fMinuit->mnexcm("SET STR",arglist,1,ierr);
806 }
807
808 if (Precision() > 0 && fMinuit->fEpsma2 != Precision() ) {
809 arglist[0] = Precision();
810 fMinuit->mnexcm("SET EPS",arglist,1,ierr);
811 }
812
813
814 // syntax of MINOS command is: MINOS [maxcalls] [parno]
815 // if parno = 0 all parameters are done
816 arglist[0] = MaxFunctionCalls();
817 arglist[1] = i+1; // par number starts from 1 in TMInuit
818
819 int nargs = 2;
820 fMinuit->mnexcm("MINOS",arglist,nargs,ierr);
821 bool isValid = (ierr == 0);
822 // check also the status from fCstatu
823 if (isValid && fMinuit->fCstatu != "SUCCESSFUL") {
824 if (fMinuit->fCstatu == "FAILURE" ) {
825 // in this case MINOS failed on all parameters, so it is not valid !
826 ierr = 5;
827 isValid = false;
828 }
829 if (fMinuit->fCstatu == "PROBLEMS") ierr = 6;
830 ierr = 7; // this should be the case UNCHANGED
831 }
832
833 fStatus += 10*ierr;
834 fMinosStatus = ierr;
835
836 fMinosRun = true;
837
838 // retrieve parameters in case a new minimum has been found
839 if (fMinuit->fCstatu == "SUCCESSFUL")
841
842 double errParab = 0;
843 double gcor = 0;
844 // what returns if parameter fixed or constant or at limit ?
845 fMinuit->mnerrs(i,errUp,errLow, errParab, gcor);
846
847 // do not flag errors case of PROBLEMS or UNCHANGED (
848 return isValid;
849
850}
851
853 // reset TMinuit
854
855 fMinuit->mncler();
856
857 //reset the internal Minuit random generator to its initial state
858 double val = 3;
859 int inseed = 12345;
860 fMinuit->mnrn15(val,inseed);
861
862 fUsed = false;
863 fgUsed = false;
864
865}
866
868 // check if a parameter is defined and in case it was fixed released
869 // TMinuit is not able to release free parameters by redefining them
870 // so we need to force the release
871 if (fMinuit == 0) return;
872 if (fMinuit->GetNumFixedPars() == 0) return;
873 // check if parameter has already been defined
874 if (int(ivar) >= fMinuit->GetNumPars() ) return;
875
876 // check if parameter is fixed
877 for (int i = 0; i < fMinuit->fNpfix; ++i) {
878 if (fMinuit->fIpfix[i] == ivar+1 ) {
879 // parameter is fixed
880 fMinuit->Release(ivar);
881 return;
882 }
883 }
884
885}
886
887
889 // print-out results using classic Minuit format (mnprin)
890 if (fMinuit == 0) return;
891
892 // print minimizer result
893 if (PrintLevel() > 2)
895 else
897}
898
900 // suppress Minuit2 warnings
901 double arglist = 0;
902 int ierr = 0;
903 if (nowarn)
904 fMinuit->mnexcm("SET NOW",&arglist,0,ierr);
905 else
906 fMinuit->mnexcm("SET WAR",&arglist,0,ierr);
907}
908
909
910bool TMinuitMinimizer::Contour(unsigned int ipar, unsigned int jpar, unsigned int &npoints, double * x, double * y) {
911 // contour plot for parameter i and j
912 // need a valid FunctionMinimum otherwise exits
913 if (fMinuit == 0) {
914 Error("TMinuitMinimizer::Contour"," invalid TMinuit instance");
915 return false;
916 }
917
918 // set error and print level
919 double arglist[1];
920 int ierr = 0;
921 arglist[0] = ErrorDef();
922 fMinuit->mnexcm("SET Err",arglist,1,ierr);
923
924 arglist[0] = PrintLevel()-1;
925 fMinuit->mnexcm("SET PRINT",arglist,1,ierr);
926
927 // suppress warning in case Printlevel() == 0
928 if (PrintLevel() == 0) fMinuit->mnexcm("SET NOW",arglist,0,ierr);
929
930 // set precision if needed
931 if (Precision() > 0) {
932 arglist[0] = Precision();
933 fMinuit->mnexcm("SET EPS",arglist,1,ierr);
934 }
935
936
937 if (npoints < 4) {
938 Error("TMinuitMinimizer::Contour","Cannot make contour with so few points");
939 return false;
940 }
941 int npfound = 0;
942 // parameter numbers in mncont start from zero
943 fMinuit->mncont( ipar,jpar,npoints, x, y,npfound);
944 if (npfound<4) {
945 // mncont did go wrong
946 Error("TMinuitMinimizer::Contour","Cannot find more than 4 points");
947 return false;
948 }
949 if (npfound!=(int)npoints) {
950 // mncont did go wrong
951 Warning("TMinuitMinimizer::Contour","Returning only %d points ",npfound);
952 npoints = npfound;
953 }
954 return true;
955
956}
957
958bool TMinuitMinimizer::Scan(unsigned int ipar, unsigned int & nstep, double * x, double * y, double xmin, double xmax) {
959 // scan a parameter (variable) around the minimum value
960 // the parameters must have been set before
961 // if xmin=0 && xmax == 0 by default scan around 2 sigma of the error
962 // if the errors are also zero then scan from min and max of parameter range
963 // (if parameters are limited Minuit scan from min and max instead of 2 sigma by default)
964 // (force in that case to use errors)
965
966 // scan is not implemented for TMinuit, the way to return the array is only via the graph
967 if (fMinuit == 0) {
968 Error("TMinuitMinimizer::Scan"," invalid TMinuit instance");
969 return false;
970 }
971
972 // case of default xmin and xmax
973 if (xmin >= xmax && (int) ipar < fMinuit->GetNumPars() ) {
974 double val = 0; double err = 0;
976 double xlow = 0; double xup = 0 ;
977 int iuint = 0;
978 // in mnpout index starts from ze
979 fMinuit->mnpout( ipar, name, val, err, xlow, xup, iuint);
980 // redefine 2 sigma for all parameters by default (TMinuit does 1 sigma and range if limited)
981 if (iuint > 0 && err > 0) {
982 xmin = val - 2.*err;
983 xmax = val + 2 * err;
984 }
985 }
986
987 double arglist[4];
988 int ierr = 0;
989
990 arglist[0] = PrintLevel()-1;
991 fMinuit->mnexcm("SET PRINT",arglist,1,ierr);
992 // suppress warning in case Printlevel() == 0
993 if (PrintLevel() == 0) fMinuit->mnexcm("SET NOW",arglist,0,ierr);
994
995 // set precision if needed
996 if (Precision() > 0) {
997 arglist[0] = Precision();
998 fMinuit->mnexcm("SET EPS",arglist,1,ierr);
999 }
1000
1001 if (nstep == 0) return false;
1002 arglist[0] = ipar+1; // TMinuit starts from 1
1003 arglist[1] = nstep; //+2; // TMinuit deletes two points
1004 int nargs = 2;
1005 if (xmax > xmin ) {
1006 arglist[2] = xmin;
1007 arglist[3] = xmax;
1008 nargs = 4;
1009 }
1010 fMinuit->mnexcm("SCAN",arglist,nargs,ierr);
1011 if (ierr) {
1012 Error("TMinuitMinimizer::Scan"," Error executing command SCAN");
1013 return false;
1014 }
1015 // get TGraph object
1016 TGraph * gr = dynamic_cast<TGraph *>(fMinuit->GetPlot() );
1017 if (!gr) {
1018 Error("TMinuitMinimizer::Scan"," Error in returned graph object");
1019 return false;
1020 }
1021 nstep = std::min(gr->GetN(), (int) nstep);
1022
1023
1024 std::copy(gr->GetX(), gr->GetX()+nstep, x);
1025 std::copy(gr->GetY(), gr->GetY()+nstep, y);
1026 nstep = gr->GetN();
1027 return true;
1028}
1029
1031 // perform calculation of Hessian
1032
1033 if (fMinuit == 0) {
1034 Error("TMinuitMinimizer::Hesse","invalid TMinuit pointer. Need to call first SetFunction and SetVariable");
1035 return false;
1036 }
1037
1038
1039 double arglist[10];
1040 int ierr = 0;
1041
1042 // set error and print level
1043 arglist[0] = ErrorDef();
1044 fMinuit->mnexcm("SET ERR",arglist,1,ierr);
1045
1046 int printlevel = PrintLevel();
1047 arglist[0] = printlevel - 1;
1048 fMinuit->mnexcm("SET PRINT",arglist,1,ierr);
1049
1050 // suppress warning in case Printlevel() == 0
1051 if (printlevel == 0) fMinuit->mnexcm("SET NOW",arglist,0,ierr);
1052
1053 // set precision if needed
1054 if (Precision() > 0) {
1055 arglist[0] = Precision();
1056 fMinuit->mnexcm("SET EPS",arglist,1,ierr);
1057 }
1058
1059 arglist[0] = MaxFunctionCalls();
1060
1061 fMinuit->mnexcm("HESSE",arglist,1,ierr);
1062 fStatus += 100*ierr;
1063
1064 if (ierr != 0) return false;
1065
1066 // retrieve results (parameter and error matrix)
1067 // only if result is OK
1068
1071
1072 return true;
1073}
1074
1076 // set debug mode
1077
1078 if (fMinuit == 0) {
1079 Error("TMinuitMinimizer::SetDebug","invalid TMinuit pointer. Need to call first SetFunction and SetVariable");
1080 return false;
1081 }
1082 int ierr = 0;
1083 double arglist[1];
1084 arglist[0] = 1;
1085 if (on)
1086 fMinuit->mnexcm("SET DEBUG",arglist,1,ierr);
1087 else
1088 fMinuit->mnexcm("SET NODEBUG",arglist,1,ierr);
1089
1090 return (ierr == 0);
1091}
1092// } // end namespace Fit
1093
1094// } // end namespace ROOT
#define f(i)
Definition RSha256.hxx:104
#define g(i)
Definition RSha256.hxx:105
PrintLevel
Definition RooMinuit.h:6
Strategy
Definition RooMinuit.h:5
#define ClassImp(name)
Definition Rtypes.h:377
#define R__ASSERT(e)
Definition TError.h:117
void Info(const char *location, const char *msgfmt,...)
Use this function for informational messages.
Definition TError.cxx:230
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:197
void Warning(const char *location, const char *msgfmt,...)
Use this function in warning situations.
Definition TError.cxx:241
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:405
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:62
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:343
virtual void Gradient(const T *x, T *grad) const
Evaluate all the vector of function derivatives (gradient) at a point x.
Definition IFunction.h:358
virtual unsigned int NDim() const=0
Retrieve the dimension of the function.
double Tolerance() const
absolute tolerance
Definition Minimizer.h:421
unsigned int MaxFunctionCalls() const
max number of function calls
Definition Minimizer.h:415
double Precision() const
precision of minimizer in the evaluation of the objective function ( a value <=0 corresponds to the l...
Definition Minimizer.h:425
int fStatus
status of minimizer
Definition Minimizer.h:497
int Strategy() const
strategy
Definition Minimizer.h:428
double ErrorDef() const
return the statistical scale used for calculate the error is typically 1 for Chi2 and 0....
Definition Minimizer.h:438
bool IsValidError() const
return true if Minimizer has performed a detailed error validation (e.g. run Hesse for Minuit)
Definition Minimizer.h:441
int PrintLevel() const
minimizer configuration parameters
Definition Minimizer.h:412
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:137
Int_t GetN() const
Definition TGraph.h:129
Double_t * GetX() const
Definition TGraph.h:136
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
TMinuitMinimizer & operator=(const TMinuitMinimizer &rhs)
Assignment operator.
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 successfull
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 ...
Implementation in C++ of the Minuit package written by Fred James.
Definition TMinuit.h:27
virtual Int_t GetParameter(Int_t parNo, Double_t &currentValue, Double_t &currentError) const
return parameter value and error
Definition TMinuit.cxx:840
virtual Int_t FixParameter(Int_t parNo)
fix a parameter
Definition TMinuit.cxx:826
virtual Int_t GetNumPars() const
returns the total number of parameters that have been defined as fixed or free.
Definition TMinuit.cxx:871
virtual Int_t GetNumFixedPars() const
returns the number of currently fixed parameters
Definition TMinuit.cxx:854
virtual Int_t Release(Int_t parNo)
release a parameter
Definition TMinuit.cxx:893
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:1102
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:919
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:1394
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:6248
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:2500
virtual void mnrn15(Double_t &val, Int_t &inseed)
This is a super-portable random number generator.
Definition TMinuit.cxx:6618
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:2577
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:2663
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:694
virtual void mnprin(Int_t inkode, Double_t fval)
Prints the values of the parameters at the time of the call.
Definition TMinuit.cxx:6305
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:380
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