Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RFitPanel.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
13#include <ROOT/RFitPanel.hxx>
14
15#include <ROOT/RLogger.hxx>
16
17#include "Fit/BinData.h"
18#include "Fit/Fitter.h"
19
20#include "TString.h"
21#include "TGraph.h"
22#include "TGraphErrors.h"
23#include "TGraph2D.h"
24#include "TGraph2DErrors.h"
25#include "TMultiGraph.h"
26#include "THStack.h"
27#include "TROOT.h"
28#include "TH1.h"
29#include "TH2.h"
30#include "TF1.h"
31#include "TF2.h"
32#include "TF3.h"
33#include "TFitResult.h"
34#include "TList.h"
35#include "TVirtualPad.h"
36#include "TCanvas.h"
37#include "TDirectory.h"
38#include "TBufferJSON.h"
39#include "Math/Minimizer.h"
40#include "HFitInterface.h"
41#include "TColor.h"
42
43#include <iostream>
44#include <iomanip>
45#include <memory>
46#include <sstream>
47
48using namespace std::string_literals;
49
50using namespace ROOT::Experimental;
51
52/** \class RFitPanel
53\ingroup webwidgets
54
55web-based FitPanel prototype.
56*/
57
58////////////////////////////////////////////////////////////////////////////////
59/// Constructor
60
61RFitPanel::FitRes::FitRes(const std::string &_objid, std::unique_ptr<TF1> &_func, TFitResultPtr &_res)
62 : objid(_objid), res(_res)
63{
64 std::swap(func, _func);
65}
66
67////////////////////////////////////////////////////////////////////////////////
68/// Destructor
69
71{
72 // to avoid dependency from TF1
73
74 // if TF1 object deleted before - prevent second delete
75 if (func && func->IsZombie())
76 func.release();
77
78}
79
80////////////////////////////////////////////////////////////////////////////////
81/// Constructor
82
83RFitPanel::RFitPanel(const std::string &title)
84{
85 model().fTitle = title;
86
88}
89
90////////////////////////////////////////////////////////////////////////////////
91/// Destructor
92
94{
95 // to avoid dependency from TF1
96}
97
98////////////////////////////////////////////////////////////////////////////////
99/// Returns RWebWindow instance, used to display FitPanel
100
101std::shared_ptr<ROOT::RWebWindow> RFitPanel::GetWindow()
102{
103 if (!fWindow) {
105
106 fWindow->SetPanelName("rootui5.fitpanel.view.FitPanel");
107
108 fWindow->SetCallBacks(
109 [this](unsigned connid)
110 {
111 fConnId = connid;
112 fWindow->Send(fConnId, "INITDONE");
113 if (!model().fInitialized)
114 SelectObject("$$$");
115 SendModel();
116 },
117 [this](unsigned connid, const std::string &arg) { ProcessData(connid, arg); },
118 [this](unsigned) { fConnId = 0; });
119
120 fWindow->SetGeometry(400, 650); // configure predefined geometry
121 }
122
123 return fWindow;
124}
125
126////////////////////////////////////////////////////////////////////////////////
127/// Update list of available data
128
130{
131 auto &m = model();
132
133 m.fDataSet.clear();
134
135 for (auto &obj : fObjects)
136 m.fDataSet.emplace_back("Panel", "panel::"s + obj->GetName(), Form("%s::%s", obj->ClassName(), obj->GetName()));
137
138 if (gDirectory) {
139 TIter iter(gDirectory->GetList());
140 TObject *obj = nullptr;
141 while ((obj = iter()) != nullptr) {
143 m.fDataSet.emplace_back("gDirectory", "gdir::"s + obj->GetName(), Form("%s::%s", obj->ClassName(), obj->GetName()));
144 }
145 }
146}
147
148////////////////////////////////////////////////////////////////////////////////
149/// Select object for fitting
150
151void RFitPanel::SelectObject(const std::string &objid)
152{
154
155 auto &m = model();
156
157 std::string id = objid;
158 if (id.compare("$$$") == 0) {
159 if (m.fDataSet.size() > 0)
160 id = m.fDataSet[0].id;
161 else
162 id.clear();
163 }
164
165 TObject *obj = GetSelectedObject(id);
166 auto kind = GetFitObjectType(obj);
167
168 m.SetObjectKind(kind);
169
170 TH1 *hist = nullptr;
171 switch (kind) {
173 hist = (TH1*)obj;
174 break;
175
177 hist = ((TGraph*)obj)->GetHistogram();
178 break;
179
181 hist = ((TMultiGraph*)obj)->GetHistogram();
182 break;
183
185 hist = ((TGraph2D*)obj)->GetHistogram("empty");
186 break;
187
189 hist = (TH1 *)((THStack *)obj)->GetHists()->First();
190 break;
191
192 //case RFitPanelModel::kObjectTree:
193 // m.fFitMethods = {{ kFP_MUBIN, "Unbinned Likelihood" }};
194 // m.fFitMethod = kFP_MUBIN;
195 // break;
196
197 default:
198 break;
199 }
200
201 if (!obj)
202 m.fSelectedData = "";
203 else
204 m.fSelectedData = id;
205
206 m.fInitialized = true;
207
208 // update list of data
209
210 m.UpdateRange(hist);
211
213
214 std::string selfunc = m.fSelectedFunc;
215
216 if (!m.HasFunction(selfunc)) {
217 if (m.fFuncList.size() > 0)
218 selfunc = m.fFuncList[0].id;
219 else
220 selfunc.clear();
221 }
222
224}
225
226////////////////////////////////////////////////////////////////////////////////
227/// Returns object based on it string id
228/// Searches either in gDirectory or in internal panel list
229
230TObject *RFitPanel::GetSelectedObject(const std::string &objid)
231{
232 if (objid.compare(0,6,"gdir::") == 0) {
233 std::string name = objid.substr(6);
234 if (gDirectory)
235 return gDirectory->GetList()->FindObject(name.c_str());
236 } else if (objid.compare(0,7,"panel::") == 0) {
237 std::string name = objid.substr(7);
238 for (auto &item : fObjects)
239 if (name.compare(item->GetName()) == 0)
240 return item;
241 }
242
243 return nullptr;
244}
245
246////////////////////////////////////////////////////////////////////////////////
247/// Returns kind of object
248
271
272////////////////////////////////////////////////////////////////////////////////
273/// Update list of available functions
274
276{
277 auto &m = model();
278
279 m.fFuncList.clear();
280
281 if (m.fDim == 1) {
282 m.fFuncList = { {"gaus"}, {"gausn"}, {"expo"}, {"landau"}, {"landaun"},
283 {"pol0"},{"pol1"},{"pol2"},{"pol3"},{"pol4"},{"pol5"},{"pol6"},{"pol7"},{"pol8"},{"pol9"},
284 {"cheb0"}, {"cheb1"}, {"cheb2"}, {"cheb3"}, {"cheb4"}, {"cheb5"}, {"cheb6"}, {"cheb7"}, {"cheb8"}, {"cheb9"} };
285 } else if (m.fDim == 2) {
286 m.fFuncList = { {"xygaus"}, {"bigaus"}, {"xyexpo"}, {"xylandau"}, {"xylandaun"} };
287 } else if (m.fDim == 3) {
288 m.fFuncList = {{"xyzgaus"}};
289 }
290
291 for (auto &func : fSystemFuncs) {
292 m.fFuncList.emplace_back("System", "system::"s + func->GetName(), func->GetName());
293 }
294
295 for (auto &entry : fPrevRes) {
296 if (entry.objid == m.fSelectedData)
297 m.fFuncList.emplace_back("Previous", "previous::"s + entry.func->GetName(), entry.func->GetName());
298 }
299}
300
301////////////////////////////////////////////////////////////////////////////////
302/// Select fit function
303
310
311////////////////////////////////////////////////////////////////////////////////
312/// Assign histogram to use with fit panel - without ownership
313
315{
316 fObjects.emplace_back(hist);
317 SelectObject("panel::"s + hist->GetName());
318 SendModel();
319}
320
321/// Assign histogram name to use with fit panel - it should be available in gDirectory
322
323void RFitPanel::AssignHistogram(const std::string &hname)
324{
325 SelectObject("gdir::" + hname);
326 SendModel();
327}
328
329////////////////////////////////////////////////////////////////////////////////
330/// assign canvas to use for drawing results of fitting or showing fitpanel itself
331
332void RFitPanel::AssignCanvas(std::shared_ptr<RCanvas> &canv)
333{
334 if (!fCanvas) {
335 fCanvas = canv;
336 } else {
337 R__LOG_ERROR(FitPanelLog()) << "FitPanel already bound to the canvas - change is not yet supported";
338 }
339}
340
341////////////////////////////////////////////////////////////////////////////////
342/// Show FitPanel
343
345{
347}
348
349////////////////////////////////////////////////////////////////////////////////
350/// Hide FitPanel
351
353{
354 if (fWindow)
355 fWindow->CloseConnections();
356}
357
358////////////////////////////////////////////////////////////////////////////////
359/// Return reference on model object
360/// Model created if was not exists before
361
362
364{
365 if (!fModel) {
366 fModel = std::make_unique<RFitPanelModel>();
367 fModel->Initialize();
368 }
369
370 return *fModel.get();
371}
372
373////////////////////////////////////////////////////////////////////////////////
374/// Send model object to the client
375
377{
378 if (fWindow && (fConnId > 0)) {
380 fWindow->Send(fConnId, "MODEL:"s + json.Data());
381 }
382}
383
384////////////////////////////////////////////////////////////////////////////////
385/// Process data from FitPanel
386/// OpenUI5-based FitPanel sends commands or status changes
387
388void RFitPanel::ProcessData(unsigned, const std::string &arg)
389{
390 if (arg == "RELOAD") {
391
393
396
397 SendModel();
398
399 } else if (arg.compare(0, 7, "UPDATE:") == 0) {
400
401 if (UpdateModel(arg.substr(7)) > 0)
402 SendModel();
403
404 } else if (arg.compare(0, 6, "DOFIT:") == 0) {
405
406 if (UpdateModel(arg.substr(6)) >= 0)
407 if (DoFit())
408 SendModel();
409
410 } else if (arg.compare(0, 7, "DODRAW:") == 0) {
411
412 if (UpdateModel(arg.substr(7)) >= 0)
413 if (DoDraw())
414 SendModel();
415
416 } else if (arg.compare(0, 8, "SETPARS:") == 0) {
417
418 auto info = TBufferJSON::FromJSON<RFitPanelModel::RFuncParsList>(arg.substr(8));
419
420 if (info) {
421 TF1 *func = FindFunction(info->id);
422
423 // copy all parameters back to the function
424 if (func)
425 info->SetParameters(func);
426 }
427
428 }
429}
430
431////////////////////////////////////////////////////////////////////////////////
432/// Search for existing functions, ownership still belongs to FitPanel or global lists
433
434TF1 *RFitPanel::FindFunction(const std::string &id)
435{
436 if (id.compare(0,8,"system::") == 0) {
437 std::string name = id.substr(8);
438
439 for (auto &item : fSystemFuncs)
440 if (name.compare(item->GetName()) == 0)
441 return item.get();
442 }
443
444 if (id.compare(0,10,"previous::") == 0) {
445 std::string name = id.substr(10);
446
447 for (auto &entry : fPrevRes)
448 if (name.compare(entry.func->GetName()) == 0)
449 return entry.func.get();
450 }
451
452 return nullptr;
453}
454
455////////////////////////////////////////////////////////////////////////////////
456/// Creates new instance to make fitting
457
459{
460 if (id.compare(0,10,"previous::") == 0) {
461 std::string name = id.substr(10);
462
463 for (auto &entry : fPrevRes)
464 if (name.compare(entry.func->GetName()) == 0)
465 return entry.res.Get();
466 }
467
468 return nullptr;
469}
470
471////////////////////////////////////////////////////////////////////////////////
472/// Creates new instance to make fitting
473
474std::unique_ptr<TF1> RFitPanel::GetFitFunction(const std::string &funcname)
475{
476 std::unique_ptr<TF1> res;
477
478 TF1 *func = FindFunction(funcname);
479
480 if (func) {
481 // Now we make a copy.
482 res.reset((TF1*) func->IsA()->New());
483 func->Copy(*res);
484 } else if (funcname.compare(0,6,"dflt::") == 0) {
485
486 std::string formula = funcname.substr(6);
487
489
490 double xmin, xmax, ymin, ymax, zmin, zmax;
491 drange.GetRange(xmin, xmax, ymin, ymax, zmin, zmax);
492
493 if ( model().fDim == 1 || model().fDim == 0 ) {
494 res.reset(new TF1(formula.c_str(), formula.c_str(), xmin, xmax));
495 } else if ( model().fDim == 2 ) {
496 res.reset(new TF2(formula.c_str(), formula.c_str(), xmin, xmax, ymin, ymax));
497 } else if ( model().fDim == 3 ) {
498 res.reset(new TF3(formula.c_str(), formula.c_str(), xmin, xmax, ymin, ymax, zmin, zmax));
499 }
500 }
501
502 return res;
503}
504
505////////////////////////////////////////////////////////////////////////////////
506/// Update fit model
507/// returns -1 if JSON fails
508/// return 0 if nothing large changed
509/// return 1 if important selection are changed and client need to be updated
510
511int RFitPanel::UpdateModel(const std::string &json)
512{
513 auto m = TBufferJSON::FromJSON<RFitPanelModel>(json);
514
515 if (!m) {
516 R__LOG_ERROR(FitPanelLog()) << "Fail to parse JSON for RFitPanelModel";
517 return -1;
518 }
519
520 m->fInitialized = true;
521
522 int res = 0; // nothing changed
523
524 if (model().fSelectedData != m->fSelectedData) {
525 res |= 1;
526 }
527
528 if (model().fSelectedFunc != m->fSelectedFunc) {
529 res |= 2;
530 }
531
532 std::swap(fModel, m); // finally replace model
533
534 if (res & 1)
535 SelectObject(model().fSelectedData);
536
537 if (res != 0)
538 SelectFunction(model().fSelectedFunc);
539
540 return res;
541}
542
543////////////////////////////////////////////////////////////////////////////////
544/// Copies f into a new TF1 to be stored in the fitpanel with it's
545/// own ownership. This is taken from Fit::StoreAndDrawFitFunction in
546/// HFitImpl.cxx
547
549{
550 double xmin = 0, xmax = 0, ymin = 0, ymax = 0, zmin = 0, zmax = 0;
551
552 // no need to use kNotGlobal bit. TF1::Copy does not add in the list by default
553 if ( dynamic_cast<TF3*>(f)) {
554 TF3* fnew = (TF3*)f->IsA()->New();
555 f->Copy(*fnew);
556 f->GetRange(xmin,ymin,zmin,xmax,ymax,zmax);
557 fnew->SetRange(xmin,ymin,zmin,xmax,ymax,zmax);
558 fnew->SetParent( nullptr );
559 fnew->AddToGlobalList(false);
560 return fnew;
561 } else if ( dynamic_cast<TF2*>(f) != 0 ) {
562 TF2* fnew = (TF2*)f->IsA()->New();
563 f->Copy(*fnew);
564 f->GetRange(xmin,ymin,xmax,ymax);
565 fnew->SetRange(xmin,ymin,xmax,ymax);
566 fnew->Save(xmin,xmax,ymin,ymax,0,0);
567 fnew->SetParent( nullptr );
568 fnew->AddToGlobalList(false);
569 return fnew;
570 }
571
572 TF1* fnew = (TF1*)f->IsA()->New();
573 f->Copy(*fnew);
574 f->GetRange(xmin,xmax);
575 fnew->SetRange(xmin,xmax);
576 // This next line is added, as fnew-Save fails with gausND! As
577 // the number of dimensions is unknown...
578 if ( '\0' != fnew->GetExpFormula()[0] )
579 fnew->Save(xmin,xmax,0,0,0,0);
580 fnew->SetParent( nullptr );
581 fnew->AddToGlobalList(false);
582 return fnew;
583}
584
585////////////////////////////////////////////////////////////////////////////////
586/// Looks for all the functions registered in the current ROOT
587/// session.
588
590{
591
592 fSystemFuncs.clear();
593
594 // Be carefull not to store functions that will be in the
595 // predefined section
596 std::vector<std::string> fnames = { "gaus" , "gausn", "expo", "landau",
597 "landaun", "pol0", "pol1", "pol2",
598 "pol3", "pol4", "pol5", "pol6",
599 "pol7", "pol8", "pol9", "user" };
600
601 // No go through all the objects registered in gROOT
602 TIter functionsIter(gROOT->GetListOfFunctions());
603 TObject* obj;
604 while( (obj = functionsIter()) != nullptr ) {
605 // And if they are TF1s
606 if ( TF1* func = dynamic_cast<TF1*>(obj) ) {
607 bool addFunction = true;
608 // And they are not already registered in fSystemFunc
609 for ( auto &name : fnames) {
610 if ( name.compare(func->GetName()) == 0 ) {
611 addFunction = false;
612 break;
613 }
614 }
615 // Add them.
616 if ( addFunction )
617 fSystemFuncs.emplace_back( copyTF1(func) );
618 }
619 }
620}
621
622////////////////////////////////////////////////////////////////////////////////
623/// Returns pad where histogram is drawn
624/// If canvas not exists, create new one
625
627{
628 if (!obj || (!force && (model().fNoDrawing || model().fNoStoreDraw)))
629 return nullptr;
630
631 std::function<TPad*(TPad*)> check = [&](TPad *pad) {
632 TPad *res = nullptr;
633 if (!pad) return res;
634 if (!fPadName.empty() && (fPadName.compare(pad->GetName()) == 0)) return pad;
635 TIter next(pad->GetListOfPrimitives());
636 TObject *prim = nullptr;
637 while (!res && (prim = next())) {
638 res = (prim == obj) ? pad : check(dynamic_cast<TPad *>(prim));
639 }
640 return res;
641 };
642
643 if (!fCanvName.empty()) {
644 auto drawcanv = dynamic_cast<TCanvas *> (gROOT->GetListOfCanvases()->FindObject(fCanvName.c_str()));
645 auto drawpad = check(drawcanv);
646 if (drawpad) {
647 drawpad->cd();
648 return drawpad;
649 }
650 if (drawcanv) {
651 drawcanv->Clear();
652 drawcanv->cd();
653 obj->Draw();
654 return drawcanv;
655 }
656 fCanvName.clear();
657 fPadName.clear();
658 }
659
660 TObject *c = nullptr;
661 TIter nextc(gROOT->GetListOfCanvases());
662 while ((c = nextc())) {
663 auto drawpad = check(dynamic_cast<TCanvas* >(c));
664 if (drawpad) {
665 drawpad->cd();
666 fCanvName = c->GetName();
667 fPadName = drawpad->GetName();
668 return drawpad;
669 }
670 }
671
672 auto canv = gROOT->MakeDefCanvas();
673 canv->SetName("fpc");
674 canv->SetTitle("Fit panel drawings");
675 fPadName = fCanvName = canv->GetName();
676
677 canv->cd();
678 obj->Draw();
679
680 return canv;
681}
682
683////////////////////////////////////////////////////////////////////////////////
684/// Perform fitting using current model settings
685/// Returns true if any action was done
686
688{
689 auto &m = model();
690
691 TObject *obj = GetSelectedObject(m.fSelectedData);
692 if (!obj) return false;
693 auto kind = GetFitObjectType(obj);
694
695 auto f1 = GetFitFunction(m.fSelectedFunc);
696 if (!f1) return false;
697
698 auto drange = m.GetRanges();
699 auto minOption = m.GetMinimizerOptions();
700 auto fitOpts = m.GetFitOptions();
701 auto drawOpts = m.GetDrawOption();
702
703 fitOpts.StoreResult = 1;
704
706
707 auto pad = GetDrawPad(obj);
708
709 TFitResultPtr res;
710
711 switch (kind) {
713
714 TH1 *hist = dynamic_cast<TH1*>(obj);
715 if (hist)
717
718 break;
719 }
721
722 TGraph *gr = dynamic_cast<TGraph*>(obj);
723 if (gr)
725 break;
726 }
728
729 TMultiGraph *mg = dynamic_cast<TMultiGraph*>(obj);
730 if (mg)
732
733 break;
734 }
736
737 TGraph2D *g2d = dynamic_cast<TGraph2D*>(obj);
738 if (g2d)
740
741 break;
742 }
744 // N/A
745 break;
746 }
747
748 default: {
749 // N/A
750 break;
751 }
752 }
753
754 // After fitting function appears in global list
755 if (f1 && gROOT->GetListOfFunctions()->FindObject(f1.get()))
756 gROOT->GetListOfFunctions()->Remove(f1.get());
757
758 if (m.fSame && f1 && pad) {
759 TF1 *copy = copyTF1(f1.get());
760 copy->SetBit(kCanDelete);
761 copy->Draw("same");
762 }
763
765
766 std::string funcname = f1->GetName();
767 if ((funcname.compare(0,4,"prev") == 0) && (funcname.find("-") > 4))
768 funcname.erase(0, funcname.find("-") + 1);
769 funcname = "prev"s + std::to_string(fPrevRes.size() + 1) + "-"s + funcname;
770 f1->SetName(funcname.c_str());
771
772 fPrevRes.emplace_back(m.fSelectedData, f1, res);
773
775
776 SelectFunction("previous::"s + funcname);
777
778 return true; // provide client with latest changes
779}
780
781////////////////////////////////////////////////////////////////////////////////
782/// Extract color from string
783/// Should be coded as #%ff00ff string
784
786{
787 if ((colorid.length() != 7) || (colorid.compare(0,1,"#") != 0)) return 0;
788
789 return TColor::GetColor(colorid.c_str());
790}
791
792////////////////////////////////////////////////////////////////////////////////
793/// Create confidence levels drawing
794/// tab. Then it call Virtual Fitter to perform it.
795
797{
798 if (!result)
799 return nullptr;
800
801 // try to use provided method
802 // auto conf = result->GetConfidenceIntervals();
803 // printf("GET INTERVALS %d\n", (int) conf.size());
804
805 const auto *function = result->FittedFunction();
806 if (!function) {
807 R__LOG_ERROR(FitPanelLog()) << "Fit Function does not exist!";
808 return nullptr;
809 }
810
811 const auto *data = result->FittedBinData();
812 if (!data) {
813 R__LOG_ERROR(FitPanelLog()) << "Unbinned data set cannot draw confidence levels.";
814 return nullptr;
815 }
816
817 std::vector<Double_t> ci(data->Size());
818 result->GetConfidenceIntervals(*data, &ci[0], model().fConfidenceLevel);
819
820 if (model().fDim == 1) {
821 TGraphErrors *g = new TGraphErrors(ci.size());
822 for (unsigned int i = 0; i < ci.size(); ++i) {
823 const Double_t *x = data->Coords(i);
824 const Double_t y = (*function)(x);
825 g->SetPoint(i, *x, y);
826 g->SetPointError(i, 0, ci[i]);
827 }
828 // std::ostringstream os;
829 // os << "Confidence Intervals with " << conflvl << " conf. band.";
830 // g->SetTitle(os.str().c_str());
831 g->SetTitle("Confidence Intervals with");
832
833 auto icol = GetColor(model().fConfidenceColor);
834 g->SetLineColor(icol);
835 g->SetFillColor(icol);
836 g->SetFillStyle(3001);
837 return g;
838 } else if (model().fDim == 2) {
839 TGraph2DErrors *g = new TGraph2DErrors(ci.size());
840 for (unsigned int i = 0; i < ci.size(); ++i) {
841 const Double_t *x = data->Coords(i);
842 const Double_t y = (*function)(x);
843 g->SetPoint(i, x[0], x[1], y);
844 g->SetPointError(i, 0, 0, ci[i]);
845 }
846 // std::ostringstream os;
847 // os << "Confidence Intervals with " << fConfLevel->GetNumber() << " conf. band.";
848 // g->SetTitle(os.str().c_str());
849
850 g->SetTitle("Confidence Intervals with");
851
852 auto icol = GetColor(model().fConfidenceColor);
853 g->SetLineColor(icol);
854 g->SetFillColor(icol);
855 g->SetFillStyle(3001);
856 return g;
857 }
858
859 return nullptr;
860}
861
862////////////////////////////////////////////////////////////////////////////////
863/// Perform drawing using current model settings
864/// Returns true if any action was done
865
867{
868 auto &m = model();
869
870 TObject *obj = GetSelectedObject(m.fSelectedData);
871 if (!obj) return false;
872
873 TObject *drawobj = nullptr;
874 std::string drawopt;
875 bool superimpose = true, objowner = true;
876
877 if (m.fHasAdvanced && (m.fSelectedTab == "Advanced")) {
878
879 TFitResult *res = FindFitResult(m.fSelectedFunc);
880 if (!res) return false;
881
882 if (m.fAdvancedTab == "Contour") {
883
884 superimpose = m.fContourSuperImpose;
885 int par1 = std::stoi(m.fContourPar1Id);
886 int par2 = std::stoi(m.fContourPar2Id);
887
888 TGraph *graph = new TGraph(m.fContourPoints);
889
890 if (!res->Contour(par1, par2, graph, m.fConfidenceLevel)) {
891 delete graph;
892 return false;
893 }
894
895 auto fillcolor = GetColor(m.fContourColor);
896 graph->SetFillColor(fillcolor);
897 graph->GetXaxis()->SetTitle( res->ParName(par1).c_str() );
898 graph->GetYaxis()->SetTitle( res->ParName(par2).c_str() );
899
900 drawobj = graph;
901 drawopt = superimpose ? "LF" : "ALF";
902
903 } else if (m.fAdvancedTab == "Scan") {
904
905 int par = std::stoi(m.fScanId);
906 TGraph *graph = new TGraph( m.fScanPoints);
907 if (!res->Scan( par, graph, m.fScanMin, m.fScanMax)) {
908 delete graph;
909 return false;
910 }
911 auto linecolor = GetColor(m.fScanColor);
912 if (!linecolor) linecolor = kBlue;
913 graph->SetLineColor(linecolor);
914 graph->SetLineWidth(2);
915 graph->GetXaxis()->SetTitle(res->ParName(par).c_str());
916 graph->GetYaxis()->SetTitle("FCN" );
917
918 superimpose = false;
919 drawobj = graph;
920 drawopt = "ALF";
921
922 } else if (m.fAdvancedTab == "Confidence") {
923
925 drawopt = "C3same";
926
927 } else {
928 return false;
929 }
930 } else {
931
932 // find already existing functions, not try to create something new
933 TF1 *func = FindFunction(m.fSelectedFunc);
934
935 // when "Pars" tab is selected, automatically update function parameters
936 if (func && (m.fSelectedTab.compare("Pars") == 0) && (m.fSelectedFunc == m.fFuncPars.id))
937 m.fFuncPars.SetParameters(func);
938
939 drawobj = func;
940 drawopt = "same";
941 objowner = true;
942 }
943
944 if (!drawobj)
945 return false;
946
947 auto pad = GetDrawPad(obj, true);
948 if (!pad) {
949 if (objowner)
950 delete drawobj;
951 return false;
952 }
953
954 if (!superimpose)
955 pad->Clear();
956
957 if (objowner)
958 drawobj->SetBit(kCanDelete);
959
960 drawobj->Draw(drawopt.c_str());
961
963
964 return true;
965}
966
967
968//////////////////////////////////////////////////////////////////////////////////////////////
969/// Mark pad modified and do update
970/// For web canvas set async mode first to avoid blocking here
971
973{
974 if (!pad) return;
975
976 pad->Modified();
977 pad->UpdateAsync();
978}
979
980
981//////////////////////////////////////////////////////////////////////////////////////////////
982/// Set handle which will be cleared when connection is closed
983
984void RFitPanel::ClearOnClose(const std::shared_ptr<void> &handle)
985{
986 GetWindow()->SetClearOnClose(handle);
987}
nlohmann::json json
#define R__LOG_ERROR(...)
Definition RLogger.hxx:357
#define f(i)
Definition RSha256.hxx:104
#define c(i)
Definition RSha256.hxx:101
#define g(i)
Definition RSha256.hxx:105
short Color_t
Color number (short)
Definition RtypesCore.h:99
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
@ kBlue
Definition Rtypes.h:67
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
#define gDirectory
Definition TDirectory.h:385
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
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 result
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize id
char name[80]
Definition TGX11.cxx:110
float xmin
float ymin
float xmax
float ymax
@ kCanDelete
Definition TObject.h:370
#define gROOT
Definition TROOT.h:411
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2495
RFitPanel(const std::string &title="Fit panel")
Constructor.
Definition RFitPanel.cxx:83
TPad * GetDrawPad(TObject *obj, bool force=false)
Returns pad where histogram is drawn If canvas not exists, create new one.
void UpdateDataSet()
Update list of available data.
void UpdateFunctionsList()
Update list of available functions.
void SelectFunction(const std::string &funcid)
Select fit function.
void GetFunctionsFromSystem()
Looks for all the functions registered in the current ROOT session.
void ClearOnClose(const std::shared_ptr< void > &handle)
Set handle which will be cleared when connection is closed.
RFitPanelModel::EFitObjectType GetFitObjectType(TObject *obj)
Returns kind of object.
std::shared_ptr< RCanvas > fCanvas
! v7 canvas used to display results
Definition RFitPanel.hxx:43
TF1 * FindFunction(const std::string &funcid)
Search for existing functions, ownership still belongs to FitPanel or global lists.
std::string fCanvName
! v6 canvas name used to display fit, will be created if not exists
Definition RFitPanel.hxx:40
std::unique_ptr< RFitPanelModel > fModel
Definition RFitPanel.hxx:37
Color_t GetColor(const std::string &colorid)
Extract color from string Should be coded as #ff00ff string.
unsigned fConnId
! client connection id
Definition RFitPanel.hxx:46
std::string fPadName
! v6 pad name in the canvas, where object is (was) drawn
Definition RFitPanel.hxx:41
bool DoFit()
Perform fitting using current model settings Returns true if any action was done.
std::shared_ptr< ROOT::RWebWindow > fWindow
! configured display
Definition RFitPanel.hxx:45
TObject * MakeConfidenceLevels(TFitResult *res)
Create confidence levels drawing tab.
std::vector< TObject * > fObjects
! objects provided directly to panel for fitting
Definition RFitPanel.hxx:39
std::unique_ptr< TF1 > GetFitFunction(const std::string &funcid)
Creates new instance to make fitting.
void AssignHistogram(TH1 *hist)
Assign histogram to use with fit panel - without ownership.
TFitResult * FindFitResult(const std::string &funcid)
Creates new instance to make fitting.
TObject * GetSelectedObject(const std::string &objid)
Returns object based on it string id Searches either in gDirectory or in internal panel list.
void SendModel()
Send model object to the client.
void Hide()
hide FitPanel
int UpdateModel(const std::string &json)
Update fit model returns -1 if JSON fails return 0 if nothing large changed return 1 if important sel...
void DoPadUpdate(TPad *pad)
Mark pad modified and do update For web canvas set async mode first to avoid blocking here.
void Show(const RWebDisplayArgs &args="")
show FitPanel in specified place
RFitPanelModel & model()
Return reference on model object Model created if was not exists before.
TF1 * copyTF1(TF1 *f)
Copies f into a new TF1 to be stored in the fitpanel with it's own ownership.
std::vector< std::unique_ptr< TF1 > > fSystemFuncs
! local copy of all internal system funcs
Definition RFitPanel.hxx:48
std::list< FitRes > fPrevRes
! all previous functions used for fitting
Definition RFitPanel.hxx:59
void ProcessData(unsigned connid, const std::string &arg)
Process data from FitPanel OpenUI5-based FitPanel sends commands or status changes.
void AssignCanvas(const std::string &cname)
bool DoDraw()
Perform drawing using current model settings Returns true if any action was done.
void SelectObject(const std::string &objid)
Select object for fitting.
std::shared_ptr< ROOT::RWebWindow > GetWindow()
Returns RWebWindow instance, used to display FitPanel.
class describing the range in the coordinates it supports multiple range in a coordinate.
Definition DataRange.h:35
std::string ParName(unsigned int i) const
name of the parameter
Holds different arguments for starting browser with RWebDisplayHandle::Display() method.
static std::shared_ptr< RWebWindow > Create()
Create new RWebWindow Using default RWebWindowsManager.
static unsigned ShowWindow(std::shared_ptr< RWebWindow > window, const RWebDisplayArgs &args="")
Static method to show web window Has to be used instead of RWebWindow::Show() when window potentially...
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:38
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:45
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:42
static TString ToJSON(const T *obj, Int_t compact=0, const char *member_name=nullptr)
Definition TBufferJSON.h:77
The Canvas class.
Definition TCanvas.h:23
void * New(ENewType defConstructor=kClassNew, Bool_t quiet=kFALSE) const
Return a pointer to a newly allocated object of this class.
Definition TClass.cxx:5017
static Int_t GetColor(const char *hexcolor)
Static method returning color number for color specified by hex color string of form: "#rrggbb",...
Definition TColor.cxx:1926
1-Dim function class
Definition TF1.h:182
void Copy(TObject &f1) const override
Copy this F1 to a new F1.
Definition TF1.cxx:1006
void Draw(Option_t *option="") override
Draw this function with its current attributes.
Definition TF1.cxx:1340
TClass * IsA() const override
Definition TF1.h:711
A 2-Dim function with parameters.
Definition TF2.h:29
TF3 defines a 3D Function with Parameters.
Definition TF3.h:28
Provides an indirection to the TFitResult class and with a semantics identical to a TFitResult pointe...
Extends the ROOT::Fit::Result class with a TNamed inheritance providing easy possibility for I/O.
Definition TFitResult.h:34
bool Contour(unsigned int ipar, unsigned int jpar, TGraph *gr, double confLevel=0.683)
Create a 2D contour around the minimum for the parameter ipar and jpar if a minimum does not exist or...
bool Scan(unsigned int ipar, TGraph *gr, double xmin=0, double xmax=0)
Scan parameter ipar between value of xmin and xmax A graph must be given which will be on return fill...
Graph 2D class with errors.
Graphics object made of three arrays X, Y and Z with the same number of points each.
Definition TGraph2D.h:41
static TClass * Class()
A TGraphErrors is a TGraph with error bars.
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
static TClass * Class()
TAxis * GetXaxis() const
Get x axis of the graph.
Definition TGraph.cxx:1593
TAxis * GetYaxis() const
Get y axis of the graph.
Definition TGraph.cxx:1602
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:109
static TClass * Class()
The Histogram stack class.
Definition THStack.h:40
static TClass * Class()
A TMultiGraph is a collection of TGraph (or derived) objects.
Definition TMultiGraph.h:34
static TClass * Class()
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:173
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition TNamed.cxx:149
Mother of all ROOT objects.
Definition TObject.h:41
virtual const char * GetName() const
Returns name of object.
Definition TObject.cxx:457
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:226
virtual Option_t * GetDrawOption() const
Get option used by the graphics system to draw this object.
Definition TObject.cxx:441
virtual TObject * FindObject(const char *name) const
Must be redefined in derived classes.
Definition TObject.cxx:421
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:881
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:543
virtual void Draw(Option_t *option="")
Default Draw method for all objects.
Definition TObject.cxx:293
The most important graphics class in the ROOT system.
Definition TPad.h:28
Basic string class.
Definition TString.h:138
small helper class to store/restore gPad context in TPad methods
Definition TVirtualPad.h:61
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
TGraphErrors * gr
Definition legend1.C:25
TF1 * f1
Definition legend1.C:11
ROOT::RLogChannel & FitPanelLog()
Log channel for FitPanel diagnostics.
TFitResultPtr FitObject(TH1 *h1, TF1 *f1, Foption_t &option, const ROOT::Math::MinimizerOptions &moption, const char *goption, ROOT::Fit::DataRange &range)
fitting function for a TH1 (called from TH1::Fit)
Definition HFitImpl.cxx:977
Data structure for the fit panel.
void UpdateAdvanced(TFitResult *res)
Update advanced parameters associated with fit function.
void SelectedFunc(const std::string &name, TF1 *func)
Select function.
std::string fTitle
title of the fit panel
TMarker m
Definition textangle.C:8