70 #if (__GNUC__==3&&__GNUC_MINOR__==2&&__GNUC_PATCHLEVEL__==3)
71 char* operator+( streampos&,
char* );
79 TVirtualFitter *RooMinuit::_theFitter = 0 ;
87 void RooMinuit::cleanup()
109 RooMinuit::RooMinuit(RooAbsReal&
function)
111 RooSentinel::activate() ;
121 _handleLocalErrors = kTRUE ;
123 _printEvalErrors = 10 ;
126 _doEvalErrorWall = kTRUE ;
129 RooArgSet* paramSet =
function.getParameters(RooArgSet()) ;
130 RooArgList paramList(*paramSet) ;
133 _floatParamList = (RooArgList*) paramList.selectByAttrib(
"Constant",kFALSE) ;
134 if (_floatParamList->getSize()>1) {
135 _floatParamList->sort() ;
137 _floatParamList->setName(
"floatParamList") ;
139 _constParamList = (RooArgList*) paramList.selectByAttrib(
"Constant",kTRUE) ;
140 if (_constParamList->getSize()>1) {
141 _constParamList->sort() ;
143 _constParamList->setName(
"constParamList") ;
146 TIterator* pIter = _floatParamList->createIterator() ;
148 while((arg=(RooAbsArg*)pIter->Next())) {
149 if (!arg->IsA()->InheritsFrom(RooAbsRealLValue::Class())) {
150 coutW(Minimization) <<
"RooMinuit::RooMinuit: removing parameter " << arg->GetName()
151 <<
" from list because it is not of type RooRealVar" << endl ;
152 _floatParamList->remove(*arg) ;
155 _nPar = _floatParamList->getSize() ;
161 _initFloatParamList = (RooArgList*) _floatParamList->snapshot(kFALSE) ;
162 _initConstParamList = (RooArgList*) _constParamList->snapshot(kFALSE) ;
165 Int_t nPar= _floatParamList->getSize() + _constParamList->getSize() ;
166 if (_theFitter)
delete _theFitter ;
167 _theFitter =
new TFitter(nPar*2+1) ;
168 _theFitter->SetObjectFit(
this) ;
175 _theFitter->SetFCN(RooMinuitGlue);
178 setErrorLevel(
function.defaultErrorLevel()) ;
181 synchronize(kFALSE) ;
188 if (RooMsgService::instance().silentMode()) {
202 RooMinuit::~RooMinuit()
204 delete _floatParamList ;
205 delete _initFloatParamList ;
206 delete _constParamList ;
207 delete _initConstParamList ;
221 void RooMinuit::setStrategy(Int_t istrat)
223 Double_t stratArg(istrat) ;
224 _theFitter->ExecuteCommand(
"SET STR",&stratArg,1) ;
235 void RooMinuit::setErrorLevel(Double_t level)
237 _theFitter->ExecuteCommand(
"SET ERR",&level,1);
245 void RooMinuit::setEps(Double_t eps)
247 _theFitter->ExecuteCommand(
"SET EPS",&eps,1) ;
255 void RooMinuit::setOffsetting(Bool_t flag)
257 _func->enableOffsetting(flag) ;
273 RooFitResult* RooMinuit::fit(
const char* options)
275 if (_floatParamList->getSize()==0) {
279 _theFitter->SetObjectFit(
this) ;
281 TString opts(options) ;
285 if (opts.Contains(
"v")) setVerbose(1) ;
286 if (opts.Contains(
"t")) setProfile(1) ;
287 if (opts.Contains(
"l")) setLogFile(Form(
"%s.log",_func->GetName())) ;
288 if (opts.Contains(
"c")) optimizeConst(1) ;
291 if (opts.Contains(
"s")) hesse() ;
292 if (opts.Contains(
"0")) setStrategy(0) ;
294 if (opts.Contains(
"0")) setStrategy(1) ;
295 if (opts.Contains(
"h")||!opts.Contains(
"m")) hesse() ;
296 if (!opts.Contains(
"m")) minos() ;
298 return (opts.Contains(
"r")) ? save() : 0 ;
309 Int_t RooMinuit::migrad()
311 if (_floatParamList->getSize()==0) {
315 _theFitter->SetObjectFit(
this) ;
318 arglist[0]= _maxEvalMult*_nPar;
321 synchronize(_verbose) ;
323 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::CollectErrors) ;
324 RooAbsReal::clearEvalErrorLog() ;
325 _status= _theFitter->ExecuteCommand(
"MIGRAD",arglist,2);
326 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::PrintErrors) ;
330 saveStatus(
"MIGRAD",_status) ;
343 Int_t RooMinuit::hesse()
345 if (_floatParamList->getSize()==0) {
349 _theFitter->SetObjectFit(
this) ;
352 arglist[0]= _maxEvalMult*_nPar;
354 synchronize(_verbose) ;
356 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::CollectErrors) ;
357 RooAbsReal::clearEvalErrorLog() ;
358 _status= _theFitter->ExecuteCommand(
"HESSE",arglist,1);
359 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::PrintErrors) ;
363 saveStatus(
"HESSE",_status) ;
376 Int_t RooMinuit::minos()
378 if (_floatParamList->getSize()==0) {
382 _theFitter->SetObjectFit(
this) ;
385 arglist[0]= _maxEvalMult*_nPar;
387 synchronize(_verbose) ;
389 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::CollectErrors) ;
390 RooAbsReal::clearEvalErrorLog() ;
391 _status= _theFitter->ExecuteCommand(
"MINOS",arglist,1);
393 if (_status == 0 && gMinuit->fCstatu !=
"SUCCESSFUL") {
394 if (gMinuit->fCstatu ==
"FAILURE" ||
395 gMinuit->fCstatu ==
"PROBLEMS") _status = 5;
399 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::PrintErrors) ;
403 saveStatus(
"MINOS",_status) ;
416 Int_t RooMinuit::minos(
const RooArgSet& minosParamList)
418 if (_floatParamList->getSize()==0) {
422 _theFitter->SetObjectFit(
this) ;
425 Double_t* arglist =
new Double_t[_nPar+1];
427 if (minosParamList.getSize()>0) {
428 TIterator* aIter = minosParamList.createIterator() ;
430 while((arg=(RooAbsArg*)aIter->Next())) {
431 RooAbsArg* par = _floatParamList->find(arg->GetName());
432 if (par && !par->isConstant()) {
433 Int_t index = _floatParamList->index(par);
435 arglist[nMinosPar]=index+1;
440 arglist[0]= _maxEvalMult*_nPar;
442 synchronize(_verbose) ;
444 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::CollectErrors) ;
445 RooAbsReal::clearEvalErrorLog() ;
446 _status= _theFitter->ExecuteCommand(
"MINOS",arglist,1+nMinosPar);
448 if (_status == 0 && gMinuit->fCstatu !=
"SUCCESSFUL") {
449 if (gMinuit->fCstatu ==
"FAILURE" ||
450 gMinuit->fCstatu ==
"PROBLEMS") _status = 5;
453 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::PrintErrors) ;
459 saveStatus(
"MINOS",_status) ;
472 Int_t RooMinuit::seek()
474 if (_floatParamList->getSize()==0) {
478 _theFitter->SetObjectFit(
this) ;
481 arglist[0]= _maxEvalMult*_nPar;
483 synchronize(_verbose) ;
485 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::CollectErrors) ;
486 RooAbsReal::clearEvalErrorLog() ;
487 _status= _theFitter->ExecuteCommand(
"SEEK",arglist,1);
488 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::PrintErrors) ;
492 saveStatus(
"SEEK",_status) ;
505 Int_t RooMinuit::simplex()
507 if (_floatParamList->getSize()==0) {
511 _theFitter->SetObjectFit(
this) ;
514 arglist[0]= _maxEvalMult*_nPar;
517 synchronize(_verbose) ;
519 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::CollectErrors) ;
520 RooAbsReal::clearEvalErrorLog() ;
521 _status= _theFitter->ExecuteCommand(
"SIMPLEX",arglist,2);
522 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::PrintErrors) ;
526 saveStatus(
"SIMPLEX",_status) ;
539 Int_t RooMinuit::improve()
541 if (_floatParamList->getSize()==0) {
545 _theFitter->SetObjectFit(
this) ;
548 arglist[0]= _maxEvalMult*_nPar;
550 synchronize(_verbose) ;
552 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::CollectErrors) ;
553 RooAbsReal::clearEvalErrorLog() ;
554 _status= _theFitter->ExecuteCommand(
"IMPROVE",arglist,1);
555 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::PrintErrors) ;
559 saveStatus(
"IMPROVE",_status) ;
569 Int_t RooMinuit::setPrintLevel(Int_t newLevel)
571 Int_t ret = _printLevel ;
572 Double_t arg(newLevel) ;
573 _theFitter->ExecuteCommand(
"SET PRINT",&arg,1);
574 _printLevel = newLevel ;
583 void RooMinuit::setNoWarn()
586 _theFitter->ExecuteCommand(
"SET NOWARNINGS",&arg,1);
595 Int_t RooMinuit::setWarnLevel(Int_t newLevel)
597 if (newLevel==_warnLevel) {
601 Int_t ret = _warnLevel ;
602 Double_t arg(newLevel) ;
605 _theFitter->ExecuteCommand(
"SET WARNINGS",&arg,1);
608 _theFitter->ExecuteCommand(
"SET NOWARNINGS",&arg2,1);
610 _warnLevel = newLevel ;
621 Bool_t RooMinuit::synchronize(Bool_t verbose)
623 Int_t oldPrint = setPrintLevel(-1) ;
624 gMinuit->fNwrmes[0] = 0;
625 Int_t oldWarn = setWarnLevel(-1) ;
627 Bool_t constValChange(kFALSE) ;
628 Bool_t constStatChange(kFALSE) ;
633 for(index= 0; index < _constParamList->getSize() ; index++) {
634 RooRealVar *par=
dynamic_cast<RooRealVar*
>(_constParamList->at(index)) ;
637 RooRealVar *oldpar=
dynamic_cast<RooRealVar*
>(_initConstParamList->at(index)) ;
638 if (!oldpar) continue ;
641 if (!par->isConstant()) {
644 _constParamList->remove(*par) ;
645 _floatParamList->add(*par) ;
646 _initFloatParamList->addClone(*oldpar) ;
647 _initConstParamList->remove(*oldpar) ;
648 constStatChange=kTRUE ;
652 coutI(Minimization) <<
"RooMinuit::synchronize: parameter " << par->GetName() <<
" is now floating." << endl ;
657 if (par->getVal()!= oldpar->getVal()) {
658 constValChange=kTRUE ;
660 coutI(Minimization) <<
"RooMinuit::synchronize: value of constant parameter " << par->GetName()
661 <<
" changed from " << oldpar->getVal() <<
" to " << par->getVal() << endl ;
668 *_initConstParamList = *_constParamList ;
672 for(index= 0; index < _nPar; index++) {
673 RooRealVar *par=
dynamic_cast<RooRealVar*
>(_floatParamList->at(index)) ;
680 if(!par->isConstant()) {
683 if (!par->IsA()->InheritsFrom(RooRealVar::Class())) {
684 coutW(Minimization) <<
"RooMinuit::fit: Error, non-constant parameter " << par->GetName()
685 <<
" is not of type RooRealVar, skipping" << endl ;
690 if (par->hasMin() && par->hasMax()) {
691 pmin = par->getMin();
692 pmax = par->getMax();
696 pstep= par->getError();
699 if (par->hasMin() && par->hasMax()) {
700 pstep= 0.1*(pmax-pmin);
703 if (pmax - par->getVal() < 2*pstep) {
704 pstep = (pmax - par->getVal())/2 ;
705 }
else if (par->getVal() - pmin < 2*pstep) {
706 pstep = (par->getVal() - pmin )/2 ;
711 pstep= 0.1*(pmax-pmin);
718 coutW(Minimization) <<
"RooMinuit::synchronize: WARNING: no initial error estimate available for "
719 << par->GetName() <<
": using " << pstep << endl;
723 pmin = par->getVal() ;
724 pmax = par->getVal() ;
728 Double_t oldVar,oldVerr,oldVlo,oldVhi ;
729 char oldParname[100] ;
730 Int_t ierr = _theFitter->GetParameter(index,oldParname,oldVar,oldVerr,oldVlo,oldVhi) ;
735 Bool_t oldFixed(kFALSE) ;
737 for (ix = 1; ix <= gMinuit->fNpfix; ++ix) {
738 if (gMinuit->fIpfix[ix-1] == index+1) oldFixed=kTRUE ;
742 if (par->isConstant() && !oldFixed) {
745 if (oldVar!=par->getVal()) {
746 Double_t arglist[2] ;
747 arglist[0] = index+1 ;
748 arglist[1] = par->getVal() ;
749 _theFitter->ExecuteCommand(
"SET PAR",arglist,2) ;
751 coutI(Minimization) <<
"RooMinuit::synchronize: value of parameter " << par->GetName() <<
" changed from " << oldVar <<
" to " << par->getVal() << endl ;
755 _theFitter->FixParameter(index) ;
756 constStatChange=kTRUE ;
758 coutI(Minimization) <<
"RooMinuit::synchronize: parameter " << par->GetName() <<
" is now fixed." << endl ;
761 }
else if (par->isConstant() && oldFixed) {
764 if (oldVar!=par->getVal()) {
765 Double_t arglist[2] ;
766 arglist[0] = index+1 ;
767 arglist[1] = par->getVal() ;
768 _theFitter->ExecuteCommand(
"SET PAR",arglist,2) ;
769 constValChange=kTRUE ;
772 coutI(Minimization) <<
"RooMinuit::synchronize: value of fixed parameter " << par->GetName() <<
" changed from " << oldVar <<
" to " << par->getVal() << endl ;
778 if (!par->isConstant() && oldFixed) {
779 _theFitter->ReleaseParameter(index) ;
780 constStatChange=kTRUE ;
783 coutI(Minimization) <<
"RooMinuit::synchronize: parameter " << par->GetName() <<
" is now floating." << endl ;
788 if (oldVar!=par->getVal() || oldVlo!=pmin || oldVhi != pmax || oldVerr!=pstep) {
789 _theFitter->SetParameter(index, par->GetName(), par->getVal(), pstep, pmin, pmax);
793 if (verbose && ierr>=0) {
796 if (oldVar!=par->getVal()) {
797 coutI(Minimization) <<
"RooMinuit::synchronize: value of parameter " << par->GetName() <<
" changed from " << oldVar <<
" to " << par->getVal() << endl ;
799 if (oldVlo!=pmin || oldVhi!=pmax) {
800 coutI(Minimization) <<
"RooMinuit::synchronize: limits of parameter " << par->GetName() <<
" changed from [" << oldVlo <<
"," << oldVhi
801 <<
"] to [" << pmin <<
"," << pmax <<
"]" << endl ;
805 if (oldVerr!=pstep && oldVerr!=0) {
806 coutI(Minimization) <<
"RooMinuit::synchronize: error/step size of parameter " << par->GetName() <<
" changed from " << oldVerr <<
" to " << pstep << endl ;
813 gMinuit->fNwrmes[0] = 0;
814 oldWarn = setWarnLevel(oldWarn) ;
815 oldPrint = setPrintLevel(oldPrint) ;
818 if (constStatChange) {
820 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::CollectErrors) ;
822 coutI(Minimization) <<
"RooMinuit::synchronize: set of constant parameters changed, rerunning const optimizer" << endl ;
823 _func->constOptimizeTestStatistic(RooAbsArg::ConfigChange) ;
824 }
else if (constValChange) {
825 coutI(Minimization) <<
"RooMinuit::synchronize: constant parameter values changed, rerunning const optimizer" << endl ;
826 _func->constOptimizeTestStatistic(RooAbsArg::ValueChange) ;
829 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::PrintErrors) ;
845 void RooMinuit::optimizeConst(Int_t flag)
847 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::CollectErrors) ;
849 if (_optConst && !flag){
850 if (_printLevel>-1) coutI(Minimization) <<
"RooMinuit::optimizeConst: deactivating const optimization" << endl ;
851 _func->constOptimizeTestStatistic(RooAbsArg::DeActivate,flag>1) ;
853 }
else if (!_optConst && flag) {
854 if (_printLevel>-1) coutI(Minimization) <<
"RooMinuit::optimizeConst: activating const optimization" << endl ;
855 _func->constOptimizeTestStatistic(RooAbsArg::Activate,flag>1) ;
857 }
else if (_optConst && flag) {
858 if (_printLevel>-1) coutI(Minimization) <<
"RooMinuit::optimizeConst: const optimization already active" << endl ;
860 if (_printLevel>-1) coutI(Minimization) <<
"RooMinuit::optimizeConst: const optimization wasn't active" << endl ;
863 RooAbsReal::setEvalErrorLoggingMode(RooAbsReal::PrintErrors) ;
877 RooFitResult* RooMinuit::save(
const char* userName,
const char* userTitle)
880 name = userName ? userName : Form(
"%s", _func->GetName()) ;
881 title = userTitle ? userTitle : Form(
"%s", _func->GetTitle()) ;
883 if (_floatParamList->getSize()==0) {
884 RooFitResult* fitRes =
new RooFitResult(name,title) ;
885 fitRes->setConstParList(*_constParamList) ;
886 fitRes->setInitParList(RooArgList()) ;
887 fitRes->setFinalParList(RooArgList()) ;
888 fitRes->setStatus(-999) ;
889 fitRes->setCovQual(-999) ;
890 fitRes->setMinNLL(_func->getVal()) ;
891 fitRes->setNumInvalidNLL(0) ;
892 fitRes->setEDM(-999) ;
896 RooFitResult* fitRes =
new RooFitResult(name,title) ;
900 RooArgList saveConstList(*_constParamList) ;
901 RooArgList saveFloatInitList(*_initFloatParamList) ;
902 RooArgList saveFloatFinalList(*_floatParamList) ;
903 for (i=0 ; i<_floatParamList->getSize() ; i++) {
904 RooAbsArg* par = _floatParamList->at(i) ;
905 if (par->isConstant()) {
906 saveFloatInitList.remove(*saveFloatInitList.find(par->GetName()),kTRUE) ;
907 saveFloatFinalList.remove(*par) ;
908 saveConstList.add(*par) ;
911 saveConstList.sort() ;
913 fitRes->setConstParList(saveConstList) ;
914 fitRes->setInitParList(saveFloatInitList) ;
916 Double_t edm, errdef, minVal;
918 Int_t icode = _theFitter->GetStats(minVal, edm, errdef, nvpar, nparx);
919 fitRes->setStatus(_status) ;
920 fitRes->setCovQual(icode) ;
921 fitRes->setMinNLL(minVal) ;
922 fitRes->setNumInvalidNLL(_numBadNLL) ;
923 fitRes->setEDM(edm) ;
924 fitRes->setFinalParList(saveFloatFinalList) ;
926 fitRes->fillCorrMatrix() ;
928 fitRes->setCovarianceMatrix(*_extV) ;
931 fitRes->setStatusHistory(_statusHistory) ;
943 RooPlot* RooMinuit::contour(RooRealVar& var1, RooRealVar& var2, Double_t n1, Double_t n2, Double_t n3, Double_t n4, Double_t n5, Double_t n6)
946 _theFitter->SetObjectFit(
this) ;
948 RooArgList* paramSave = (RooArgList*) _floatParamList->snapshot() ;
951 Int_t index1= _floatParamList->index(&var1);
953 coutE(Minimization) <<
"RooMinuit::contour(" << GetName()
954 <<
") ERROR: " << var1.GetName() <<
" is not a floating parameter of " << _func->GetName() << endl ;
958 Int_t index2= _floatParamList->index(&var2);
960 coutE(Minimization) <<
"RooMinuit::contour(" << GetName()
961 <<
") ERROR: " << var2.GetName() <<
" is not a floating parameter of PDF " << _func->GetName() << endl ;
966 RooPlot* frame =
new RooPlot(var1,var2) ;
969 TMarker *point=
new TMarker(var1.getVal(), var2.getVal(), 8);
970 frame->addObject(point) ;
973 Double_t errdef= gMinuit->fUp;
976 n[0] = n1 ; n[1] = n2 ; n[2] = n3 ; n[3] = n4 ; n[4] = n5 ; n[5] = n6 ;
979 for (Int_t ic = 0 ; ic<6 ; ic++) {
982 gMinuit->SetErrorDef(n[ic]*n[ic]*errdef);
984 TGraph* graph= (TGraph*)gMinuit->Contour(50, index1, index2);
986 coutE(Minimization) <<
"RooMinuit::contour(" << GetName() <<
") ERROR: MINUIT did not return a contour graph for n=" << n[ic] << endl ;
988 graph->SetName(Form(
"contour_%s_n%f",_func->GetName(),n[ic])) ;
989 graph->SetLineStyle(ic+1) ;
990 graph->SetLineWidth(2) ;
991 graph->SetLineColor(kBlue) ;
992 frame->addObject(graph,
"L") ;
998 gMinuit->SetErrorDef(errdef);
1001 *_floatParamList = *paramSave ;
1015 Bool_t RooMinuit::setLogFile(
const char* inLogfile)
1018 coutI(Minimization) <<
"RooMinuit::setLogFile: closing previous log file" << endl ;
1023 _logfile =
new ofstream(inLogfile) ;
1024 if (!_logfile->good()) {
1025 coutI(Minimization) <<
"RooMinuit::setLogFile: cannot open file " << inLogfile << endl ;
1038 Double_t RooMinuit::getPdfParamVal(Int_t index)
1040 return ((RooRealVar*)_floatParamList->at(index))->getVal() ;
1048 Double_t RooMinuit::getPdfParamErr(Int_t index)
1050 return ((RooRealVar*)_floatParamList->at(index))->getError() ;
1058 Bool_t RooMinuit::setPdfParamVal(Int_t index, Double_t value, Bool_t verbose)
1061 RooRealVar* par = (RooRealVar*)_floatParamVec[index] ;
1063 if (par->getVal()!=value) {
1064 if (verbose) cout << par->GetName() <<
"=" << value <<
", " ;
1065 par->setVal(value) ;
1077 void RooMinuit::setPdfParamErr(Int_t index, Double_t value)
1079 ((RooRealVar*)_floatParamList->at(index))->setError(value) ;
1087 void RooMinuit::clearPdfParamAsymErr(Int_t index)
1089 ((RooRealVar*)_floatParamList->at(index))->removeAsymError() ;
1096 void RooMinuit::setPdfParamErr(Int_t index, Double_t loVal, Double_t hiVal)
1098 ((RooRealVar*)_floatParamList->at(index))->setAsymError(loVal,hiVal) ;
1106 void RooMinuit::profileStart()
1110 _cumulTimer.Start(kFALSE) ;
1120 void RooMinuit::profileStop()
1124 _cumulTimer.Stop() ;
1125 coutI(Minimization) <<
"Command timer: " ; _timer.Print() ;
1126 coutI(Minimization) <<
"Session timer: " ; _cumulTimer.Print() ;
1137 void RooMinuit::backProp()
1139 Double_t val,err,vlo,vhi, eplus, eminus, eparab, globcc;
1142 for(index= 0; index < _nPar; index++) {
1143 _theFitter->GetParameter(index, buffer, val, err, vlo, vhi);
1144 setPdfParamVal(index, val);
1145 _theFitter->GetErrors(index, eplus, eminus, eparab, globcc);
1148 setPdfParamErr(index, err);
1150 if(eplus > 0 || eminus < 0) {
1152 setPdfParamErr(index, eminus,eplus);
1155 clearPdfParamAsymErr(index) ;
1163 void RooMinuit::updateFloatVec()
1165 _floatParamVec.clear() ;
1166 RooFIter iter = _floatParamList->fwdIterator() ;
1168 _floatParamVec.resize(_floatParamList->getSize()) ;
1170 while((arg=iter.next())) {
1171 _floatParamVec[i++] = arg ;
1182 void RooMinuit::applyCovarianceMatrix(TMatrixDSym& V)
1184 _extV = (TMatrixDSym*) V.Clone() ;
1186 for (Int_t i=0 ; i<getNPar() ; i++) {
1188 if (_floatParamList->at(i)->isConstant()) {
1191 RooMinuit* context = (RooMinuit*) RooMinuit::_theFitter->GetObjectFit() ;
1192 if (context && context->_verbose)
1193 cout <<
"setting parameter " << i <<
" error to " << sqrt((*_extV)(i,i)) << endl ;
1194 setPdfParamErr(i, sqrt((*_extV)(i,i))) ;
1202 void RooMinuitGlue(Int_t& , Double_t* ,
1203 Double_t &f, Double_t *par, Int_t )
1208 RooMinuit* context = (RooMinuit*) RooMinuit::_theFitter->GetObjectFit() ;
1209 ofstream* logf = context->logfile() ;
1210 Double_t& maxFCN = context->maxFCN() ;
1211 Bool_t verbose = context->_verbose ;
1214 Int_t nPar= context->getNPar();
1215 for(Int_t index= 0; index < nPar; index++) {
1216 if (logf) (*logf) << par[index] <<
" " ;
1217 context->setPdfParamVal(index, par[index],verbose);
1221 RooAbsReal::setHideOffset(kFALSE) ;
1222 f= context->_func->getVal() ;
1223 RooAbsReal::setHideOffset(kTRUE) ;
1224 context->_evalCounter++ ;
1225 if (RooAbsReal::numEvalErrors()>0 || f>1e30) {
1227 if (context->_printEvalErrors>=0) {
1229 if (context->_doEvalErrorWall) {
1230 oocoutW(context,Minimization) <<
"RooMinuitGlue: Minimized function has error status." << endl
1231 <<
"Returning maximum FCN so far (" << maxFCN
1232 <<
") to force MIGRAD to back out of this region. Error log follows" << endl ;
1234 oocoutW(context,Minimization) <<
"RooMinuitGlue: Minimized function has error status but is ignored" << endl ;
1237 TIterator* iter = context->_floatParamList->createIterator() ;
1239 Bool_t first(kTRUE) ;
1240 ooccoutW(context,Minimization) <<
"Parameter values: " ;
1241 while((var=(RooRealVar*)iter->Next())) {
1242 if (first) { first = kFALSE ; }
else ooccoutW(context,Minimization) <<
", " ;
1243 ooccoutW(context,Minimization) << var->GetName() <<
"=" << var->getVal() ;
1246 ooccoutW(context,Minimization) << endl ;
1248 RooAbsReal::printEvalErrors(ooccoutW(context,Minimization),context->_printEvalErrors) ;
1249 ooccoutW(context,Minimization) << endl ;
1252 if (context->_doEvalErrorWall) {
1256 RooAbsReal::clearEvalErrorLog() ;
1257 context->_numBadNLL++ ;
1258 }
else if (f>maxFCN) {
1263 if (logf) (*logf) << setprecision(15) << f << setprecision(4) << endl;
1265 cout <<
"\nprevFCN" << (context->_func->isOffsetting()?
"-offset":
"") <<
" = " << setprecision(10) << f << setprecision(4) <<
" " ;