Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
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 }
288
289 for (auto &func : fSystemFuncs) {
290 m.fFuncList.emplace_back("System", "system::"s + func->GetName(), func->GetName());
291 }
292
293 for (auto &entry : fPrevRes) {
294 if (entry.objid == m.fSelectedData)
295 m.fFuncList.emplace_back("Previous", "previous::"s + entry.func->GetName(), entry.func->GetName());
296 }
297}
298
299////////////////////////////////////////////////////////////////////////////////
300/// Select fit function
301
308
309////////////////////////////////////////////////////////////////////////////////
310/// Assign histogram to use with fit panel - without ownership
311
313{
314 fObjects.emplace_back(hist);
315 SelectObject("panel::"s + hist->GetName());
316 SendModel();
317}
318
319/// Assign histogram name to use with fit panel - it should be available in gDirectory
320
321void RFitPanel::AssignHistogram(const std::string &hname)
322{
323 SelectObject("gdir::" + hname);
324 SendModel();
325}
326
327////////////////////////////////////////////////////////////////////////////////
328/// assign canvas to use for drawing results of fitting or showing fitpanel itself
329
330void RFitPanel::AssignCanvas(std::shared_ptr<RCanvas> &canv)
331{
332 if (!fCanvas) {
333 fCanvas = canv;
334 } else {
335 R__LOG_ERROR(FitPanelLog()) << "FitPanel already bound to the canvas - change is not yet supported";
336 }
337}
338
339////////////////////////////////////////////////////////////////////////////////
340/// Show FitPanel
341
343{
345}
346
347////////////////////////////////////////////////////////////////////////////////
348/// Hide FitPanel
349
351{
352 if (fWindow)
353 fWindow->CloseConnections();
354}
355
356////////////////////////////////////////////////////////////////////////////////
357/// Return reference on model object
358/// Model created if was not exists before
359
360
362{
363 if (!fModel) {
364 fModel = std::make_unique<RFitPanelModel>();
365 fModel->Initialize();
366 }
367
368 return *fModel.get();
369}
370
371////////////////////////////////////////////////////////////////////////////////
372/// Send model object to the client
373
375{
376 if (fWindow && (fConnId > 0)) {
378 fWindow->Send(fConnId, "MODEL:"s + json.Data());
379 }
380}
381
382////////////////////////////////////////////////////////////////////////////////
383/// Process data from FitPanel
384/// OpenUI5-based FitPanel sends commands or status changes
385
386void RFitPanel::ProcessData(unsigned, const std::string &arg)
387{
388 if (arg == "RELOAD") {
389
391
394
395 SendModel();
396
397 } else if (arg.compare(0, 7, "UPDATE:") == 0) {
398
399 if (UpdateModel(arg.substr(7)) > 0)
400 SendModel();
401
402 } else if (arg.compare(0, 6, "DOFIT:") == 0) {
403
404 if (UpdateModel(arg.substr(6)) >= 0)
405 if (DoFit())
406 SendModel();
407
408 } else if (arg.compare(0, 7, "DODRAW:") == 0) {
409
410 if (UpdateModel(arg.substr(7)) >= 0)
411 if (DoDraw())
412 SendModel();
413
414 } else if (arg.compare(0, 8, "SETPARS:") == 0) {
415
416 auto info = TBufferJSON::FromJSON<RFitPanelModel::RFuncParsList>(arg.substr(8));
417
418 if (info) {
419 TF1 *func = FindFunction(info->id);
420
421 // copy all parameters back to the function
422 if (func)
423 info->SetParameters(func);
424 }
425
426 }
427}
428
429////////////////////////////////////////////////////////////////////////////////
430/// Search for existing functions, ownership still belongs to FitPanel or global lists
431
432TF1 *RFitPanel::FindFunction(const std::string &id)
433{
434 if (id.compare(0,8,"system::") == 0) {
435 std::string name = id.substr(8);
436
437 for (auto &item : fSystemFuncs)
438 if (name.compare(item->GetName()) == 0)
439 return item.get();
440 }
441
442 if (id.compare(0,10,"previous::") == 0) {
443 std::string name = id.substr(10);
444
445 for (auto &entry : fPrevRes)
446 if (name.compare(entry.func->GetName()) == 0)
447 return entry.func.get();
448 }
449
450 return nullptr;
451}
452
453////////////////////////////////////////////////////////////////////////////////
454/// Creates new instance to make fitting
455
457{
458 if (id.compare(0,10,"previous::") == 0) {
459 std::string name = id.substr(10);
460
461 for (auto &entry : fPrevRes)
462 if (name.compare(entry.func->GetName()) == 0)
463 return entry.res.Get();
464 }
465
466 return nullptr;
467}
468
469////////////////////////////////////////////////////////////////////////////////
470/// Creates new instance to make fitting
471
472std::unique_ptr<TF1> RFitPanel::GetFitFunction(const std::string &funcname)
473{
474 std::unique_ptr<TF1> res;
475
476 TF1 *func = FindFunction(funcname);
477
478 if (func) {
479 // Now we make a copy.
480 res.reset((TF1*) func->IsA()->New());
481 func->Copy(*res);
482 } else if (funcname.compare(0,6,"dflt::") == 0) {
483
484 std::string formula = funcname.substr(6);
485
487
488 double xmin, xmax, ymin, ymax, zmin, zmax;
489 drange.GetRange(xmin, xmax, ymin, ymax, zmin, zmax);
490
491 if ( model().fDim == 1 || model().fDim == 0 ) {
492 res.reset(new TF1(formula.c_str(), formula.c_str(), xmin, xmax));
493 } else if ( model().fDim == 2 ) {
494 res.reset(new TF2(formula.c_str(), formula.c_str(), xmin, xmax, ymin, ymax));
495 } else if ( model().fDim == 3 ) {
496 res.reset(new TF3(formula.c_str(), formula.c_str(), xmin, xmax, ymin, ymax, zmin, zmax));
497 }
498 }
499
500 return res;
501}
502
503////////////////////////////////////////////////////////////////////////////////
504/// Update fit model
505/// returns -1 if JSON fails
506/// return 0 if nothing large changed
507/// return 1 if important selection are changed and client need to be updated
508
509int RFitPanel::UpdateModel(const std::string &json)
510{
511 auto m = TBufferJSON::FromJSON<RFitPanelModel>(json);
512
513 if (!m) {
514 R__LOG_ERROR(FitPanelLog()) << "Fail to parse JSON for RFitPanelModel";
515 return -1;
516 }
517
518 m->fInitialized = true;
519
520 int res = 0; // nothing changed
521
522 if (model().fSelectedData != m->fSelectedData) {
523 res |= 1;
524 }
525
526 if (model().fSelectedFunc != m->fSelectedFunc) {
527 res |= 2;
528 }
529
530 std::swap(fModel, m); // finally replace model
531
532 if (res & 1)
533 SelectObject(model().fSelectedData);
534
535 if (res != 0)
536 SelectFunction(model().fSelectedFunc);
537
538 return res;
539}
540
541////////////////////////////////////////////////////////////////////////////////
542/// Copies f into a new TF1 to be stored in the fitpanel with it's
543/// own ownership. This is taken from Fit::StoreAndDrawFitFunction in
544/// HFitImpl.cxx
545
547{
548 double xmin = 0, xmax = 0, ymin = 0, ymax = 0, zmin = 0, zmax = 0;
549
550 // no need to use kNotGlobal bit. TF1::Copy does not add in the list by default
551 if ( dynamic_cast<TF3*>(f)) {
552 TF3* fnew = (TF3*)f->IsA()->New();
553 f->Copy(*fnew);
554 f->GetRange(xmin,ymin,zmin,xmax,ymax,zmax);
555 fnew->SetRange(xmin,ymin,zmin,xmax,ymax,zmax);
556 fnew->SetParent( nullptr );
557 fnew->AddToGlobalList(false);
558 return fnew;
559 } else if ( dynamic_cast<TF2*>(f) != 0 ) {
560 TF2* fnew = (TF2*)f->IsA()->New();
561 f->Copy(*fnew);
562 f->GetRange(xmin,ymin,xmax,ymax);
563 fnew->SetRange(xmin,ymin,xmax,ymax);
564 fnew->Save(xmin,xmax,ymin,ymax,0,0);
565 fnew->SetParent( nullptr );
566 fnew->AddToGlobalList(false);
567 return fnew;
568 }
569
570 TF1* fnew = (TF1*)f->IsA()->New();
571 f->Copy(*fnew);
572 f->GetRange(xmin,xmax);
573 fnew->SetRange(xmin,xmax);
574 // This next line is added, as fnew-Save fails with gausND! As
575 // the number of dimensions is unknown...
576 if ( '\0' != fnew->GetExpFormula()[0] )
577 fnew->Save(xmin,xmax,0,0,0,0);
578 fnew->SetParent( nullptr );
579 fnew->AddToGlobalList(false);
580 return fnew;
581}
582
583////////////////////////////////////////////////////////////////////////////////
584/// Looks for all the functions registered in the current ROOT
585/// session.
586
588{
589
590 fSystemFuncs.clear();
591
592 // Be carefull not to store functions that will be in the
593 // predefined section
594 std::vector<std::string> fnames = { "gaus" , "gausn", "expo", "landau",
595 "landaun", "pol0", "pol1", "pol2",
596 "pol3", "pol4", "pol5", "pol6",
597 "pol7", "pol8", "pol9", "user" };
598
599 // No go through all the objects registered in gROOT
600 TIter functionsIter(gROOT->GetListOfFunctions());
601 TObject* obj;
602 while( (obj = functionsIter()) != nullptr ) {
603 // And if they are TF1s
604 if ( TF1* func = dynamic_cast<TF1*>(obj) ) {
605 bool addFunction = true;
606 // And they are not already registered in fSystemFunc
607 for ( auto &name : fnames) {
608 if ( name.compare(func->GetName()) == 0 ) {
609 addFunction = false;
610 break;
611 }
612 }
613 // Add them.
614 if ( addFunction )
615 fSystemFuncs.emplace_back( copyTF1(func) );
616 }
617 }
618}
619
620////////////////////////////////////////////////////////////////////////////////
621/// Returns pad where histogram is drawn
622/// If canvas not exists, create new one
623
625{
626 if (!obj || (!force && (model().fNoDrawing || model().fNoStoreDraw)))
627 return nullptr;
628
629 std::function<TPad*(TPad*)> check = [&](TPad *pad) {
630 TPad *res = nullptr;
631 if (!pad) return res;
632 if (!fPadName.empty() && (fPadName.compare(pad->GetName()) == 0)) return pad;
633 TIter next(pad->GetListOfPrimitives());
634 TObject *prim = nullptr;
635 while (!res && (prim = next())) {
636 res = (prim == obj) ? pad : check(dynamic_cast<TPad *>(prim));
637 }
638 return res;
639 };
640
641 if (!fCanvName.empty()) {
642 auto drawcanv = dynamic_cast<TCanvas *> (gROOT->GetListOfCanvases()->FindObject(fCanvName.c_str()));
643 auto drawpad = check(drawcanv);
644 if (drawpad) {
645 drawpad->cd();
646 return drawpad;
647 }
648 if (drawcanv) {
649 drawcanv->Clear();
650 drawcanv->cd();
651 obj->Draw();
652 return drawcanv;
653 }
654 fCanvName.clear();
655 fPadName.clear();
656 }
657
658 TObject *c = nullptr;
659 TIter nextc(gROOT->GetListOfCanvases());
660 while ((c = nextc())) {
661 auto drawpad = check(dynamic_cast<TCanvas* >(c));
662 if (drawpad) {
663 drawpad->cd();
664 fCanvName = c->GetName();
665 fPadName = drawpad->GetName();
666 return drawpad;
667 }
668 }
669
670 auto canv = gROOT->MakeDefCanvas();
671 canv->SetName("fpc");
672 canv->SetTitle("Fit panel drawings");
673 fPadName = fCanvName = canv->GetName();
674
675 canv->cd();
676 obj->Draw();
677
678 return canv;
679}
680
681////////////////////////////////////////////////////////////////////////////////
682/// Perform fitting using current model settings
683/// Returns true if any action was done
684
686{
687 auto &m = model();
688
689 TObject *obj = GetSelectedObject(m.fSelectedData);
690 if (!obj) return false;
691 auto kind = GetFitObjectType(obj);
692
693 auto f1 = GetFitFunction(m.fSelectedFunc);
694 if (!f1) return false;
695
696 auto drange = m.GetRanges();
697 auto minOption = m.GetMinimizerOptions();
698 auto fitOpts = m.GetFitOptions();
699 auto drawOpts = m.GetDrawOption();
700
701 fitOpts.StoreResult = 1;
702
704
705 auto pad = GetDrawPad(obj);
706
707 TFitResultPtr res;
708
709 switch (kind) {
711
712 TH1 *hist = dynamic_cast<TH1*>(obj);
713 if (hist)
715
716 break;
717 }
719
720 TGraph *gr = dynamic_cast<TGraph*>(obj);
721 if (gr)
723 break;
724 }
726
727 TMultiGraph *mg = dynamic_cast<TMultiGraph*>(obj);
728 if (mg)
730
731 break;
732 }
734
735 TGraph2D *g2d = dynamic_cast<TGraph2D*>(obj);
736 if (g2d)
738
739 break;
740 }
742 // N/A
743 break;
744 }
745
746 default: {
747 // N/A
748 break;
749 }
750 }
751
752 // After fitting function appears in global list
753 if (f1 && gROOT->GetListOfFunctions()->FindObject(f1.get()))
754 gROOT->GetListOfFunctions()->Remove(f1.get());
755
756 if (m.fSame && f1 && pad) {
757 TF1 *copy = copyTF1(f1.get());
758 copy->SetBit(kCanDelete);
759 copy->Draw("same");
760 }
761
763
764 std::string funcname = f1->GetName();
765 if ((funcname.compare(0,4,"prev") == 0) && (funcname.find("-") > 4))
766 funcname.erase(0, funcname.find("-") + 1);
767 funcname = "prev"s + std::to_string(fPrevRes.size() + 1) + "-"s + funcname;
768 f1->SetName(funcname.c_str());
769
770 fPrevRes.emplace_back(m.fSelectedData, f1, res);
771
773
774 SelectFunction("previous::"s + funcname);
775
776 return true; // provide client with latest changes
777}
778
779////////////////////////////////////////////////////////////////////////////////
780/// Extract color from string
781/// Should be coded as #%ff00ff string
782
784{
785 if ((colorid.length() != 7) || (colorid.compare(0,1,"#") != 0)) return 0;
786
787 return TColor::GetColor(colorid.c_str());
788}
789
790////////////////////////////////////////////////////////////////////////////////
791/// Create confidence levels drawing
792/// tab. Then it call Virtual Fitter to perform it.
793
795{
796 if (!result)
797 return nullptr;
798
799 // try to use provided method
800 // auto conf = result->GetConfidenceIntervals();
801 // printf("GET INTERVALS %d\n", (int) conf.size());
802
803 const auto *function = result->FittedFunction();
804 if (!function) {
805 R__LOG_ERROR(FitPanelLog()) << "Fit Function does not exist!";
806 return nullptr;
807 }
808
809 const auto *data = result->FittedBinData();
810 if (!data) {
811 R__LOG_ERROR(FitPanelLog()) << "Unbinned data set cannot draw confidence levels.";
812 return nullptr;
813 }
814
815 std::vector<Double_t> ci(data->Size());
816 result->GetConfidenceIntervals(*data, &ci[0], model().fConfidenceLevel);
817
818 if (model().fDim == 1) {
819 TGraphErrors *g = new TGraphErrors(ci.size());
820 for (unsigned int i = 0; i < ci.size(); ++i) {
821 const Double_t *x = data->Coords(i);
822 const Double_t y = (*function)(x);
823 g->SetPoint(i, *x, y);
824 g->SetPointError(i, 0, ci[i]);
825 }
826 // std::ostringstream os;
827 // os << "Confidence Intervals with " << conflvl << " conf. band.";
828 // g->SetTitle(os.str().c_str());
829 g->SetTitle("Confidence Intervals with");
830
831 auto icol = GetColor(model().fConfidenceColor);
832 g->SetLineColor(icol);
833 g->SetFillColor(icol);
834 g->SetFillStyle(3001);
835 return g;
836 } else if (model().fDim == 2) {
837 TGraph2DErrors *g = new TGraph2DErrors(ci.size());
838 for (unsigned int i = 0; i < ci.size(); ++i) {
839 const Double_t *x = data->Coords(i);
840 const Double_t y = (*function)(x);
841 g->SetPoint(i, x[0], x[1], y);
842 g->SetPointError(i, 0, 0, ci[i]);
843 }
844 // std::ostringstream os;
845 // os << "Confidence Intervals with " << fConfLevel->GetNumber() << " conf. band.";
846 // g->SetTitle(os.str().c_str());
847
848 g->SetTitle("Confidence Intervals with");
849
850 auto icol = GetColor(model().fConfidenceColor);
851 g->SetLineColor(icol);
852 g->SetFillColor(icol);
853 g->SetFillStyle(3001);
854 return g;
855 }
856
857 return nullptr;
858}
859
860////////////////////////////////////////////////////////////////////////////////
861/// Perform drawing using current model settings
862/// Returns true if any action was done
863
865{
866 auto &m = model();
867
868 TObject *obj = GetSelectedObject(m.fSelectedData);
869 if (!obj) return false;
870
871 TObject *drawobj = nullptr;
872 std::string drawopt;
873 bool superimpose = true, objowner = true;
874
875 if (m.fHasAdvanced && (m.fSelectedTab == "Advanced")) {
876
877 TFitResult *res = FindFitResult(m.fSelectedFunc);
878 if (!res) return false;
879
880 if (m.fAdvancedTab == "Contour") {
881
882 superimpose = m.fContourSuperImpose;
883 int par1 = std::stoi(m.fContourPar1Id);
884 int par2 = std::stoi(m.fContourPar2Id);
885
886 TGraph *graph = new TGraph(m.fContourPoints);
887
888 if (!res->Contour(par1, par2, graph, m.fConfidenceLevel)) {
889 delete graph;
890 return false;
891 }
892
893 auto fillcolor = GetColor(m.fContourColor);
894 graph->SetFillColor(fillcolor);
895 graph->GetXaxis()->SetTitle( res->ParName(par1).c_str() );
896 graph->GetYaxis()->SetTitle( res->ParName(par2).c_str() );
897
898 drawobj = graph;
899 drawopt = superimpose ? "LF" : "ALF";
900
901 } else if (m.fAdvancedTab == "Scan") {
902
903 int par = std::stoi(m.fScanId);
904 TGraph *graph = new TGraph( m.fScanPoints);
905 if (!res->Scan( par, graph, m.fScanMin, m.fScanMax)) {
906 delete graph;
907 return false;
908 }
909 auto linecolor = GetColor(m.fScanColor);
910 if (!linecolor) linecolor = kBlue;
911 graph->SetLineColor(linecolor);
912 graph->SetLineWidth(2);
913 graph->GetXaxis()->SetTitle(res->ParName(par).c_str());
914 graph->GetYaxis()->SetTitle("FCN" );
915
916 superimpose = false;
917 drawobj = graph;
918 drawopt = "ALF";
919
920 } else if (m.fAdvancedTab == "Confidence") {
921
923 drawopt = "C3same";
924
925 } else {
926 return false;
927 }
928 } else {
929
930 // find already existing functions, not try to create something new
931 TF1 *func = FindFunction(m.fSelectedFunc);
932
933 // when "Pars" tab is selected, automatically update function parameters
934 if (func && (m.fSelectedTab.compare("Pars") == 0) && (m.fSelectedFunc == m.fFuncPars.id))
935 m.fFuncPars.SetParameters(func);
936
937 drawobj = func;
938 drawopt = "same";
939 objowner = true;
940 }
941
942 if (!drawobj)
943 return false;
944
945 auto pad = GetDrawPad(obj, true);
946 if (!pad) {
947 if (objowner)
948 delete drawobj;
949 return false;
950 }
951
952 if (!superimpose)
953 pad->Clear();
954
955 if (objowner)
956 drawobj->SetBit(kCanDelete);
957
958 drawobj->Draw(drawopt.c_str());
959
961
962 return true;
963}
964
965
966//////////////////////////////////////////////////////////////////////////////////////////////
967/// Mark pad modified and do update
968/// For web canvas set async mode first to avoid blocking here
969
971{
972 if (!pad) return;
973
974 pad->Modified();
975 pad->UpdateAsync();
976}
977
978
979//////////////////////////////////////////////////////////////////////////////////////////////
980/// Set handle which will be cleared when connection is closed
981
982void RFitPanel::ClearOnClose(const std::shared_ptr<void> &handle)
983{
984 GetWindow()->SetClearOnClose(handle);
985}
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
Definition RtypesCore.h:85
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
@ kBlue
Definition Rtypes.h:66
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
#define gDirectory
Definition TDirectory.h:384
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:373
#define gROOT
Definition TROOT.h:406
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2489
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:75
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:5100
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:1927
1-Dim function class
Definition TF1.h:234
void Copy(TObject &f1) const override
Copy this F1 to a new F1.
Definition TF1.cxx:1007
void Draw(Option_t *option="") override
Draw this function with its current attributes.
Definition TF1.cxx:1335
TClass * IsA() const override
Definition TF1.h:763
A 2-Dim function with parameters.
Definition TF2.h:29
A 3-Dim 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:1567
TAxis * GetYaxis() const
Get y axis of the graph.
Definition TGraph.cxx:1576
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
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:174
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:150
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:864
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:139
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