Logo ROOT   6.12/07
Reference Guide
TFitPanel.cxx
Go to the documentation of this file.
1 /// \file ROOT/TFitPanel.cxx
2 /// \ingroup WebGui ROOT7
3 /// \author Sergey Linev <S.Linev@gsi.de>
4 /// \date 2017-10-24
5 /// \warning This is part of the ROOT 7 prototype! It will change without notice. It might trigger earthquakes. Feedback
6 /// is welcome!
7 
8 /*************************************************************************
9  * Copyright (C) 1995-2017, Rene Brun and Fons Rademakers. *
10  * All rights reserved. *
11  * *
12  * For the licensing terms see $ROOTSYS/LICENSE. *
13  * For the list of contributors see $ROOTSYS/README/CREDITS. *
14  *************************************************************************/
15 
16 #include "ROOT/TFitPanel.hxx"
17 
18 #include <ROOT/TWebWindowsManager.hxx>
19 #include <ROOT/TLogger.hxx>
20 #include "ROOT/TDirectory.hxx"
21 
22 #include "TString.h"
23 #include "TROOT.h"
24 #include "TBufferJSON.h"
25 
26 /** \class ROOT::Experimental::TFitPanel
27 \ingroup webdisplay
28 
29 web-based FitPanel prototype.
30 */
31 
32 ///////////////////////////////////////////////////////////////////////////////////////////////////////
33 /// Returns TWebWindow instance, used to display FitPanel
34 
35 std::shared_ptr<ROOT::Experimental::TWebWindow> ROOT::Experimental::TFitPanel::GetWindow()
36 {
37  if (!fWindow) {
38  fWindow = TWebWindowsManager::Instance()->CreateWindow(gROOT->IsBatch());
39 
40  fWindow->SetPanelName("FitPanel");
41 
42  fWindow->SetDataCallBack([this](unsigned connid, const std::string &arg) { ProcessData(connid, arg); });
43 
44  fWindow->SetGeometry(300, 500); // configure predefined geometry
45  }
46 
47  return fWindow;
48 }
49 
50 ///////////////////////////////////////////////////////////////////////////////////////////////////////
51 /// Show FitPanel
52 
53 void ROOT::Experimental::TFitPanel::Show(const std::string &where)
54 {
55  GetWindow()->Show(where);
56 }
57 
58 ///////////////////////////////////////////////////////////////////////////////////////////////////////
59 /// Hide FitPanel
60 
62 {
63  if (!fWindow)
64  return;
65 
66  fWindow->CloseConnections();
67 }
68 
69 ///////////////////////////////////////////////////////////////////////////////////////////////////////
70 /// Process data from FitPanel
71 /// OpenUI5-based FitPanel sends commands or status changes
72 
73 void ROOT::Experimental::TFitPanel::ProcessData(unsigned connid, const std::string &arg)
74 {
75  if (arg == "CONN_READY") {
76  fConnId = connid;
77  printf("FitPanel connection established %u\n", fConnId);
78  fWindow->Send("INITDONE", fConnId);
79 
81  model.fDataNames.push_back(ComboBoxItem("1", "RootData1"));
82  model.fDataNames.push_back(ComboBoxItem("2", "RootData2"));
83  model.fDataNames.push_back(ComboBoxItem("3", "RootData3"));
84  model.fSelectDataId = "1";
85 
86  model.fModelNames.push_back(ComboBoxItem("1", "RootModel1"));
87  model.fModelNames.push_back(ComboBoxItem("2", "RootModel2"));
88  model.fModelNames.push_back(ComboBoxItem("3", "RootModel3"));
89  model.fSelectModelId = "3";
90 
91  TString json = TBufferJSON::ConvertToJSON(&model, gROOT->GetClass("ROOT::Experimental::TFitPanelModel"));
92 
93  fWindow->Send(std::string("MODEL:") + json.Data(), fConnId);
94 
95  return;
96  }
97 
98  if (arg == "CONN_CLOSED") {
99  printf("FitPanel connection closed\n");
100  fConnId = 0;
101  return;
102  }
103 
104  if (arg.find("DOFIT:") == 0) {
105  TString exec;
106  exec.Form("((ROOT::Experimental::TFitPanel *) %p)->DoFit(%s);", this, arg.c_str() + 6);
107  printf("Execute %s\n", exec.Data());
108  gROOT->ProcessLine(exec.Data());
109  return;
110  }
111 }
112 
113 ///////////////////////////////////////////////////////////////////////////////////////////////////////
114 /// Let use canvas to display fit results
115 
116 void ROOT::Experimental::TFitPanel::UseCanvas(std::shared_ptr<TCanvas> &canv)
117 {
118  if (!fCanvas) {
119  fCanvas = canv;
120  } else {
121  R__ERROR_HERE("ShowIn") << "FitPanel already bound to the canvas - change is not yet supported";
122  }
123 }
124 
125 ///////////////////////////////////////////////////////////////////////////////////////////////////////
126 /// Dummy function, called when "Fit" button pressed in UI
127 
128 void ROOT::Experimental::TFitPanel::DoFit(const std::string &dname, const std::string &mname)
129 {
130  printf("DoFit %s %s\n", dname.c_str(), mname.c_str());
131 
132  bool first_time = false;
133 
134  if (!fCanvas) {
135  fCanvas = Experimental::TCanvas::Create("FitPanel Canvas");
136  first_time = true;
137  }
138 
139  if (!fFitHist) {
140 
141  // Create the histogram.
142  auto xaxis = std::make_shared<ROOT::Experimental::TAxisConfig>(10, 0., 10.);
143 
144  fFitHist = std::make_shared<ROOT::Experimental::TH1D>(*xaxis.get());
145 
146  // Fill a few points.
147  fFitHist->Fill(5);
148  fFitHist->Fill(6);
149  fFitHist->Fill(6);
150  fFitHist->Fill(7);
151 
152  fCanvas->Draw(fFitHist).SetLineColor(Experimental::TColor::kBlue);
153 
154  // workaround to keep histogram in the lists
155  ROOT::Experimental::TDirectory::Heap().Add("fitaxis", xaxis);
156 
157  if (first_time) {
158  fCanvas->Show();
159  // fCanvas->Update();
160  } else {
161  fCanvas->Modified();
162  // fCanvas->Update();
163  }
164  }
165 }
std::vector< ComboBoxItem > fModelNames
Definition: TFitPanel.hxx:50
void Hide()
hide FitPanel
Definition: TFitPanel.cxx:61
std::vector< ComboBoxItem > fDataNames
Definition: TFitPanel.hxx:48
void Add(std::string_view name, const std::shared_ptr< T > &ptr)
Add an existing object (rather a shared_ptr to it) to the TDirectory.
Definition: TDirectory.hxx:172
#define gROOT
Definition: TROOT.h:402
Basic string class.
Definition: TString.h:125
std::shared_ptr< TH1D > fFitHist
!< canvas used to display results
Definition: TFitPanel.hxx:64
static TDirectory & Heap()
Dedicated, process-wide TDirectory.
Definition: TFile.cxx:23
unsigned fConnId
! connection id
Definition: TFitPanel.hxx:58
static constexpr PredefinedRGB kBlue
Definition: TColor.hxx:187
struct ROOT::Experimental::TFitPanelModelModel, used to initialized openui5 FitPanel ...
Definition: TFitPanel.hxx:47
std::shared_ptr< TWebWindow > GetWindow()
Returns TWebWindow instance, used to display FitPanel.
Definition: TFitPanel.cxx:35
void Form(const char *fmt,...)
Formats a string using a printf style format descriptor.
Definition: TString.cxx:2343
std::shared_ptr< TWebWindow > fWindow
Definition: TFitPanel.hxx:60
void Show(const std::string &where="")
show FitPanel in specified place
Definition: TFitPanel.cxx:53
void ProcessData(unsigned connid, const std::string &arg)
process data from UI
Definition: TFitPanel.cxx:73
static TString ConvertToJSON(const TObject *obj, Int_t compact=0, const char *member_name=0)
Converts object, inherited from TObject class, to JSON string Lower digit of compact parameter define...
void DoFit(const std::string &dname, const std::string &mname)
Dummy function, called when "Fit" button pressed in UI.
Definition: TFitPanel.cxx:128
struct ROOT::Experimental::ComboBoxItemDescriptor for the openui5 ComboBox, used in FitPanel ...
Definition: TFitPanel.hxx:35
void UseCanvas(std::shared_ptr< TCanvas > &canv)
let use canvas to display fit results
Definition: TFitPanel.cxx:116
static std::shared_ptr< TCanvas > Create(const std::string &title)
Definition: TCanvas.cxx:63
#define R__ERROR_HERE(GROUP)
Definition: TLogger.hxx:126
std::shared_ptr< TCanvas > fCanvas
!< configured display
Definition: TFitPanel.hxx:62
const char * Data() const
Definition: TString.h:345