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
27
28
39
40using namespace std::string_literals;
41
42
43/////////////////////////////////////////////////////////
44/// Clear list of patameters
45
47{
48 haspars = false;
49 id.clear();
50 name.clear();
51 pars.clear();
52}
53
54/////////////////////////////////////////////////////////
55/// Extract TF1 parameters, used in editor on client sides
56
58{
59 pars.clear();
60 haspars = true;
61 name = func->GetName();
62
63 for (int n = 0; n < func->GetNpar(); ++n) {
64 pars.emplace_back(n, func->GetParName(n));
65 auto &par = pars.back();
66
67 par.value = std::to_string(func->GetParameter(n));
68 par.error = std::to_string(func->GetParError(n));
69 double min, max;
70 func->GetParLimits(n, min, max);
71 par.min = std::to_string(min);
72 par.max = std::to_string(max);
73
74 if ((min >= max) && ((min != 0) || (max != 0)))
75 par.fixed = true;
76 }
77}
78
79/////////////////////////////////////////////////////////
80/// Set parameters to TF1
81
83{
84 if (func->GetNpar() != (int) pars.size()) {
85 ::Error("RFitFuncParsList::SetParameters", "Mismatch in parameters numbers");
86 return;
87 }
88
89 for (int n = 0; n < func->GetNpar(); ++n) {
90 if (pars[n].name.compare(func->GetParName(n)) != 0) {
91 ::Error("RFitFuncParsList::SetParameters", "Mismatch in parameter %d name %s %s", n, pars[n].name.c_str(), func->GetParName(n));
92 return;
93 }
94
95 func->SetParameter(n, std::stod(pars[n].value));
96 func->SetParError(n, std::stod(pars[n].error));
97 if (pars[n].fixed) {
98 func->FixParameter(n, std::stod(pars[n].value));
99 } else {
100 func->ReleaseParameter(n);
101 double min = std::stod(pars[n].min);
102 double max = std::stod(pars[n].max);
103 if (min < max)
104 func->SetParLimits(n, min, max);
105 }
106 }
107}
108
109/////////////////////////////////////////////////////////
110/// Update range values
111
113{
114 fDim = hist ? hist->GetDimension() : 0;
115
116 fMinRangeX = 0.;
117 fMaxRangeX = 100.;
118 fMinRangeY = 0.;
119 fMaxRangeY = 100.;
120
121 if (hist && (fDim > 0)) {
122 fMinRangeX = hist->GetXaxis()->GetXmin();
123 fMaxRangeX = hist->GetXaxis()->GetXmax();
124 }
125 if (hist && (fDim > 1)) {
126 fMinRangeY = hist->GetYaxis()->GetXmin();
127 fMaxRangeY = hist->GetYaxis()->GetXmax();
128 }
129
130 // defined values
131 fStepX = (fMaxRangeX - fMinRangeX) / 100;
132 fRangeX[0] = fMinRangeX;
133 fRangeX[1] = fMaxRangeX;
134
135 fStepY = (fMaxRangeY - fMinRangeY) / 100;
136 fRangeY[0] = fMinRangeY;
137 fRangeY[1] = fMaxRangeY;
138}
139
140/////////////////////////////////////////////////////////
141/// Check if function id is exists
142
144{
145 if (id.empty())
146 return false;
147
148 for (auto &item : fFuncList)
149 if (item.id == id)
150 return true;
151
152 return false;
153}
154
155/////////////////////////////////////////////////////////
156/// Select function
157
159{
160 if (HasFunction(id))
162 else
163 fSelectedFunc.clear();
164
166 if (func) {
167 fFuncPars.id = id;
169 }
170}
171
172/////////////////////////////////////////////////////////
173/// Initialize model - set some meaningful default values
174
176{
177 // build list of available histograms, as id use name from gdir
178 fSelectedData = "";
179
180 // build list of available functions
181
182 // Sub ComboBox for Type Function
183 fSelectedFunc = "";
184 fDim = 1;
185
186 fSelectedTab = "General";
187
188 // corresponds when Type == User Func (fSelectedTypeID == 1)
189
190 // ComboBox for General Tab --- Method
191 fFitMethods.clear();
192 fFitMethod = 0;
193
194 fLinearFit = false;
195 fRobust = false;
196 fRobustLevel = 0.95;
197
198 fIntegral = false;
199 fAllWeights1 = false;
200 fAddToList = false;
201 fEmptyBins1 = false;
202 fUseGradient = false;
203
204 fSame = false;
205 fNoDrawing = false;
206 fNoStoreDraw = false;
207
208 // Minimization method
209 fLibrary = 0;
210 // corresponds to library == 0
211 fMethodMinAll = {
212 {0, kFP_MIGRAD, "MIGRAD"}, {0, kFP_SIMPLX, "SIMPLEX"}, {0, kFP_SCAN, "SCAN"}, {0, kFP_COMBINATION, "Combination"},
213 {1, kFP_MIGRAD, "MIGRAD"}, {1, kFP_SIMPLX, "SIMPLEX"}, {1, kFP_FUMILI2, "FUMILI"}, {1, kFP_SCAN, "SCAN"}, {1, kFP_COMBINATION, "Combination"},
214 {2, kFP_FUMILI, "FUMILI"},
215
216 {3, kFP_GSLFR, "Fletcher-Reeves conjugate gradient"},
217 {3, kFP_GSLPR, "Polak-Ribiere conjugate gradient"},
218 {3, kFP_BFGS, "BFGS conjugate gradient"},
219 {3, kFP_BFGS2, "BFGS conjugate gradient (Version 2)"},
220 {3, kFP_GSLLM, "Levenberg-Marquardt"},
221 {3, kFP_GSLSA, "Simulated Annealing"}
222 };
223
224 fHasGenetics = false;
225 if ( gPluginMgr->FindHandler("ROOT::Math::Minimizer","GAlibMin") ) {
226 fMethodMinAll.emplace_back(4, kFP_GALIB, "GA Lib Genetic Algorithm");
227 fHasGenetics = true;
228 }
229 if (gPluginMgr->FindHandler("ROOT::Math::Minimizer","Genetic")) {
230 fMethodMinAll.emplace_back(4, kFP_TMVAGA, "TMVA Genetic Algorithm");
231 fHasGenetics = true;
232 }
233
235
236 fPrint = 0;
237
238 fAdvancedTab = "Contour";
239
240 fContourPoints = 40;
241 fScanPoints = 40;
242 fConfidenceLevel = 0.683;
243 fContourColor = "#c0b6ac"; // #11
244 fScanColor = "#0000ff"; // kBlue
245 fConfidenceColor = "#c0b6ac"; // #11
246}
247
248
249////////////////////////////////////////////////////////////
250/// Update setting dependent from object type
251
253{
254 fDataType = kind;
255
256 fFitMethods.clear();
257 fFitMethod = 0;
258 fRobust = false;
259
260 switch (kind) {
261 case kObjectHisto:
262 fFitMethods = {{kFP_MCHIS, "Chi-square"}, {kFP_MBINL, "Binned Likelihood"}};
264 break;
265
266 case kObjectGraph:
267 fFitMethods = {{kFP_MCHIS, "Chi-square"}};
269 fRobust = true;
270 break;
271
273 fFitMethods = {{kFP_MCHIS, "Chi-square"}};
275 fRobust = true;
276 break;
277
278 case kObjectGraph2D:
279 fFitMethods = {{kFP_MCHIS, "Chi-square"}};
281 break;
282
283 case kObjectHStack:
284 fFitMethods = {{kFP_MCHIS, "Chi-square"}};
286 break;
287
288 //case kObjectTree:
289 // fFitMethods = {{ kFP_MUBIN, "Unbinned Likelihood" }};
290 // fFitMethod = kFP_MUBIN;
291 // break;
292
293 default:
294 break;
295 }
296
297}
298
299
300////////////////////////////////////////////////////////////
301/// Update advanced parameters associated with fit function
302
304{
305 fAdvancedPars.clear();
306
307 fHasAdvanced = (res!=nullptr);
308
309 auto checkid = [&](std::string &id, const std::string &dflt) {
310 if (!res) { id.clear(); return; }
311 for (auto &item : fAdvancedPars)
312 if (item.key == id) return;
313 id = dflt;
314 };
315
316 if (res)
317 for (unsigned n = 0; n < res->NPar(); ++n)
318 fAdvancedPars.emplace_back(std::to_string(n), res->ParName(n));
319
322 checkid(fScanId, "0");
323}
324
325
327{
329
330 if (fDim > 0)
331 drange.AddRange(0, fRangeX[0], fRangeX[1]);
332
333 if ( fDim > 1 )
334 drange.AddRange(1, fRangeY[0], fRangeY[1]);
335
336 return drange;
337}
338
339/////////////////////////////////////////////////////////
340/// Provide initialized Foption_t instance
341
343{
345 fitOpts.Range = fUseRange;
346 fitOpts.Integral = fIntegral;
348 fitOpts.Errors = fBestErrors;
349 fitOpts.Like = fFitMethod != kFP_MCHIS;
350
351 if (fEmptyBins1)
352 fitOpts.W1 = 2;
353 else if (fAllWeights1)
354 fitOpts.W1 = 1;
355
356 // TODO: fEnteredFunc->GetText();
357 TString tmpStr = ""; // fEnteredFunc->GetText();
358 if ( !fLinearFit && (tmpStr.Contains("pol") || tmpStr.Contains("++")) )
359 fitOpts.Minuit = 1;
360
361 // TODO: fChangedParams
362 bool fChangedParams = false;
363 if (fChangedParams) {
364 fitOpts.Bound = 1;
365 fChangedParams = false; // reset
366 }
367
368 //fitOpts.Nochisq = (fNoChi2->GetState() == kButtonDown);
369 fitOpts.Nostore = fNoStoreDraw;
370 fitOpts.Nograph = fNoDrawing;
371 fitOpts.Plus = false; // TODO: (fAdd2FuncList->GetState() == kButtonDown);
372 fitOpts.Gradient = fUseGradient;
373 fitOpts.Quiet = fPrint == 2;
374 fitOpts.Verbose = fPrint == 1;
375
376 if (fRobust) {
377 fitOpts.Robust = 1;
378 fitOpts.hRobust = fRobustLevel;
379 }
380
381 return fitOpts;
382}
383
384/////////////////////////////////////////////////////////
385/// Provide initialized ROOT::Math::MinimizerOptions instance
386
388{
390
391 switch (fLibrary) {
392 case 0: minOpts.SetMinimizerType("Minuit"); break;
393 case 1: minOpts.SetMinimizerType("Minuit2"); break;
394 case 2: minOpts.SetMinimizerType("Fumili"); break;
395 case 3: minOpts.SetMinimizerType("GSLMultiMin"); break;
396 case 4: minOpts.SetMinimizerType("Geneti2c"); break;
397 }
398
399 switch(fSelectMethodMin) {
400 case kFP_MIGRAD: minOpts.SetMinimizerAlgorithm( "Migrad" ); break;
401 case kFP_FUMILI: minOpts.SetMinimizerAlgorithm("Fumili"); break;
402 case kFP_FUMILI2: minOpts.SetMinimizerAlgorithm("Fumili2"); break;
403 case kFP_SIMPLX: minOpts.SetMinimizerAlgorithm("Simplex"); break;
404 case kFP_SCAN: minOpts.SetMinimizerAlgorithm("Scan"); break;
405 case kFP_COMBINATION: minOpts.SetMinimizerAlgorithm("Minimize"); break;
406 case kFP_GSLFR: minOpts.SetMinimizerAlgorithm("conjugatefr"); break;
407 case kFP_GSLPR: minOpts.SetMinimizerAlgorithm("conjugatepr"); break;
408 case kFP_BFGS: minOpts.SetMinimizerAlgorithm("bfgs"); break;
409 case kFP_BFGS2: minOpts.SetMinimizerAlgorithm("bfgs2"); break;
410
411 case kFP_GSLLM:
412 minOpts.SetMinimizerType("GSLMultiFit");
413 minOpts.SetMinimizerAlgorithm("");
414 break;
415 case kFP_GSLSA:
416 minOpts.SetMinimizerType("GSLSimAn");
417 minOpts.SetMinimizerAlgorithm("");
418 break;
419 case kFP_TMVAGA:
420 minOpts.SetMinimizerType("Geneti2c");
421 minOpts.SetMinimizerAlgorithm("");
422 break;
423 case kFP_GALIB:
424 minOpts.SetMinimizerType("GAlibMin");
425 minOpts.SetMinimizerAlgorithm("");
426 break;
427 default: minOpts.SetMinimizerAlgorithm(""); break;
428 }
429
430 minOpts.SetErrorDef(fErrorDef);
431 minOpts.SetTolerance(fMaxTolerance);
432 minOpts.SetMaxIterations(fMaxIterations);
433 minOpts.SetMaxFunctionCalls(fMaxIterations);
434
435 return minOpts;
436}
437
438/////////////////////////////////////////////////////////
439/// Retrun draw option - dummy now
440
@ 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
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:185
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
unsigned int NPar() const
total number of parameters (abbreviation)
Definition FitResult.h:122
std::string ParName(unsigned int i) const
name of the parameter
Double_t GetXmax() const
Definition TAxis.h:140
Double_t GetXmin() const
Definition TAxis.h:139
1-Dim function class
Definition TF1.h:233
virtual void ReleaseParameter(Int_t ipar)
Release parameter number ipar during a fit operation.
Definition TF1.cxx:3151
virtual void SetParError(Int_t ipar, Double_t error)
Set error for parameter number ipar.
Definition TF1.cxx:3479
virtual void GetParLimits(Int_t ipar, Double_t &parmin, Double_t &parmax) const
Return limits for parameter ipar.
Definition TF1.cxx:1942
virtual Double_t GetParError(Int_t ipar) const
Return value of parameter number ipar.
Definition TF1.cxx:1932
virtual Int_t GetNpar() const
Definition TF1.h:507
virtual const char * GetParName(Int_t ipar) const
Definition TF1.h:555
virtual void SetParLimits(Int_t ipar, Double_t parmin, Double_t parmax)
Set lower and upper limits for parameter ipar.
Definition TF1.cxx:3507
virtual void SetParameter(Int_t param, Double_t value)
Definition TF1.h:660
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:1559
virtual Double_t GetParameter(Int_t ipar) const
Definition TF1.h:538
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:59
virtual Int_t GetDimension() const
Definition TH1.h:283
TAxis * GetXaxis()
Definition TH1.h:324
TAxis * GetYaxis()
Definition TH1.h:325
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
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.
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