Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RFitPanelModel.cxx
Go to the documentation of this file.
1// Authors: Sergey Linev <S.Linev@gsi.de> Iliana Betsou <Iliana.Betsou@cern.ch>
2// Date: 2019-04-11
3// Warning: This is part of the ROOT 7 prototype! It will change without notice. It might trigger earthquakes. Feedback is welcome!
4
5/*************************************************************************
6 * Copyright (C) 1995-2019, Rene Brun and Fons Rademakers. *
7 * All rights reserved. *
8 * *
9 * For the licensing terms see $ROOTSYS/LICENSE. *
10 * For the list of contributors see $ROOTSYS/README/CREDITS. *
11 *************************************************************************/
12
14
15#include <ROOT/RLogger.hxx>
16
17#include "TH1.h"
18#include "TPluginManager.h"
19#include "TFitResult.h"
20#include "TF1.h"
21
23 static RLogChannel sLog("ROOT.FitPanel");
24 return sLog;
25}
26
27
34
35 kFP_MCHIS, kFP_MBINL, kFP_MUBIN // Fit methods
36
37};
38
39using namespace std::string_literals;
40
41
42/////////////////////////////////////////////////////////
43/// Clear list of patameters
44
46{
47 haspars = false;
48 id.clear();
49 name.clear();
50 pars.clear();
51}
52
53/////////////////////////////////////////////////////////
54/// Extract TF1 parameters, used in editor on client sides
55
57{
58 pars.clear();
59 haspars = true;
60 name = func->GetName();
61
62 for (int n = 0; n < func->GetNpar(); ++n) {
63 pars.emplace_back(n, func->GetParName(n));
64 auto &par = pars.back();
65
66 par.value = std::to_string(func->GetParameter(n));
67 par.error = std::to_string(func->GetParError(n));
68 double min, max;
69 func->GetParLimits(n, min, max);
70 par.min = std::to_string(min);
71 par.max = std::to_string(max);
72
73 if ((min >= max) && ((min != 0) || (max != 0)))
74 par.fixed = true;
75 }
76}
77
78/////////////////////////////////////////////////////////
79/// Set parameters to TF1
80
82{
83 if (func->GetNpar() != (int) pars.size()) {
84 ::Error("RFitFuncParsList::SetParameters", "Mismatch in parameters numbers");
85 return;
86 }
87
88 for (int n = 0; n < func->GetNpar(); ++n) {
89 if (pars[n].name.compare(func->GetParName(n)) != 0) {
90 ::Error("RFitFuncParsList::SetParameters", "Mismatch in parameter %d name %s %s", n, pars[n].name.c_str(), func->GetParName(n));
91 return;
92 }
93
94 func->SetParameter(n, std::stod(pars[n].value));
95 func->SetParError(n, std::stod(pars[n].error));
96 if (pars[n].fixed) {
97 func->FixParameter(n, std::stod(pars[n].value));
98 } else {
99 func->ReleaseParameter(n);
100 double min = std::stod(pars[n].min);
101 double max = std::stod(pars[n].max);
102 if (min < max)
103 func->SetParLimits(n, min, max);
104 }
105 }
106}
107
108/////////////////////////////////////////////////////////
109/// Update range values
110
112{
113 fDim = hist ? hist->GetDimension() : 0;
114
115 fMinRangeX = 0.;
116 fMaxRangeX = 100.;
117 fMinRangeY = 0.;
118 fMaxRangeY = 100.;
119
120 if (hist && (fDim > 0)) {
121 fMinRangeX = hist->GetXaxis()->GetXmin();
122 fMaxRangeX = hist->GetXaxis()->GetXmax();
123 }
124 if (hist && (fDim > 1)) {
125 fMinRangeY = hist->GetYaxis()->GetXmin();
126 fMaxRangeY = hist->GetYaxis()->GetXmax();
127 }
128
129 // defined values
130 fStepX = (fMaxRangeX - fMinRangeX) / 100;
131 fRangeX[0] = fMinRangeX;
132 fRangeX[1] = fMaxRangeX;
133
134 fStepY = (fMaxRangeY - fMinRangeY) / 100;
135 fRangeY[0] = fMinRangeY;
136 fRangeY[1] = fMaxRangeY;
137}
138
139/////////////////////////////////////////////////////////
140/// Check if function id is exists
141
143{
144 if (id.empty())
145 return false;
146
147 for (auto &item : fFuncList)
148 if (item.id == id)
149 return true;
150
151 return false;
152}
153
154/////////////////////////////////////////////////////////
155/// Select function
156
158{
159 if (HasFunction(id))
161 else
162 fSelectedFunc.clear();
163
165 if (func) {
166 fFuncPars.id = id;
168 }
169}
170
171/////////////////////////////////////////////////////////
172/// Initialize model - set some meaningful default values
173
175{
176 // build list of available histograms, as id use name from gdir
177 fSelectedData = "";
178
179 // build list of available functions
180
181 // Sub ComboBox for Type Function
182 fSelectedFunc = "";
183 fDim = 1;
184
185 fSelectedTab = "General";
186
187 // corresponds when Type == User Func (fSelectedTypeID == 1)
188
189 // ComboBox for General Tab --- Method
190 fFitMethods.clear();
191 fFitMethod = 0;
192
193 fLinearFit = false;
194 fRobust = false;
195 fRobustLevel = 0.95;
196
197 fIntegral = false;
198 fAllWeights1 = false;
199 fAddToList = false;
200 fEmptyBins1 = false;
201 fUseGradient = false;
202
203 fSame = false;
204 fNoDrawing = false;
205 fNoStoreDraw = false;
206
207 // Minimization method
208 fLibrary = 0;
209 // corresponds to library == 0
210 fMethodMinAll = {
211 {0, kFP_MIGRAD, "MIGRAD"}, {0, kFP_SIMPLX, "SIMPLEX"}, {0, kFP_SCAN, "SCAN"}, {0, kFP_COMBINATION, "Combination"},
212 {1, kFP_MIGRAD, "MIGRAD"}, {1, kFP_SIMPLX, "SIMPLEX"}, {1, kFP_FUMILI2, "FUMILI"}, {1, kFP_SCAN, "SCAN"}, {1, kFP_COMBINATION, "Combination"},
213 {2, kFP_FUMILI, "FUMILI"},
214
215 {3, kFP_GSLFR, "Fletcher-Reeves conjugate gradient"},
216 {3, kFP_GSLPR, "Polak-Ribiere conjugate gradient"},
217 {3, kFP_BFGS, "BFGS conjugate gradient"},
218 {3, kFP_BFGS2, "BFGS conjugate gradient (Version 2)"},
219 {3, kFP_GSLLM, "Levenberg-Marquardt"},
220 {3, kFP_GSLSA, "Simulated Annealing"}
221 };
222
223 fHasGenetics = false;
224 if ( gPluginMgr->FindHandler("ROOT::Math::Minimizer","GAlibMin") ) {
225 fMethodMinAll.emplace_back(4, kFP_GALIB, "GA Lib Genetic Algorithm");
226 fHasGenetics = true;
227 }
228 if (gPluginMgr->FindHandler("ROOT::Math::Minimizer","Genetic")) {
229 fMethodMinAll.emplace_back(4, kFP_TMVAGA, "TMVA Genetic Algorithm");
230 fHasGenetics = true;
231 }
232
234
235 fPrint = 0;
236
237 fAdvancedTab = "Contour";
238
239 fContourPoints = 40;
240 fScanPoints = 40;
241 fConfidenceLevel = 0.683;
242 fContourColor = "#c0b6ac"; // #11
243 fScanColor = "#0000ff"; // kBlue
244 fConfidenceColor = "#c0b6ac"; // #11
245}
246
247
248////////////////////////////////////////////////////////////
249/// Update setting dependent from object type
250
252{
253 fDataType = kind;
254
255 fFitMethods.clear();
256 fFitMethod = 0;
257 fRobust = false;
258
259 switch (kind) {
260 case kObjectHisto:
261 fFitMethods = {{kFP_MCHIS, "Chi-square"}, {kFP_MBINL, "Binned Likelihood"}};
263 break;
264
265 case kObjectGraph:
266 fFitMethods = {{kFP_MCHIS, "Chi-square"}};
268 fRobust = true;
269 break;
270
272 fFitMethods = {{kFP_MCHIS, "Chi-square"}};
274 fRobust = true;
275 break;
276
277 case kObjectGraph2D:
278 fFitMethods = {{kFP_MCHIS, "Chi-square"}};
280 break;
281
282 case kObjectHStack:
283 fFitMethods = {{kFP_MCHIS, "Chi-square"}};
285 break;
286
287 //case kObjectTree:
288 // fFitMethods = {{ kFP_MUBIN, "Unbinned Likelihood" }};
289 // fFitMethod = kFP_MUBIN;
290 // break;
291
292 default:
293 break;
294 }
295
296}
297
298
299////////////////////////////////////////////////////////////
300/// Update advanced parameters associated with fit function
301
303{
304 fAdvancedPars.clear();
305
306 fHasAdvanced = (res!=nullptr);
307
308 auto checkid = [&](std::string &id, const std::string &dflt) {
309 if (!res) { id.clear(); return; }
310 for (auto &item : fAdvancedPars)
311 if (item.key == id) return;
312 id = dflt;
313 };
314
315 if (res)
316 for (unsigned n = 0; n < res->NPar(); ++n)
317 fAdvancedPars.emplace_back(std::to_string(n), res->ParName(n));
318
319 checkid(fContourPar1Id, "0");
320 checkid(fContourPar2Id, "1");
321 checkid(fScanId, "0");
322}
323
324
326{
328
329 if (fDim > 0)
330 drange.AddRange(0, fRangeX[0], fRangeX[1]);
331
332 if ( fDim > 1 )
333 drange.AddRange(1, fRangeY[0], fRangeY[1]);
334
335 return drange;
336}
337
338/////////////////////////////////////////////////////////
339/// Provide initialized Foption_t instance
340
342{
343 Foption_t fitOpts;
344 fitOpts.Range = fUseRange;
345 fitOpts.Integral = fIntegral;
346 fitOpts.More = fImproveFitResults;
347 fitOpts.Errors = fBestErrors;
348 fitOpts.Like = fFitMethod != kFP_MCHIS;
349
350 if (fEmptyBins1)
351 fitOpts.W1 = 2;
352 else if (fAllWeights1)
353 fitOpts.W1 = 1;
354
355 // TODO: fEnteredFunc->GetText();
356 TString tmpStr = ""; // fEnteredFunc->GetText();
357 if ( !fLinearFit && (tmpStr.Contains("pol") || tmpStr.Contains("++")) )
358 fitOpts.Minuit = 1;
359
360 // TODO: fChangedParams
361 bool fChangedParams = false;
362 if (fChangedParams) {
363 fitOpts.Bound = 1;
364 fChangedParams = false; // reset
365 }
366
367 //fitOpts.Nochisq = (fNoChi2->GetState() == kButtonDown);
368 fitOpts.Nostore = fNoStoreDraw;
369 fitOpts.Nograph = fNoDrawing;
370 fitOpts.Plus = false; // TODO: (fAdd2FuncList->GetState() == kButtonDown);
371 fitOpts.Gradient = fUseGradient;
372 fitOpts.Quiet = fPrint == 2;
373 fitOpts.Verbose = fPrint == 1;
374
375 if (fRobust) {
376 fitOpts.Robust = 1;
377 fitOpts.hRobust = fRobustLevel;
378 }
379
380 return fitOpts;
381}
382
383/////////////////////////////////////////////////////////
384/// Provide initialized ROOT::Math::MinimizerOptions instance
385
387{
389
390 switch (fLibrary) {
391 case 0: minOpts.SetMinimizerType("Minuit"); break;
392 case 1: minOpts.SetMinimizerType("Minuit2"); break;
393 case 2: minOpts.SetMinimizerType("Fumili"); break;
394 case 3: minOpts.SetMinimizerType("GSLMultiMin"); break;
395 case 4: minOpts.SetMinimizerType("Geneti2c"); break;
396 }
397
398 switch(fSelectMethodMin) {
399 case kFP_MIGRAD: minOpts.SetMinimizerAlgorithm( "Migrad" ); break;
400 case kFP_FUMILI: minOpts.SetMinimizerAlgorithm("Fumili"); break;
401 case kFP_FUMILI2: minOpts.SetMinimizerAlgorithm("Fumili2"); break;
402 case kFP_SIMPLX: minOpts.SetMinimizerAlgorithm("Simplex"); break;
403 case kFP_SCAN: minOpts.SetMinimizerAlgorithm("Scan"); break;
404 case kFP_COMBINATION: minOpts.SetMinimizerAlgorithm("Minimize"); break;
405 case kFP_GSLFR: minOpts.SetMinimizerAlgorithm("conjugatefr"); break;
406 case kFP_GSLPR: minOpts.SetMinimizerAlgorithm("conjugatepr"); break;
407 case kFP_BFGS: minOpts.SetMinimizerAlgorithm("bfgs"); break;
408 case kFP_BFGS2: minOpts.SetMinimizerAlgorithm("bfgs2"); break;
409
410 case kFP_GSLLM:
411 minOpts.SetMinimizerType("GSLMultiFit");
412 minOpts.SetMinimizerAlgorithm("");
413 break;
414 case kFP_GSLSA:
415 minOpts.SetMinimizerType("GSLSimAn");
416 minOpts.SetMinimizerAlgorithm("");
417 break;
418 case kFP_TMVAGA:
419 minOpts.SetMinimizerType("Geneti2c");
420 minOpts.SetMinimizerAlgorithm("");
421 break;
422 case kFP_GALIB:
423 minOpts.SetMinimizerType("GAlibMin");
424 minOpts.SetMinimizerAlgorithm("");
425 break;
426 default: minOpts.SetMinimizerAlgorithm(""); break;
427 }
428
429 minOpts.SetErrorDef(fErrorDef);
433
434 return minOpts;
435}
436
437/////////////////////////////////////////////////////////
438/// Retrun draw option - dummy now
439
441{
442 TString res;
443 return res;
444}
@ kFP_MCHIS
@ kFP_NONE
@ kFP_GSLSA
@ kFP_FUMILI2
@ kFP_FUMILI
@ kFP_SIMPLX
@ kFP_MBINL
@ kFP_GSLLM
@ kFP_BFGS
@ kFP_MUBIN
@ kFP_GALIB
@ kFP_GSLPR
@ kFP_SCAN
@ kFP_COMBINATION
@ kFP_GSLFR
@ kFP_BFGS2
@ kFP_MIGRAD
@ kFP_TMVAGA
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:197
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize id
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
char name[80]
Definition TGX11.cxx:110
R__EXTERN TPluginManager * gPluginMgr
A log configuration for a channel, e.g.
Definition RLogger.hxx:101
class describing the range in the coordinates it supports multiple range in a coordinate.
Definition DataRange.h:35
void AddRange(unsigned int icoord, double xmin, double xmax)
add a range [xmin,xmax] for the new coordinate icoord Adding a range does not delete existing one,...
Definition DataRange.cxx:94
unsigned int NPar() const
total number of parameters (abbreviation)
Definition FitResult.h:122
std::string ParName(unsigned int i) const
name of the parameter
void SetMaxFunctionCalls(unsigned int maxfcn)
set maximum of function calls
void SetMaxIterations(unsigned int maxiter)
set maximum iterations (one iteration can have many function calls)
void SetMinimizerType(const char *type)
set minimizer type
void SetErrorDef(double err)
set error def
void SetMinimizerAlgorithm(const char *type)
set minimizer algorithm
void SetTolerance(double tol)
set the tolerance
Double_t GetXmax() const
Definition TAxis.h:135
Double_t GetXmin() const
Definition TAxis.h:134
1-Dim function class
Definition TF1.h:213
virtual void ReleaseParameter(Int_t ipar)
Release parameter number ipar during a fit operation.
Definition TF1.cxx:3153
virtual void SetParError(Int_t ipar, Double_t error)
Set error for parameter number ipar.
Definition TF1.cxx:3471
virtual void GetParLimits(Int_t ipar, Double_t &parmin, Double_t &parmax) const
Return limits for parameter ipar.
Definition TF1.cxx:1941
virtual Double_t GetParError(Int_t ipar) const
Return value of parameter number ipar.
Definition TF1.cxx:1931
virtual Int_t GetNpar() const
Definition TF1.h:486
virtual const char * GetParName(Int_t ipar) const
Definition TF1.h:534
virtual void SetParLimits(Int_t ipar, Double_t parmin, Double_t parmax)
Set lower and upper limits for parameter ipar.
Definition TF1.cxx:3499
virtual void SetParameter(Int_t param, Double_t value)
Definition TF1.h:639
virtual void FixParameter(Int_t ipar, Double_t value)
Fix the value of a parameter for a fit operation The specified value will be used in the fit and the ...
Definition TF1.cxx:1558
virtual Double_t GetParameter(Int_t ipar) const
Definition TF1.h:517
Extends the ROOT::Fit::Result class with a TNamed inheritance providing easy possibility for I/O.
Definition TFitResult.h:34
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:58
virtual Int_t GetDimension() const
Definition TH1.h:281
TAxis * GetXaxis()
Definition TH1.h:322
TAxis * GetYaxis()
Definition TH1.h:323
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
TPluginHandler * FindHandler(const char *base, const char *uri=nullptr)
Returns the handler if there exists a handler for the specified URI.
Basic string class.
Definition TString.h:139
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:636
EFitPanel
Definition CommonDefs.h:7
@ kFP_MCHIS
Definition CommonDefs.h:18
@ kFP_GSLSA
Definition CommonDefs.h:24
@ kFP_FUMILI
Definition CommonDefs.h:23
@ kFP_SIMPLX
Definition CommonDefs.h:22
@ kFP_MBINL
Definition CommonDefs.h:18
@ kFP_GSLLM
Definition CommonDefs.h:24
@ kFP_BFGS
Definition CommonDefs.h:24
@ kFP_GALIB
Definition CommonDefs.h:25
@ kFP_GSLPR
Definition CommonDefs.h:24
@ kFP_SCAN
Definition CommonDefs.h:25
@ kFP_COMBINATION
Definition CommonDefs.h:23
@ kFP_GSLFR
Definition CommonDefs.h:24
@ kFP_BFGS2
Definition CommonDefs.h:24
@ kFP_MIGRAD
Definition CommonDefs.h:22
@ kFP_TMVAGA
Definition CommonDefs.h:25
const Int_t n
Definition legend1.C:16
RLogChannel & FitPanelLog()
Log channel for FitPanel diagnostics.
int Range
Definition Foption.h:39
int Nograph
Definition Foption.h:42
int Quiet
Definition Foption.h:29
int Like
Definition Foption.h:34
int W1
Definition Foption.h:36
int Gradient
Definition Foption.h:40
int Robust
Definition Foption.h:48
double hRobust
Definition Foption.h:51
int Plus
Definition Foption.h:43
int Integral
Definition Foption.h:44
int Bound
Definition Foption.h:31
int Nostore
Definition Foption.h:41
int More
Definition Foption.h:38
int Minuit
Definition Foption.h:46
int Errors
Definition Foption.h:37
int Verbose
Definition Foption.h:30
void GetParameters(TF1 *f1)
Extract TF1 parameters, used in editor on client sides.
std::string id
function id in the FitPanel
void SetParameters(TF1 *f1)
Set parameters to TF1.
std::vector< RItemInfo > fFuncList
all available fit functions
EFitObjectType fDataType
selected object type, provided by server
RFuncParsList fFuncPars
Parameters.
int fLibrary
selected minimization library
TString GetDrawOption()
Retrun draw option - dummy now.
ROOT::Math::MinimizerOptions GetMinimizerOptions()
Provide initialized ROOT::Math::MinimizerOptions instance.
void Initialize()
Initialize model - set some meaningful default values.
std::vector< RComboBoxItem > fAdvancedPars
std::string fSelectedData
selected data
bool HasFunction(const std::string &id)
Check if function id is exists.
std::string fSelectedFunc
id of selected fit function like dflt::gaus
void UpdateAdvanced(TFitResult *res)
Update advanced parameters associated with fit function.
void SelectedFunc(const std::string &name, TF1 *func)
Select function.
std::string fSelectedTab
key of selected tab, useful for drawing
void UpdateRange(TH1 *hist)
Update range values.
void SetObjectKind(EFitObjectType kind)
Update setting dependent from object type.
Foption_t GetFitOptions()
Provide initialized Foption_t instance.
std::vector< RMinimezerAlgorithm > fMethodMinAll
all items for all methods
std::string fConfidenceColor
Confidence sub-tab.
int fDim
number of dimensions in selected data object
std::vector< RMethodInfo > fFitMethods
all supported for selected data
bool fHasGenetics
is genetics available