Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RFitPanel.hxx
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#ifndef ROOT_RFitPanel
14#define ROOT_RFitPanel
15
16#include <ROOT/RWebWindow.hxx>
17
19
20#include <ROOT/RCanvas.hxx>
21
22#include "ROOT/RHist.hxx"
23
24#include "TFitResultPtr.h"
25
26#include <memory>
27#include <vector>
28#include <list>
29
30class TPad;
31class TH1;
32class TF1;
33
34namespace ROOT {
35namespace Experimental {
36
37class RFitPanel {
38
39 std::unique_ptr<RFitPanelModel> fModel;
40
41 std::vector<TObject*> fObjects; ///<! objects provided directly to panel for fitting
42 std::string fCanvName; ///<! v6 canvas name used to display fit, will be created if not exists
43 std::string fPadName; ///<! v6 pad name in the canvas, where object is (was) drawn
44
45 std::shared_ptr<RCanvas> fCanvas; ///<! v7 canvas used to display results
46 std::shared_ptr<RH1D> fFitHist; ///<! v7 histogram for fitting
47
48 std::shared_ptr<RWebWindow> fWindow; ///<! configured display
49 unsigned fConnId{0}; ///<! client connection id
50
51 std::vector<std::unique_ptr<TF1>> fSystemFuncs; ///<! local copy of all internal system funcs
52
53 struct FitRes {
54 std::string objid; // object used for fitting
55 std::unique_ptr<TF1> func; // function
56 TFitResultPtr res; // result
57 FitRes() = default;
58 FitRes(const std::string &_objid, std::unique_ptr<TF1> &_func, TFitResultPtr &_res);
59 ~FitRes();
60 };
61
62 std::list<FitRes> fPrevRes; ///<! all previous functions used for fitting
63
64 TF1 *copyTF1(TF1 *f);
65
67
68 void ProcessData(unsigned connid, const std::string &arg);
69
71
72 int UpdateModel(const std::string &json);
73
74 void SelectObject(const std::string &objid);
75 TObject *GetSelectedObject(const std::string &objid);
77 void UpdateDataSet();
78
80 TF1 *FindFunction(const std::string &funcid);
81 TFitResult *FindFitResult(const std::string &funcid);
82 std::unique_ptr<TF1> GetFitFunction(const std::string &funcid);
83 void SelectFunction(const std::string &funcid);
84
86
87 Color_t GetColor(const std::string &colorid);
88
89 TPad *GetDrawPad(TObject *obj, bool force = false);
90
91 void SendModel();
92
93 bool DoFit();
94 bool DoDraw();
95
96public:
97 RFitPanel(const std::string &title = "Fit panel");
98 ~RFitPanel();
99
100 // method required when any panel want to be inserted into the RCanvas
101 std::shared_ptr<RWebWindow> GetWindow();
102
103 void AssignHistogram(TH1 *hist);
104
105 void AssignHistogram(const std::string &hname);
106
107 void AssignCanvas(const std::string &cname) { fCanvName = cname; }
108
109 void AssignCanvas(std::shared_ptr<RCanvas> &canv);
110
111 void AssignHistogram(std::shared_ptr<RH1D> &hist);
112
113 /// show FitPanel in specified place
114 void Show(const std::string &where = "");
115
116 /// hide FitPanel
117 void Hide();
118};
119
120} // namespace Experimental
121} // namespace ROOT
122
123#endif
#define f(i)
Definition RSha256.hxx:104
short Color_t
Definition RtypesCore.h:92
web-based FitPanel prototype.
Definition RFitPanel.hxx:37
void Show(const std::string &where="")
show FitPanel in specified place
bool DoFit()
Perform fitting using current model settings Returns true if any action was done.
void AssignHistogram(TH1 *hist)
Assign histogram to use with fit panel - without ownership.
std::unique_ptr< TF1 > GetFitFunction(const std::string &funcid)
Creates new instance to make fitting.
void SelectFunction(const std::string &funcid)
Select fit function.
RFitPanelModel::EFitObjectType GetFitObjectType(TObject *obj)
Returns kind of object.
std::shared_ptr< RCanvas > fCanvas
! v7 canvas used to display results
Definition RFitPanel.hxx:45
void SendModel()
Send model object to the client.
std::shared_ptr< RWebWindow > fWindow
! configured display
Definition RFitPanel.hxx:48
std::string fCanvName
! v6 canvas name used to display fit, will be created if not exists
Definition RFitPanel.hxx:42
void GetFunctionsFromSystem()
Looks for all the functions registered in the current ROOT session.
std::unique_ptr< RFitPanelModel > fModel
Definition RFitPanel.hxx:39
void Hide()
hide FitPanel
unsigned fConnId
! client connection id
Definition RFitPanel.hxx:49
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...
TF1 * copyTF1(TF1 *f)
Copies f into a new TF1 to be stored in the fitpanel with it's own ownership.
TPad * GetDrawPad(TObject *obj, bool force=false)
Returns pad where histogram is drawn If canvas not exists, create new one.
std::shared_ptr< RWebWindow > GetWindow()
Returns RWebWindow instance, used to display FitPanel.
Definition RFitPanel.cxx:94
void SelectObject(const std::string &objid)
Select object for fitting.
std::string fPadName
! v6 pad name in the canvas, where object is (was) drawn
Definition RFitPanel.hxx:43
TObject * GetSelectedObject(const std::string &objid)
Returns object based on it string id Searches either in gDirectory or in internal panel list.
std::vector< TObject * > fObjects
! objects provided directly to panel for fitting
Definition RFitPanel.hxx:41
bool DoDraw()
Perform drawing using current model settings Returns true if any action was done.
Color_t GetColor(const std::string &colorid)
Extract color from string Should be coded as #ff00ff string.
void UpdateFunctionsList()
Update list of available functions.
RFitPanelModel & model()
Return reference on model object Model created if was not exists before.
std::shared_ptr< RH1D > fFitHist
! v7 histogram for fitting
Definition RFitPanel.hxx:46
TObject * MakeConfidenceLevels(TFitResult *res)
Create confidence levels drawing tab.
TF1 * FindFunction(const std::string &funcid)
Search for existing functions, ownership still belongs to FitPanel or global lists.
void ProcessData(unsigned connid, const std::string &arg)
Process data from FitPanel OpenUI5-based FitPanel sends commands or status changes.
TFitResult * FindFitResult(const std::string &funcid)
Creates new instance to make fitting.
void UpdateDataSet()
Update list of available data.
std::vector< std::unique_ptr< TF1 > > fSystemFuncs
! local copy of all internal system funcs
Definition RFitPanel.hxx:51
std::list< FitRes > fPrevRes
! all previous functions used for fitting
Definition RFitPanel.hxx:62
void AssignCanvas(const std::string &cname)
1-Dim function class
Definition TF1.h:213
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
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:58
Mother of all ROOT objects.
Definition TObject.h:41
The most important graphics class in the ROOT system.
Definition TPad.h:26
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...
Data structure for the fit panel.