Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TWebCanvas.cxx
Go to the documentation of this file.
1// Author: Sergey Linev, GSI 7/12/2016
2
3/*************************************************************************
4 * Copyright (C) 1995-2023, Rene Brun and Fons Rademakers. *
5 * All rights reserved. *
6 * *
7 * For the licensing terms see $ROOTSYS/LICENSE. *
8 * For the list of contributors see $ROOTSYS/README/CREDITS. *
9 *************************************************************************/
10
11#include "TWebCanvas.h"
12
13#include "TWebSnapshot.h"
14#include "TWebPadPainter.h"
15#include "TWebPS.h"
16#include "TWebMenuItem.h"
17
18#include "TSystem.h"
19#include "TStyle.h"
20#include "TCanvas.h"
21#include "TButton.h"
22#include "TFrame.h"
23#include "TPaveText.h"
24#include "TPaveStats.h"
25#include "TText.h"
26#include "TROOT.h"
27#include "TClass.h"
28#include "TColor.h"
29#include "TObjArray.h"
30#include "TArrayI.h"
31#include "TList.h"
32#include "TF1.h"
33#include "TF2.h"
34#include "TH1.h"
35#include "TH2.h"
36#include "TH1K.h"
37#include "TH2.h"
38#include "THStack.h"
39#include "TMultiGraph.h"
40#include "TEnv.h"
41#include "TError.h"
42#include "TGraph.h"
43#include "TGraph2D.h"
44#include "TGaxis.h"
45#include "TScatter.h"
46#include "TCutG.h"
47#include "TBufferJSON.h"
48#include "TBase64.h"
49#include "TAtt3D.h"
50#include "TView.h"
51#include "TExec.h"
52#include "TVirtualX.h"
53#include "TMath.h"
54#include "TTimer.h"
55
56#include <cstdio>
57#include <cstring>
58#include <fstream>
59#include <iostream>
60#include <memory>
61#include <sstream>
62
63class TWebCanvasTimer : public TTimer {
68public:
70
71 Bool_t IsSlow() const { return fSlow; }
72 void SetSlow(Bool_t slow = kTRUE)
73 {
74 fSlow = slow;
75 fSlowCnt = 0;
76 SetTime(slow ? 1000 : 10);
77 }
78
79 /// used to send control messages to clients
80 void Timeout() override
81 {
82 if (fProcessing || fCanv.fProcessingData) return;
86 if (res) {
87 fSlowCnt = 0;
88 } else if (++fSlowCnt > 10 && !IsSlow()) {
90 }
91 }
92};
93
94
95/** \class TWebCanvas
96\ingroup webgui6
97
98Basic TCanvasImp ABI implementation for Web-based Graphics
99Provides painting of main ROOT classes in web browsers using [JSROOT](https://root.cern/js/)
100
101Following settings parameters can be useful for TWebCanvas:
102
103 WebGui.FullCanvas: 1 read-only mode (0), full-functional canvas (1) (default - 1)
104 WebGui.StyleDelivery: 1 provide gStyle object to JSROOT client (default - 1)
105 WebGui.PaletteDelivery: 1 provide color palette to JSROOT client (default - 1)
106 WebGui.TF1UseSave: 0 used saved values for function drawing (1) or calculate function on the client side (0) (default - 0)
107
108TWebCanvas is used by default in interactive ROOT session. To use web-based canvas in batch mode for image
109generation, one should explicitly specify `--web` option when starting ROOT:
110
111 [shell] root -b --web tutorials/hsimple.root -e 'hpxpy->Draw("colz"); c1->SaveAs("image.png");'
112
113If for any reasons TWebCanvas does not provide required functionality, one always can disable it.
114Either by specifying `root --web=off` when starting ROOT or by setting `Canvas.Name: TRootCanvas` in rootrc file.
115
116*/
117
118using namespace std::string_literals;
119
120////////////////////////////////////////////////////////////////////////////////
121/// Constructor
122
124 : TCanvasImp(c, name, x, y, width, height)
125{
126 fTimer = new TWebCanvasTimer(*this);
127
128 fReadOnly = readonly;
129 fStyleDelivery = gEnv->GetValue("WebGui.StyleDelivery", 1);
130 fPaletteDelivery = gEnv->GetValue("WebGui.PaletteDelivery", 1);
131 fPrimitivesMerge = gEnv->GetValue("WebGui.PrimitivesMerge", 100);
132 fTF1UseSave = gEnv->GetValue("WebGui.TF1UseSave", (Int_t) 1) > 0;
134
135 fWebConn.emplace_back(0); // add special connection which only used to perform updates
136
137 fTimer->TurnOn();
138
139 // fAsyncMode = kTRUE;
140}
141
142
143////////////////////////////////////////////////////////////////////////////////
144/// Destructor
145
147{
148 delete fTimer;
149}
150
151////////////////////////////////////////////////////////////////////////////////
152/// Initialize window for the web canvas
153/// At this place canvas is not yet register to the list of canvases - one cannot call RWebWindow::Show()
154
156{
157 return 111222333; // should not be used at all
158}
159
160////////////////////////////////////////////////////////////////////////////////
161/// Creates web-based pad painter
162
164{
165 return new TWebPadPainter();
166}
167
168////////////////////////////////////////////////////////////////////////////////
169/// Returns kTRUE when object is fully supported on JSROOT side
170/// In ROOT7 Paint function will just return appropriate flag that object can be displayed on JSROOT side
171
173{
174 if (!obj)
175 return kTRUE;
176
177 static const struct {
178 const char *name{nullptr};
179 bool with_derived{false};
180 bool reduse_by_many{false};
181 } supported_classes[] = {{"TH1", true},
182 {"TF1", true},
183 {"TGraph", true},
184 {"TScatter"},
185 {"TFrame"},
186 {"THStack"},
187 {"TMultiGraph"},
188 {"TGraphPolargram", true},
189 {"TPave", true},
190 {"TGaxis"},
191 {"TPave", true},
192 {"TArrow"},
193 {"TBox", false, true}, // can be handled via TWebPainter, disable for large number of primitives (like in greyscale.C)
194 {"TWbox"}, // some extra calls which cannot be handled via TWebPainter
195 {"TLine", false, true}, // can be handler via TWebPainter, disable for large number of primitives (like in greyscale.C)
196 {"TEllipse", true, true}, // can be handled via TWebPainter, disable for large number of primitives (like in greyscale.C)
197 {"TText"},
198 {"TLatex"},
199 {"TAnnotation"},
200 {"TMathText"},
201 {"TMarker"},
202 {"TPolyMarker"},
203 {"TPolyLine", true, true}, // can be handled via TWebPainter, simplify colors handling
204 {"TPolyMarker3D"},
205 {"TPolyLine3D"},
206 {"TGraphTime"},
207 {"TGraph2D"},
208 {"TGraph2DErrors"},
209 {"TGraphTime"},
210 {"TASImage"},
211 {"TRatioPlot"},
212 {"TSpline"},
213 {"TSpline3"},
214 {"TSpline5"},
215 {"TGeoManager"},
216 {"TGeoVolume"},
217 {}};
218
219 // fast check of class name
220 for (int i = 0; supported_classes[i].name != nullptr; ++i)
221 if ((!many_primitives || !supported_classes[i].reduse_by_many) && (strcmp(supported_classes[i].name, obj->ClassName()) == 0))
222 return kTRUE;
223
224 // now check inheritance only for configured classes
225 for (int i = 0; supported_classes[i].name != nullptr; ++i)
226 if (supported_classes[i].with_derived && (!many_primitives || !supported_classes[i].reduse_by_many))
227 if (obj->InheritsFrom(supported_classes[i].name))
228 return kTRUE;
229
230 return IsCustomClass(obj->IsA());
231}
232
233//////////////////////////////////////////////////////////////////////////////////////////////////
234/// Configures custom script for canvas.
235/// If started from "load:" or "assert:" prefix will be loaded with JSROOT.AssertPrerequisites function
236/// Script should implement custom user classes, which transferred as is to client
237/// In the script draw handler for appropriate classes would be assigned
238
239void TWebCanvas::SetCustomScripts(const std::string &src)
240{
242}
243
244//////////////////////////////////////////////////////////////////////////////////////////////////
245/// Assign custom class
246
247void TWebCanvas::AddCustomClass(const std::string &clname, bool with_derived)
248{
249 if (with_derived)
250 fCustomClasses.emplace_back("+"s + clname);
251 else
252 fCustomClasses.emplace_back(clname);
253}
254
255//////////////////////////////////////////////////////////////////////////////////////////////////
256/// Checks if class belongs to custom
257
259{
260 for (auto &name : fCustomClasses) {
261 if (name[0] == '+') {
262 if (cl->InheritsFrom(name.substr(1).c_str()))
263 return true;
264 } else if (name.compare(cl->GetName()) == 0) {
265 return true;
266 }
267 }
268 return false;
269}
270
271//////////////////////////////////////////////////////////////////////////////////////////////////
272/// Creates representation of the object for painting in web browser
273
274void TWebCanvas::CreateObjectSnapshot(TPadWebSnapshot &master, TPad *pad, TObject *obj, const char *opt, TWebPS *masterps)
275{
276 if (IsJSSupportedClass(obj, masterps != nullptr)) {
277 master.NewPrimitive(obj, opt).SetSnapshot(TWebSnapshot::kObject, obj);
278 return;
279 }
280
281 // painter is not necessary for batch canvas, but keep configuring it for a while
282 auto *painter = dynamic_cast<TWebPadPainter *>(Canvas()->GetCanvasPainter());
283
284 TView *view = nullptr;
285
287
288 gPad = pad;
289
290 if (obj->InheritsFrom(TAtt3D::Class()) && !pad->GetView()) {
291 pad->GetViewer3D("pad");
292 view = TView::CreateView(1, 0, 0); // Cartesian view by default
293 pad->SetView(view);
294
295 // Set view to perform first auto-range (scaling) pass
296 view->SetAutoRange(kTRUE);
297 }
298
299 TVirtualPS *saveps = gVirtualPS;
300
301 TWebPS ps;
302 gVirtualPS = masterps ? masterps : &ps;
303 if (painter)
304 painter->SetPainting(ps.GetPainting());
305
306 // calling Paint function for the object
307 obj->Paint(opt);
308
309 if (view) {
310 view->SetAutoRange(kFALSE);
311 // call 3D paint once again to make real drawing
312 obj->Paint(opt);
313 pad->SetView(nullptr);
314 }
315
316 if (painter)
317 painter->SetPainting(nullptr);
318
319 gVirtualPS = saveps;
320
321 fPadsStatus[pad]._has_specials = true;
322
323 // if there are master PS, do not create separate entries
324 if (!masterps && !ps.IsEmptyPainting())
326}
327
328//////////////////////////////////////////////////////////////////////////////////////////////////
329/// Calculate hash function for all colors and palette
330
332{
333 UInt_t hash = 0;
334
335 TObjArray *colors = (TObjArray *)gROOT->GetListOfColors();
336
337 if (colors) {
338 for (Int_t n = 0; n <= colors->GetLast(); ++n)
339 if (colors->At(n))
340 hash += TString::Hash(colors->At(n), TColor::Class()->Size());
341 }
342
344
345 hash += TString::Hash(pal.GetArray(), pal.GetSize() * sizeof(Int_t));
346
347 return hash;
348}
349
350
351//////////////////////////////////////////////////////////////////////////////////////////////////
352/// Add special canvas objects like colors list at selected palette
353
355{
356 TObjArray *colors = (TObjArray *)gROOT->GetListOfColors();
357
358 if (!colors)
359 return;
360
361 //Int_t cnt = 0;
362 //for (Int_t n = 0; n <= colors->GetLast(); ++n)
363 // if (colors->At(n))
364 // cnt++;
365 //if (cnt <= 598)
366 // return; // normally there are 598 colors defined
367
369
370 auto *listofcols = new TWebPainting;
371 for (Int_t n = 0; n <= colors->GetLast(); ++n)
372 if (colors->At(n))
373 listofcols->AddColor(n, (TColor *)colors->At(n));
374
375 // store palette in the buffer
376 auto *tgt = listofcols->Reserve(pal.GetSize());
377 for (Int_t i = 0; i < pal.GetSize(); i++)
378 tgt[i] = pal[i];
379 listofcols->FixSize();
380
381 master.NewSpecials().SetSnapshot(TWebSnapshot::kColors, listofcols, kTRUE);
382}
383
384//////////////////////////////////////////////////////////////////////////////////////////////////
385/// Create snapshot for pad and all primitives
386/// Callback function is used to create JSON in the middle of data processing -
387/// when all misc objects removed from canvas list of primitives or histogram list of functions
388/// After that objects are moved back to their places
389
391{
392 auto &pad_status = fPadsStatus[pad];
393
394 // send primitives if version 0 or actual pad version grater than already send version
395 bool process_primitives = (version == 0) || (pad_status.fVersion > version);
396
397 if (paddata.IsSetObjectIds()) {
398 paddata.SetActive(pad == gPad);
399 paddata.SetObjectIDAsPtr(pad);
400 }
401 paddata.SetSnapshot(TWebSnapshot::kSubPad, pad); // add ref to the pad
402 paddata.SetWithoutPrimitives(!process_primitives);
403 paddata.SetHasExecs(pad->GetListOfExecs()); // if pad execs are there provide more events from client
404
405 // check style changes every time when creating canvas snapshot
406 if (resfunc && (GetStyleDelivery() > 0)) {
407
409 auto hash = TString::Hash(gStyle, TStyle::Class()->Size());
410 if ((hash != fStyleHash) || (fStyleVersion == 0)) {
411 fStyleHash = hash;
413 }
414 }
415
416 if (fStyleVersion > version)
418 }
419
420 TList *primitives = pad->GetListOfPrimitives();
421
422 if (primitives) fPrimitivesLists.Add(primitives); // add list of primitives
423
424 TWebPS masterps;
425 bool usemaster = primitives ? (primitives->GetSize() > fPrimitivesMerge) : false;
426
427 TIter iter(primitives);
428 TObject *obj = nullptr;
429 TFrame *frame = nullptr;
430 TPaveText *title = nullptr;
431 bool need_frame = false, has_histo = false, need_palette = false;
432 std::string need_title;
433
434 while (process_primitives && ((obj = iter()) != nullptr)) {
435 TString opt = iter.GetOption();
436 opt.ToUpper();
437
438 if (obj->InheritsFrom(THStack::Class())) {
439 // workaround for THStack, create extra components before sending to client
440 auto hs = static_cast<THStack *>(obj);
441 TVirtualPad::TContext ctxt(pad, kFALSE);
442 hs->BuildPrimitives(iter.GetOption());
443 has_histo = true;
444 } else if (obj->InheritsFrom(TMultiGraph::Class())) {
445 // workaround for TMultiGraph
446 if (opt.Contains("A")) {
447 auto mg = static_cast<TMultiGraph *>(obj);
449 mg->GetHistogram(); // force creation of histogram without any drawings
450 has_histo = true;
451 if (strlen(obj->GetTitle()) > 0)
452 need_title = obj->GetTitle();
453 }
454 } else if (obj->InheritsFrom(TFrame::Class())) {
455 if (!frame)
456 frame = static_cast<TFrame *>(obj);
457 } else if (obj->InheritsFrom(TH1::Class())) {
458 need_frame = true;
459 has_histo = true;
460 if (!obj->TestBit(TH1::kNoTitle) && !opt.Contains("SAME") && !opt.Contains("AXIS") && !opt.Contains("AXIG") && (strlen(obj->GetTitle()) > 0))
461 need_title = obj->GetTitle();
462 if (obj->InheritsFrom(TH2::Class()) && (opt.Contains("COLZ") || opt.Contains("LEGO2Z") || opt.Contains("LEGO4Z") || opt.Contains("SURF2Z")))
463 need_palette = true;
464 } else if (obj->InheritsFrom(TGraph::Class())) {
465 if (opt.Contains("A")) {
466 need_frame = true;
467 if (!has_histo && (strlen(obj->GetTitle()) > 0))
468 need_title = obj->GetTitle();
469 }
470 } else if (obj->InheritsFrom(TGraph2D::Class())) {
471 if (!has_histo && (strlen(obj->GetTitle()) > 0))
472 need_title = obj->GetTitle();
473 } else if (obj->InheritsFrom(TScatter::Class())) {
474 need_frame = need_palette = true;
475 if (strlen(obj->GetTitle()) > 0)
476 need_title = obj->GetTitle();
477 } else if (obj->InheritsFrom(TF1::Class())) {
478 need_frame = !obj->InheritsFrom(TF2::Class());
479 if (!has_histo && (strlen(obj->GetTitle()) > 0))
480 need_title = obj->GetTitle();
481 } else if (obj->InheritsFrom(TPaveText::Class())) {
482 if (strcmp(obj->GetName(), "title") == 0)
483 title = static_cast<TPaveText *>(obj);
484 }
485 }
486
487 if (need_frame && !frame && primitives && CanCreateObject("TFrame")) {
488 if (!IsReadOnly() && need_palette && (pad->GetRightMargin() < 0.12) && (pad->GetRightMargin() == gStyle->GetPadRightMargin()))
489 pad->SetRightMargin(0.12);
490
491 frame = pad->GetFrame();
492 if(frame)
493 primitives->AddFirst(frame);
494 }
495
496 if (!need_title.empty() && gStyle->GetOptTitle()) {
497 if (title) {
498 auto line0 = title->GetLine(0);
499 if (line0 && !IsReadOnly()) line0->SetTitle(need_title.c_str());
500 } else if (primitives && CanCreateObject("TPaveText")) {
501 title = new TPaveText(0, 0, 0, 0, "blNDC");
504 title->SetName("title");
507 title->SetTextFont(gStyle->GetTitleFont(""));
508 if (gStyle->GetTitleFont("") % 10 > 2)
510 title->AddText(need_title.c_str());
511 title->SetBit(kCanDelete);
512 primitives->Add(title);
513 }
514 }
515
516 auto flush_master = [&]() {
517 if (!usemaster || masterps.IsEmptyPainting()) return;
518
520 masterps.CreatePainting(); // create for next operations
521 };
522
523 auto check_save_tf1 = [&](TObject *fobj, bool ignore_nodraw = false) {
524 if (!paddata.IsBatchMode() && !fTF1UseSave)
525 return;
526 if (!ignore_nodraw && fobj->TestBit(TF1::kNotDraw))
527 return;
528
529 auto f1 = static_cast<TF1 *>(fobj);
530 Double_t xmin = 0., ymin = 0., zmin = 0., xmax = 0., ymax = 0., zmax = 0.;
531 f1->GetRange(xmin, ymin, zmin, xmax, ymax, zmax);
532 f1->Save(xmin, xmax, ymin, ymax, zmin, zmax);
533 };
534
535 iter.Reset();
536
537 bool first_obj = true;
538
539 if (process_primitives)
540 pad_status._has_specials = false;
541
542 while ((obj = iter()) != nullptr) {
543 if (obj->InheritsFrom(TPad::Class())) {
544 flush_master();
545 CreatePadSnapshot(paddata.NewSubPad(), (TPad *)obj, version, nullptr);
546 } else if (!process_primitives) {
547 continue;
548 } else if (obj->InheritsFrom(TH1K::Class())) {
549 flush_master();
550 TH1K *hist = (TH1K *)obj;
551
552 Int_t nbins = hist->GetXaxis()->GetNbins();
553
554 TH1D *h1 = new TH1D("__dummy_name__", hist->GetTitle(), nbins, hist->GetXaxis()->GetXmin(), hist->GetXaxis()->GetXmax());
555 h1->SetDirectory(nullptr);
556 h1->SetName(hist->GetName());
557 hist->TAttLine::Copy(*h1);
558 hist->TAttFill::Copy(*h1);
559 hist->TAttMarker::Copy(*h1);
560 for (Int_t n = 1; n <= nbins; ++n)
561 h1->SetBinContent(n, hist->GetBinContent(n));
562
563 TIter fiter(hist->GetListOfFunctions());
564 while (auto fobj = fiter())
565 h1->GetListOfFunctions()->Add(fobj->Clone());
566
567 TString hopt = iter.GetOption();
568 if (title && first_obj) hopt.Append(";;use_pad_title");
569
571
572 } else if (obj->InheritsFrom(TH1::Class())) {
573 flush_master();
574
575 TH1 *hist = static_cast<TH1 *>(obj);
576 TIter fiter(hist->GetListOfFunctions());
577 TObject *fobj = nullptr;
578 TPaveStats *stats = nullptr;
579 TObject *palette = nullptr;
580
581 hist->BufferEmpty();
582
583 while ((fobj = fiter()) != nullptr) {
584 if (fobj->InheritsFrom(TPaveStats::Class()))
585 stats = dynamic_cast<TPaveStats *> (fobj);
586 else if (fobj->InheritsFrom("TPaletteAxis"))
587 palette = fobj;
588 else if (fobj->InheritsFrom(TF1::Class()))
589 check_save_tf1(fobj);
590 }
591
592 TString hopt = iter.GetOption();
593 TString o = hopt;
594 o.ToUpper();
595
596 if (!stats && (first_obj || o.Contains("SAMES")) && (gStyle->GetOptStat() > 0) && CanCreateObject("TPaveStats")) {
597 stats = new TPaveStats(
600 gStyle->GetStatX(),
601 gStyle->GetStatY(), "brNDC");
602
603 stats->SetParent(hist);
604 // do not set optfit and optstat, they calling pad->Update,
605 // values correctly set already in TPaveStats constructor
606 // stats->SetOptFit(gStyle->GetOptFit());
607 // stats->SetOptStat(gStyle->GetOptStat());
611 stats->SetTextFont(gStyle->GetStatFont());
612 if (gStyle->GetStatFont()%10 > 2)
616 stats->SetName("stats");
617
619 stats->SetTextAlign(12);
620 stats->SetBit(kCanDelete);
621 stats->SetBit(kMustCleanup);
622
623 hist->GetListOfFunctions()->Add(stats);
624 }
625
626 if (!palette && CanCreateObject("TPaletteAxis") && (hist->GetDimension() > 1) &&
627 (o.Contains("COLZ") || o.Contains("LEGO2Z") || o.Contains("LEGO4Z") || o.Contains("SURF2Z"))) {
628 std::stringstream exec;
629 exec << "new TPaletteAxis(0,0,0,0, (TH1*)" << std::hex << std::showbase << (size_t)hist << ");";
630 palette = (TObject *)gROOT->ProcessLine(exec.str().c_str());
631 if (palette)
632 hist->GetListOfFunctions()->AddFirst(palette);
633 }
634
635 if (title && first_obj) hopt.Append(";;use_pad_title");
636
637 // if (stats) hopt.Append(";;use_pad_stats");
638
639 if (palette) hopt.Append(";;use_pad_palette");
640
641 paddata.NewPrimitive(obj, hopt.Data()).SetSnapshot(TWebSnapshot::kObject, obj);
642
643 if (hist->GetDimension() == 2) {
644 TString opt = iter.GetOption();
645 auto p1 = opt.Index("["), p2 = opt.Index("]");
646 if ((p1 != kNPOS) && (p2 != kNPOS) && p2 > p1 + 1) {
647 TString cutname = opt(p1 + 1, p2 - p1 - 1);
648 TObject *cutg = primitives->FindObject(cutname.Data());
649 if (!cutg || (cutg->IsA() != TCutG::Class())) {
650 cutg = gROOT->GetListOfSpecials()->FindObject(cutname.Data());
651 if (cutg && cutg->IsA() == TCutG::Class())
652 paddata.NewPrimitive(cutg, "__ignore_drawing__").SetSnapshot(TWebSnapshot::kObject, cutg);
653 }
654 }
655 }
656
657 // do not extract objects from list of functions - stats and func need to be handled together with hist
658 //
659 // fiter.Reset();
660 // while ((fobj = fiter()) != nullptr)
661 // CreateObjectSnapshot(paddata, pad, fobj, fiter.GetOption());
662
663 // fPrimitivesLists.Add(hist->GetListOfFunctions());
664
665 first_obj = false;
666 } else if (obj->InheritsFrom(TGraph::Class())) {
667 flush_master();
668
669 TGraph *gr = static_cast<TGraph *>(obj);
670 auto funcs = gr->GetListOfFunctions();
671
672 TIter fiter(funcs);
673 TPaveStats *stats = nullptr;
674
675 while (auto fobj = fiter()) {
676 if (fobj->InheritsFrom(TPaveStats::Class()))
677 stats = dynamic_cast<TPaveStats *> (fobj);
678 else if (fobj->InheritsFrom(TF1::Class()))
679 check_save_tf1(fobj);
680 }
681
682 TString gropt = iter.GetOption();
683
684 // ensure histogram exists on server to draw it properly on clients side
685 if (!IsReadOnly() && (first_obj || gropt.Index("A", 0, TString::kIgnoreCase) != kNPOS ||
686 (gropt.Index("X+", 0, TString::kIgnoreCase) != kNPOS) || (gropt.Index("X+", 0, TString::kIgnoreCase) != kNPOS)))
687 gr->GetHistogram();
688
689 if (title && first_obj) gropt.Append(";;use_pad_title");
690 if (stats) gropt.Append(";;use_pad_stats");
691
692 paddata.NewPrimitive(obj, gropt.Data()).SetSnapshot(TWebSnapshot::kObject, obj);
693
694 fiter.Reset();
695 while (auto fobj = fiter())
696 CreateObjectSnapshot(paddata, pad, fobj, fiter.GetOption());
697
698 if (funcs)
700 first_obj = false;
701 } else if (obj->InheritsFrom(TGraph2D::Class())) {
702 flush_master();
703
704 TGraph2D *gr2d = static_cast<TGraph2D *>(obj);
705
706 // ensure correct range of histogram
707 if (!IsReadOnly() && first_obj) {
708 TString gropt = iter.GetOption();
709 gropt.ToUpper();
710 Bool_t zscale = gropt.Contains("TRI1") || gropt.Contains("TRI2") || gropt.Contains("COL");
711 Bool_t real_draw = gropt.Contains("TRI") || gropt.Contains("LINE") || gropt.Contains("ERR") || gropt.Contains("P0");
712
713 TString hopt = !real_draw ? iter.GetOption() : (zscale ? "lego2z" : "lego2");
714 if (title) hopt.Append(";;use_pad_title");
715
716 // if gr2d not draw - let create histogram with correspondent content
717 auto hist = gr2d->GetHistogram(real_draw ? "empty" : "");
718
719 paddata.NewPrimitive(gr2d, hopt.Data(), "#hist").SetSnapshot(TWebSnapshot::kObject, hist);
720 }
721
722 paddata.NewPrimitive(obj, iter.GetOption()).SetSnapshot(TWebSnapshot::kObject, obj);
723 first_obj = false;
724 } else if (obj->InheritsFrom(TMultiGraph::Class())) {
725 flush_master();
726
727 TMultiGraph *mgr = static_cast<TMultiGraph *>(obj);
728 auto funcs = mgr->GetListOfFunctions();
729
730 TIter fiter(funcs);
731 while (auto fobj = fiter()) {
732 if (fobj->InheritsFrom(TF1::Class()))
733 check_save_tf1(fobj);
734 }
735
736 paddata.NewPrimitive(obj, iter.GetOption()).SetSnapshot(TWebSnapshot::kObject, obj);
737
738 fiter.Reset();
739 while (auto fobj = fiter())
740 CreateObjectSnapshot(paddata, pad, fobj, fiter.GetOption());
741
742 if (funcs)
744
745 first_obj = false;
746 } else if (obj->InheritsFrom(TScatter::Class())) {
747 flush_master();
748
749 TScatter *scatter = (TScatter *)obj;
750 auto funcs = scatter->GetGraph()->GetListOfFunctions();
751
752 TIter fiter(funcs);
753 TObject *fobj = nullptr, *palette = nullptr;
754
755 while ((fobj = fiter()) != nullptr) {
756 if (fobj->InheritsFrom("TPaletteAxis"))
757 palette = fobj;
758 }
759
760 // ensure histogram exists on server to draw it properly on clients side
761 if (!IsReadOnly() && first_obj)
762 scatter->GetHistogram();
763
764 if (!palette && CanCreateObject("TPaletteAxis")) {
765 std::stringstream exec;
766 exec << "new TPaletteAxis(0,0,0,0,0,0);";
767 palette = (TObject *)gROOT->ProcessLine(exec.str().c_str());
768 if (palette)
769 funcs->AddFirst(palette);
770 }
771
772 TString scopt = iter.GetOption();
773 if (title && first_obj) scopt.Append(";;use_pad_title");
774 if (palette) scopt.Append(";;use_pad_palette");
775
776 paddata.NewPrimitive(obj, scopt.Data()).SetSnapshot(TWebSnapshot::kObject, obj);
777
778 fiter.Reset();
779 while ((fobj = fiter()) != nullptr)
780 CreateObjectSnapshot(paddata, pad, fobj, fiter.GetOption());
781
782 if (funcs)
784
785 first_obj = false;
786 } else if (obj->InheritsFrom(TF1::Class())) {
787 flush_master();
788 auto f1 = static_cast<TF1 *> (obj);
789
790 TString f1opt = iter.GetOption();
791
792 check_save_tf1(obj, true);
793 // if (fTF1UseSave)
794 // f1opt.Append(";force_saved");
795
796 if (first_obj) {
797 auto hist = f1->GetHistogram();
798 paddata.NewPrimitive(hist, "__ignore_drawing__").SetSnapshot(TWebSnapshot::kObject, hist);
799 f1opt.Append(";webcanv_hist");
800 }
801
803
804 first_obj = false;
805
806 } else if (obj->InheritsFrom(TGaxis::Class())) {
807 flush_master();
808 auto gaxis = static_cast<TGaxis *> (obj);
809 auto func = gaxis->GetFunction();
810 if (func)
811 paddata.NewPrimitive(func, "__ignore_drawing__").SetSnapshot(TWebSnapshot::kObject, func);
812
813 paddata.NewPrimitive(obj, iter.GetOption()).SetSnapshot(TWebSnapshot::kObject, obj);
814 } else if (obj->InheritsFrom(TFrame::Class())) {
815 flush_master();
816 if (frame && (obj == frame)) {
817 paddata.NewPrimitive(obj, iter.GetOption()).SetSnapshot(TWebSnapshot::kObject, obj);
818 frame = nullptr; // add frame only once
819 }
820 } else if (IsJSSupportedClass(obj, usemaster)) {
821 flush_master();
822 paddata.NewPrimitive(obj, iter.GetOption()).SetSnapshot(TWebSnapshot::kObject, obj);
823 } else {
824 CreateObjectSnapshot(paddata, pad, obj, iter.GetOption(), usemaster ? &masterps : nullptr);
825 }
826 }
827
828 flush_master();
829
830 bool provide_colors = false;
831
832 if ((GetPaletteDelivery() > 2) || ((GetPaletteDelivery() == 2) && resfunc)) {
833 // provide colors: either for each subpad (> 2) or only for canvas (== 2)
834 provide_colors = process_primitives;
835 } else if ((GetPaletteDelivery() == 1) && resfunc) {
836 // check that colors really changing, using hash
837
839 auto hash = CalculateColorsHash();
840 if ((hash != fColorsHash) || (fColorsVersion == 0)) {
841 fColorsHash = hash;
843 }
844 }
845
846 provide_colors = fColorsVersion > version;
847 }
848
849 // add colors after painting is performed - new colors may be generated only during painting
850 if (provide_colors)
851 AddColorsPalette(paddata);
852
853 if (!resfunc)
854 return;
855
856 // now move all primitives and functions into separate list to perform I/O
857
858 TList save_lst;
859 TIter diter(&fPrimitivesLists);
860 TList *dlst = nullptr;
861 while ((dlst = (TList *)diter()) != nullptr) {
862 TIter fiter(dlst);
863 while ((obj = fiter()) != nullptr)
864 save_lst.Add(obj, fiter.GetOption());
865 save_lst.Add(dlst); // add list itself to have marker
866 dlst->Clear("nodelete");
867 }
868
869 // execute function to prevent storing of colors with custom TCanvas streamer
870 // TODO: Olivier - we need to change logic here!
872
873 // invoke callback for master painting
874 resfunc(&paddata);
875
876 TIter siter(&save_lst);
877 diter.Reset();
878 while ((dlst = (TList *)diter()) != nullptr) {
879 while ((obj = siter()) != nullptr) {
880 if (obj == dlst)
881 break;
882 dlst->Add(obj, siter.GetOption());
883 }
884 }
885
886 save_lst.Clear("nodelete");
887
888 fPrimitivesLists.Clear("nodelete");
889}
890
891//////////////////////////////////////////////////////////////////////////////////////////////////
892/// Add control message for specified connection
893/// Same control message can be overwritten many time before it really sends to the client
894/// If connid == 0, message will be add to all connections
895/// After ctrl message is add to the output, short timer is activated and message send afterwards
896
897void TWebCanvas::AddCtrlMsg(unsigned connid, const std::string &key, const std::string &value)
898{
899 Bool_t new_ctrl = kFALSE;
900
901 for (auto &conn : fWebConn) {
902 if (conn.match(connid)) {
903 conn.fCtrl[key] = value;
904 new_ctrl = kTRUE;
905 }
906 }
907
908 if (new_ctrl && fTimer->IsSlow())
910}
911
912
913//////////////////////////////////////////////////////////////////////////////////////////////////
914/// Add message to send queue for specified connection
915/// If connid == 0, message will be add to all connections
916
917void TWebCanvas::AddSendQueue(unsigned connid, const std::string &msg)
918{
919 for (auto &conn : fWebConn) {
920 if (conn.match(connid))
921 conn.fSend.emplace(msg);
922 }
923}
924
925
926//////////////////////////////////////////////////////////////////////////////////////////////////
927/// Check if any data should be send to client
928/// If connid != 0, only selected connection will be checked
929
931{
932 if (!Canvas())
933 return kFALSE;
934
935 bool isMoreData = false, isAnySend = false;
936
937 for (auto &conn : fWebConn) {
938
939 bool isConnData = !conn.fCtrl.empty() || !conn.fSend.empty() ||
940 ((conn.fCheckedVersion < fCanvVersion) && (conn.fSendVersion == conn.fDrawVersion));
941
942 while ((conn.is_batch() && !connid) || (conn.match(connid) && fWindow && fWindow->CanSend(conn.fConnId, true))) {
943 // check if any control messages still there to keep timer running
944
945 std::string buf;
946
947 if (!conn.fCtrl.empty()) {
949 conn.fCtrl.clear();
950 } else if (!conn.fSend.empty()) {
951 std::swap(buf, conn.fSend.front());
952 conn.fSend.pop();
953 } else if ((conn.fCheckedVersion < fCanvVersion) && (conn.fSendVersion == conn.fDrawVersion)) {
954
955 buf = "SNAP6:"s + std::to_string(fCanvVersion) + ":"s;
956
957 TCanvasWebSnapshot holder(IsReadOnly(), true, false); // readonly, set ids, batchmode
958
959 holder.SetFixedSize(fFixedSize); // set fixed size flag
960
961 // scripts send only when canvas drawn for the first time
962 if (!conn.fSendVersion)
964
965 holder.SetHighlightConnect(Canvas()->HasConnection("Highlighted(TVirtualPad*,TObject*,Int_t,Int_t)"));
966
967 CreatePadSnapshot(holder, Canvas(), conn.fSendVersion, [&buf, &conn, this](TPadWebSnapshot *snap) {
968 if (conn.is_batch()) {
969 // for batch connection only calling of CreatePadSnapshot is important
970 buf.clear();
971 return;
972 }
973
974 auto json = TBufferJSON::ToJSON(snap, fJsonComp);
975 auto hash = json.Hash();
976 if (conn.fLastSendHash && (conn.fLastSendHash == hash) && conn.fSendVersion) {
977 // prevent looping when same data send many times
978 buf.clear();
979 } else {
980 buf.append(json.Data());
981 conn.fLastSendHash = hash;
982 }
983 });
984
985 conn.fCheckedVersion = fCanvVersion;
986
987 conn.fSendVersion = fCanvVersion;
988
989 if (buf.empty())
990 conn.fDrawVersion = fCanvVersion;
991 } else {
992 isConnData = false;
993 break;
994 }
995
996 if (!buf.empty() && !conn.is_batch()) {
997 fWindow->Send(conn.fConnId, buf);
998 isAnySend = true;
999 }
1000 }
1001
1002 if (isConnData)
1003 isMoreData = true;
1004 }
1005
1006 if (fTimer->IsSlow() && isMoreData)
1007 fTimer->SetSlow(kFALSE);
1008
1009 return isAnySend;
1010}
1011
1012//////////////////////////////////////////////////////////////////////////////////////////
1013/// Close web canvas - not implemented
1014
1016{
1017}
1018
1019//////////////////////////////////////////////////////////////////////////////////////////
1020/// Show canvas in specified place.
1021/// If parameter args not specified, default ROOT web display will be used
1022
1024{
1025 if (!fWindow) {
1027
1028 fWindow->SetConnLimit(0); // configure connections limit
1029
1030 fWindow->SetDefaultPage("file:rootui5sys/canv/canvas6.html");
1031
1032 fWindow->SetCallBacks(
1033 // connection
1034 [this](unsigned connid) {
1035 fWebConn.emplace_back(connid);
1036 CheckDataToSend(connid);
1037 },
1038 // data
1039 [this](unsigned connid, const std::string &arg) {
1040 ProcessData(connid, arg);
1042 },
1043 // disconnect
1044 [this](unsigned connid) {
1045 unsigned indx = 0;
1046 for (auto &c : fWebConn) {
1047 if (c.fConnId == connid) {
1048 fWebConn.erase(fWebConn.begin() + indx);
1049 break;
1050 }
1051 indx++;
1052 }
1053 });
1054 }
1055
1056 auto w = Canvas()->GetWindowWidth(), h = Canvas()->GetWindowHeight();
1057 if ((w > 0) && (w < 50000) && (h > 0) && (h < 30000))
1058 fWindow->SetGeometry(w, h);
1059
1064
1065 fWindow->Show(args);
1066}
1067
1068//////////////////////////////////////////////////////////////////////////////////////////
1069/// Show canvas in browser window
1070
1072{
1073 if (gROOT->IsWebDisplayBatch())
1074 return;
1075
1077 args.SetWidgetKind("TCanvas");
1078 args.SetSize(Canvas()->GetWindowWidth(), Canvas()->GetWindowHeight());
1079 args.SetPos(Canvas()->GetWindowTopX(), Canvas()->GetWindowTopY());
1080
1081 ShowWebWindow(args);
1082}
1083
1084//////////////////////////////////////////////////////////////////////////////////////////
1085/// Function used to send command to browser to toggle menu, toolbar, editors, ...
1086
1087void TWebCanvas::ShowCmd(const std::string &arg, Bool_t show)
1088{
1089 AddCtrlMsg(0, arg, show ? "1"s : "0"s);
1090}
1091
1092//////////////////////////////////////////////////////////////////////////////////////////
1093/// Activate object in editor in web browser
1094
1096{
1097 if (!pad || !obj) return;
1098
1099 UInt_t hash = TString::Hash(&obj, sizeof(obj));
1100
1101 AddCtrlMsg(0, "edit"s, std::to_string(hash));
1102}
1103
1104//////////////////////////////////////////////////////////////////////////////////////////
1105/// Returns kTRUE if web canvas has graphical editor
1106
1108{
1109 return (fClientBits & TCanvas::kShowEditor) != 0;
1110}
1111
1112//////////////////////////////////////////////////////////////////////////////////////////
1113/// Returns kTRUE if web canvas has menu bar
1114
1116{
1117 return (fClientBits & TCanvas::kMenuBar) != 0;
1118}
1119
1120//////////////////////////////////////////////////////////////////////////////////////////
1121/// Returns kTRUE if web canvas has status bar
1122
1124{
1125 return (fClientBits & TCanvas::kShowEventStatus) != 0;
1126}
1127
1128//////////////////////////////////////////////////////////////////////////////////////////
1129/// Returns kTRUE if tooltips are activated in web canvas
1130
1132{
1133 return (fClientBits & TCanvas::kShowToolTips) != 0;
1134}
1135
1136//////////////////////////////////////////////////////////////////////////////////////////
1137/// Set window position of web canvas
1138
1140{
1141 AddCtrlMsg(0, "x"s, std::to_string(x));
1142 AddCtrlMsg(0, "y"s, std::to_string(y));
1143}
1144
1145//////////////////////////////////////////////////////////////////////////////////////////
1146/// Set window size of web canvas
1147
1149{
1150 AddCtrlMsg(0, "w"s, std::to_string(w));
1151 AddCtrlMsg(0, "h"s, std::to_string(h));
1152}
1153
1154//////////////////////////////////////////////////////////////////////////////////////////
1155/// Set window title of web canvas
1156
1157void TWebCanvas::SetWindowTitle(const char *newTitle)
1158{
1159 AddCtrlMsg(0, "title"s, newTitle);
1160}
1161
1162//////////////////////////////////////////////////////////////////////////////////////////
1163/// Set canvas size of web canvas
1164
1166{
1167 fFixedSize = kTRUE;
1168 AddCtrlMsg(0, "cw"s, std::to_string(cw));
1169 AddCtrlMsg(0, "ch"s, std::to_string(ch));
1170 if ((cw > 0) && (ch > 0)) {
1171 Canvas()->fCw = cw;
1172 Canvas()->fCh = ch;
1173 } else {
1174 // temporary value, will be reported back from client
1175 Canvas()->fCw = Canvas()->fWindowWidth;
1177 }
1178}
1179
1180//////////////////////////////////////////////////////////////////////////////////////////
1181/// Iconify browser window
1182
1184{
1185 AddCtrlMsg(0, "winstate"s, "iconify"s);
1186}
1187
1188//////////////////////////////////////////////////////////////////////////////////////////
1189/// Raise browser window
1190
1192{
1193 AddCtrlMsg(0, "winstate"s, "raise"s);
1194}
1195
1196//////////////////////////////////////////////////////////////////////////////////////////
1197/// Assign clients bits
1198
1200{
1201 fClientBits = bits;
1206}
1207
1208//////////////////////////////////////////////////////////////////////////////////////////////////
1209/// Decode all pad options, which includes ranges plus objects options
1210
1211Bool_t TWebCanvas::DecodePadOptions(const std::string &msg, bool process_execs)
1212{
1213 if (IsReadOnly() || msg.empty())
1214 return kFALSE;
1215
1216 auto arr = TBufferJSON::FromJSON<std::vector<TWebPadOptions>>(msg);
1217
1218 if (!arr)
1219 return kFALSE;
1220
1221 Bool_t need_update = kFALSE;
1222
1223 TPad *pad_with_execs = nullptr;
1224 TExec *hist_exec = nullptr;
1225
1226 for (unsigned n = 0; n < arr->size(); ++n) {
1227 auto &r = arr->at(n);
1228 TPad *pad = dynamic_cast<TPad *>(FindPrimitive(r.snapid));
1229
1230 if (!pad)
1231 continue;
1232
1233 if (pad == Canvas()) {
1234 AssignStatusBits(r.bits);
1235 Canvas()->fCw = r.cw;
1236 Canvas()->fCh = r.ch;
1237 if (r.w.size() == 4) {
1238 fWindowGeometry = r.w;
1243 }
1244 }
1245
1246 if (r.active && (pad != gPad)) gPad = pad;
1247
1248 if ((pad->GetTickx() != r.tickx) || (pad->GetTicky() != r.ticky))
1249 pad->SetTicks(r.tickx, r.ticky);
1250 if ((pad->GetGridx() != (r.gridx > 0)) || (pad->GetGridy() != (r.gridy > 0)))
1251 pad->SetGrid(r.gridx, r.gridy);
1252 pad->fLogx = r.logx;
1253 pad->fLogy = r.logy;
1254 pad->fLogz = r.logz;
1255
1256 pad->SetLeftMargin(r.mleft);
1257 pad->SetRightMargin(r.mright);
1258 pad->SetTopMargin(r.mtop);
1259 pad->SetBottomMargin(r.mbottom);
1260
1261 if (r.ranges) {
1262 // avoid call of original methods, set members directly
1263 // pad->Range(r.px1, r.py1, r.px2, r.py2);
1264 // pad->RangeAxis(r.ux1, r.uy1, r.ux2, r.uy2);
1265
1266 pad->fX1 = r.px1;
1267 pad->fX2 = r.px2;
1268 pad->fY1 = r.py1;
1269 pad->fY2 = r.py2;
1270
1271 pad->fUxmin = r.ux1;
1272 pad->fUxmax = r.ux2;
1273 pad->fUymin = r.uy1;
1274 pad->fUymax = r.uy2;
1275 }
1276
1277 // pad->SetPad(r.mleft, r.mbottom, 1-r.mright, 1-r.mtop);
1278
1279 pad->fAbsXlowNDC = r.xlow;
1280 pad->fAbsYlowNDC = r.ylow;
1281 pad->fAbsWNDC = r.xup - r.xlow;
1282 pad->fAbsHNDC = r.yup - r.ylow;
1283
1284 if (pad == Canvas()) {
1285 pad->fXlowNDC = r.xlow;
1286 pad->fYlowNDC = r.ylow;
1287 pad->fXUpNDC = r.xup;
1288 pad->fYUpNDC = r.yup;
1289 pad->fWNDC = r.xup - r.xlow;
1290 pad->fHNDC = r.yup - r.ylow;
1291 } else {
1292 auto mother = pad->GetMother();
1293 if (mother->GetAbsWNDC() > 0. && mother->GetAbsHNDC() > 0.) {
1294 pad->fXlowNDC = (r.xlow - mother->GetAbsXlowNDC()) / mother->GetAbsWNDC();
1295 pad->fYlowNDC = (r.ylow - mother->GetAbsYlowNDC()) / mother->GetAbsHNDC();
1296 pad->fXUpNDC = (r.xup - mother->GetAbsXlowNDC()) / mother->GetAbsWNDC();
1297 pad->fYUpNDC = (r.yup - mother->GetAbsYlowNDC()) / mother->GetAbsHNDC();
1298 pad->fWNDC = (r.xup - r.xlow) / mother->GetAbsWNDC();
1299 pad->fHNDC = (r.yup - r.ylow) / mother->GetAbsHNDC();
1300 }
1301 }
1302
1303 // copy of code from TPad::ResizePad()
1304
1305 Double_t pxlow = r.xlow * r.cw;
1306 Double_t pylow = (1-r.ylow) * r.ch;
1307 Double_t pxrange = (r.xup - r.xlow) * r.cw;
1308 Double_t pyrange = -1*(r.yup - r.ylow) * r.ch;
1309
1310 Double_t rounding = 0.00005;
1311 Double_t xrange = r.px2 - r.px1;
1312 Double_t yrange = r.py2 - r.py1;
1313
1314 if ((xrange != 0.) && (pxrange != 0)) {
1315 // Linear X axis
1316 pad->fXtoAbsPixelk = rounding + pxlow - pxrange*r.px1/xrange; //origin at left
1317 pad->fXtoPixelk = rounding + -pxrange*r.px1/xrange;
1318 pad->fXtoPixel = pxrange/xrange;
1319 pad->fAbsPixeltoXk = r.px1 - pxlow*xrange/pxrange;
1320 pad->fPixeltoXk = r.px1;
1321 pad->fPixeltoX = xrange/pxrange;
1322 }
1323
1324 if ((yrange != 0.) && (pyrange != 0.)) {
1325 // Linear Y axis
1326 pad->fYtoAbsPixelk = rounding + pylow - pyrange*r.py1/yrange; //origin at top
1327 pad->fYtoPixelk = rounding + -pyrange - pyrange*r.py1/yrange;
1328 pad->fYtoPixel = pyrange/yrange;
1329 pad->fAbsPixeltoYk = r.py1 - pylow*yrange/pyrange;
1330 pad->fPixeltoYk = r.py1;
1331 pad->fPixeltoY = yrange/pyrange;
1332 }
1333
1335
1336 TObjLink *objlnk = nullptr;
1337
1338 TH1 *hist = static_cast<TH1 *>(FindPrimitive("histogram", 1, pad, &objlnk));
1339
1340 if (hist) {
1341
1342 TObject *hist_holder = objlnk ? objlnk->GetObject() : nullptr;
1343 if (hist_holder == hist)
1344 hist_holder = nullptr;
1345
1346 Bool_t no_entries = hist->GetEntries();
1347
1348 Double_t hmin = 0., hmax = 0.;
1349
1350 if (r.zx1 == r.zx2)
1351 hist->GetXaxis()->SetRange(0,0);
1352 else
1353 hist->GetXaxis()->SetRangeUser(r.zx1, r.zx2);
1354
1355 if (hist->GetDimension() == 1) {
1356 hmin = r.zy1;
1357 hmax = r.zy2;
1358 if ((hmin == hmax) && !no_entries) {
1359 // if there are no zooming on Y and histogram has no entries, hmin/hmax should be set to full range
1360 hmin = pad->fLogy ? TMath::Power(pad->fLogy < 2 ? 10 : pad->fLogy, r.uy1) : r.uy1;
1361 hmax = pad->fLogy ? TMath::Power(pad->fLogy < 2 ? 10 : pad->fLogy, r.uy2) : r.uy2;
1362 }
1363 } else if (r.zy1 == r.zy2) {
1364 hist->GetYaxis()->SetRange(0., 0.);
1365 } else {
1366 hist->GetYaxis()->SetRangeUser(r.zy1, r.zy2);
1367 }
1368
1369 if (hist->GetDimension() == 2) {
1370 hmin = r.zz1;
1371 hmax = r.zz2;
1372 if ((hmin == hmax) && !no_entries) {
1373 // z scale is not transformed
1374 hmin = r.uz1;
1375 hmax = r.uz2;
1376 }
1377 } else if (hist->GetDimension() == 3) {
1378 if (r.zz1 == r.zz2) {
1379 hist->GetZaxis()->SetRange(0., 0.);
1380 } else {
1381 hist->GetZaxis()->SetRangeUser(r.zz1, r.zz2);
1382 }
1383 }
1384
1385 if (hmin == hmax)
1386 hmin = hmax = -1111;
1387
1388 if (!hist_holder || (hist_holder->IsA() == TScatter::Class())) {
1389 hist->SetMinimum(hmin);
1390 hist->SetMaximum(hmax);
1391 } else {
1392 auto SetMember = [hist_holder](const char *name, Double_t value) {
1393 auto offset = hist_holder->IsA()->GetDataMemberOffset(name);
1394 if (offset > 0)
1395 *((Double_t *)((char*) hist_holder + offset)) = value;
1396 else
1397 ::Error("SetMember", "Cannot find %s data member in %s", name, hist_holder->ClassName());
1398 };
1399
1400 // directly set min/max in classes like THStack, TGraph, TMultiGraph
1401 SetMember("fMinimum", hmin);
1402 SetMember("fMaximum", hmax);
1403 }
1404
1405 TIter next(hist->GetListOfFunctions());
1406 while (auto fobj = next())
1407 if (!hist_exec && fobj->InheritsFrom(TExec::Class())) {
1408 hist_exec = (TExec *) fobj;
1409 need_update = kTRUE;
1410 }
1411 }
1412
1413 std::map<std::string, int> idmap;
1414
1415 for (auto &item : r.primitives) {
1416 auto iter = idmap.find(item.snapid);
1417 int idcnt = 1;
1418 if (iter == idmap.end())
1419 idmap[item.snapid] = 1;
1420 else
1421 idcnt = ++iter->second;
1422
1423 ProcessObjectOptions(item, pad, idcnt);
1424 }
1425
1426 // without special objects no need for explicit update of the pad
1427 if (fPadsStatus[pad]._has_specials) {
1428 pad->Modified(kTRUE);
1429 need_update = kTRUE;
1430 }
1431
1432 if (process_execs && (gPad == pad))
1433 pad_with_execs = pad;
1434 }
1435
1436 ProcessExecs(pad_with_execs, hist_exec);
1437
1438 if (fUpdatedSignal) fUpdatedSignal(); // invoke signal
1439
1440 return need_update;
1441}
1442
1443//////////////////////////////////////////////////////////////////////////////////////////////////
1444/// Process TExec objects in the pad
1445
1447{
1448 auto execs = pad ? pad->GetListOfExecs() : nullptr;
1449
1450 if ((!execs || !execs->GetSize()) && !extra)
1451 return;
1452
1453 auto saveps = gVirtualPS;
1454 TWebPS ps;
1455 gVirtualPS = &ps;
1456
1457 auto savex = gVirtualX;
1458 TVirtualX x;
1459 gVirtualX = &x;
1460
1461 TIter next(execs);
1462 while (auto obj = next()) {
1463 auto exec = dynamic_cast<TExec *>(obj);
1464 if (exec)
1465 exec->Exec();
1466 }
1467
1468 if (extra)
1469 extra->Exec();
1470
1471 gVirtualPS = saveps;
1472 gVirtualX = savex;
1473}
1474
1475//////////////////////////////////////////////////////////////////////////////////////////
1476/// Execute one or several methods for selected object
1477/// String can be separated by ";;" to let execute several methods at once
1478void TWebCanvas::ProcessLinesForObject(TObject *obj, const std::string &lines)
1479{
1480 std::string buf = lines;
1481
1482 Int_t indx = 0;
1483
1484 while (obj && !buf.empty()) {
1485 std::string sub = buf;
1486 auto pos = buf.find(";;");
1487 if (pos == std::string::npos) {
1488 sub = buf;
1489 buf.clear();
1490 } else {
1491 sub = buf.substr(0,pos);
1492 buf = buf.substr(pos+2);
1493 }
1494 if (sub.empty()) continue;
1495
1496 std::stringstream exec;
1497 exec << "((" << obj->ClassName() << " *) " << std::hex << std::showbase << (size_t)obj << ")->" << sub << ";";
1498 if (indx < 3 || gDebug > 0)
1499 Info("ProcessLinesForObject", "Obj %s Execute %s", obj->GetName(), exec.str().c_str());
1500 gROOT->ProcessLine(exec.str().c_str());
1501 indx++;
1502 }
1503}
1504
1505//////////////////////////////////////////////////////////////////////////////////////////
1506/// Handle data from web browser
1507/// Returns kFALSE if message was not processed
1508
1509Bool_t TWebCanvas::ProcessData(unsigned connid, const std::string &arg)
1510{
1511 if (arg.empty())
1512 return kTRUE;
1513
1514 // try to identify connection for given WS request
1515 unsigned indx = 0; // first connection is batch and excluded
1516 while(++indx < fWebConn.size()) {
1517 if (fWebConn[indx].fConnId == connid)
1518 break;
1519 }
1520 if (indx >= fWebConn.size())
1521 return kTRUE;
1522
1523 Bool_t is_main_connection = indx == 1; // first connection allow to make changes
1524
1525 struct FlagGuard {
1526 Bool_t &flag;
1527 FlagGuard(Bool_t &_flag) : flag(_flag) { flag = true; }
1528 ~FlagGuard() { flag = false; }
1529 };
1530
1531 FlagGuard guard(fProcessingData);
1532
1533 const char *cdata = arg.c_str();
1534
1535 if (arg == "KEEPALIVE") {
1536 // do nothing
1537
1538 } else if (arg == "QUIT") {
1539
1540 // use window manager to correctly terminate http server
1541 fWindow->TerminateROOT();
1542
1543 } else if (arg.compare(0, 7, "READY6:") == 0) {
1544
1545 // this is reply on drawing of ROOT6 snapshot
1546 // it confirms when drawing of specific canvas version is completed
1547
1548 cdata += 7;
1549
1550 const char *separ = strchr(cdata, ':');
1551 if (!separ) {
1552 fWebConn[indx].fDrawVersion = std::stoll(cdata);
1553 } else {
1554 fWebConn[indx].fDrawVersion = std::stoll(std::string(cdata, separ - cdata));
1555 if (is_main_connection && !IsReadOnly())
1556 if (DecodePadOptions(separ+1, false))
1558 }
1559
1560 } else if (arg == "RELOAD") {
1561
1562 // trigger reload of canvas data
1563 fWebConn[indx].reset();
1564
1565 } else if (arg.compare(0, 5, "SAVE:") == 0) {
1566
1567 // save image produced by the client side - like png or svg
1568 const char *img = cdata + 5;
1569
1570 const char *separ = strchr(img, ':');
1571 if (separ) {
1572 TString filename(img, separ - img);
1573 img = separ + 1;
1574
1575 std::ofstream ofs(filename.Data());
1576
1577 if (filename.Index(".svg") != kNPOS) {
1578 // ofs << "<?xml version=\"1.0\" standalone=\"no\"?>";
1579 ofs << img;
1580 } else {
1581 TString binary = TBase64::Decode(img);
1582 ofs.write(binary.Data(), binary.Length());
1583 }
1584 ofs.close();
1585
1586 Info("ProcessData", "File %s has been created", filename.Data());
1587 }
1588
1589 } else if (arg.compare(0, 8, "PRODUCE:") == 0) {
1590
1591 // create ROOT, PDF, ... files using native ROOT functionality
1592 Canvas()->Print(arg.c_str() + 8);
1593
1594 } else if (arg.compare(0, 9, "OPTIONS6:") == 0) {
1595
1596 if (is_main_connection && !IsReadOnly())
1597 if (DecodePadOptions(arg.substr(9), true))
1599
1600 } else if (arg.compare(0, 11, "STATUSBITS:") == 0) {
1601
1602 if (is_main_connection) {
1603 AssignStatusBits(std::stoul(arg.substr(11)));
1604 if (fUpdatedSignal) fUpdatedSignal(); // invoke signal
1605 }
1606 } else if (arg.compare(0, 10, "HIGHLIGHT:") == 0) {
1607 if (is_main_connection) {
1608 auto arr = TBufferJSON::FromJSON<std::vector<std::string>>(arg.substr(10));
1609 if (!arr || (arr->size() != 4)) {
1610 Error("ProcessData", "Wrong arguments count %d in highlight message", (int)(arr ? arr->size() : -1));
1611 } else {
1612 auto pad = dynamic_cast<TVirtualPad *>(FindPrimitive(arr->at(0)));
1613 auto obj = FindPrimitive(arr->at(1));
1614 int argx = std::stoi(arr->at(2));
1615 int argy = std::stoi(arr->at(3));
1616 if (pad && obj) {
1617 Canvas()->Highlighted(pad, obj, argx, argy);
1619 }
1620 }
1621 }
1622 } else if (ROOT::RWebWindow::IsFileDialogMessage(arg)) {
1623
1625
1626 } else if (arg == "FITPANEL"s) {
1627
1628 TH1 *hist = nullptr;
1629 TIter iter(Canvas()->GetListOfPrimitives());
1630 while (auto obj = iter()) {
1631 hist = dynamic_cast<TH1 *>(obj);
1632 if (hist) break;
1633 }
1634
1635 TString cmd = TString::Format("auto panel = std::make_shared<ROOT::Experimental::RFitPanel>(\"FitPanel\");"
1636 "panel->AssignCanvas(\"%s\");"
1637 "panel->AssignHistogram((TH1 *)0x%zx);"
1638 "panel->Show();"
1639 "panel->ClearOnClose(panel);", Canvas()->GetName(), (size_t) hist);
1640
1641 gROOT->ProcessLine(cmd.Data());
1642
1643 } else if (arg == "START_BROWSER"s) {
1644
1645 gROOT->ProcessLine("new TBrowser;");
1646
1647 } else if (IsReadOnly()) {
1648
1649 // all following messages are not allowed in readonly mode
1650 return kFALSE;
1651
1652 } else if (arg.compare(0, 6, "EVENT:") == 0) {
1653 auto arr = TBufferJSON::FromJSON<std::vector<std::string>>(arg.substr(6));
1654 if (!arr || (arr->size() != 5)) {
1655 Error("ProcessData", "Wrong arguments count %d in event message", (int)(arr ? arr->size() : -1));
1656 } else {
1657 auto pad = dynamic_cast<TPad *>(FindPrimitive(arr->at(0)));
1658 std::string kind = arr->at(1);
1659 int event = -1;
1660 if (kind == "move"s) event = kMouseMotion;
1661 int argx = std::stoi(arr->at(2));
1662 int argy = std::stoi(arr->at(3));
1663 auto selobj = FindPrimitive(arr->at(4));
1664
1665 if ((event >= 0) && pad && (pad == gPad)) {
1666 Canvas()->fEvent = event;
1667 Canvas()->fEventX = argx;
1668 Canvas()->fEventY = argy;
1669
1670 Canvas()->fSelected = selobj;
1671
1672 ProcessExecs(pad);
1673 }
1674 }
1675
1676 } else if (arg.compare(0, 8, "GETMENU:") == 0) {
1677
1678 TObject *obj = FindPrimitive(arg.substr(8));
1679 if (!obj)
1680 obj = Canvas();
1681
1682 TWebMenuItems items(arg.c_str() + 8);
1683 items.PopulateObjectMenu(obj, obj->IsA());
1684 std::string buf = "MENU:";
1685 buf.append(TBufferJSON::ToJSON(&items, 103).Data());
1686
1687 AddSendQueue(connid, buf);
1688
1689 } else if (arg.compare(0, 8, "PRIMIT6:") == 0) {
1690
1691 if (IsFirstConn(connid) && !IsReadOnly()) { // only first connection can modify object
1692
1693 auto opt = TBufferJSON::FromJSON<TWebObjectOptions>(arg.c_str() + 8);
1694
1695 if (opt) {
1696 TPad *modpad = ProcessObjectOptions(*opt, nullptr);
1697
1698 // indicate that pad was modified
1699 if (modpad)
1700 modpad->Modified();
1701 }
1702 }
1703
1704 } else if (arg.compare(0, 11, "PADCLICKED:") == 0) {
1705
1706 auto click = TBufferJSON::FromJSON<TWebPadClick>(arg.c_str() + 11);
1707
1708 if (click && IsFirstConn(connid) && !IsReadOnly()) {
1709
1710 TPad *pad = dynamic_cast<TPad *>(FindPrimitive(click->padid));
1711
1712 if (pad && pad->InheritsFrom(TButton::Class())) {
1713 auto btn = (TButton *) pad;
1714 const char *mthd = btn->GetMethod();
1715 if (mthd && *mthd) {
1716 TVirtualPad::TContext ctxt(gROOT->GetSelectedPad(), kTRUE, kTRUE);
1717 gROOT->ProcessLine(mthd);
1718 }
1719 return kTRUE;
1720 }
1721
1722 if (pad && (pad != gPad)) {
1723 gPad = pad;
1727 }
1728
1729 if (!click->objid.empty()) {
1730 auto selobj = FindPrimitive(click->objid);
1731 Canvas()->SetClickSelected(selobj);
1732 Canvas()->fSelected = selobj;
1733 if (pad && selobj && fObjSelectSignal)
1734 fObjSelectSignal(pad, selobj);
1735 }
1736
1737 if ((click->x >= 0) && (click->y >= 0)) {
1738 Canvas()->fEvent = click->dbl ? kButton1Double : kButton1Up;
1739 Canvas()->fEventX = click->x;
1740 Canvas()->fEventY = click->y;
1741 if (click->dbl && fPadDblClickedSignal)
1742 fPadDblClickedSignal(pad, click->x, click->y);
1743 else if (!click->dbl && fPadClickedSignal)
1744 fPadClickedSignal(pad, click->x, click->y);
1745 }
1746
1747 ProcessExecs(pad);
1748 }
1749
1750 } else if (arg.compare(0, 8, "OBJEXEC:") == 0) {
1751
1752 auto buf = arg.substr(8);
1753 auto pos = buf.find(":");
1754
1755 if ((pos > 0) && IsFirstConn(connid) && !IsReadOnly()) { // only first client can execute commands
1756 auto sid = buf.substr(0, pos);
1757 buf.erase(0, pos + 1);
1758
1759 TObjLink *lnk = nullptr;
1760 TPad *objpad = nullptr;
1761
1762 TObject *obj = FindPrimitive(sid, 1, nullptr, &lnk, &objpad);
1763
1764 if (obj && !buf.empty()) {
1765
1766 ProcessLinesForObject(obj, buf);
1767
1768 if (objpad)
1769 objpad->Modified();
1770 else
1771 Canvas()->Modified();
1772
1774 }
1775 }
1776
1777 } else if (arg.compare(0, 12, "EXECANDSEND:") == 0) {
1778
1779 // execute method and send data, used by drawing projections
1780
1781 std::string buf = arg.substr(12);
1782 std::string reply;
1783 TObject *obj = nullptr;
1784
1785 auto pos = buf.find(":");
1786
1787 if ((pos > 0) && IsFirstConn(connid) && !IsReadOnly()) {
1788 // only first client can execute commands
1789 reply = buf.substr(0, pos);
1790 buf.erase(0, pos + 1);
1791 pos = buf.find(":");
1792 if (pos > 0) {
1793 auto sid = buf.substr(0, pos);
1794 buf.erase(0, pos + 1);
1795 obj = FindPrimitive(sid);
1796 }
1797 }
1798
1799 if (obj && !buf.empty() && !reply.empty()) {
1800 std::stringstream exec;
1801 exec << "((" << obj->ClassName() << " *) " << std::hex << std::showbase << (size_t)obj
1802 << ")->" << buf << ";";
1803 if (gDebug > 1)
1804 Info("ProcessData", "Obj %s Exec %s", obj->GetName(), exec.str().c_str());
1805
1806 auto res = gROOT->ProcessLine(exec.str().c_str());
1807 TObject *resobj = (TObject *)(res);
1808 if (resobj) {
1809 std::string send = reply;
1810 send.append(":");
1811 send.append(TBufferJSON::ToJSON(resobj, 23).Data());
1812 AddSendQueue(connid, send);
1813 if (reply[0] == 'D')
1814 delete resobj; // delete object if first symbol in reply is D
1815 }
1816 }
1817
1818 } else if (arg.compare(0, 6, "CLEAR:") == 0) {
1819 std::string snapid = arg.substr(6);
1820
1821 TPad *pad = dynamic_cast<TPad *>(FindPrimitive(snapid));
1822
1823 if (pad) {
1824 pad->Clear();
1825 pad->Modified();
1827 } else {
1828 Error("ProcessData", "Not found pad with id %s to clear\n", snapid.c_str());
1829 }
1830 } else if (arg.compare(0, 7, "DIVIDE:") == 0) {
1831 auto arr = TBufferJSON::FromJSON<std::vector<std::string>>(arg.substr(7));
1832 if (arr && arr->size() == 2) {
1833 TPad *pad = dynamic_cast<TPad *>(FindPrimitive(arr->at(0)));
1834 int nn = 0, n1 = 0, n2 = 0;
1835
1836 std::string divide = arr->at(1);
1837 auto p = divide.find('x');
1838 if (p == std::string::npos)
1839 p = divide.find('X');
1840
1841 if (p != std::string::npos) {
1842 n1 = std::stoi(divide.substr(0,p));
1843 n2 = std::stoi(divide.substr(p+1));
1844 } else {
1845 nn = std::stoi(divide);
1846 }
1847
1848 if (pad && ((nn > 1) || (n1*n2 > 1))) {
1849 pad->Clear();
1850 pad->Modified();
1851 if (nn > 1)
1852 pad->DivideSquare(nn);
1853 else
1854 pad->Divide(n1, n2);
1855 pad->cd(1);
1857 }
1858 }
1859
1860 } else if (arg.compare(0, 8, "DRAWOPT:") == 0) {
1861 auto arr = TBufferJSON::FromJSON<std::vector<std::string>>(arg.substr(8));
1862 if (arr && arr->size() == 2) {
1863 TObjLink *objlnk = nullptr;
1864 FindPrimitive(arr->at(0), 1, nullptr, &objlnk);
1865 if (objlnk)
1866 objlnk->SetOption(arr->at(1).c_str());
1867 }
1868 } else if (arg.compare(0, 8, "RESIZED:") == 0) {
1869 auto arr = TBufferJSON::FromJSON<std::vector<int>>(arg.substr(8));
1870 if (arr && arr->size() == 7) {
1871 // set members directly to avoid redrawing of the client again
1872 Canvas()->fCw = arr->at(4);
1873 Canvas()->fCh = arr->at(5);
1874 fFixedSize = arr->at(6) > 0;
1875 arr->resize(4);
1876 fWindowGeometry = *arr;
1881 }
1882 } else if (arg.compare(0, 7, "POPOBJ:") == 0) {
1883 auto arr = TBufferJSON::FromJSON<std::vector<std::string>>(arg.substr(7));
1884 if (arr && arr->size() == 2) {
1885 TPad *pad = dynamic_cast<TPad *>(FindPrimitive(arr->at(0)));
1886 TObject *obj = FindPrimitive(arr->at(1), 0, pad);
1887 if (pad && obj && (obj != pad->GetListOfPrimitives()->Last())) {
1888 TIter next(pad->GetListOfPrimitives());
1889 while (auto o = next())
1890 if (obj == o) {
1891 TString opt = next.GetOption();
1892 pad->GetListOfPrimitives()->Remove(obj);
1893 pad->GetListOfPrimitives()->AddLast(obj, opt.Data());
1894 pad->Modified();
1895 break;
1896 }
1897 }
1898 }
1899 } else if (arg == "INTERRUPT"s) {
1900 gROOT->SetInterrupt();
1901 } else {
1902 // unknown message, probably should be processed by other implementation
1903 return kFALSE;
1904 }
1905
1906 return kTRUE;
1907}
1908
1909//////////////////////////////////////////////////////////////////////////////////////////
1910/// Returns true if any pad in the canvas were modified
1911/// Reset modified flags, increment canvas version (if inc_version is true)
1912
1914{
1915 if (fPadsStatus.find(pad) == fPadsStatus.end())
1916 fPadsStatus[pad] = PadStatus{0, true, true};
1917
1918 auto &entry = fPadsStatus[pad];
1919 entry._detected = true;
1920 if (pad->IsModified()) {
1921 pad->Modified(kFALSE);
1922 entry._modified = true;
1923 }
1924
1925 TIter iter(pad->GetListOfPrimitives());
1926 while (auto obj = iter()) {
1927 if (obj->InheritsFrom(TPad::Class()))
1928 CheckPadModified(static_cast<TPad *>(obj));
1929 }
1930}
1931
1932//////////////////////////////////////////////////////////////////////////////////////////
1933/// Check if any pad on the canvas was modified
1934/// If yes, increment version of correspondent pad
1935/// Returns true when canvas really modified
1936
1938{
1939 // clear temporary flags
1940 for (auto &entry : fPadsStatus) {
1941 entry.second._detected = false;
1942 entry.second._modified = force_modified;
1943 }
1944
1945 // scan sub-pads
1947
1948 // remove no-longer existing pads
1949 bool is_any_modified = false;
1950 for(auto iter = fPadsStatus.begin(); iter != fPadsStatus.end(); ) {
1951 if (iter->second._modified)
1952 is_any_modified = true;
1953 if (!iter->second._detected)
1954 fPadsStatus.erase(iter++);
1955 else
1956 iter++;
1957 }
1958
1959 // if any pad modified, increment canvas version and set version of modified pads
1960 if (is_any_modified) {
1961 fCanvVersion++;
1962 for(auto &entry : fPadsStatus)
1963 if (entry.second._modified)
1964 entry.second.fVersion = fCanvVersion;
1965 }
1966
1967 return is_any_modified;
1968}
1969
1970
1971//////////////////////////////////////////////////////////////////////////////////////////
1972/// Returns window geometry including borders and menus
1973
1975{
1976 if (fWindowGeometry.size() == 4) {
1977 x = fWindowGeometry[0];
1978 y = fWindowGeometry[1];
1979 w = fWindowGeometry[2];
1980 h = fWindowGeometry[3];
1981 } else {
1982 x = Canvas()->fWindowTopX;
1983 y = Canvas()->fWindowTopY;
1984 w = Canvas()->fWindowWidth;
1985 h = Canvas()->fWindowHeight;
1986 }
1987 return 0;
1988}
1989
1990
1991//////////////////////////////////////////////////////////////////////////////////////////
1992/// if canvas or any subpad was modified,
1993/// scan all primitives in the TCanvas and subpads and convert them into
1994/// the structure which will be delivered to JSROOT client
1995
1997{
1998 if (CheckCanvasModified()) {
1999 // configure selected pad that method like TPad::WaitPrimitive() can correctly work
2000 // can be removed again once WaitPrimitive implemented differently
2001 if (gPad && (gPad->GetCanvas() == Canvas()))
2002 gROOT->SetSelectedPad(gPad);
2003 }
2004
2006
2007 if (!fProcessingData && !IsAsyncMode() && !async)
2009
2010 return kTRUE;
2011}
2012
2013//////////////////////////////////////////////////////////////////////////////////////////
2014/// Increment canvas version and force sending data to client - do not wait for reply
2015
2017{
2018 CheckCanvasModified(true);
2019
2020 if (!fWindow) {
2021 TCanvasWebSnapshot holder(IsReadOnly(), false, true); // readonly, set ids, batchmode
2022 CreatePadSnapshot(holder, Canvas(), 0, nullptr);
2023 } else {
2025 }
2026}
2027
2028//////////////////////////////////////////////////////////////////////////////////////////
2029/// Wait when specified version of canvas was painted and confirmed by browser
2030
2032{
2033 if (!fWindow)
2034 return kTRUE;
2035
2036 // simple polling loop until specified version delivered to the clients
2037 // first 500 loops done without sleep, then with 1ms sleep and last 500 with 100 ms sleep
2038
2039 long cnt = 0, cnt_limit = GetLongerPolling() ? 5500 : 1500;
2040
2041 if (gDebug > 2)
2042 Info("WaitWhenCanvasPainted", "version %ld", (long)ver);
2043
2044 while (cnt++ < cnt_limit) {
2045
2046 if (!fWindow->HasConnection(0, false)) {
2047 if (gDebug > 2)
2048 Info("WaitWhenCanvasPainted", "no connections - abort");
2049 return kFALSE; // wait ~1 min if no new connection established
2050 }
2051
2052 if ((fWebConn.size() > 1) && (fWebConn[1].fDrawVersion >= ver)) {
2053 if (gDebug > 2)
2054 Info("WaitWhenCanvasPainted", "ver %ld got painted", (long)ver);
2055 return kTRUE;
2056 }
2057
2059 if (cnt > 500)
2060 gSystem->Sleep((cnt < cnt_limit - 500) ? 1 : 100); // increase sleep interval when do very often
2061 }
2062
2063 if (gDebug > 2)
2064 Info("WaitWhenCanvasPainted", "timeout");
2065
2066 return kFALSE;
2067}
2068
2069//////////////////////////////////////////////////////////////////////////////////////////
2070/// Create JSON painting output for given pad
2071/// Produce JSON can be used for offline drawing with JSROOT
2072
2073TString TWebCanvas::CreatePadJSON(TPad *pad, Int_t json_compression, Bool_t batchmode)
2074{
2075 TString res;
2076 if (!pad)
2077 return res;
2078
2079 TCanvas *c = dynamic_cast<TCanvas *>(pad);
2080 if (c) {
2081 res = CreateCanvasJSON(c, json_compression, batchmode);
2082 } else {
2083 auto imp = std::make_unique<TWebCanvas>(pad->GetCanvas(), pad->GetName(), 0, 0, pad->GetWw(), pad->GetWh(), kTRUE);
2084
2085 TPadWebSnapshot holder(true, false, batchmode); // readonly, no ids, batchmode
2086
2087 imp->CreatePadSnapshot(holder, pad, 0, [&res, json_compression](TPadWebSnapshot *snap) {
2088 res = TBufferJSON::ToJSON(snap, json_compression);
2089 });
2090 }
2091
2092 return res;
2093}
2094
2095//////////////////////////////////////////////////////////////////////////////////////////
2096/// Create JSON painting output for given canvas
2097/// Produce JSON can be used for offline drawing with JSROOT
2098
2100{
2101 TString res;
2102
2103 if (!c)
2104 return res;
2105
2106 {
2107 auto imp = std::make_unique<TWebCanvas>(c, c->GetName(), 0, 0, c->GetWw(), c->GetWh(), kTRUE);
2108
2109 TCanvasWebSnapshot holder(true, false, batchmode); // readonly, no ids, batchmode
2110
2111 imp->CreatePadSnapshot(holder, c, 0, [&res, json_compression](TPadWebSnapshot *snap) {
2112 res = TBufferJSON::ToJSON(snap, json_compression);
2113 });
2114 }
2115
2116 return res;
2117}
2118
2119//////////////////////////////////////////////////////////////////////////////////////////
2120/// Create JSON painting output for given canvas and store into the file
2121/// See TBufferJSON::ExportToFile() method for more details about option
2122/// If option string starts with symbol 'b', JSON for batch mode will be generated
2123
2125{
2126 Int_t res = 0;
2127 Bool_t batchmode = kFALSE;
2128 if (option && *option == 'b') {
2129 batchmode = kTRUE;
2130 ++option;
2131 }
2132
2133 if (!c)
2134 return res;
2135
2136 {
2137 auto imp = std::make_unique<TWebCanvas>(c, c->GetName(), 0, 0, c->GetWw(), c->GetWh(), kTRUE);
2138
2139 TCanvasWebSnapshot holder(true, false, batchmode); // readonly, no ids, batchmode
2140
2141 imp->CreatePadSnapshot(holder, c, 0, [&res, filename, option](TPadWebSnapshot *snap) {
2143 });
2144 }
2145
2146 return res;
2147}
2148
2149//////////////////////////////////////////////////////////////////////////////////////////
2150/// Create image using batch (headless) capability of Chrome or Firefox browsers
2151/// Supported png, jpeg, svg, pdf formats
2152
2153bool TWebCanvas::ProduceImage(TPad *pad, const char *fileName, Int_t width, Int_t height)
2154{
2155 if (!pad)
2156 return false;
2157
2159 if (!json.Length())
2160 return false;
2161
2162 if (!width && !height) {
2163 if ((pad->GetCanvas() == pad) || (pad->IsA() == TCanvas::Class())) {
2164 width = pad->GetWw();
2165 height = pad->GetWh();
2166 } else {
2167 width = (Int_t) (pad->GetAbsWNDC() * pad->GetCanvas()->GetWw());
2168 height = (Int_t) (pad->GetAbsHNDC() * pad->GetCanvas()->GetWh());
2169 }
2170 }
2171
2172 return ROOT::RWebDisplayHandle::ProduceImage(fileName, json.Data(), width, height);
2173}
2174
2175//////////////////////////////////////////////////////////////////////////////////////////
2176/// Create images for several pads using batch (headless) capability of Chrome or Firefox browsers
2177/// Supported png, jpeg, svg, pdf, webp formats
2178/// One can include %d qualifier which will be replaced by image index using printf functionality.
2179/// If for pdf format %d qualifier not specified, all images will be stored in single PDF file.
2180/// For all other formats %d qualifier will be add before extension automatically
2181
2182bool TWebCanvas::ProduceImages(std::vector<TPad *> pads, const char *filename, Int_t width, Int_t height)
2183{
2184 if (pads.empty())
2185 return false;
2186
2187 std::vector<std::string> jsons;
2188 std::vector<Int_t> widths, heights;
2189
2190 bool isMultiPdf = (strstr(filename, ".pdf") || strstr(filename, ".PDF")) && strstr(filename, "%");
2191 bool is_multipdf_ok = true;
2192
2193 for (unsigned n = 0; n < pads.size(); ++n) {
2194 auto pad = pads[n];
2195
2197 if (!json.Length())
2198 continue;
2199
2200 Int_t w = width, h = height;
2201
2202 if (!w && !h) {
2203 if ((pad->GetCanvas() == pad) || (pad->IsA() == TCanvas::Class())) {
2204 w = pad->GetWw();
2205 h = pad->GetWh();
2206 } else {
2207 w = (Int_t) (pad->GetAbsWNDC() * pad->GetCanvas()->GetWw());
2208 h = (Int_t) (pad->GetAbsHNDC() * pad->GetCanvas()->GetWh());
2209 }
2210 }
2211
2212 if (isMultiPdf) {
2213 TString pdfname = TString::Format(filename, (int)n);
2214 if (!ROOT::RWebDisplayHandle::ProduceImage(pdfname.Data(), json.Data(), w, h))
2215 is_multipdf_ok = false;
2216 } else {
2217 jsons.emplace_back(json.Data());
2218 widths.emplace_back(w);
2219 heights.emplace_back(h);
2220 }
2221 }
2222
2223 if (isMultiPdf)
2224 return is_multipdf_ok;
2225
2226 return ROOT::RWebDisplayHandle::ProduceImages(filename, jsons, widths, heights);
2227}
2228
2229
2230//////////////////////////////////////////////////////////////////////////////////////////
2231/// Process data for single primitive
2232/// Returns object pad if object was modified
2233
2235{
2236 TObjLink *lnk = nullptr;
2237 TPad *objpad = nullptr;
2238 TObject *obj = FindPrimitive(item.snapid, idcnt, pad, &lnk, &objpad);
2239
2240 if (item.fcust.compare("exec") == 0) {
2241 auto pos = item.opt.find("(");
2242 if (obj && (pos != std::string::npos) && obj->IsA()->GetMethodAllAny(item.opt.substr(0,pos).c_str())) {
2243 std::stringstream exec;
2244 exec << "((" << obj->ClassName() << " *) " << std::hex << std::showbase
2245 << (size_t)obj << ")->" << item.opt << ";";
2246 Info("ProcessObjectOptions", "Obj %s Execute %s", obj->GetName(), exec.str().c_str());
2247 gROOT->ProcessLine(exec.str().c_str());
2248 } else {
2249 Error("ProcessObjectOptions", "Fail to execute %s for object %p %s", item.opt.c_str(), obj, obj ? obj->ClassName() : "---");
2250 objpad = nullptr;
2251 }
2252 return objpad;
2253 }
2254
2255 bool modified = false;
2256
2257 if (obj && lnk) {
2258 auto pos = item.opt.find(";;use_"); // special coding of extra options
2259 if (pos != std::string::npos) item.opt.resize(pos);
2260
2261 if (gDebug > 1)
2262 Info("ProcessObjectOptions", "Set draw option %s for object %s %s", item.opt.c_str(),
2263 obj->ClassName(), obj->GetName());
2264
2265 lnk->SetOption(item.opt.c_str());
2266
2267 modified = true;
2268 }
2269
2270 if (item.fcust.compare(0,10,"auto_exec:") == 0) {
2271 ProcessLinesForObject(obj, item.fcust.substr(10));
2272 } else if (item.fcust.compare("frame") == 0) {
2273 if (obj && obj->InheritsFrom(TFrame::Class())) {
2274 TFrame *frame = static_cast<TFrame *>(obj);
2275 if (item.fopt.size() >= 4) {
2276 frame->SetX1(item.fopt[0]);
2277 frame->SetY1(item.fopt[1]);
2278 frame->SetX2(item.fopt[2]);
2279 frame->SetY2(item.fopt[3]);
2280 modified = true;
2281 }
2282 }
2283 } else if (item.fcust.compare(0,4,"pave") == 0) {
2284 if (obj && obj->InheritsFrom(TPave::Class())) {
2285 TPave *pave = static_cast<TPave *>(obj);
2286 if ((item.fopt.size() >= 4) && objpad) {
2287 TVirtualPad::TContext ctxt(objpad, kFALSE);
2288
2289 // first time need to overcome init problem
2290 pave->ConvertNDCtoPad();
2291
2292 pave->SetX1NDC(item.fopt[0]);
2293 pave->SetY1NDC(item.fopt[1]);
2294 pave->SetX2NDC(item.fopt[2]);
2295 pave->SetY2NDC(item.fopt[3]);
2296 modified = true;
2297
2298 pave->ConvertNDCtoPad();
2299 }
2300 if ((item.fcust.length() > 4) && pave->InheritsFrom(TPaveStats::Class())) {
2301 // add text lines for statsbox
2302 auto stats = static_cast<TPaveStats *>(pave);
2303 stats->Clear();
2304 size_t pos_start = 6, pos_end;
2305 while ((pos_end = item.fcust.find(";;", pos_start)) != std::string::npos) {
2306 stats->AddText(item.fcust.substr(pos_start, pos_end - pos_start).c_str());
2307 pos_start = pos_end + 2;
2308 }
2309 stats->AddText(item.fcust.substr(pos_start).c_str());
2310 }
2311 }
2312 } else if (item.fcust.compare(0,9,"func_fail") == 0) {
2313 if (!fTF1UseSave) {
2315 modified = true;
2316 } else {
2317 Error("ProcessObjectOptions", "Client fails to calculate function %s cl %s but it should not try!", obj ? obj->GetName() : "---", obj ? obj->ClassName() : "---");
2318 }
2319 }
2320
2321 return modified ? objpad : nullptr;
2322}
2323
2324//////////////////////////////////////////////////////////////////////////////////////////////////
2325/// Search of object with given id in list of primitives
2326/// One could specify pad where search could be start
2327/// Also if object is in list of primitives, one could ask for entry link for such object,
2328/// This can allow to change draw option
2329
2330TObject *TWebCanvas::FindPrimitive(const std::string &sid, int idcnt, TPad *pad, TObjLink **objlnk, TPad **objpad)
2331{
2332 if (sid.empty() || (sid == "0"s))
2333 return nullptr;
2334
2335 if (!pad)
2336 pad = Canvas();
2337
2338 std::string kind;
2339 auto separ = sid.find("#");
2340 long unsigned id = 0;
2341 bool search_hist = false;
2342
2343 if (sid == "histogram") {
2344 search_hist = true;
2345 } else if (separ == std::string::npos) {
2346 id = std::stoul(sid);
2347 } else {
2348 kind = sid.substr(separ + 1);
2349 id = std::stoul(sid.substr(0, separ));
2350 }
2351
2352 if (!search_hist && TString::Hash(&pad, sizeof(pad)) == id)
2353 return pad;
2354
2355 auto getHistogram = [](TObject *obj) -> TH1* {
2356 auto offset = obj->IsA()->GetDataMemberOffset("fHistogram");
2357 if (offset > 0)
2358 return *((TH1 **)((char*) obj + offset));
2359 ::Error("getHistogram", "Cannot access fHistogram data member in %s", obj->ClassName());
2360 return nullptr;
2361 };
2362
2363 for (auto lnk = pad->GetListOfPrimitives()->FirstLink(); lnk != nullptr; lnk = lnk->Next()) {
2364 TObject *obj = lnk->GetObject();
2365 if (!obj) continue;
2366
2367 TString opt = lnk->GetOption();
2368 opt.ToUpper();
2369
2370 TH1 *h1 = obj->InheritsFrom(TH1::Class()) ? static_cast<TH1 *>(obj) : nullptr;
2371 TGraph *gr = obj->InheritsFrom(TGraph::Class()) ? static_cast<TGraph *>(obj) : nullptr;
2372 TGraph2D *gr2d = obj->InheritsFrom(TGraph2D::Class()) ? static_cast<TGraph2D *>(obj) : nullptr;
2373 TScatter *scatter = obj->InheritsFrom(TScatter::Class()) ? static_cast<TScatter *>(obj) : nullptr;
2374 TMultiGraph *mg = obj->InheritsFrom(TMultiGraph::Class()) ? static_cast<TMultiGraph *>(obj) : nullptr;
2375 THStack *hs = obj->InheritsFrom(THStack::Class()) ? static_cast<THStack *>(obj) : nullptr;
2376 TF1 *f1 = obj->InheritsFrom(TF1::Class()) ? static_cast<TF1 *>(obj) : nullptr;
2377
2378 if (search_hist) {
2379 if (objlnk)
2380 *objlnk = lnk;
2381
2382 if (h1)
2383 return h1;
2384 if (gr)
2385 return getHistogram(gr);
2386 if (gr2d)
2387 return getHistogram(gr2d);
2388 if (scatter)
2389 return getHistogram(scatter);
2390 if (mg && opt.Contains("A"))
2391 return getHistogram(mg);
2392 if (hs)
2393 return getHistogram(hs);
2394 if (f1)
2395 return getHistogram(f1);
2396
2397 if (objlnk)
2398 *objlnk = nullptr;
2399
2400 continue;
2401 }
2402
2403 if ((TString::Hash(&obj, sizeof(obj)) == id) && (--idcnt <= 0)) {
2404 if (objpad)
2405 *objpad = pad;
2406
2407 if (kind.compare(0, 4, "hist") == 0) {
2408 if (gr)
2409 obj = h1 = getHistogram(gr);
2410 else if (mg)
2411 obj = h1 = getHistogram(mg);
2412 else if (hs)
2413 obj = h1 = getHistogram(hs);
2414 else if (scatter)
2415 obj = h1 = getHistogram(scatter);
2416 else if (f1)
2417 obj = h1 = getHistogram(f1);
2418 else if (gr2d)
2419 obj = h1 = getHistogram(gr2d);
2420
2421 kind.erase(0,4);
2422 if (!kind.empty() && (kind[0]=='#')) kind.erase(0,1);
2423 objlnk = nullptr;
2424 }
2425
2426 if (h1 && (kind == "x"))
2427 return h1->GetXaxis();
2428 if (h1 && (kind == "y"))
2429 return h1->GetYaxis();
2430 if (h1 && (kind == "z"))
2431 return h1->GetZaxis();
2432
2433 if ((h1 || gr || scatter) && !kind.empty() && (kind.compare(0,5,"func_") == 0)) {
2434 auto funcname = kind.substr(5);
2436 return col ? col->FindObject(funcname.c_str()) : nullptr;
2437 }
2438
2439 if ((h1 || gr) && !kind.empty() && (kind.compare(0,5,"indx_") == 0)) {
2440 auto col = h1 ? h1->GetListOfFunctions() : gr->GetListOfFunctions();
2441 return col ? col->At(std::stoi(kind.substr(5))) : nullptr;
2442 }
2443
2444 if (!kind.empty() && (kind.compare(0,7,"member_") == 0)) {
2445 auto member = kind.substr(7);
2446 auto offset = obj->IsA() ? obj->IsA()->GetDataMemberOffset(member.c_str()) : 0;
2447 if (offset > 0) {
2448 TObject **mobj = (TObject **)((char*) obj + offset);
2449 return *mobj;
2450 }
2451 return nullptr;
2452 }
2453
2454 if (objlnk)
2455 *objlnk = lnk;
2456 return obj;
2457 }
2458
2459 if (h1 || gr || mg || scatter) {
2460 TList *funcs = nullptr;
2461 if (h1)
2463 else if (gr)
2465 else if (mg)
2466 funcs = mg->GetListOfFunctions();
2467 else if (scatter && scatter->GetGraph())
2468 funcs = scatter->GetGraph()->GetListOfFunctions();
2469 TIter fiter(funcs);
2470 while (auto fobj = fiter())
2471 if (TString::Hash(&fobj, sizeof(fobj)) == id) {
2472 if (objpad)
2473 *objpad = pad;
2474 return fobj;
2475 }
2476 } else if (obj->InheritsFrom(TPad::Class())) {
2477 obj = FindPrimitive(sid, idcnt, (TPad *)obj, objlnk, objpad);
2478 if (objpad && !*objpad)
2479 *objpad = pad;
2480 if (obj)
2481 return obj;
2482 }
2483 }
2484
2485 return nullptr;
2486}
2487
2488//////////////////////////////////////////////////////////////////////////////////////////////////
2489/// Static method to create TWebCanvas instance
2490/// Used by plugin manager
2491
2493{
2494 Bool_t readonly = gEnv->GetValue("WebGui.FullCanvas", (Int_t) 1) == 0;
2495
2496 auto imp = new TWebCanvas(c, name, x, y, width, height, readonly);
2497
2498 c->fWindowTopX = x;
2499 c->fWindowTopY = y;
2500 c->fWindowWidth = width;
2501 c->fWindowHeight = height;
2502 if (!gROOT->IsBatch() && (height > 25))
2503 height -= 25;
2504 c->fCw = width;
2505 c->fCh = height;
2506
2507 return imp;
2508}
2509
@ kMouseMotion
Definition Buttons.h:23
@ kButton1Double
Definition Buttons.h:24
@ kButton1Up
Definition Buttons.h:19
nlohmann::json json
#define c(i)
Definition RSha256.hxx:101
#define h(i)
Definition RSha256.hxx:106
int Int_t
Definition RtypesCore.h:45
constexpr Bool_t kFALSE
Definition RtypesCore.h:101
constexpr Ssiz_t kNPOS
Definition RtypesCore.h:124
long long Long64_t
Definition RtypesCore.h:80
constexpr Bool_t kTRUE
Definition RtypesCore.h:100
R__EXTERN TEnv * gEnv
Definition TEnv.h:170
void Info(const char *location, const char *msgfmt,...)
Use this function for informational messages.
Definition TError.cxx:218
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:185
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t hmin
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 Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char filename
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t hmax
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 Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h offset
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 r
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void funcs
Option_t Option_t width
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t src
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t height
char name[80]
Definition TGX11.cxx:110
float xmin
float ymin
float xmax
float ymax
@ kCanDelete
Definition TObject.h:367
@ kMustCleanup
Definition TObject.h:368
Int_t gDebug
Definition TROOT.cxx:597
#define gROOT
Definition TROOT.h:407
R__EXTERN TStyle * gStyle
Definition TStyle.h:433
R__EXTERN TSystem * gSystem
Definition TSystem.h:560
R__EXTERN TVirtualPS * gVirtualPS
Definition TVirtualPS.h:81
#define gPad
#define gVirtualX
Definition TVirtualX.h:338
Color * colors
Definition X3DBuffer.c:21
Holds different arguments for starting browser with RWebDisplayHandle::Display() method.
EBrowserKind GetBrowserKind() const
returns configured browser kind, see EBrowserKind for supported values
RWebDisplayArgs & SetWidgetKind(const std::string &kind)
set widget kind
RWebDisplayArgs & SetSize(int w, int h)
set preferable web window width and height
RWebDisplayArgs & SetPos(int x=-1, int y=-1)
set preferable web window x and y position, negative is default
@ kCEF
Chromium Embedded Framework - local display with CEF libs.
@ kQt5
Qt5 QWebEngine libraries - Chromium code packed in qt5.
@ kQt6
Qt6 QWebEngine libraries - Chromium code packed in qt6.
static bool ProduceImages(const std::string &fname, const std::vector< std::string > &jsons, const std::vector< int > &widths, const std::vector< int > &heights, const char *batch_file=nullptr)
Produce image file(s) using JSON data as source Invokes JSROOT drawing functionality in headless brow...
static bool ProduceImage(const std::string &fname, const std::string &json, int width=800, int height=600, const char *batch_file=nullptr)
Produce image file using JSON data as source Invokes JSROOT drawing functionality in headless browser...
static std::shared_ptr< RWebWindow > Create()
Create new RWebWindow Using default RWebWindowsManager.
static bool EmbedFileDialog(const std::shared_ptr< RWebWindow > &window, unsigned connid, const std::string &args)
Create dialog instance to use as embedded dialog inside provided widget Loads libROOTBrowserv7 and tr...
static bool IsFileDialogMessage(const std::string &msg)
Check if this could be the message send by client to start new file dialog If returns true,...
Array of integers (32 bits per element).
Definition TArrayI.h:27
const Int_t * GetArray() const
Definition TArrayI.h:43
Int_t GetSize() const
Definition TArray.h:47
static TClass * Class()
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:37
virtual void SetFillStyle(Style_t fstyle)
Set the fill area style.
Definition TAttFill.h:39
virtual void SetBottomMargin(Float_t bottommargin)
Set Pad bottom margin in fraction of the pad height.
Definition TAttPad.cxx:99
virtual void SetLeftMargin(Float_t leftmargin)
Set Pad left margin in fraction of the pad width.
Definition TAttPad.cxx:109
virtual void SetRightMargin(Float_t rightmargin)
Set Pad right margin in fraction of the pad width.
Definition TAttPad.cxx:119
Float_t GetRightMargin() const
Definition TAttPad.h:45
virtual void SetTopMargin(Float_t topmargin)
Set Pad top margin in fraction of the pad height.
Definition TAttPad.cxx:129
virtual void SetTextAlign(Short_t align=11)
Set the text alignment.
Definition TAttText.h:42
virtual void SetTextColor(Color_t tcolor=1)
Set the text color.
Definition TAttText.h:44
virtual void SetTextFont(Font_t tfont=62)
Set the text font.
Definition TAttText.h:46
virtual void SetTextSize(Float_t tsize=1)
Set the text size.
Definition TAttText.h:47
Double_t GetXmax() const
Definition TAxis.h:140
Double_t GetXmin() const
Definition TAxis.h:139
Int_t GetNbins() const
Definition TAxis.h:125
virtual void SetRangeUser(Double_t ufirst, Double_t ulast)
Set the viewing range for the axis from ufirst to ulast (in user coordinates, that is,...
Definition TAxis.cxx:1078
virtual void SetRange(Int_t first=0, Int_t last=0)
Set the viewing range for the axis using bin numbers.
Definition TAxis.cxx:1052
static TString Decode(const char *data)
Decode a base64 string date into a generic TString.
Definition TBase64.cxx:131
virtual void SetY2(Double_t y2)
Definition TBox.h:65
virtual void SetX1(Double_t x1)
Definition TBox.h:62
virtual void SetX2(Double_t x2)
Definition TBox.h:63
virtual void SetY1(Double_t y1)
Definition TBox.h:64
static Int_t ExportToFile(const char *filename, const TObject *obj, const char *option=nullptr)
Convert object into JSON and store in text file Returns size of the produce file Used in TObject::Sav...
static TString ToJSON(const T *obj, Int_t compact=0, const char *member_name=nullptr)
Definition TBufferJSON.h:75
@ kNoSpaces
no new lines plus remove all spaces around "," and ":" symbols
Definition TBufferJSON.h:39
@ kMapAsObject
store std::map, std::unordered_map as JSON object
Definition TBufferJSON.h:41
@ kSameSuppression
zero suppression plus compress many similar values together
Definition TBufferJSON.h:45
A TButton object is a user interface object.
Definition TButton.h:18
static TClass * Class()
ABC describing GUI independent main window (with menubar, scrollbars and a drawing area).
Definition TCanvasImp.h:30
TCanvas * Canvas() const
Definition TCanvasImp.h:58
void SetScripts(const std::string &src)
void SetFixedSize(bool on=true)
void SetHighlightConnect(bool on=true)
The Canvas class.
Definition TCanvas.h:23
UInt_t fCw
Width of the canvas along X (pixels)
Definition TCanvas.h:43
UInt_t GetWindowHeight() const
Definition TCanvas.h:162
void SetClickSelectedPad(TPad *pad)
Definition TCanvas.h:211
Int_t fWindowTopX
Top X position of window (in pixels)
Definition TCanvas.h:39
Int_t fEventX
! Last X mouse position in canvas
Definition TCanvas.h:46
TVirtualPadPainter * GetCanvasPainter()
Access and (probably) creation of pad painter.
Definition TCanvas.cxx:2602
UInt_t fWindowWidth
Width of window (including borders, etc.)
Definition TCanvas.h:41
Int_t fEventY
! Last Y mouse position in canvas
Definition TCanvas.h:47
UInt_t fWindowHeight
Height of window (including menubar, borders, etc.)
Definition TCanvas.h:42
TObject * fSelected
! Currently selected object
Definition TCanvas.h:49
UInt_t fCh
Height of the canvas along Y (pixels)
Definition TCanvas.h:44
UInt_t GetWindowWidth() const
Definition TCanvas.h:161
Int_t fWindowTopY
Top Y position of window (in pixels)
Definition TCanvas.h:40
void SetClickSelected(TObject *obj)
Definition TCanvas.h:209
@ kShowToolTips
Definition TCanvas.h:97
@ kShowEventStatus
Definition TCanvas.h:89
@ kMenuBar
Definition TCanvas.h:91
@ kShowEditor
Definition TCanvas.h:93
UInt_t GetWw() const override
Definition TCanvas.h:163
UInt_t GetWh() const override
Definition TCanvas.h:164
virtual void Highlighted(TVirtualPad *pad, TObject *obj, Int_t x, Int_t y)
Emit Highlighted() signal.
Definition TCanvas.cxx:1610
static TClass * Class()
Int_t fEvent
! Type of current or last handled event
Definition TCanvas.h:45
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition TClass.h:81
Longptr_t GetDataMemberOffset(const char *membername) const
return offset for member name.
Definition TClass.cxx:3477
Int_t Size() const
Return size of object of this class.
Definition TClass.cxx:5704
Bool_t InheritsFrom(const char *cl) const override
Return kTRUE if this class inherits from a class with name "classname".
Definition TClass.cxx:4874
TMethod * GetMethodAllAny(const char *method)
Return pointer to method without looking at parameters.
Definition TClass.cxx:4384
Collection abstract base class.
Definition TCollection.h:65
TObject * FindObject(const char *name) const override
Find an object in this collection using its name.
virtual Int_t GetSize() const
Return the capacity of the collection, i.e.
The color creation and management class.
Definition TColor.h:19
static const TArrayI & GetPalette()
Static function returning the current active palette.
Definition TColor.cxx:1462
static TClass * Class()
static Bool_t DefinedColors()
Static function returning kTRUE if some new colors have been defined after initialisation or since th...
Definition TColor.cxx:1480
static TClass * Class()
virtual Int_t GetValue(const char *name, Int_t dflt) const
Returns the integer value for a resource.
Definition TEnv.cxx:491
TExec is a utility class that can be used to execute a C++ command when some event happens in a pad.
Definition TExec.h:26
virtual void Exec(const char *command="")
Execute the command referenced by this object.
Definition TExec.cxx:143
static TClass * Class()
1-Dim function class
Definition TF1.h:214
virtual TH1 * GetHistogram() const
Return a pointer to the histogram used to visualise the function Note that this histogram is managed ...
Definition TF1.cxx:1586
static TClass * Class()
@ kNotDraw
Definition TF1.h:327
virtual void GetRange(Double_t *xmin, Double_t *xmax) const
Return range of a generic N-D function.
Definition TF1.cxx:2281
virtual void Save(Double_t xmin, Double_t xmax, Double_t ymin, Double_t ymax, Double_t zmin, Double_t zmax)
Save values of function in array fSave.
Definition TF1.cxx:3161
static TClass * Class()
Define a Frame.
Definition TFrame.h:19
static TClass * Class()
The axis painter class.
Definition TGaxis.h:24
static TClass * Class()
TF1 * GetFunction() const
Definition TGaxis.h:77
Graphics object made of three arrays X, Y and Z with the same number of points each.
Definition TGraph2D.h:41
static TClass * Class()
TH2D * GetHistogram(Option_t *option="")
By default returns a pointer to the Delaunay histogram.
Definition TGraph2D.cxx:979
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
static TClass * Class()
TList * GetListOfFunctions() const
Definition TGraph.h:124
virtual TH1F * GetHistogram() const
Returns a pointer to the histogram used to draw the axis Takes into account the two following cases.
Definition TGraph.cxx:1402
1-D histogram with a double per channel (see TH1 documentation)}
Definition TH1.h:620
TH1K class supports the nearest K Neighbours method, widely used in cluster analysis.
Definition TH1K.h:26
Double_t GetBinContent(Int_t bin) const override
Return content of global bin number bin.
Definition TH1K.cxx:116
static TClass * Class()
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:58
virtual void SetDirectory(TDirectory *dir)
By default, when a histogram is created, it is added to the list of histogram objects in the current ...
Definition TH1.cxx:8854
TAxis * GetZaxis()
Definition TH1.h:324
static TClass * Class()
virtual Int_t GetDimension() const
Definition TH1.h:281
@ kNoTitle
Don't draw the histogram title.
Definition TH1.h:168
TAxis * GetXaxis()
Definition TH1.h:322
virtual void SetMaximum(Double_t maximum=-1111)
Definition TH1.h:400
TAxis * GetYaxis()
Definition TH1.h:323
virtual void SetMinimum(Double_t minimum=-1111)
Definition TH1.h:401
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content see convention for numbering bins in TH1::GetBin In case the bin number is greater th...
Definition TH1.cxx:9139
virtual Double_t GetEntries() const
Return the current number of entries.
Definition TH1.cxx:4426
TList * GetListOfFunctions() const
Definition TH1.h:242
void SetName(const char *name) override
Change the name of this histogram.
Definition TH1.cxx:8877
virtual Int_t BufferEmpty(Int_t action=0)
Fill histogram with all entries in the buffer.
Definition TH1.cxx:1387
static TClass * Class()
The Histogram stack class.
Definition THStack.h:38
static TClass * Class()
Option_t * GetOption() const
void Reset()
A doubly linked list.
Definition TList.h:38
void Clear(Option_t *option="") override
Remove all objects from the list.
Definition TList.cxx:402
TObject * FindObject(const char *name) const override
Find an object in this list using its name.
Definition TList.cxx:578
void Add(TObject *obj) override
Definition TList.h:81
TObject * Remove(TObject *obj) override
Remove object from the list.
Definition TList.cxx:822
void AddLast(TObject *obj) override
Add object at the end of the list.
Definition TList.cxx:152
TObject * Last() const override
Return the last object in the list. Returns 0 when list is empty.
Definition TList.cxx:693
virtual TObjLink * FirstLink() const
Definition TList.h:102
TObject * At(Int_t idx) const override
Returns the object at position idx. Returns 0 if idx is out of range.
Definition TList.cxx:357
void AddFirst(TObject *obj) override
Add object at the beginning of the list.
Definition TList.cxx:100
A TMultiGraph is a collection of TGraph (or derived) objects.
Definition TMultiGraph.h:34
static TClass * Class()
TList * GetListOfFunctions()
Return pointer to list of functions.
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:164
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
const char * GetTitle() const override
Returns title of object.
Definition TNamed.h:48
An array of TObjects.
Definition TObjArray.h:31
Mother of all ROOT objects.
Definition TObject.h:41
virtual const char * GetName() const
Returns name of object.
Definition TObject.cxx:439
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:199
virtual TObject * Clone(const char *newname="") const
Make a clone of an object using the Streamer facility.
Definition TObject.cxx:223
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:207
virtual TObject * FindObject(const char *name) const
Must be redefined in derived classes.
Definition TObject.cxx:403
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:780
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:525
virtual const char * GetTitle() const
Returns title of object.
Definition TObject.cxx:483
virtual TClass * IsA() const
Definition TObject.h:243
virtual void Paint(Option_t *option="")
This method must be overridden if a class wants to paint itself.
Definition TObject.cxx:607
TPadWebSnapshot & NewSubPad()
Create new entry for subpad.
TWebSnapshot & NewPrimitive(TObject *obj=nullptr, const std::string &opt="", const std::string &suffix="")
Create new entry in list of primitives.
TWebSnapshot & NewSpecials()
Create new entry in list of primitives in the front.
void SetHasExecs(bool on=true)
void SetWithoutPrimitives(bool on=true)
void SetActive(bool on=true)
bool IsSetObjectIds() const
bool IsBatchMode() const
The most important graphics class in the ROOT system.
Definition TPad.h:28
Int_t GetTicky() const override
Definition TPad.h:237
Double_t fAbsYlowNDC
Absolute Y top left corner of pad in NDC [0,1].
Definition TPad.h:70
Double_t fXtoAbsPixelk
Conversion coefficient for X World to absolute pixel.
Definition TPad.h:41
virtual void DivideSquare(Int_t n, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
"n" is the total number of sub-pads.
Definition TPad.cxx:1245
static TClass * Class()
Double_t fWNDC
Width of pad along X in Normalized Coordinates (NDC)
Definition TPad.h:66
void SetView(TView *view=nullptr) override
Set the current TView. Delete previous view if view=0.
Definition TPad.cxx:6070
TVirtualViewer3D * GetViewer3D(Option_t *type="") override
Create/obtain handle to 3D viewer.
Definition TPad.cxx:7025
Double_t fPixeltoYk
Conversion coefficient for pixel to Y World.
Definition TPad.h:59
void SetGrid(Int_t valuex=1, Int_t valuey=1) override
Definition TPad.h:331
Double_t fPixeltoY
yworld = fPixeltoYk + fPixeltoY*ypixel
Definition TPad.h:60
Double_t fAbsXlowNDC
Absolute X top left corner of pad in NDC [0,1].
Definition TPad.h:69
TList * GetListOfExecs() const override
Definition TPad.h:244
void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0) override
Automatic pad generation by division.
Definition TPad.cxx:1153
Double_t fXtoPixel
xpixel = fXtoPixelk + fXtoPixel*xworld
Definition TPad.h:43
Bool_t GetGridx() const override
Definition TPad.h:233
Double_t fX2
X of upper X coordinate.
Definition TPad.h:38
Double_t fPixeltoX
xworld = fPixeltoXk + fPixeltoX*xpixel
Definition TPad.h:57
Double_t fYtoPixel
ypixel = fYtoPixelk + fYtoPixel*yworld
Definition TPad.h:46
Double_t fAbsWNDC
Absolute Width of pad along X in NDC.
Definition TPad.h:71
UInt_t GetWw() const override
Get Ww.
Definition TPad.cxx:2741
Double_t fX1
X of lower X coordinate.
Definition TPad.h:36
TList * GetListOfPrimitives() const override
Definition TPad.h:243
Double_t fUymin
Minimum value on the Y axis.
Definition TPad.h:75
Int_t fLogz
(=0 if Z linear scale, =1 if log scale)
Definition TPad.h:93
Double_t fYtoPixelk
Conversion coefficient for Y World to pixel.
Definition TPad.h:45
Double_t fPixeltoXk
Conversion coefficient for pixel to X World.
Definition TPad.h:56
Bool_t IsModified() const override
Definition TPad.h:272
Double_t fY1
Y of lower Y coordinate.
Definition TPad.h:37
Double_t fYlowNDC
Y bottom left corner of pad in NDC [0,1].
Definition TPad.h:63
Double_t fAbsPixeltoXk
Conversion coefficient for absolute pixel to X World.
Definition TPad.h:55
void Clear(Option_t *option="") override
Delete all pad primitives.
Definition TPad.cxx:626
Int_t GetTickx() const override
Definition TPad.h:236
Double_t fUymax
Maximum value on the Y axis.
Definition TPad.h:77
void Modified(Bool_t flag=1) override
Definition TPad.h:420
TVirtualPad * GetMother() const override
Definition TPad.h:257
TView * GetView() const override
Definition TPad.h:252
TClass * IsA() const override
Definition TPad.h:415
Bool_t GetGridy() const override
Definition TPad.h:234
Double_t fAbsHNDC
Absolute Height of pad along Y in NDC.
Definition TPad.h:72
void SetFixedAspectRatio(Bool_t fixed=kTRUE) override
Fix pad aspect ratio to current value if fixed is true.
Definition TPad.cxx:5906
Int_t fLogx
(=0 if X linear scale, =1 if log scale)
Definition TPad.h:91
Double_t GetAbsWNDC() const override
Definition TPad.h:220
UInt_t GetWh() const override
Get Wh.
Definition TPad.cxx:2733
TCanvas * GetCanvas() const override
Definition TPad.h:260
Double_t fXUpNDC
Definition TPad.h:64
TVirtualPad * cd(Int_t subpadnumber=0) override
Set Current pad.
Definition TPad.cxx:597
void Print(const char *filename="") const override
This method is equivalent to SaveAs("filename"). See TPad::SaveAs for details.
Definition TPad.cxx:4688
TFrame * GetFrame() override
Get frame.
Definition TPad.cxx:2859
Double_t fYUpNDC
Definition TPad.h:65
Double_t fYtoAbsPixelk
Conversion coefficient for Y World to absolute pixel.
Definition TPad.h:44
Double_t fXtoPixelk
Conversion coefficient for X World to pixel.
Definition TPad.h:42
Int_t fLogy
(=0 if Y linear scale, =1 if log scale)
Definition TPad.h:92
Double_t fHNDC
Height of pad along Y in Normalized Coordinates (NDC)
Definition TPad.h:67
Double_t fXlowNDC
X bottom left corner of pad in NDC [0,1].
Definition TPad.h:62
Double_t fUxmin
Minimum value on the X axis.
Definition TPad.h:74
Double_t GetAbsHNDC() const override
Definition TPad.h:221
void SetTicks(Int_t valuex=1, Int_t valuey=1) override
Definition TPad.h:351
Double_t fUxmax
Maximum value on the X axis.
Definition TPad.h:76
Double_t fY2
Y of upper Y coordinate.
Definition TPad.h:39
Double_t fAbsPixeltoYk
Conversion coefficient for absolute pixel to Y World.
Definition TPad.h:58
const char * GetName() const override
Returns name of object.
Definition TPad.h:258
The histogram statistics painter class.
Definition TPaveStats.h:18
virtual void SetStatFormat(const char *format="6.4g")
Change (i.e. set) the format for printing statistics.
virtual void SetFitFormat(const char *format="5.4g")
Change (i.e. set) the format for printing fit parameters in statistics box.
void SetParent(TObject *obj) override
Definition TPaveStats.h:52
static TClass * Class()
A Pave (see TPave) with text, lines or/and boxes inside.
Definition TPaveText.h:21
virtual TText * AddText(Double_t x1, Double_t y1, const char *label)
Add a new Text line to this pavetext at given coordinates.
static TClass * Class()
void Clear(Option_t *option="") override
Clear all lines in this pavetext.
virtual TText * GetLine(Int_t number) const
Get Pointer to line number in this pavetext.
A TBox with a bordersize and a shadow option.
Definition TPave.h:19
virtual void SetY1NDC(Double_t y1)
Definition TPave.h:80
virtual void ConvertNDCtoPad()
Convert pave coordinates from NDC to Pad coordinates.
Definition TPave.cxx:139
virtual void SetName(const char *name="")
Definition TPave.h:75
virtual void SetBorderSize(Int_t bordersize=4)
Definition TPave.h:73
virtual void SetY2NDC(Double_t y2)
Definition TPave.h:81
static TClass * Class()
virtual void SetX1NDC(Double_t x1)
Definition TPave.h:78
virtual void SetX2NDC(Double_t x2)
Definition TPave.h:79
A TScatter is able to draw four variables scatter plot on a single plot.
Definition TScatter.h:32
TGraph * GetGraph() const
Get the graph holding X and Y positions.
Definition TScatter.h:58
TH2F * GetHistogram() const
Get the graph histogram used for drawing axis.
Definition TScatter.cxx:159
static TClass * Class()
Basic string class.
Definition TString.h:139
Ssiz_t Length() const
Definition TString.h:421
const char * Data() const
Definition TString.h:380
@ kIgnoreCase
Definition TString.h:279
void ToUpper()
Change string to upper case.
Definition TString.cxx:1183
UInt_t Hash(ECaseCompare cmp=kExact) const
Return hash value.
Definition TString.cxx:670
TString & Append(const char *cs)
Definition TString.h:576
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2356
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:636
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition TString.h:651
Int_t GetOptStat() const
Definition TStyle.h:243
Color_t GetStatTextColor() const
Definition TStyle.h:256
Int_t GetOptTitle() const
Definition TStyle.h:244
Float_t GetStatFontSize() const
Definition TStyle.h:259
Float_t GetStatX() const
Definition TStyle.h:262
Float_t GetPadRightMargin() const
Definition TStyle.h:212
Style_t GetTitleFont(Option_t *axis="X") const
Return title font.
Definition TStyle.cxx:1212
Float_t GetStatY() const
Definition TStyle.h:263
Color_t GetTitleFillColor() const
Definition TStyle.h:269
Style_t GetTitleStyle() const
Definition TStyle.h:271
Color_t GetStatColor() const
Definition TStyle.h:255
Float_t GetStatH() const
Definition TStyle.h:265
static TClass * Class()
Width_t GetTitleBorderSize() const
Definition TStyle.h:273
Width_t GetStatBorderSize() const
Definition TStyle.h:257
Color_t GetTitleTextColor() const
Definition TStyle.h:270
Style_t GetStatStyle() const
Definition TStyle.h:260
Float_t GetStatW() const
Definition TStyle.h:264
const char * GetFitFormat() const
Definition TStyle.h:198
const char * GetStatFormat() const
Definition TStyle.h:261
Style_t GetStatFont() const
Definition TStyle.h:258
Float_t GetTitleFontSize() const
Definition TStyle.h:272
virtual void Sleep(UInt_t milliSec)
Sleep milliSec milli seconds.
Definition TSystem.cxx:424
virtual Bool_t ProcessEvents()
Process pending events (GUI, timers, sockets).
Definition TSystem.cxx:403
Handles synchronous and a-synchronous timer events.
Definition TTimer.h:51
virtual void TurnOn()
Add the timer to the system timer list.
Definition TTimer.cxx:243
void SetTime(Long_t milliSec)
Definition TTimer.h:91
See TView3D.
Definition TView.h:25
static TView * CreateView(Int_t system=1, const Double_t *rmin=nullptr, const Double_t *rmax=nullptr)
Create a concrete default 3-d view via the plug-in manager.
Definition TView.cxx:27
virtual void SetAutoRange(Bool_t autorange=kTRUE)=0
TVirtualPS is an abstract interface to Postscript, PDF, SVG.
Definition TVirtualPS.h:30
To make it possible to use GL for 2D graphic in a TPad/TCanvas.
small helper class to store/restore gPad context in TPad methods
Definition TVirtualPad.h:61
TVirtualPad is an abstract base class for the Pad and Canvas classes.
Definition TVirtualPad.h:51
Semi-Abstract base class defining a generic interface to the underlying, low level,...
Definition TVirtualX.h:46
void SetSlow(Bool_t slow=kTRUE)
TWebCanvasTimer(TWebCanvas &canv)
Bool_t IsSlow() const
void Timeout() override
used to send control messages to clients
TWebCanvas & fCanv
Basic TCanvasImp ABI implementation for Web-based Graphics Provides painting of main ROOT classes in ...
Definition TWebCanvas.h:35
TVirtualPadPainter * CreatePadPainter() override
Creates web-based pad painter.
void ForceUpdate() override
Increment canvas version and force sending data to client - do not wait for reply.
void AddCustomClass(const std::string &clname, bool with_derived=false)
Assign custom class.
static TString CreatePadJSON(TPad *pad, Int_t json_compression=0, Bool_t batchmode=kFALSE)
Create JSON painting output for given pad Produce JSON can be used for offline drawing with JSROOT.
Bool_t fTF1UseSave
! use save buffer for TF1/TF2, need when evaluation failed on client side
Definition TWebCanvas.h:108
void SetCanvasSize(UInt_t w, UInt_t h) override
Set canvas size of web canvas.
UInt_t fColorsHash
! last hash of colors/palette
Definition TWebCanvas.h:107
void ShowCmd(const std::string &arg, Bool_t show)
Function used to send command to browser to toggle menu, toolbar, editors, ...
Long64_t fColorsVersion
! current colors/palette version, checked every time when new snapshot created
Definition TWebCanvas.h:106
std::string fCustomScripts
! custom JavaScript code or URL on JavaScript files to load before start drawing
Definition TWebCanvas.h:98
virtual Bool_t IsReadOnly() const
Definition TWebCanvas.h:176
std::shared_ptr< ROOT::RWebWindow > fWindow
Definition TWebCanvas.h:88
virtual Bool_t IsJSSupportedClass(TObject *obj, Bool_t many_primitives=kFALSE)
Returns kTRUE when object is fully supported on JSROOT side In ROOT7 Paint function will just return ...
void AddCtrlMsg(unsigned connid, const std::string &key, const std::string &value)
Add control message for specified connection Same control message can be overwritten many time before...
void SetCustomScripts(const std::string &src)
Configures custom script for canvas.
ObjectSelectSignal_t fObjSelectSignal
! signal emitted when new object selected in the pad
Definition TWebCanvas.h:116
PadClickedSignal_t fPadClickedSignal
! signal emitted when simple mouse click performed on the pad
Definition TWebCanvas.h:114
void SetLongerPolling(Bool_t on)
Definition TWebCanvas.h:232
UInt_t fStyleHash
! last hash of gStyle
Definition TWebCanvas.h:105
virtual Bool_t CanCreateObject(const std::string &)
Definition TWebCanvas.h:158
void ShowWebWindow(const ROOT::RWebDisplayArgs &user_args="")
Show canvas in specified place.
Int_t fPrimitivesMerge
! number of PS primitives, which will be merged together
Definition TWebCanvas.h:96
void Show() override
Show canvas in browser window.
Bool_t WaitWhenCanvasPainted(Long64_t ver)
Wait when specified version of canvas was painted and confirmed by browser.
std::vector< std::string > fCustomClasses
! list of custom classes, which can be delivered as is to client
Definition TWebCanvas.h:99
Bool_t IsAsyncMode() const
Definition TWebCanvas.h:241
UInt_t CalculateColorsHash()
Calculate hash function for all colors and palette.
Bool_t HasStatusBar() const override
Returns kTRUE if web canvas has status bar.
void Close() override
Close web canvas - not implemented.
static bool ProduceImages(std::vector< TPad * > pads, const char *filename, Int_t width=0, Int_t height=0)
Create images for several pads using batch (headless) capability of Chrome or Firefox browsers Suppor...
Bool_t HasMenuBar() const override
Returns kTRUE if web canvas has menu bar.
Int_t InitWindow() override
Initialize window for the web canvas At this place canvas is not yet register to the list of canvases...
void CheckPadModified(TPad *pad)
Returns true if any pad in the canvas were modified Reset modified flags, increment canvas version (i...
void RaiseWindow() override
Raise browser window.
static bool ProduceImage(TPad *pad, const char *filename, Int_t width=0, Int_t height=0)
Create image using batch (headless) capability of Chrome or Firefox browsers Supported png,...
void ActivateInEditor(TPad *pad, TObject *obj)
Activate object in editor in web browser.
std::vector< WebConn > fWebConn
! connections
Definition TWebCanvas.h:83
PadSignal_t fActivePadChangedSignal
! signal emitted when active pad changed in the canvas
Definition TWebCanvas.h:113
Bool_t GetLongerPolling() const
Definition TWebCanvas.h:233
UInt_t fClientBits
! latest status bits from client like editor visible or not
Definition TWebCanvas.h:92
std::function< void(TPadWebSnapshot *)> PadPaintingReady_t
Function called when pad painting produced.
Definition TWebCanvas.h:55
Int_t fPaletteDelivery
! colors palette delivery 0:never, 1:once, 2:always, 3:per subpad
Definition TWebCanvas.h:95
Bool_t fProcessingData
! flag used to prevent blocking methods when process data is invoked
Definition TWebCanvas.h:102
Bool_t HasToolTips() const override
Returns kTRUE if tooltips are activated in web canvas.
friend class TWebCanvasTimer
Definition TWebCanvas.h:37
TWebCanvasTimer * fTimer
! timer to submit control messages
Definition TWebCanvas.h:84
Long64_t fCanvVersion
! actual canvas version, changed with every new Modified() call
Definition TWebCanvas.h:91
bool IsCustomClass(const TClass *cl) const
Checks if class belongs to custom.
std::vector< int > fWindowGeometry
! last received window geometry
Definition TWebCanvas.h:109
TPad * ProcessObjectOptions(TWebObjectOptions &item, TPad *pad, int idcnt=1)
Process data for single primitive Returns object pad if object was modified.
void CreateObjectSnapshot(TPadWebSnapshot &master, TPad *pad, TObject *obj, const char *opt, TWebPS *masterps=nullptr)
Creates representation of the object for painting in web browser.
void AddColorsPalette(TPadWebSnapshot &master)
Add special canvas objects like colors list at selected palette.
Long64_t fStyleVersion
! current gStyle object version, checked every time when new snapshot created
Definition TWebCanvas.h:104
void AddSendQueue(unsigned connid, const std::string &msg)
Add message to send queue for specified connection If connid == 0, message will be add to all connect...
void SetWindowPosition(Int_t x, Int_t y) override
Set window position of web canvas.
UpdatedSignal_t fUpdatedSignal
! signal emitted when canvas updated or state is changed
Definition TWebCanvas.h:112
Int_t fJsonComp
! compression factor for messages send to the client
Definition TWebCanvas.h:97
Bool_t IsFirstConn(unsigned connid) const
Definition TWebCanvas.h:146
~TWebCanvas() override
Destructor.
Bool_t fReadOnly
!< configured display
Definition TWebCanvas.h:90
std::map< TPad *, PadStatus > fPadsStatus
! map of pads in canvas and their status flags
Definition TWebCanvas.h:86
Bool_t CheckDataToSend(unsigned connid=0)
Check if any data should be send to client If connid != 0, only selected connection will be checked.
Bool_t PerformUpdate(Bool_t async) override
if canvas or any subpad was modified, scan all primitives in the TCanvas and subpads and convert them...
void AssignStatusBits(UInt_t bits)
Assign clients bits.
virtual Bool_t ProcessData(unsigned connid, const std::string &arg)
Handle data from web browser Returns kFALSE if message was not processed.
void ProcessLinesForObject(TObject *obj, const std::string &lines)
Execute one or several methods for selected object String can be separated by ";;" to let execute sev...
TWebCanvas(TCanvas *c, const char *name, Int_t x, Int_t y, UInt_t width, UInt_t height, Bool_t readonly=kTRUE)
Constructor.
static Int_t StoreCanvasJSON(TCanvas *c, const char *filename, const char *option="")
Create JSON painting output for given canvas and store into the file See TBufferJSON::ExportToFile() ...
void Iconify() override
Iconify browser window.
void SetWindowTitle(const char *newTitle) override
Set window title of web canvas.
UInt_t GetWindowGeometry(Int_t &x, Int_t &y, UInt_t &w, UInt_t &h) override
Returns window geometry including borders and menus.
static TCanvasImp * NewCanvas(TCanvas *c, const char *name, Int_t x, Int_t y, UInt_t width, UInt_t height)
Static method to create TWebCanvas instance Used by plugin manager.
static TString CreateCanvasJSON(TCanvas *c, Int_t json_compression=0, Bool_t batchmode=kFALSE)
Create JSON painting output for given canvas Produce JSON can be used for offline drawing with JSROOT...
Bool_t HasEditor() const override
Returns kTRUE if web canvas has graphical editor.
Int_t fStyleDelivery
! gStyle delivery to clients: 0:never, 1:once, 2:always
Definition TWebCanvas.h:94
PadClickedSignal_t fPadDblClickedSignal
! signal emitted when simple mouse click performed on the pad
Definition TWebCanvas.h:115
void ProcessExecs(TPad *pad, TExec *extra=nullptr)
Process TExec objects in the pad.
TList fPrimitivesLists
! list of lists of primitives, temporary collected during painting
Definition TWebCanvas.h:93
void CreatePadSnapshot(TPadWebSnapshot &paddata, TPad *pad, Long64_t version, PadPaintingReady_t func)
Create snapshot for pad and all primitives Callback function is used to create JSON in the middle of ...
Bool_t CheckCanvasModified(bool force_modified=false)
Check if any pad on the canvas was modified If yes, increment version of correspondent pad Returns tr...
virtual Bool_t DecodePadOptions(const std::string &, bool process_execs=false)
Decode all pad options, which includes ranges plus objects options.
Int_t GetPaletteDelivery() const
Definition TWebCanvas.h:227
void SetWindowSize(UInt_t w, UInt_t h) override
Set window size of web canvas.
TObject * FindPrimitive(const std::string &id, int idcnt=1, TPad *pad=nullptr, TObjLink **objlnk=nullptr, TPad **objpad=nullptr)
Search of object with given id in list of primitives One could specify pad where search could be star...
Bool_t fFixedSize
! is canvas size fixed
Definition TWebCanvas.h:110
Int_t GetStyleDelivery() const
Definition TWebCanvas.h:224
void PopulateObjectMenu(void *obj, TClass *cl)
Class used to transport drawing options from the client.
std::vector< double > fopt
custom float array
std::string fcust
custom string
std::string snapid
id of the object
std::string opt
drawing options
void CreatePainting()
Definition TWebPS.cxx:26
Bool_t IsEmptyPainting() const
Definition TWebPS.h:32
TWebPainting * TakePainting()
Definition TWebPS.h:34
TWebPainting * GetPainting()
Definition TWebPS.h:33
Implement TVirtualPadPainter which abstracts painting operations.
Object used to store paint operations and deliver them to JSROOT.
void AddColor(Int_t indx, TColor *col)
Add custom color to operations.
@ kStyle
gStyle object
@ kObject
object itself
@ kSVG
list of SVG primitives
@ kSubPad
subpad
@ kColors
list of ROOT colors + palette
void SetSnapshot(Int_t kind, TObject *snapshot, Bool_t owner=kFALSE)
SetUse pointer to assign object id - TString::Hash.
void SetObjectIDAsPtr(void *ptr, const std::string &suffix="")
Use pointer to assign object id - TString::Hash.
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
TGraphErrors * gr
Definition legend1.C:25
TH1F * h1
Definition legend1.C:5
TF1 * f1
Definition legend1.C:11
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Definition TMath.h:721
bool _detected
! if pad was detected during last scan
Definition TWebCanvas.h:78