Logo ROOT   6.18/05
Reference Guide
TTreeViewer.cxx
Go to the documentation of this file.
1// @(#)root/treeviewer:$Id: c8e226dde2f9b6f39946bfe90cabcb778d63dc4f $
2//Author : Andrei Gheata 16/08/00
3
4/*************************************************************************
5 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12/** \class TTreeViewer
13A graphic user interface designed to handle ROOT trees and to take advantage of
14TTree class features.
15
16It uses ROOT native GUI widgets adapted for "drag and drop" functionality.
17in the same session.
18
19### The following capabilities are making the viewer a helpful tool for analysis:
20
21 - several trees may be opened in the same session;
22 - branches and leaves can be easily browsed or scanned;
23 - fast drawing of branch expressions by double-clicking;
24 - new variables/selections easy to compose with the built-in editor;
25 - histograms can be composed by dragging leaves or user-defined expressions
26 to X, Y and Z axis items;
27 - the tree entries to be processed can be selected with a double slider;
28 - selections can be defined and activated by dragging them to the 'Cut' item;
29 - all expressions can be aliased and aliases can be used in composing others;
30 - input/output event lists easy to handle;
31 - menu with histogram drawing options;
32 - user commands may be executed within the viewer and the current command
33 can be echoed;
34 - current 'Draw' event loop is reflected by a progress bar and may be
35 interrupted by the user;
36 - all widgets have self-explaining tool tips and/or context menus;
37 - expressions/leaves can be dragged to a 'scan box' and scanned by
38 double-clicking this item. The result can be redirected to an ASCII file;
39
40### The layout has the following items:
41
42 - a menu bar with entries : File, Edit, Run, Options and Help;
43 - a toolbar in the upper part where you can issue user commands, change
44 the drawing option and the histogram name, three check buttons Hist, Rec
45 and Scan.HIST toggles histogram drawing mode, REC enables recording of the
46 last command issued and SCAN enables redirecting of TTree::Scan command in
47 an ASCII file (see -Scanning expressions-);
48 - a button bar in the lower part with : buttons DRAW/STOP that issue histogram
49 drawing and stop the current command respectively, two text widgets where
50 input and output event lists can be specified, a message box and a RESET
51 button on the right that clear edited expression content (see Editing...)
52 - a tree-type list on the main left panel where you can select among trees or
53 branches. The tree/branch will be detailed in the right panel.
54 Mapped trees are provided with context menus, activated by right-clicking;
55 - a view-type list on the right panel. The first column contain X, Y and
56 Z expression items, an optional cut and ten optional editable expressions.
57 Expressions and leaf-type items can be dragged or deleted. A right click on
58 the list-box or item activates context menus.
59
60### Opening a new tree and saving a session :
61
62 To open a new tree in the viewer use `<File/Open tree file>` menu
63The content of the file (keys) will be listed. Use `<SetTreeName>` function
64from the context menu of the right panel, entering a tree name among those
65listed.
66
67 To save the current session, use `<File/Save>` menu or the `<SaveSource>`
68function from the context menu of the right panel (to specify the name of the
69file - name.C)
70
71 To open a previously saved session for the tree MyTree, first open MyTree
72in the browser, then use `<File/Open session>` menu.
73
74### Dragging items:
75
76Items that can be dragged from the list in the right : expressions and
77leaves. Dragging an item and dropping to another will copy the content of first
78to the last (leaf->expression, expression->expression). Items far to the right
79side of the list can be easily dragged to the left (where expressions are
80placed) by dragging them to the left at least 10 pixels.
81
82### Editing expressions:
83
84 Any editable expression from the right panel has two components : a
85true name (that will be used when TTree::Draw() commands are issued) and an
86alias. The visible name is the alias. Aliases of user defined expressions have
87a leading ~ and may be used in new expressions. Expressions containing boolean
88operators have a specific icon and may be dragged to the active cut (scissors
89item) position.
90
91 The expression editor can be activated by double-clicking empty expression,
92using `<EditExpression>` from the selected expression context menu or using
93`<Edit/Expression>` menu.
94
95 The editor will pop-up in the left part, but it can be moved.
96The editor usage is the following :
97
98 - you can write C expressions made of leaf names by hand or you can insert
99 any item from the right panel by clicking on it (recommandable);
100 - you can click on other expressions/leaves to paste them in the editor;
101 - you should write the item alias by hand since it not only make the
102 expression meaningful, but it also highly improve the layout for big
103 expressions
104 - you may redefine an old alias - the other expressions depending on it will
105 be modified accordingly. An alias must not be the leading string of other
106 aliases. When Draw commands are issued, the name of the corresponding
107 histogram axes will become the aliases of the expressions.
108
109User commands can be issued directly from the textbox labeled "Command"
110from the upper-left toolbar by typing and pressing Enter at the end.
111
112 Another way is from the right panel context menu : ExecuteCommand.
113All commands can be interrupted at any time by pressing the STOP button
114from the bottom-left
115You can toggle recording of the current command in the history file by
116checking the Rec button from the top-right
117
118### Context menus
119
120 You can activate context menus by right-clicking on items or inside the
121right panel.
122
123Context menus for mapped items from the left tree-type list :
124 The items from the left that are provided with context menus are tree and
125branch items. You can directly activate the *MENU* marked methods of TTree
126from this menu.
127
128Context menu for the right panel:
129
130 A general context menu is activated if the user right-clicks the right panel.
131
132 Commands are :
133 - EmptyAll : clears the content of all expressions;
134 - ExecuteCommand : execute a ROOT command;
135 - MakeSelector : equivalent of TTree::MakeSelector();
136 - NewExpression : add an expression item in the right panel;
137 - Process : equivalent of TTree::Process();
138 - SaveSource : save the current session as a C++ macro;
139 - SetScanFileName : define a name for the file where TTree::Scan command
140 is redirected when the `<Scan>` button is checked;
141 - SetTreeName : open a new tree with this name in the viewer;
142
143 A specific context menu is activated if expressions/leaves are right-clicked.
144
145 Commands are :
146 - Draw : draw a histogram for this item;
147 - EditExpression : pops-up the expression editor;
148 - Empty : empty the name and alias of this item;
149 - RemoveItem : removes clicked item from the list;
150 - Scan : scan this expression;
151 - SetExpression : edit name and alias for this item by hand;
152
153Starting the viewer
154
155 1. From the TBrowser: Select a tree in the TBrowser, then call the
156 StartViewer() method from its context menu (right-click on the tree).
157 2. From the command line: Start a ROOT session in the directory where you have
158 your tree. You will need first to load the library for TTreeViewer and
159 optionally other libraries for user defined classes (you can do this later in
160 the session) :
161~~~ {.cpp}
162 root [0] gSystem->Load(\"TTreeViewer\");
163~~~
164Supposing you have the tree MyTree in the file MyFile, you can do :
165~~~ {.cpp}
166 root [1] TFile file("Myfile");
167 root [2] new TTreeViewer("Mytree");
168~~~
169or :
170~~~ {.cpp}
171 root [2] TTreeViewer *tv = new TTreeViewer();
172 root [3] tv->SetTreeName("Mytree");
173~~~
174\image html ttree_treeview.png
175*/
176
177#include "RConfigure.h"
178
179#include "Riostream.h"
180#include "TTreeViewer.h"
181#include "HelpText.h"
182#include "HelpTextTV.h"
183#include "TTVLVContainer.h"
184#include "TTVSession.h"
185
186#include "TROOT.h"
187#include "TError.h"
188#include "TGMsgBox.h"
189#include "TTreePlayer.h"
190#include "TContextMenu.h"
191#include "TInterpreter.h"
192#include "TLeaf.h"
193#include "TRootHelpDialog.h"
194#include "TSystem.h"
195#include "TApplication.h"
196#include "TVirtualX.h"
197#include "TGClient.h"
198#include "TKey.h"
199#include "TFile.h"
200#include "TGMenu.h"
201#include "TGFrame.h"
202#include "TCanvas.h"
203#include "TH1.h"
204#include "TTree.h"
205#include "TFriendElement.h"
206#include "TObjArray.h"
207#include "TObjString.h"
208#include "TGButton.h"
209#include "TGButtonGroup.h"
210#include "TGTextEntry.h"
211#include "TGComboBox.h"
212#include "TGLabel.h"
213#include "TGListView.h"
214#include "TGListTree.h"
215#include "TGMimeTypes.h"
216#include "TGSplitter.h"
217#include "TGDoubleSlider.h"
218#include "TGToolBar.h"
219#include "TGStatusBar.h"
220#include "Getline.h"
221#include "TTimer.h"
222#include "TG3DLine.h"
223#include "TGFileDialog.h"
224#include "TGProgressBar.h"
225#include "TClonesArray.h"
226#include "TSpider.h"
227
228#ifdef WIN32
229#include "TWin32SplashThread.h"
230#endif
231
232// drawing options
233static const char* gOptgen[16] =
234{
235 "","AXIS","HIST","SAME","CYL","POL","SPH","PSR","LEGO","LEGO1","LEGO2",
236 "SURF","SURF1","SURF2","SURF3","SURF4"
237};
238static const char* gOpt1D[12] =
239{
240 "","AH","B","C","E","E1","E2","E3","E4","L","P","*H"
241};
242static const char* gOpt2D[14] =
243{
244 "","ARR","BOX","COL","COL2","CONT","CONT0","CONT1","CONT2","CONT3",
245 "FB","BB","SCAT","PROF"
246};
247
248static const char* gOpenTypes[] = {"Root files", "*.root",
249 0, 0 };
250
251static const char* gMacroTypes[] = {"C++ macros", "*.C",
252 0, 0 };
253
254// Menu command id's
264
269
272
277
290
294 kAxis
296
297// button Id's
308 kBGLast
310
312
313////////////////////////////////////////////////////////////////////////////////
314/// TTreeViewer default constructor
315
316TTreeViewer::TTreeViewer(const char* treeName) :
318 fDimension(0), fVarDraw(0), fScanMode(0),
319 fTreeIndex(0), fDefaultCursor(0), fWatchCursor(0),
320 fCounting(0), fStopMapping(0), fEnableCut(0),fNexpressions(0)
321{
322 fTree = 0;
323 if (!gClient) return;
324 char command[128];
325 gROOT->ProcessLine("#ifndef GTV_DEFINED\n\
326 TTreeViewer *gTV = 0;\n\
327 TTree *tv__tree = 0;\n\
328 TList *tv__tree_list = 0;\n\
329 TFile *tv__tree_file = 0;\n\
330 #define GTV_DEFINED\n\
331 #endif");
332 snprintf(command,128, "gTV = (TTreeViewer*)0x%lx", (ULong_t)this);
333 gROOT->ProcessLine(command);
334 fTreeList = new TList;
335 gROOT->ProcessLine("tv__tree_list = new TList;");
336 fFilename = "";
337 gInterpreter->SaveContext();
339 SetTreeName(treeName);
340}
341
342////////////////////////////////////////////////////////////////////////////////
343
345 TGMainFrame(0, 10, 10, kVerticalFrame),
346 fDimension(0), fVarDraw(0), fScanMode(0),
347 fTreeIndex(0), fDefaultCursor(0), fWatchCursor(0),
348 fCounting(0), fStopMapping(0), fEnableCut(0),fNexpressions(0)
349
350{
351 // TTreeViewer constructor with a pointer to a Tree
352
353 fTree = 0;
354 char command[128];
355 gROOT->ProcessLine("#ifndef GTV_DEFINED\n\
356 TTreeViewer *gTV = 0;\n\
357 TTree *tv__tree = 0;\n\
358 TList *tv__tree_list = 0;\n\
359 TFile *tv__tree_file = 0;\n\
360 #define GTV_DEFINED\n\
361 #endif");
362 snprintf(command,128, "gTV = (TTreeViewer*)0x%lx", (ULong_t)this);
363 gROOT->ProcessLine(command);
364 if (!tree) return;
365 fTreeList = new TList;
366 gROOT->ProcessLine("tv__tree_list = new TList;");
367 fFilename = "";
368 gInterpreter->SaveContext();
370 TDirectory *dirsav = gDirectory;
371 TDirectory *cdir = tree->GetDirectory();
372 if (cdir) cdir->cd();
373
374 SetTree((TTree *)tree);
375 // If the tree is a chain, the tree directory will be changed by SwitchTree
376 // (called by SetTreeName)
377 cdir = tree->GetDirectory();
378 if (cdir) {
379 if (cdir->GetFile()) fFilename = cdir->GetFile()->GetName();
380 }
381 if (dirsav) dirsav->cd();
382}
383////////////////////////////////////////////////////////////////////////////////
384/// Allow geting the tree from the context menu.
385
387{
388 if (!tree) return;
389 TTree *ftree;
390 if (fTreeList) {
391 if (fTreeList->FindObject(tree)) {
392 printf("Tree found\n");
393 TIter next(fTreeList);
394 Int_t index = 0;
395 while ((ftree = (TTree*)next())) {
396 if (ftree==tree) {printf("found at index %i\n", index);break;}
397 index++;
398 }
399 SwitchTree(index);
400 if (fTree != fMappedTree) {
401 // switch also the global "tree" variable
403 // map it on the right panel
404 MapTree(fTree);
405 fListView->Layout();
406 TGListTreeItem *base = 0;
407 TGListTreeItem *parent = fLt->FindChildByName(base, "TreeList");
408 TGListTreeItem *item = fLt->FindChildByName(parent, fTree->GetName());
410 fLt->HighlightItem(item);
412 }
413 return;
414 }
415 }
416 if (fTree != tree) {
417 fTree = tree;
418 // load the tree via the interpreter
419 char command[100];
420 command[0] = 0;
421 // define a global "tree" variable for the same tree
422 snprintf(command,100, "tv__tree = (TTree *)0x%lx;", (ULong_t)tree);
423 ExecuteCommand(command);
424 }
425 //--- add the tree to the list if it is not already in
427 ExecuteCommand("tv__tree_list->Add(tv__tree);");
428 //--- map this tree
429 TGListTreeItem *base = 0;
430 TGListTreeItem *parent = fLt->FindChildByName(base, "TreeList");
431 if (!parent) parent = fLt->AddItem(base, "TreeList", new ULong_t(kLTNoType));
432 ULong_t *itemType = new ULong_t((fTreeIndex << 8) | kLTTreeType);
433 fTreeIndex++;
434 TGListTreeItem *lTreeItem = fLt->AddItem(parent, tree->GetName(), itemType,
435 gClient->GetPicture("tree_t.xpm"), gClient->GetPicture("tree_t.xpm"));
436 MapTree(fTree, lTreeItem, kFALSE);
437 fLt->OpenItem(parent);
438 fLt->HighlightItem(lTreeItem);
440
441 //--- map slider and list view
444 MapTree(fTree);
445 fListView->Layout();
446 SetFile();
447}
448////////////////////////////////////////////////////////////////////////////////
449/// Change the number of expression widgets.
450
452{
453 Int_t diff = expr - fNexpressions;
454 if (diff <= 0) return;
455 if (!fLVContainer) return;
456 for (Int_t i=0; i<TMath::Abs(diff); i++) NewExpression();
457}
458////////////////////////////////////////////////////////////////////////////////
459/// Set the name of the file where to redirect `<Scan>` output.
460
461void TTreeViewer::SetScanFileName(const char *name)
462{
464}
465////////////////////////////////////////////////////////////////////////////////
466/// Set the state of Scan check button.
467
469{
470 if (mode)
472 else
474}
475////////////////////////////////////////////////////////////////////////////////
476/// Assign the fTree member from existing tree, e.g. when calling
477/// tree->StartViewer() from the browser, or even from the command line.
478
480{
481 if (!tree) return;
482 if (fTree != tree) {
483 fTree = tree;
484 // load the tree via the interpreter
485 // define a global "tree" variable for the same tree
486 TString command = TString::Format("tv__tree = (TTree *)0x%lx;", (ULong_t)tree);
487 ExecuteCommand(command.Data());
488 }
489 //--- add the tree to the list if it is not already in
491 ExecuteCommand("tv__tree_list->Add(tv__tree);");
492 //--- map this tree
493 TGListTreeItem *base = 0;
494 TGListTreeItem *parent = fLt->FindChildByName(base, "TreeList");
495 if (!parent) parent = fLt->AddItem(base, "TreeList", new ULong_t(kLTNoType));
496 ULong_t *itemType = new ULong_t((fTreeIndex << 8) | kLTTreeType);
497 fTreeIndex++;
498 TGListTreeItem *lTreeItem = fLt->AddItem(parent, tree->GetName(), itemType,
499 gClient->GetPicture("tree_t.xpm"), gClient->GetPicture("tree_t.xpm"));
500 MapTree(fTree, lTreeItem, kFALSE);
501 fLt->OpenItem(parent);
502 fLt->HighlightItem(lTreeItem);
504
505 //--- map slider and list view
508 MapTree(fTree);
509 fListView->Layout();
510 SetFile();
511}
512////////////////////////////////////////////////////////////////////////////////
513/// Allow geting the tree from the context menu.
514
515void TTreeViewer::SetTreeName(const char* treeName)
516{
517 if (!treeName) return;
518 TTree *tree = (TTree *) gROOT->FindObject(treeName);
519 if (fTreeList) {
520 if (fTreeList->FindObject(treeName)) {
521 printf("Tree found\n");
522 TIter next(fTreeList);
523 Int_t index = 0;
524 while ((tree = (TTree*)next())) {
525 if (!strcmp(treeName, tree->GetName())) {printf("found at index %i\n", index);break;}
526 index++;
527 }
528 SwitchTree(index);
529 if (fTree != fMappedTree) {
530 // switch also the global "tree" variable
532 // map it on the right panel
533 MapTree(fTree);
534 fListView->Layout();
535 TGListTreeItem *base = 0;
536 TGListTreeItem *parent = fLt->FindChildByName(base, "TreeList");
537 TGListTreeItem *item = fLt->FindChildByName(parent, fTree->GetName());
539 fLt->HighlightItem(item);
541 }
542 return;
543 }
544 }
545 if (!tree) return;
546// ((TTreePlayer *)tree->GetPlayer())->SetViewer(this);
547 if (fTree != tree) {
548 fTree = tree;
549 // load the tree via the interpreter
550 // define a global "tree" variable for the same tree
551 TString command = TString::Format("tv__tree = (TTree *) gROOT->FindObject(\"%s\");", treeName);
552 ExecuteCommand(command.Data());
553 }
554 //--- add the tree to the list if it is not already in
556 ExecuteCommand("tv__tree_list->Add(tv__tree);");
557 //--- map this tree
558 TGListTreeItem *base = 0;
559 TGListTreeItem *parent = fLt->FindChildByName(base, "TreeList");
560 if (!parent) parent = fLt->AddItem(base, "TreeList", new ULong_t(kLTNoType));
561 ULong_t *itemType = new ULong_t((fTreeIndex << 8) | kLTTreeType);
562 fTreeIndex++;
563 TGListTreeItem *lTreeItem = fLt->AddItem(parent, treeName, itemType,
564 gClient->GetPicture("tree_t.xpm"), gClient->GetPicture("tree_t.xpm"));
565 MapTree(fTree, lTreeItem, kFALSE);
566 fLt->OpenItem(parent);
567 fLt->HighlightItem(lTreeItem);
569
570 //--- map slider and list view
573 MapTree(fTree);
574 fListView->Layout();
575 SetFile();
576}
577////////////////////////////////////////////////////////////////////////////////
578/// Set file name containing the tree.
579
581{
582 if (!fTree) return;
583 TSeqCollection *list = gROOT->GetListOfFiles();
584 TTree *tree;
585 TIter next(list);
586 TObject *obj;
587 TFile *file;
588 while ((obj=next())) {
589 file = (TFile*)obj;
590 if (file) {
591 tree = (TTree*)file->Get(fTree->GetName());
592 if (tree) {
593 fFilename = file->GetName();
594 std::cout << "File name : "<< fFilename << std::endl;
595 return;
596 } else {
597 fFilename = "";
598 }
599 }
600 }
601 fFilename = "";
602}
603////////////////////////////////////////////////////////////////////////////////
604/// Create all viewer widgets.
605
607{
608 //--- timer & misc
612 fTimer = new TTimer(this, 20, kTRUE);
613 fLastOption = "";
614 fSession = new TTVSession(this);
615 //--- cursors
616 fDefaultCursor = gVirtualX->CreateCursor(kPointer);
617 fWatchCursor = gVirtualX->CreateCursor(kWatch);
618 //--- colours
619 ULong_t color;
620 gClient->GetColorByName("blue",color);
621 //--- pictures for X, Y and Z expression items
622 fPicX = gClient->GetPicture("x_pic.xpm");
623 fPicY = gClient->GetPicture("y_pic.xpm");
624 fPicZ = gClient->GetPicture("z_pic.xpm");
625
626 //--- general context menu
627 fContextMenu = new TContextMenu("TreeViewer context menu","");
628 fMappedTree = 0;
629 fMappedBranch = 0;
630 fDialogBox = 0;
631 fDimension = 0;
634// fFilename = "";
635 fSourceFile = "treeviewer.C";
636 //--- lists : trees and widgets to be removed
637// fTreeList = 0;
638 fTreeIndex = 0;
639 fWidgets = new TList();
640 //--- create menus --------------------------------------------------------
641 //--- File menu
643 fFileMenu->AddEntry("&New canvas", kFileCanvas);
644 fFileMenu->AddEntry("Open &tree file...", kFileBrowse);
645 fFileMenu->AddEntry("&Load Library...", kFileLoadLibrary);
646 fFileMenu->AddEntry("&Open session", kFileOpenSession);
647 fFileMenu->AddEntry("&Save source...", kFileSaveMacro);
649 fFileMenu->AddEntry("&Print", kFilePrint);
650 fFileMenu->AddEntry("&Close", kFileClose);
652 fFileMenu->AddEntry("&Quit ROOT", kFileQuit);
653
655
656 //--- Edit menu
657 fEditMenu = new TGPopupMenu(gClient->GetRoot());
658 fEditMenu->AddEntry("&Expression...", kEditExpression);
659 fEditMenu->AddEntry("&Cut...", kEditCut);
660 fEditMenu->AddEntry("&Macro...", kEditMacro);
661 fEditMenu->AddEntry("E&Vent...", kEditEvent);
662
665 //---Run menu
666 fRunMenu = new TGPopupMenu(gClient->GetRoot());
667 fRunMenu->AddEntry("&Macro...", kRunMacro);
669 //--- Options menu
670 //--- General options
671 fOptionsGen = new TGPopupMenu(gClient->GetRoot());
674 fOptionsGen->AddEntry("Axis only", kOptionsGeneral+1); // "AXIS"
675 fOptionsGen->AddEntry("Contour only", kOptionsGeneral+2); // "HIST"
676 fOptionsGen->AddEntry("Superimpose", kOptionsGeneral+3); //"SAME"
677 fOptionsGen->AddEntry("Cylindrical", kOptionsGeneral+4); //"CYL"
678 fOptionsGen->AddEntry("Polar", kOptionsGeneral+5); //"POL"
679 fOptionsGen->AddEntry("Spherical", kOptionsGeneral+6); //"SPH"
680 fOptionsGen->AddEntry("PsRap/Phi", kOptionsGeneral+7); //"PSR"
681 fOptionsGen->AddEntry("Lego HLR", kOptionsGeneral+8); //"LEGO"
682 fOptionsGen->AddEntry("Lego HSR", kOptionsGeneral+9); //"LEGO1"
683 fOptionsGen->AddEntry("Lego Color", kOptionsGeneral+10); //"LEGO2"
684 fOptionsGen->AddEntry("Surface HLR", kOptionsGeneral+11); //"SURF"
685 fOptionsGen->AddEntry("Surface HSR", kOptionsGeneral+12); //"SURF1"
686 fOptionsGen->AddEntry("Surface Col", kOptionsGeneral+13); //"SURF2"
687 fOptionsGen->AddEntry("Surf+Cont", kOptionsGeneral+14); //"SURF3"
688 fOptionsGen->AddEntry("Gouraud", kOptionsGeneral+15); //"SURF4"
689 fOptionsGen->Associate(this);
690 //--- 1D options
691 fOptions1D = new TGPopupMenu(gClient->GetRoot());
692 fOptions1D->AddEntry("Default", kOptions1D);
694 fOptions1D->AddEntry("No labels/ticks", kOptions1D+1); // "AH"
695 fOptions1D->AddEntry("Bar chart", kOptions1D+2); // "B"
696 fOptions1D->AddEntry("Smooth curve", kOptions1D+3); // "C"
697 fOptions1D->AddEntry("Errors", kOptions1D+4); // "E"
698 fOptions1D->AddEntry("Errors 1", kOptions1D+5); // "E1"
699 fOptions1D->AddEntry("Errors 2", kOptions1D+6); // "E2"
700 fOptions1D->AddEntry("Errors 3", kOptions1D+7); // "E3"
701 fOptions1D->AddEntry("Errors 4", kOptions1D+8); // "E4"
702 fOptions1D->AddEntry("Line", kOptions1D+9); // "L"
703 fOptions1D->AddEntry("Markers", kOptions1D+10); // "P"
704 fOptions1D->AddEntry("Stars", kOptions1D+11); // "*H"
705 fOptions1D->Associate(this);
706 //--- 2D options
707 fOptions2D = new TGPopupMenu(gClient->GetRoot());
708 fOptions2D->AddEntry("Default", kOptions2D);
710 fOptions2D->AddEntry("Arrows", kOptions2D+1); // "ARR"
711 fOptions2D->AddEntry("Box/Surf", kOptions2D+2); // "BOX"
712 fOptions2D->AddEntry("Box/Color", kOptions2D+3); // "COL"
713 fOptions2D->AddEntry("Box/ColMap", kOptions2D+4); // "COLZ"
714 fOptions2D->AddEntry("Contour", kOptions2D+5); // "CONT"
715 fOptions2D->AddEntry("Contour 0", kOptions2D+6); // "CONT0"
716 fOptions2D->AddEntry("Contour 1", kOptions2D+7); // "CONT1"
717 fOptions2D->AddEntry("Contour 2", kOptions2D+8); // "CONT2"
718 fOptions2D->AddEntry("Contour 3", kOptions2D+9); // "CONT3"
719 fOptions2D->AddEntry("No front-box", kOptions2D+10); // "FB"
720 fOptions2D->AddEntry("No back-box", kOptions2D+11); // "BB"
721 fOptions2D->AddEntry("Scatter", kOptions2D+12); // "SCAT"
722 fOptions2D->AddEntry("Profile", kOptions2D+13); // "SCAT"
723 fOptions2D->Associate(this);
724
725 fOptionsMenu = new TGPopupMenu(gClient->GetRoot());
726 fOptionsMenu->AddPopup("&General Options...", fOptionsGen);
727 fOptionsMenu->AddPopup("&1D Options", fOptions1D);
728 fOptionsMenu->AddPopup("&2D Options", fOptions2D);
730 fOptionsMenu->AddEntry("&Reset options", kOptionsReset);
731 //--- Help menu
732 fHelpMenu = new TGPopupMenu(gClient->GetRoot());
733 fHelpMenu->AddEntry("&About ROOT...", kHelpAbout);
734 fHelpMenu->AddEntry("&About TreeViewer...", kHelpAboutTV);
736 fHelpMenu->AddEntry("&Starting...", kHelpStart);
737 fHelpMenu->AddEntry("&Layout...", kHelpLayout);
738 fHelpMenu->AddEntry("&Open/Save", kHelpOpenSave);
739 fHelpMenu->AddEntry("&Dragging...", kHelpDragging);
740 fHelpMenu->AddEntry("&Editing expressions...",kHelpEditing);
741 fHelpMenu->AddEntry("&Session...", kHelpSession);
742 fHelpMenu->AddEntry("&User commands...", kHelpCommands);
743 fHelpMenu->AddEntry("&Context menus...", kHelpContext);
744 fHelpMenu->AddEntry("D&rawing...", kHelpDrawing);
745 fHelpMenu->AddEntry("&Macros...", kHelpMacros);
746
747 fFileMenu->Associate(this);
748 fEditMenu->Associate(this);
749 fRunMenu->Associate(this);
750 fOptionsMenu->Associate(this);
751 fHelpMenu->Associate(this);
752
753 //--- menubar layout hints
757 //--- create menubar and add popup menus
758 fMenuBar = new TGMenuBar(this, 1, 1, kHorizontalFrame);
759
765
767 //--- toolbar ----------------------------------------------------------------
768 fToolBar = new TGToolBar(this, 10, 10, kHorizontalFrame);
770
771 TGLayoutHints *lo;
772 lo = new TGLayoutHints(kLHintsLeft | kLHintsCenterY, 4,4,0,0);
773 fWidgets->Add(lo);
774 //--- label for Command text entry
775 fBarLbl1 = new TGLabel(fToolBar,"Command");
777 //--- command text entry
779 fBarCommand->SetWidth(120);
780 fBarCommand->Associate(this);
781 fBarCommand->SetToolTipText("User commands executed via interpreter. Type <ENTER> to execute");
783 //--- first vertical separator
784 TGVertical3DLine *vSeparator = new TGVertical3DLine(fToolBar);
785 lo = new TGLayoutHints(kLHintsLeft | kLHintsExpandY, 4,4,0,0);
786 fWidgets->Add(lo);
787 fWidgets->Add(vSeparator);
788 fToolBar->AddFrame(vSeparator, lo);
789
790 lo = new TGLayoutHints(kLHintsLeft | kLHintsCenterY, 4,4,0,0);
791 fWidgets->Add(lo);
792 //--- label for Option text entry
793 fBarLbl2 = new TGLabel(fToolBar,"Option");
795 //--- drawing option text entry
797 fBarOption->SetWidth(100);
798 fBarOption->Associate(this);
799 fBarOption->SetToolTipText("Histogram graphics option. Type option here and click <Draw> (or <ENTER> to update current histogram).");
801 //--- second vertical separator
802 vSeparator = new TGVertical3DLine(fToolBar);
803 lo = new TGLayoutHints(kLHintsLeft | kLHintsExpandY, 4,4,0,0);
804 fWidgets->Add(lo);
805 fWidgets->Add(vSeparator);
806 fToolBar->AddFrame(vSeparator, lo);
807
808 lo = new TGLayoutHints(kLHintsLeft | kLHintsCenterY, 4,4,0,0);
809 fWidgets->Add(lo);
810 //--- label for Histogram text entry
811 fBarLbl3 = new TGLabel(fToolBar,"Histogram");
813
814 //--- histogram name text entry
815 lo = new TGLayoutHints(kLHintsCenterY | kLHintsExpandX, 4,4,0,0);
816 fWidgets->Add(lo);
817 fBarHist = new TGTextEntry(fToolBar, new TGTextBuffer(100));
820 fBarHist->SetText("htemp");
822
823 //--- Hist check button
824 lo = new TGLayoutHints(kLHintsLeft | kLHintsCenterY, 4,4,0,0);
825 fWidgets->Add(lo);
826 fBarH = new TGCheckButton(fToolBar, "Hist");
827 fBarH->SetToolTipText("Checked : redraw only current histogram");
829 fToolBar->AddFrame(fBarH, lo);
830 //--- Scan check button
831 fBarScan = new TGCheckButton(fToolBar, "Scan");
833 fBarScan->SetToolTipText("Check to redirect TTree::Scan command in a file");
835 //--- Rec check button
836 fBarRec = new TGCheckButton(fToolBar, "Rec");
838 fBarRec->SetToolTipText("Check to record commands in history file and be verbose");
840 //--- 1'st horizontal tool bar separator ----------------------------------------
841 TGHorizontal3DLine *toolBarSep = new TGHorizontal3DLine(this);
842 fWidgets->Add(toolBarSep);
843 AddFrame(toolBarSep, fBarLayout);
845 //--- 2'nd horizontal tool bar separator ----------------------------------------
846 toolBarSep = new TGHorizontal3DLine(this);
847 fWidgets->Add(toolBarSep);
848 AddFrame(toolBarSep, fBarLayout);
849
850 //--- Horizontal mother frame ---------------------------------------------------
851 fHf = new TGHorizontalFrame(this, 10, 10);
852 //--- Vertical frames
854// fSlider->SetBackgroundColor(color);
855 fSlider->Associate(this);
856
857 //--- fV1 -----------------------------------------------------------------------
858 fV1 = new TGVerticalFrame(fHf, 10, 10, kFixedWidth);
860
861 fLbl1 = new TGLabel(fTreeHdr, "Current Folder");
862 lo = new TGLayoutHints(kLHintsLeft | kLHintsTop | kLHintsCenterY, 3, 0, 0, 0);
863 fWidgets->Add(lo);
864 fTreeHdr->AddFrame(fLbl1, lo);
865
866 lo = new TGLayoutHints(kLHintsTop | kLHintsExpandX, 2, 0, 1, 0);
867 fWidgets->Add(lo);
868 fV1->AddFrame(fTreeHdr, lo);
869
870 //--- tree view canvas on the left
872 //--- container frame
874 GetWhitePixel());
875 fLt->Associate(this);
877
878 lo = new TGLayoutHints(kLHintsExpandX | kLHintsExpandY, 2,0,0,0);
879 fWidgets->Add(lo);
880 fV1->AddFrame(fTreeView, lo);
881
882 //--- button horizontal frame
884
885 //--- DRAW button
886 fPicDraw = gClient->GetPicture("draw_t.xpm");
888 fDRAW->SetToolTipText("Draw current selection");
889 fDRAW->Associate(this);
890
891 lo = new TGLayoutHints(kLHintsTop | kLHintsLeft, 2,2,4,2);
892 fWidgets->Add(lo);
893 fHpb->AddFrame(fDRAW, lo);
894
895 //--- SPIDER button
896 fSPIDER = new TGTextButton(fHpb,"SPIDER");
897 fSPIDER->SetToolTipText("Scan current selection using a spider plot");
898 fSPIDER->Associate(this);
899
900 lo = new TGLayoutHints(kLHintsTop | kLHintsLeft, 2,2,4,2);
901 fWidgets->Add(lo);
902 fHpb->AddFrame(fSPIDER,lo);
903 //---connect SPIDER button to ExecuteScan() method
904 fSPIDER->Connect("Clicked()","TTreeViewer",this,"ExecuteSpider()");
905
906 //--- STOP button (breaks current operation)
907// fPicStop = gClient->GetPicture("mb_stop_s.xpm");
908 fPicStop = gClient->GetPicture("stop_t.xpm");
910 fSTOP->SetToolTipText("Abort current operation");
911 fSTOP->Associate(this);
912
913 lo = new TGLayoutHints(kLHintsTop | kLHintsLeft, 2,2,4,2);
914 fWidgets->Add(lo);
915 fHpb->AddFrame(fSTOP, lo);
916
917 //--- REFR button (breaks current operation)
918 fPicRefr = gClient->GetPicture("refresh2.xpm");
920 fREFR->SetToolTipText("Update the tree viewer");
921 lo = new TGLayoutHints(kLHintsTop | kLHintsLeft, 2,2,4,2);
922 fWidgets->Add(lo);
923 fHpb->AddFrame(fREFR, lo);
924 //---connect REFR button to DoRefresh() method
925 fREFR->Connect("Clicked()", "TTreeViewer", this, "DoRefresh()");
926
927 lo = new TGLayoutHints(kLHintsTop | kLHintsLeft, 2,2,2,2);
928 fWidgets->Add(lo);
929 fV1->AddFrame(fHpb, lo);
930
931 //--- fV2
932 fV2 = new TGVerticalFrame(fHf, 10, 10);
934 fLbl2 = new TGLabel(fListHdr, "Current Tree: ");
935 lo = new TGLayoutHints(kLHintsTop | kLHintsLeft, 3, 0, 0, 0);
936 fWidgets->Add(lo);
937 fListHdr->AddFrame(fLbl2, lo);
938
939 //--- progress bar
943 lo = new TGLayoutHints(kLHintsBottom | kLHintsExpandX, 2,2,4,2);
944 fWidgets->Add(lo);
947 fWidgets->Add(lo);
948 fV2->AddFrame(fListHdr, lo);
949
952 fWidgets->Add(lo);
953 fHf->AddFrame(fSlider, lo);
955 fWidgets->Add(lo);
956 fHf->AddFrame(fV1, lo);
957
958 //--- vertical splitter
960 splitter->SetFrame(fV1,kTRUE);
963 fWidgets->Add(lo);
964 fHf->AddFrame(splitter,lo);
965
966
967
968 //-- listview for the content of the tree/branch -----------------------------
969 fListView = new TGListView(fListHdr,400,300);
970 //--- container frame
972 fLVContainer->Associate(this);
974 fLVContainer->SetViewer(this);
980 fWidgets->Add(lo);
981
983
985 fWidgets->Add(lo);
986 fHf->AddFrame(fV2,lo);
987
988 AddFrame(fHf, lo);
989 //--- 3rd horizontal tool bar separator ----------------------------------------
990 toolBarSep = new TGHorizontal3DLine(this);
991 fWidgets->Add(toolBarSep);
992 AddFrame(toolBarSep, fBarLayout);
993
994 //--- label for IList text entry
995 fBFrame = new TGHorizontalFrame(this,10,10);
996 fBLbl4 = new TGLabel(fBFrame,"IList");
997 lo = new TGLayoutHints(kLHintsLeft | kLHintsBottom, 2,2,2,2);
998 fWidgets->Add(lo);
999 fBFrame->AddFrame(fBLbl4, lo);
1000 //--- IList text entry
1001 fBarListIn = new TGTextEntry(fBFrame, new TGTextBuffer(100));
1002 fBarListIn->SetWidth(60);
1003 fBarListIn->SetToolTipText("Name of a previously created event list");
1005 //--- label for OList text entry
1006 fBLbl5 = new TGLabel(fBFrame,"OList");
1007 fBFrame->AddFrame(fBLbl5, lo);
1008 //--- OList text entry
1010 fBarListOut->SetWidth(60);
1011 fBarListOut->SetToolTipText("Output event list. Use <Draw> to generate it.");
1013 //--- Status bar
1014 fStatusBar = new TGStatusBar(fBFrame, 10, 10);
1015 fStatusBar->SetWidth(200);
1018 fWidgets->Add(lo);
1020 //--- RESET button
1021 fReset = new TGTextButton(fBFrame,"RESET",kRESET);
1022 fReset->SetToolTipText("Reset variable's fields and drawing options");
1023 fReset->Associate(this);
1024 lo = new TGLayoutHints(kLHintsTop | kLHintsRight, 2,2,2,2);
1025 fWidgets->Add(lo);
1026 fBFrame->AddFrame(fReset,lo);
1027 //--- group of buttons for session handling
1029 gClient->GetPicture("first_t.xpm"), kBGFirst);
1030 fBGFirst->SetToolTipText("First record");
1031 fBGFirst->Associate(this);
1033 gClient->GetPicture("previous_t.xpm"), kBGPrevious);
1034 fBGPrevious->SetToolTipText("Previous record");
1035 fBGPrevious->Associate(this);
1037 gClient->GetPicture("record_t.xpm"), kBGRecord);
1038 fBGRecord->SetToolTipText("Record");
1039 fBGRecord->Associate(this);
1041 gClient->GetPicture("next_t.xpm"), kBGNext);
1042 fBGNext->SetToolTipText("Next record");
1043 fBGNext->Associate(this);
1045 gClient->GetPicture("last_t.xpm"), kBGLast);
1046 fBGLast->SetToolTipText("Last record");
1047 fBGLast->Associate(this);
1048
1049 fCombo = new TGComboBox(fBFrame, 0);
1051 fCombo->SetWidth(100);
1052 fCombo->Associate(this);
1053
1054 lo = new TGLayoutHints(kLHintsCenterY | kLHintsRight, 0,0,2,0);
1055 fWidgets->Add(lo);
1056 fBFrame->AddFrame(fCombo, lo);
1057 fBFrame->AddFrame(fBGLast, lo);
1058 fBFrame->AddFrame(fBGNext, lo);
1062 lo = new TGLayoutHints(kLHintsExpandX,2,2,2,0);
1063 fWidgets->Add(lo);
1064 AddFrame(fBFrame,lo);
1065
1066 // map the window
1067 SetWindowName("TreeViewer");
1068 MapSubwindows();
1070 MapWindow();
1071
1072 // put default items in the listview on the right
1073 const TGPicture *pic, *spic;
1074
1076 TTVLVEntry* entry;
1077 Char_t symbol;
1079 symbol = 'X';
1080 entry->SetUserData(new ULong_t((symbol << 8) | kLTExpressionType | kLTTreeType));
1081 entry->SetToolTipText("X expression. Drag and drop expressions here");
1082 //--- X item
1083 fLVContainer->AddThisItem(entry);
1084 entry->Empty();
1085 entry->MapWindow();
1086
1088 symbol = 'Y';
1089 entry->SetUserData(new ULong_t((symbol << 8) | kLTExpressionType | kLTTreeType));
1090 entry->SetToolTipText("Y expression. Drag and drop expressions here");
1091 //--- Y item
1092 fLVContainer->AddThisItem(entry);
1093 entry->Empty();
1094 entry->MapWindow();
1095
1097 symbol = 'Z';
1098 entry->SetUserData(new ULong_t((symbol << 8) | kLTExpressionType | kLTTreeType));
1099 entry->SetToolTipText("Z expression. Drag and drop expressions here");
1100 //--- Z item
1101 fLVContainer->AddThisItem(entry);
1102 entry->Empty();
1103 entry->MapWindow();
1104
1105 pic = gClient->GetPicture("cut_t.xpm");
1106 spic = gClient->GetPicture("cut_t.xpm");
1107 entry = new TTVLVEntry(fLVContainer,pic,spic,new TGString(),0,kLVSmallIcons);
1109 entry->SetToolTipText("Active cut. Double-click to enable/disable");
1110 //--- Cut item (scissors icon)
1111 fLVContainer->AddThisItem(entry);
1112 entry->Empty();
1113 entry->MapWindow();
1114
1115 pic = gClient->GetPicture("pack_t.xpm");
1116 spic = gClient->GetPicture("pack-empty_t.xpm");
1117 entry = new TTVLVEntry(fLVContainer,pic,spic,new TGString("Scan box"),0,kLVSmallIcons);
1119 entry->SetToolTipText("Drag and drop expressions/leaves here. Double-click to scan. Check <Scan> to redirect on file.");
1120 //--- Scan Box
1121 fLVContainer->AddThisItem(entry);
1122 entry->MapWindow();
1123 entry->SetTrueName("");
1124
1125 //--- 10 expression items
1126 fNexpressions = 10;
1127 for (Int_t i=0; i<fNexpressions; i++) {
1128 pic = gClient->GetPicture("expression_t.xpm");
1129 spic = gClient->GetPicture("expression_t.xpm");
1130 entry = new TTVLVEntry(fLVContainer,pic,spic,new TGString(),0,kLVSmallIcons);
1132 entry->SetToolTipText("User defined expression/cut. Double-click to edit");
1133 fLVContainer->AddThisItem(entry);
1134 entry->Empty();
1135 entry->MapWindow();
1136 }
1137
1138 fListView->Layout();
1139 fListView->Resize();
1140// EmptyAll();
1141 // map the tree if it was supplied in the constructor
1142
1143 if (!fTree) {
1144 fSlider->SetRange(0,1000000);
1145 fSlider->SetPosition(0,1000000);
1146 } else {
1149 }
1150 PrintEntries();
1154
1155 // map the window
1156 ///SetWindowName("TreeViewer");
1157 MapSubwindows();
1159 MapWindow();
1160}
1161
1162////////////////////////////////////////////////////////////////////////////////
1163/// TTreeViewer destructor.
1164
1166{
1167 if (!gClient) return;
1168 gClient->FreePicture(fPicX);
1169 gClient->FreePicture(fPicY);
1170 gClient->FreePicture(fPicZ);
1171 gClient->FreePicture(fPicDraw);
1172 gClient->FreePicture(fPicStop);
1173 gClient->FreePicture(fPicRefr);
1174
1176 if (fDialogBox) delete fDialogBox;
1177
1178 delete fContextMenu;
1179
1180 delete fBarLbl1;
1181 delete fBarLbl2;
1182 delete fBarLbl3;
1183 delete fBLbl4;
1184 delete fBLbl5;
1185 delete fBarCommand;
1186 delete fBarOption;
1187 delete fBarHist;
1188 delete fBarListIn;
1189 delete fBarListOut;
1190
1191 delete fBarH;
1192 delete fBarScan;
1193 delete fBarRec;
1194
1195 delete fToolBar;
1196
1197 delete fSlider;
1198 delete fV1;
1199 delete fV2;
1200 delete fLbl1;
1201 delete fLbl2;
1202 delete fHf;
1203 delete fTreeHdr;
1204 delete fListHdr;
1205 delete fLt;
1206 delete fTreeView;
1207 delete fLVContainer;
1208 delete fListView;
1209
1210 delete fProgressBar;
1211 delete fHpb;
1212
1213 delete fDRAW;
1214 delete fSPIDER;
1215 delete fSTOP;
1216 delete fReset;
1217 delete fBGFirst;
1218 delete fBGPrevious;
1219 delete fBGRecord;
1220 delete fBGNext;
1221 delete fBGLast;
1222 delete fCombo;
1223 delete fBFrame;
1224
1225 delete fMenuBar;
1226 delete fFileMenu;
1227 delete fEditMenu;
1228
1229 delete fOptionsGen;
1230 delete fOptions1D;
1231 delete fOptions2D;
1232 delete fOptionsMenu;
1233 delete fHelpMenu;
1234 delete fMenuBarLayout;
1235 delete fMenuBarItemLayout;
1236 delete fMenuBarHelpLayout;
1237 delete fBarLayout;
1238
1239 fWidgets->Delete();
1240 delete fWidgets;
1241 if (fTreeList) {
1242 delete fTreeList;
1243 }
1244 delete fTimer;
1245 delete fSession;
1246}
1247
1248////////////////////////////////////////////////////////////////////////////////
1249/// Enable/disable session buttons.
1250
1252 Bool_t next, Bool_t last)
1253{
1256 if (previous) fBGPrevious->SetState(kButtonUp);
1258 if (next) fBGNext->SetState(kButtonUp);
1260 if (last) fBGLast->SetState(kButtonUp);
1262}
1263
1264////////////////////////////////////////////////////////////////////////////////
1265/// Apply Cut
1266
1267const char* TTreeViewer::Cut()
1268{
1269 return fLVContainer->Cut();
1270}
1271
1272////////////////////////////////////////////////////////////////////////////////
1273/// returns scanlist
1274
1275const char* TTreeViewer::ScanList()
1276{
1277 return fLVContainer->ScanList();
1278}
1279
1280////////////////////////////////////////////////////////////////////////////////
1281/// Set current session
1282
1284{
1285 if (session) {
1286 delete fSession;
1287 fSession = session;
1288 }
1289}
1290
1291////////////////////////////////////////////////////////////////////////////////
1292/// Empty the bracket content of a string.
1293
1294const char* TTreeViewer::EmptyBrackets(const char* name)
1295{
1296 TString stripped(name);
1297 if (!stripped.Contains("[")) return name;
1298 TString retstr(name);
1299 TObjString *objstr;
1300 Int_t index = 0;
1301 while (stripped.Index("[", index) != kNPOS) {
1302 Int_t start = stripped.Index("[", index);
1303 Int_t end = stripped.Index("]", index);
1304 if (end == kNPOS) {
1305 objstr = new TObjString(retstr.Data());
1306 fWidgets->Add(objstr);
1307 return (objstr->String()).Data();
1308 }
1309 index = start+2;
1310 retstr = stripped.Remove(start+1, end-start-1);
1311 stripped = retstr;
1312 }
1313 objstr = new TObjString(retstr.Data());
1314 fWidgets->Add(objstr);
1315 return (objstr->String()).Data();
1316}
1317
1318////////////////////////////////////////////////////////////////////////////////
1319/// Clear the content of all items in the list view.
1320
1322{
1324}
1325
1326////////////////////////////////////////////////////////////////////////////////
1327/// Empty the content of the selected expression.
1328
1329void TTreeViewer::Empty()
1330{
1331 void *p = 0;
1332 TTVLVEntry *item = 0;
1333 if ((item = (TTVLVEntry *) fLVContainer->GetNextSelected(&p)) == 0) {
1334 Warning("Empty", "No item selected.");
1335 return;
1336 }
1337 ULong_t *itemType = (ULong_t *) item->GetUserData();
1338 if (!(*itemType & kLTExpressionType)) {
1339 Warning("Empty", "Not expression type.");
1340 return;
1341 }
1342 if (*itemType & kLTPackType) {
1343 item->SetSmallPic(fClient->GetPicture("pack-empty_t.xpm"));
1344 item->SetTrueName("");
1345 return;
1346 }
1347 item->Empty();
1348}
1349
1350////////////////////////////////////////////////////////////////////////////////
1351/// Get the item from a specific position.
1352
1354{
1355 return fLVContainer->ExpressionItem(index);
1356}
1357
1358////////////////////////////////////////////////////////////////////////////////
1359/// Get the list of expression items.
1360
1362{
1363 return fLVContainer->ExpressionList();
1364}
1365
1366////////////////////////////////////////////////////////////////////////////////
1367/// Compute dimension of the histogram.
1368
1370{
1371 fDimension = 0;
1372 if (Ex() && strlen(Ex())) fDimension++;
1373 if (Ey() && strlen(Ey())) fDimension++;
1374 if (Ez() && strlen(Ez())) fDimension++;
1375 return fDimension;
1376}
1377
1378////////////////////////////////////////////////////////////////////////////////
1379/// Called when the DRAW button is executed.
1380
1382{
1383 TString varexp;
1384 TString command;
1385 Int_t dimension = 0;
1386 TString alias[3];
1387 TTVLVEntry *item;
1388 Int_t i;
1389 // fill in expressions
1390 if (fVarDraw) {
1391 void *p = 0;
1392 dimension = 1;
1393 if (!(item = (TTVLVEntry *) fLVContainer->GetNextSelected(&p))) return;
1394 alias[0] = item->GetAlias();
1395 if (alias[0].BeginsWith("~")) alias[0].Remove(0, 1);
1396 varexp = item->ConvertAliases();
1397 } else {
1398 if (Ez() && strlen(Ez())) {
1399 dimension++;
1400 varexp = Ez();
1401 item = ExpressionItem(2);
1402 alias[2] = item->GetAlias();
1403 if (alias[2].BeginsWith("~")) alias[2].Remove(0, 1);
1404 }
1405 if ((Ez() && strlen(Ez())) && ((Ex() &&strlen(Ex())) || (Ey() && strlen(Ey())))) varexp += ":";
1406 if (Ey() && strlen(Ey())) {
1407 dimension++;
1408 varexp += Ey();
1409 item = ExpressionItem(1);
1410 alias[1] = item->GetAlias();
1411 if (alias[1].BeginsWith("~")) alias[1].Remove(0, 1);
1412 }
1413 if (Ey() && strlen(Ey()) && Ex() && strlen(Ex())) varexp += ":";
1414 if (Ex () && strlen(Ex())) {
1415 dimension++;
1416 varexp += Ex();
1417 item = ExpressionItem(0);
1418 alias[0] = item->GetAlias();
1419 if (alias[0].BeginsWith("~")) alias[0].Remove(0, 1);
1420 }
1421 }
1422 if (!dimension && !fScanMode) {
1423 Warning("ExecuteDraw", "Nothing to draw on X,Y,Z.");
1424 return;
1425 }
1426 // find ListIn
1427 fTree->SetEventList(0);
1428 TEventList *elist = 0;
1429 if (strlen(fBarListIn->GetText())) {
1430 elist = (TEventList *) gROOT->FindObject(fBarListIn->GetText());
1431 if (elist) fTree->SetEventList(elist);
1432 }
1433 // find ListOut
1434 if (strlen(fBarListOut->GetText())) varexp = TString::Format(">>%s", fBarListOut->GetText());
1435 // find histogram name
1436 if (strcmp("htemp", fBarHist->GetText())) {
1437 varexp += ">>";
1438 varexp += fBarHist->GetText();
1439 }
1440 // find canvas/pad where to draw
1441 TPad *pad = (TPad*)gROOT->GetSelectedPad();
1442 if (pad) pad->cd();
1443 // find graphics option
1444 const char* gopt = fBarOption->GetText();
1445 // just in case a previous interrupt was posted
1446 gROOT->SetInterrupt(kFALSE);
1447 // check if cut is enabled
1448 const char *cut = "";
1449 if (fEnableCut) cut = Cut();
1450
1451 // get entries to be processed
1453 fSlider->GetMinPosition() + 1);
1454 Long64_t firstentry =(Long64_t) fSlider->GetMinPosition();
1455//printf("firstentry=%lld, nentries=%lld\n",firstentry,nentries);
1456 // check if Scan is checked and if there is something in the box
1457 if (fScanMode) {
1458// fBarScan->SetState(kButtonUp);
1459 fScanMode = kFALSE;
1460 if (ScanList() && strlen(ScanList())) varexp = ScanList();
1461 command = TString::Format("tv__tree->Scan(\"%s\",\"%s\",\"%s\", %lld, %lld);",
1462 varexp.Data(), cut, gopt, nentries, firstentry);
1463 if (fBarScan->GetState() == kButtonDown) {
1465 } else {
1467 }
1468 ExecuteCommand(command.Data(), kTRUE);
1469 return;
1470 }
1471 // check if only histogram has to be updated
1472 if (fBarH->GetState() == kButtonDown) {
1473 // reset 'Hist' mode
1475 TH1 *hist = fTree->GetHistogram();
1476 if (hist && gPad) {
1477 //hist = (TH1*)gPad->GetListOfPrimitives()->FindObject(fBarHist->GetText());
1478 if (hist) {
1479 // check if graphic option was modified
1480 TString last(fLastOption);
1481 TString current(gopt);
1482 current.ToUpper();
1483 last.ToUpper();
1484 if (current == last) {
1485 gPad->Update();
1486 return;
1487 }
1488 if (dimension == 3 && strlen(gopt)) {
1489 std::cout << "Graphics option " << gopt << " not valid for 3D histograms" << std::endl;
1490 return;
1491 }
1492 std::cout << " Graphics option for current histogram changed to " << gopt << std::endl;
1493 hist->Draw(gopt);
1495 gPad->Update();
1496 return;
1497 }
1498 }
1499 }
1500 // send draw command
1502 //if (!gopt[0] && dimension!=3) {
1503 // gopt = "hist";
1504 // fLastOption = "hist";
1505 //}
1506 if (dimension == 3 && strlen(gopt)) {
1507 std::cout << "Graphics option " << gopt << " not valid for 3D histograms" << std::endl;
1508 gopt = "";
1509 fLastOption = "";
1510 }
1511 command = TString::Format("tv__tree->Draw(\"%s\",\"%s\",\"%s\", %lld, %lld);",
1512 varexp.Data(), cut, gopt, nentries, firstentry);
1513 if (fCounting) return;
1514 fCounting = kTRUE;
1515 fTree->SetTimerInterval(200);
1516 fTimer->TurnOn();
1517 ExecuteCommand(command.Data());
1519 fTimer->TurnOff();
1521 fCounting = kFALSE;
1524 TH1 *hist = fTree->GetHistogram();
1525 if (hist) {
1526 // put expressions aliases on axes
1527 Int_t current = 0;
1528 for (i=0; i<3; i++) {
1529 if (alias[i].Length()) {
1530 if (i != current) {
1531 alias[current] = alias[i];
1532 alias[i] = "";
1533 }
1534 current++;
1535 }
1536 }
1537 //hist = (TH1*)gPad->GetListOfPrimitives()->FindObject(fBarHist->GetText());
1538 TAxis *axis[3];
1539 axis[0] = hist->GetXaxis();
1540 axis[1] = hist->GetYaxis();
1541 axis[2] = hist->GetZaxis();
1542 for (Int_t ind=0; ind<3; ind++) axis[ind]->SetTitle(alias[ind].Data());
1543 }
1544 if (gPad) gPad->Update();
1545}
1546
1547////////////////////////////////////////////////////////////////////////////////
1548/// Draw a spider plot for the selected entries.
1549
1551{
1552 TString varexp;
1553 Int_t dimension = 0;
1554 TString alias[3];
1555 TTVLVEntry *item;
1556 Bool_t previousexp = kFALSE;
1557 // fill in expressions
1558 if (Ez() && strlen(Ez())) {
1559 previousexp = kTRUE;
1560 dimension++;
1561 varexp = Ez();
1562 item = ExpressionItem(2);
1563 alias[2] = item->GetAlias();
1564 if (alias[2].BeginsWith("~")) alias[2].Remove(0, 1);
1565 }
1566 if ((Ez() && strlen(Ez())) && ((Ex() && strlen(Ex())) || (Ey() && strlen(Ey())))) varexp += ":";
1567 if (Ey() && strlen(Ey())) {
1568 previousexp = kTRUE;
1569 dimension++;
1570 varexp += Ey();
1571 item = ExpressionItem(1);
1572 alias[1] = item->GetAlias();
1573 if (alias[1].BeginsWith("~")) alias[1].Remove(0, 1);
1574 }
1575 if (Ey() && strlen(Ey()) && Ex() && strlen(Ex())) varexp += ":";
1576 if (Ex() && strlen(Ex())) {
1577 previousexp = kTRUE;
1578 dimension++;
1579 varexp += Ex();
1580 item = ExpressionItem(0);
1581 alias[0] = item->GetAlias();
1582 if (alias[0].BeginsWith("~")) alias[0].Remove(0, 1);
1583 }
1584 for(Int_t i=0;i<10;++i){
1585 if(En(i+5) && strlen(En(i+5))){
1586 ++dimension;
1587 if(previousexp){
1588 varexp += ":";
1589 varexp += En(i+5);
1590 } else varexp = En(i+5);
1591 previousexp = kTRUE;
1592 }
1593 }
1594 if (dimension<3) {
1595 Warning("ExecuteSpider", "Need at least 3 variables");
1596 return;
1597 }
1598 // find ListIn
1599 fTree->SetEventList(0);
1600 TEventList *elist = 0;
1601 if (strlen(fBarListIn->GetText())) {
1602 elist = (TEventList *) gROOT->FindObject(fBarListIn->GetText());
1603 if (elist) fTree->SetEventList(elist);
1604 }
1605 // find ListOut
1606 if (strlen(fBarListOut->GetText())) varexp = TString::Format(">>%s", fBarListOut->GetText());
1607 // find canvas/pad where to draw
1608 TPad *pad = (TPad*)gROOT->GetSelectedPad();
1609 if (pad) pad->cd();
1610 // find graphics option
1611 const char* gopt = fBarOption->GetText();
1612 // just in case a previous interrupt was posted
1613 gROOT->SetInterrupt(kFALSE);
1614 // check if cut is enabled
1615 const char *cut = "";
1616 if (fEnableCut) cut = Cut();
1617
1618 // get entries to be processed
1620 fSlider->GetMinPosition() + 1);
1621 Long64_t firstentry =(Long64_t) fSlider->GetMinPosition();
1622
1623 // create the spider plot
1624
1625 TSpider* spider = new TSpider(fTree,varexp.Data(),cut,Form("%s spider average",gopt),nentries,firstentry);
1626 spider->Draw();
1627
1628 if (gPad) gPad->Update();
1629}
1630
1631////////////////////////////////////////////////////////////////////////////////
1632/// Get the expression to be drawn on X axis.
1633
1634const char* TTreeViewer::Ex()
1635{
1636 return fLVContainer->Ex();
1637}
1638
1639////////////////////////////////////////////////////////////////////////////////
1640/// Get the expression to be drawn on Y axis.
1641
1642const char* TTreeViewer::Ey()
1643{
1644 return fLVContainer->Ey();
1645}
1646
1647////////////////////////////////////////////////////////////////////////////////
1648/// Get the expression to be drawn on Z axis.
1649
1650const char* TTreeViewer::Ez()
1651{
1652 return fLVContainer->Ez();
1653}
1654
1655////////////////////////////////////////////////////////////////////////////////
1656/// Get the n'th expression
1657
1658const char* TTreeViewer::En(Int_t n)
1659{
1661 if(e) return e->ConvertAliases();
1662 return "";
1663}
1664
1665////////////////////////////////////////////////////////////////////////////////
1666/// Start the expression editor.
1667
1669{
1670 void *p = 0;
1671 // get the selected item
1672 TTVLVEntry *item = 0;
1673 if ((item = (TTVLVEntry *) fLVContainer->GetNextSelected(&p)) == 0) {
1674 Warning("EditExpression", "No item selected.");
1675 return;
1676 }
1677 // check if it is an expression
1678 ULong_t *itemType = (ULong_t *) item->GetUserData();
1679 if (!(*itemType & kLTExpressionType)) {
1680 Warning("EditExpression", "Not expression type.");
1681 return;
1682 }
1683 // check if the editor is already active
1685 if (!fDialogBox) {
1686 fDialogBox = new TGSelectBox(fClient->GetRoot(), this, fV1->GetWidth() - 10);
1687 }
1688 // copy current item data into editor boxes
1689 fDialogBox->SetEntry(item);
1690 fDialogBox->SetWindowName("Expression editor");
1691 // check if you are editing the cut expression
1692 if (*itemType & kLTCutType || item->IsCut()) {
1693 fDialogBox->SetLabel("Selection");
1694 } else {
1695 fDialogBox->SetLabel("Expression");
1696 }
1697}
1698
1699////////////////////////////////////////////////////////////////////////////////
1700/// Get use of TTree::MakeSelector() via the context menu.
1701
1702Int_t TTreeViewer::MakeSelector(const char* selector)
1703{
1704 if (!fTree) return 0;
1705 return fTree->MakeSelector(selector);
1706}
1707
1708////////////////////////////////////////////////////////////////////////////////
1709/// Get use of TTree::Process() via the context menu.
1710
1711Long64_t TTreeViewer::Process(const char* filename, Option_t *option, Long64_t nentries, Long64_t firstentry)
1712{
1713 if (!fTree) return 0;
1714 return fTree->Process(filename, option, nentries, firstentry);
1715}
1716
1717////////////////////////////////////////////////////////////////////////////////
1718/// Get graph option
1719
1720const char *TTreeViewer::GetGrOpt()
1721{
1722 return fBarOption->GetText();
1723}
1724
1725////////////////////////////////////////////////////////////////////////////////
1726/// Set graph option
1727
1728void TTreeViewer::SetGrOpt(const char *option)
1729{
1730 fBarOption->SetText(option);
1731}
1732
1733////////////////////////////////////////////////////////////////////////////////
1734/// Return kTRUE if scan is redirected
1735
1737{
1738 return (fBarScan->GetState()==kButtonDown);
1739}
1740
1741////////////////////////////////////////////////////////////////////////////////
1742/// Remove the selected item from the list.
1743
1745{
1746 void *p = 0;
1747 TTVLVEntry *item = 0;
1748 // get the selected item
1749 if ((item = (TTVLVEntry *) fLVContainer->GetNextSelected(&p)) == 0) {
1750 Warning("RemoveItem", "No item selected.");
1751 return;
1752 }
1753 // check if it is removable
1754 ULong_t *itemType = (ULong_t *) item->GetUserData();
1755 if (!(*itemType & kLTDragType)) {
1756 Warning("RemoveItem", "Not removable type.");
1757 return;
1758 }
1759 fLVContainer->RemoveItem(item);
1760 fListView->Layout();
1761}
1762
1763////////////////////////////////////////////////////////////////////////////////
1764/// Remove the current record.
1765
1767{
1769}
1770
1771////////////////////////////////////////////////////////////////////////////////
1772/// This function is called by the fTimer object.
1773
1775{
1776 if (fCounting) {
1778 Float_t last = fSlider->GetMaxPosition();
1779 Float_t current = (Float_t)fTree->GetReadEntry();
1780 Float_t percent = (current-first+1)/(last-first+1);
1781 fProgressBar->SetPosition(100.*percent);
1783 }
1784 timer->Reset();
1785 return kFALSE;
1786}
1787
1788////////////////////////////////////////////////////////////////////////////////
1789/// Handle menu and other commands generated.
1790
1792{
1793 TRootHelpDialog *hd;
1794 TTVRecord *record;
1795
1796 switch (GET_MSG(msg)) {
1797 case kC_VSLIDER :
1798 // handle slider messages
1799 PrintEntries();
1800 break;
1801 case kC_TEXTENTRY:
1802 switch (GET_SUBMSG(msg)) {
1803 // handle enter posted by the Command text entry
1804 case kTE_ENTER:
1805 if ((ERootTreeViewerCommands)parm1 == kBarCommand) {
1807 fBarCommand->Clear();
1808 }
1809 if ((ERootTreeViewerCommands)parm1 == kBarOption) {
1810 fVarDraw = kFALSE;
1812 ExecuteDraw();
1814 }
1815 break;
1816 default:
1817 break;
1818 }
1819 break;
1820 case kC_LISTTREE:
1821 switch (GET_SUBMSG(msg)) {
1822 // handle mouse messages in the list-tree (left panel)
1823 case kCT_ITEMCLICK :
1824 // tell coverity that parm1 is a Long_t, and not an enum (even
1825 // if we compare it with an enum value) and the meaning of
1826 // parm1 depends on GET_MSG(msg) and GET_SUBMSG(msg)
1827 // coverity[mixed_enums]
1828 if (((EMouseButton)parm1==kButton1) ||
1829 ((EMouseButton)parm1==kButton3)) {
1830 TGListTreeItem *ltItem = 0;
1831 // get item that sent this
1832 if ((ltItem = fLt->GetSelected()) != 0) {
1833 // get item type
1834 ULong_t *itemType = (ULong_t *)ltItem->GetUserData();
1835 if (!itemType)
1836 break;
1837 if (*itemType & kLTTreeType) {
1838 // already mapped tree item clicked
1839 Int_t index = (Int_t)(*itemType >> 8);
1840 SwitchTree(index);
1841 if (fTree != fMappedTree) {
1842 // switch also the global "tree" variable
1844 // map it on the right panel
1845 MapTree(fTree);
1846 fListView->Layout();
1847 }
1848 // activate context menu for this tree
1849 if (parm1 == kButton3) {
1850 Int_t x = (Int_t)(parm2 &0xffff);
1851 Int_t y = (Int_t)((parm2 >> 16) & 0xffff);
1853 }
1854 }
1855
1856 if (*itemType & kLTBranchType) {
1857 // branch item clicked
1858 SetParentTree(ltItem);
1859 if (!fTree) break; // really needed ?
1860 TBranch *branch = fTree->GetBranch(ltItem->GetText());
1861 if (!branch) break;
1862 // check if it is mapped on the right panel
1863 if (branch != fMappedBranch) {
1865 MapBranch(branch);
1867 fListView->Layout();
1868 }
1869 // activate context menu for this branch (no *MENU* methods ):)
1870 if (parm1 == kButton3) {
1871 Int_t x = (Int_t)(parm2 &0xffff);
1872 Int_t y = (Int_t)((parm2 >> 16) & 0xffff);
1873 fContextMenu->Popup(x, y, branch);
1874 }
1875 }
1876
1877 if (*itemType & kLTLeafType) {
1878 // leaf item clicked
1879 SetParentTree(ltItem);
1880 if (!fTree) break;
1881 // find parent branch
1882 TBranch *branch = fTree->GetBranch(ltItem->GetParent()->GetText());
1883 if (!branch) {
1884 if (fTree != fMappedTree) {
1886 MapTree(fTree);
1887 fListView->Layout();
1888 }
1889 } else {
1890 // check if it is already mapped
1891 if (branch!=fMappedBranch) {
1893 MapBranch(branch);
1895 fListView->Layout();
1896 }
1897 }
1898 // select corresponding leaf on the right panel
1899 fLVContainer->SelectItem(ltItem->GetText());
1900 if (parm1 == kButton3) {
1901 // activate context menu for this leaf
1903 }
1904 }
1905 }
1906 }
1907 break;
1908 case kCT_ITEMDBLCLICK :
1910 if (parm1 == kButton1) {
1911 // execute double-click action for corresponding item in the right panel
1913 }
1914 break;
1915 default:
1916 break;
1917 }
1918 break;
1919 case kC_COMMAND:
1920 switch (GET_SUBMSG(msg)) {
1921 case kCM_COMBOBOX:
1922 if ((record = fSession->GetRecord((Int_t)parm2)))
1923 fSession->Show(record);
1924 break;
1925 case kCM_BUTTON:
1926 switch (parm1) {
1927 // handle button messages
1928 case kRESET:
1929 EmptyAll();
1930 break;
1931 case kDRAW:
1932 fVarDraw = kFALSE;
1933 ExecuteDraw();
1934 break;
1935 case kSTOP:
1936 if (fCounting)
1937 gROOT->SetInterrupt(kTRUE);
1938 break;
1939 case kCLOSE:
1941 break;
1942 case kBGFirst:
1943 if ((record = fSession->First()))
1944 fSession->Show(record);
1945 break;
1946 case kBGPrevious:
1947 if ((record = fSession->Previous()))
1948 fSession->Show(record);
1949 break;
1950 case kBGRecord:
1952 break;
1953 case kBGNext:
1954 if ((record = fSession->Next()))
1955 fSession->Show(record);
1956 break;
1957 case kBGLast:
1958 if ((record = fSession->Last()))
1959 fSession->Show(record);
1960 break;
1961 default:
1962 break;
1963 }
1964 break;
1965 case kCM_MENU:
1966 // handle menu messages
1967 // check if sent by Options menu
1968 if ((parm1>=kOptionsReset) && (parm1<kHelpAbout)) {
1969 Dimension();
1970 if ((fDimension==0) && (parm1>=kOptions1D)) {
1971 Warning("ProcessMessage", "Edit expressions first.");
1972 break;
1973 }
1974 if ((fDimension==1) && (parm1>=kOptions2D)) {
1975 Warning("ProcessMessage", "You have only one expression active.");
1976 break;
1977 }
1978 if ((fDimension==2) && (parm1>=kOptions1D) &&(parm1<kOptions2D)) {
1979 Warning("ProcessMessage", "1D drawing options not apply to 2D histograms.");
1980 break;
1981 }
1982 // make composed option
1983 MapOptions(parm1);
1984 break;
1985 }
1986 switch (parm1) {
1987 case kFileCanvas:
1988 gROOT->MakeDefCanvas();
1989 break;
1990 case kFileBrowse:
1991 if (1) {
1992 static TString dir(".");
1993 TGFileInfo info;
1994 info.fFileTypes = gOpenTypes;
1995 info.fIniDir = StrDup(dir);
1996 new TGFileDialog(fClient->GetRoot(), this, kFDOpen, &info);
1997 if (!info.fFilename) return kTRUE;
1998 dir = info.fIniDir;
1999 TString command = TString::Format("tv__tree_file = new TFile(\"%s\");",
2001 ExecuteCommand(command.Data());
2002 ExecuteCommand("tv__tree_file->ls();");
2003 std::cout << "Use SetTreeName() from context menu and supply a tree name" << std::endl;
2004 std::cout << "The context menu is activated by right-clicking the panel from right" << std::endl;
2005 }
2006 break;
2007 case kFileLoadLibrary:
2008 fBarCommand->SetText("gSystem->Load(\"\");");
2009 if (1) {
2010 Event_t event;
2011 event.fType = kButtonPress;
2012 event.fCode = kButton1;
2013 event.fX = event.fY = 1;
2014 fBarCommand->HandleButton(&event);
2015 }
2017 break;
2018 case kFileOpenSession:
2019 if (1) {
2020 static TString dir(".");
2021 TGFileInfo info;
2022 info.fFileTypes = gMacroTypes;
2023 info.fIniDir = StrDup(dir);
2024 new TGFileDialog(fClient->GetRoot(), this, kFDOpen, &info);
2025 if (!info.fFilename) return kTRUE;
2026 dir = info.fIniDir;
2027 gInterpreter->Reset();
2028 if (!gInterpreter->IsLoaded(info.fFilename)) gInterpreter->LoadMacro(info.fFilename);
2029 char command[1024];
2030 command[0] = 0;
2031 snprintf(command,1024,"open_session((void*)0x%lx);", (Long_t)this);
2032 ExecuteCommand(command);
2033 }
2034 break;
2035 case kFileSaveMacro:
2036 fContextMenu->Action(this,(TMethod*)IsA()->GetListOfMethods()->FindObject("SaveSource"));
2037 break;
2038 case kFilePrint:
2039 break;
2040 case kFileClose:
2042 break;
2043 case kFileQuit:
2045 break;
2046 case kEditExpression:
2048 break;
2049 case kEditCut:
2051 break;
2052 case kEditMacro:
2053 break;
2054 case kEditEvent:
2055 break;
2056 case kRunMacro:
2057 break;
2058 case kHelpAbout:
2059 {
2060#ifdef R__UNIX
2061 TString rootx = TROOT::GetBinDir() + "/root -a &";
2062 gSystem->Exec(rootx);
2063#else
2064#ifdef WIN32
2066#else
2067 char str[32];
2068 snprintf(str,32, "About ROOT %s...", gROOT->GetVersion());
2069 hd = new TRootHelpDialog(this, str, 600, 400);
2070 hd->SetText(gHelpAbout);
2071 hd->Popup();
2072#endif
2073#endif
2074 }
2075 break;
2076 case kHelpAboutTV:
2077 hd = new TRootHelpDialog(this, "About TreeViewer...", 600, 400);
2078 hd->SetText(gTVHelpAbout);
2079 hd->Resize(hd->GetDefaultSize());
2080 hd->Popup();
2081 break;
2082 case kHelpStart:
2083 hd = new TRootHelpDialog(this, "Quick start...", 600, 400);
2084 hd->SetText(gTVHelpStart);
2085 hd->Popup();
2086 break;
2087 case kHelpLayout:
2088 hd = new TRootHelpDialog(this, "Layout...", 600, 400);
2090 hd->Popup();
2091 break;
2092 case kHelpOpenSave:
2093 hd = new TRootHelpDialog(this, "Open/Save...", 600, 400);
2095 hd->Popup();
2096 break;
2097 case kHelpDragging:
2098 hd = new TRootHelpDialog(this, "Dragging items...", 600, 400);
2100 hd->Popup();
2101 break;
2102 case kHelpEditing:
2103 hd = new TRootHelpDialog(this, "Editing expressions...", 600, 400);
2105 hd->Popup();
2106 break;
2107 case kHelpSession:
2108 hd = new TRootHelpDialog(this, "Session...", 600, 400);
2110 hd->Popup();
2111 break;
2112 case kHelpCommands:
2113 hd = new TRootHelpDialog(this, "Executing user commands...", 600, 400);
2115 hd->Popup();
2116 break;
2117 case kHelpContext:
2118 hd = new TRootHelpDialog(this, "Context menus...", 600, 400);
2120 hd->Popup();
2121 break;
2122 case kHelpDrawing:
2123 hd = new TRootHelpDialog(this, "Drawing histograms...", 600, 400);
2125 hd->Popup();
2126 break;
2127 case kHelpMacros:
2128 hd = new TRootHelpDialog(this, "Using macros...", 600, 400);
2130 hd->Popup();
2131 break;
2132 default:
2133 break;
2134 }
2135 break;
2136 default:
2137 break;
2138 }
2139 break;
2140 case kC_CONTAINER:
2141 switch (GET_SUBMSG(msg)) {
2142 // handle messages sent from the listview (right panel)
2143 case kCT_SELCHANGED:
2144 break;
2145 case kCT_ITEMCLICK:
2146 // handle mouse messages
2147 switch (parm1) {
2148 case kButton1:
2149 if (fLVContainer->NumSelected()) {
2150 // get item that sent this
2151 void *p = 0;
2152 TTVLVEntry *item;
2153 if ((item = (TTVLVEntry *) fLVContainer->GetNextSelected(&p)) != 0) {
2154 const char* vname = item->GetTrueName();
2155 TString trueName(vname);
2156 if (trueName.Contains("[]")) {
2157 TIter next(fTree->GetListOfLeaves());
2158 TLeaf *leaf;
2159 while((leaf=(TLeaf*)next())) {
2160 if (!strcmp(vname, EmptyBrackets(leaf->GetName())))
2161 vname = leaf->GetName();
2162 }
2163 }
2164 char* msg2 = new char[2000];
2165 // get item type
2166 ULong_t *itemType = (ULong_t *) item->GetUserData();
2167 if (*itemType & kLTTreeType) {
2168 // X, Y or Z clicked
2169 char symbol = (char)((*itemType) >> 8);
2170 snprintf(msg2,2000, "%c expression : %s", symbol, vname);
2171 } else {
2172 if (*itemType & kLTCutType) {
2173 // scissors clicked
2174 snprintf(msg2,2000, "Cut : %s", vname);
2175 } else {
2176 if (*itemType & kLTPackType) {
2177 snprintf(msg2,2000, "Box : %s", vname);
2178 } else {
2179 if (*itemType & kLTExpressionType) {
2180 // expression clicked
2181 snprintf(msg2,2000, "Expression : %s", vname);
2182 } else {
2183 if (*itemType & kLTBranchType) {
2184 snprintf(msg2,2000, "Branch : %s", vname);
2185 } else {
2186 snprintf(msg2,2000, "Leaf : %s", vname);
2187 }
2188 }
2189 }
2190 }
2191 }
2192 // write who is responsable for this
2193 TString message = msg2;
2194 message = message(0,150);
2195 Message(msg2);
2196 delete[] msg2;
2197 // check if this should be pasted into the expression editor
2198 if ((*itemType & kLTBranchType) || (*itemType & kLTCutType)) break;
2200 if (!fDialogBox || !vname[0]) break;
2201 if (item == fDialogBox->EditedEntry()) break;
2202 // paste it
2203// char first = (char) vname[0];
2204 TString insert(item->GetAlias());
2205// if (first != '(') insert += "(";
2206// insert += item->GetAlias();
2207// if (first != '(') insert += ")";
2208
2210 fDialogBox->InsertText(insert.Data());
2211 // put the cursor at the right position
2212 }
2213 }
2214 break;
2215 case kButton2:
2216 break;
2217 case kButton3:
2218 // activate general context menu
2219 if (fLVContainer->NumSelected()) {
2220 void *p = 0;
2221 Int_t x = (Int_t)(parm2 &0xffff);
2222 Int_t y = (Int_t)((parm2 >> 16) & 0xffff);
2223 TTVLVEntry *item = 0;
2224 if ((item = (TTVLVEntry *) fLVContainer->GetNextSelected(&p)) != 0) {
2225 fContextMenu->Popup(x, y, item->GetContext());
2226 }
2227 } else { // empty click
2228 Int_t x = (Int_t)(parm2 &0xffff);
2229 Int_t y = (Int_t)((parm2 >> 16) & 0xffff);
2230 fContextMenu->Popup(x, y, this);
2231 }
2232 break;
2233 default:
2234 break;
2235 }
2236 break;
2237 case kCT_ITEMDBLCLICK:
2238 switch (parm1) {
2239 case kButton1:
2240 if (fLVContainer->NumSelected()) {
2241 // get item that sent this
2242 void *p = 0;
2243 TTVLVEntry *item;
2244 if ((item = (TTVLVEntry *) fLVContainer->GetNextSelected(&p)) != 0) {
2245 // get item type
2246 ULong_t *itemType = (ULong_t *) item->GetUserData();
2247 if (!(*itemType & kLTCutType) && !(*itemType & kLTBranchType)
2248 && !(*itemType & kLTPackType)) {
2249 if (strlen(item->GetTrueName())) {
2250 fVarDraw = kTRUE;
2251 // draw on double-click
2252 ExecuteDraw();
2253 break;
2254 } else {
2255 // open expression in editor
2257 }
2258 }
2259 if (*itemType & kLTCutType) {
2261 if (fEnableCut) {
2262 item->SetSmallPic(gClient->GetPicture("cut_t.xpm"));
2263 } else {
2264 item->SetSmallPic(gClient->GetPicture("cut-disable_t.xpm"));
2265 }
2266 }
2267 if (*itemType & kLTPackType) {
2268 fScanMode = kTRUE;
2269 ExecuteDraw();
2270 }
2271 }
2272 }
2273 break;
2274 case kButton2:
2275 break;
2276 case kButton3:
2277 break;
2278 default:
2279 break;
2280 }
2281 break;
2282 case 4:
2283// std::cout << "Dragging Item" << std::endl;
2284 default:
2285 break;
2286 }
2287 break;
2288 default:
2289 break;
2290 }
2291 return kTRUE;
2292}
2293
2294////////////////////////////////////////////////////////////////////////////////
2295/// Close the viewer.
2296
2298{
2299 DeleteWindow();
2300}
2301
2302////////////////////////////////////////////////////////////////////////////////
2303/// Execute all user commands.
2304
2305void TTreeViewer::ExecuteCommand(const char* command, Bool_t fast)
2306{
2307 // Execute the command, write it to history file and echo it to output
2308 if (fBarRec->GetState() == kButtonDown) {
2309 // show the command on the command line
2310 //printf("%s\n", command);
2311 char comm[2000];
2312 comm[0] = 0;
2313 if (strlen(command) > 1999) {
2314 Warning("ExecuteCommand", "Command too long: aborting.");
2315 return;
2316 }
2317 snprintf(comm,2000, "%s", command);
2318 // print the command to history file
2319 Gl_histadd(comm);
2320 }
2321 // execute it
2322 if (fast) {
2323 gROOT->ProcessLineFast(command);
2324 } else {
2325 gROOT->ProcessLine(command);
2326 }
2327 // make sure that 'draw on double-click' flag is reset
2328 fVarDraw = kFALSE;
2329}
2330
2331////////////////////////////////////////////////////////////////////////////////
2332/// Scan the selected options from option menu.
2333
2335{
2336 Int_t ind;
2337 if (parm1 == kOptionsReset) {
2338 for (ind=kOptionsGeneral; ind<kOptionsGeneral+16; ind++)
2340 for (ind=kOptions1D; ind<kOptions1D+12; ind++)
2342 for (ind=kOptions2D; ind<kOptions2D+14; ind++)
2344 }
2345 if ((parm1 < kOptions1D) && (parm1 != kOptionsReset)) {
2346 if (fOptionsGen->IsEntryChecked((Int_t)parm1)) {
2348 } else {
2349 fOptionsGen->CheckEntry((Int_t)parm1);
2351 }
2353 // uncheck all in this menu
2354 for (ind=kOptionsGeneral+1; ind<kOptionsGeneral+16; ind++) {
2356 }
2357 }
2358 }
2359
2360 if ((parm1 < kOptions2D) && (parm1 >= kOptions1D)) {
2361 if (fOptions1D->IsEntryChecked((Int_t)parm1)) {
2362 fOptions1D->UnCheckEntry((Int_t)parm1);
2363 } else {
2364 fOptions1D->CheckEntry((Int_t)parm1);
2366 }
2368 // uncheck all in this menu
2369 for (ind=kOptions1D+1; ind<kOptions1D+12; ind++) {
2371 }
2372 }
2373 }
2374
2375 if (parm1 >= kOptions2D) {
2376 if (fOptions2D->IsEntryChecked((Int_t)parm1)) {
2377 fOptions2D->UnCheckEntry((Int_t)parm1);
2378 } else {
2379 fOptions2D->CheckEntry((Int_t)parm1);
2381 }
2383 // uncheck all in this menu
2384 for (ind=kOptions2D+1; ind<kOptions2D+14; ind++) {
2386 }
2387 }
2388 }
2389 // concatenate options
2390 fBarOption->SetText("");
2391 for (ind=kOptionsGeneral; ind<kOptionsGeneral+16; ind++) {
2392 if (fOptionsGen->IsEntryChecked(ind))
2394 }
2395 if (Dimension() == 1) {
2396 for (ind=kOptions1D; ind<kOptions1D+12; ind++) {
2397 if (fOptions1D->IsEntryChecked(ind))
2399 }
2400 }
2401 if (Dimension() == 2) {
2402 for (ind=kOptions2D; ind<kOptions2D+14; ind++) {
2403 if (fOptions2D->IsEntryChecked(ind))
2405 }
2406 }
2407}
2408
2409////////////////////////////////////////////////////////////////////////////////
2410/// Map current tree and expand its content (including friends) in the lists.
2411
2412void TTreeViewer::MapTree(TTree *tree, TGListTreeItem *parent, Bool_t listIt)
2413{
2414 if (!tree) return;
2415 TObjArray *branches = tree->GetListOfBranches();
2416 if (!branches) return; // A Chain with no underlying trees.
2417 TBranch *branch;
2418 // loop on branches
2419 Int_t id;
2420 for (id=0; id<branches->GetEntries(); id++) {
2421 branch = (TBranch *)branches->At(id);
2422 if (branch->TestBit(kDoNotProcess)) continue;
2423 TString name = branch->GetName();
2424 if (name.Contains("fBits") || name.Contains("fUniqueID")) continue;
2425 // now map sub-branches
2426 MapBranch(branch, "", parent, listIt);
2428 }
2429 //Map branches of friend Trees (if any)
2430 //Look at tree->GetTree() to insure we see both the friendss of a chain
2431 //and the friends of the chain members
2432 TIter nextf( tree->GetTree()->GetListOfFriends() );
2433 TFriendElement *fr;
2434 while ((fr = (TFriendElement*)nextf())) {
2435 TTree * t = fr->GetTree();
2436 branches = t->GetListOfBranches();
2437 for (id=0; id<branches->GetEntries(); id++) {
2438 branch = (TBranch *)branches->At(id);
2439 if (branch->TestBit(kDoNotProcess)) continue;
2440 TString name = branch->GetName();
2441 if (name.Contains("fBits") || name.Contains("fUniqueID")) continue;
2442 // now map sub-branches
2443 MapBranch(branch, fr->GetName(), parent, listIt);
2445 }
2446 }
2447
2448 // tell who was last mapped
2449 if (listIt) {
2450 fMappedTree = tree;
2451 fMappedBranch = 0;
2452 }
2453}
2454
2455////////////////////////////////////////////////////////////////////////////////
2456/// Map current branch and expand its content in the list view.
2457
2458void TTreeViewer::MapBranch(TBranch *branch, const char *prefix, TGListTreeItem *parent, Bool_t listIt)
2459{
2460 if (!branch) return;
2461 TString name;
2462 if (prefix && strlen(prefix) > 0) {
2463 name = prefix;
2464 if (!name.EndsWith(".")) name += ".";
2465 name += branch->GetName();
2466 }
2467 else name = branch->GetName();
2468 Int_t ind;
2469 TGListTreeItem *branchItem = 0;
2470 ULong_t *itemType;
2471 // map this branch
2472 if (name.Contains("fBits") || name.Contains("fUniqueID")) return;
2473 if (parent) {
2474 // make list tree items for each branch according to the type
2475 const TGPicture *pic, *spic;
2476 if ((branch->GetListOfBranches()->GetEntries()) ||
2477 (branch->GetNleaves())) {
2478 if (branch->GetListOfBranches()->GetEntries()) {
2479 itemType = new ULong_t(kLTBranchType);
2480 if (branch->InheritsFrom("TBranchObject")) {
2481 pic = gClient->GetPicture("branch-ob_t.xpm");
2482 spic = gClient->GetPicture("branch-ob_t.xpm");
2483 } else {
2484 if (branch->InheritsFrom("TBranchClones")) {
2485 pic = gClient->GetPicture("branch-cl_t.xpm");
2486 spic = gClient->GetPicture("branch-cl_t.xpm");
2487 } else {
2488 pic = gClient->GetPicture("branch_t.xpm");
2489 spic = gClient->GetPicture("branch_t.xpm");
2490 }
2491 }
2492 branchItem = fLt->AddItem(parent, EmptyBrackets(name), itemType, pic, spic);
2493 } else {
2494 if (branch->GetNleaves() > 1) {
2495 itemType = new ULong_t(kLTBranchType);
2496 pic = gClient->GetPicture("branch_t.xpm");
2497 spic = gClient->GetPicture("branch_t.xpm");
2498 branchItem = fLt->AddItem(parent, EmptyBrackets(name), itemType,pic, spic);
2499 TObjArray *leaves = branch->GetListOfLeaves();
2500 TLeaf *leaf = 0;
2501 TString leafName;
2502 for (Int_t lf=0; lf<leaves->GetEntries(); lf++) {
2503 leaf = (TLeaf *)leaves->At(lf);
2504 leafName = name;
2505 if (!leafName.EndsWith(".")) leafName.Append(".");
2506 leafName.Append(EmptyBrackets(leaf->GetName()));
2507 itemType = new ULong_t(kLTLeafType);
2508 pic = gClient->GetPicture("leaf_t.xpm");
2509 spic = gClient->GetPicture("leaf_t.xpm");
2510 fLt->AddItem(branchItem, leafName.Data(), itemType, pic, spic);
2511 }
2512 } else {
2513 itemType = new ULong_t(kLTLeafType);
2514 pic = gClient->GetPicture("leaf_t.xpm");
2515 spic = gClient->GetPicture("leaf_t.xpm");
2516 branchItem = fLt->AddItem(parent, EmptyBrackets(name), itemType, pic, spic);
2517 }
2518 }
2519 }
2520 }
2521 // list branch in list view if necessary
2522 if (listIt) {
2523 TGString *textEntry = 0;
2524 const TGPicture *pic, *spic;
2525 TTVLVEntry *entry;
2526 // make list view items in the right frame
2527 if (!fStopMapping) {
2528 fMappedBranch = branch;
2529 fMappedTree = 0;
2531 }
2532 if ((branch->GetListOfBranches()->GetEntries()) ||
2533 (branch->GetNleaves())) {
2534 textEntry = new TGString(EmptyBrackets(name.Data()));
2535 if (branch->GetListOfBranches()->GetEntries()) {
2536 if (branch->InheritsFrom("TBranchObject")) {
2537 pic = gClient->GetPicture("branch-ob_t.xpm");
2538 spic = gClient->GetPicture("branch-ob_t.xpm");
2539 } else {
2540 if (branch->InheritsFrom("TBranchClones")) {
2541 pic = gClient->GetPicture("branch-cl_t.xpm");
2542 spic = gClient->GetPicture("branch-cl_t.xpm");
2543 } else {
2544 pic = gClient->GetPicture("branch_t.xpm");
2545 spic = gClient->GetPicture("branch_t.xpm");
2546 }
2547 }
2548 entry = new TTVLVEntry(fLVContainer,pic,spic,textEntry,0,kLVSmallIcons);
2549 entry->SetUserData(new UInt_t(kLTBranchType));
2550 entry->SetToolTipText("Branch with sub-branches. Can not be dragged");
2551 fLVContainer->AddThisItem(entry);
2552 entry->MapWindow();
2553 entry->SetAlias(textEntry->GetString());
2554 } else {
2555 if (branch->GetNleaves() > 1) {
2556 if (textEntry) delete textEntry;
2557 textEntry = new TGString(EmptyBrackets(name.Data()));
2558 pic = gClient->GetPicture("branch_t.xpm");
2559 spic = gClient->GetPicture("branch_t.xpm");
2560 entry = new TTVLVEntry(fLVContainer, pic, spic, textEntry,0,kLVSmallIcons);
2561 entry->SetUserData(new UInt_t(kLTBranchType));
2562 entry->SetToolTipText("Branch with more than one leaf. Can not be dragged");
2563 fLVContainer->AddThisItem(entry);
2564 entry->MapWindow();
2565 entry->SetAlias(textEntry->GetString());
2566
2567 TObjArray *leaves = branch->GetListOfLeaves();
2568 TLeaf *leaf = 0;
2569 TString leafName;
2570 for (Int_t lf=0; lf<leaves->GetEntries(); lf++) {
2571 leaf = (TLeaf *)leaves->At(lf);
2572 leafName = name;
2573 if (!leafName.EndsWith(".")) leafName.Append(".");
2574 leafName.Append(EmptyBrackets(leaf->GetName()));
2575 textEntry = new TGString(leafName.Data());
2576 pic = gClient->GetPicture("leaf_t.xpm");
2577 spic = gClient->GetPicture("leaf_t.xpm");
2578 entry = new TTVLVEntry(fLVContainer, pic, spic, textEntry,0,kLVSmallIcons);
2580 entry->SetToolTipText("Double-click to draw. Drag to X, Y, Z or scan box.");
2581 fLVContainer->AddThisItem(entry);
2582 entry->MapWindow();
2583 entry->SetAlias(textEntry->GetString());
2584 }
2585 } else {
2586 pic = (gClient->GetMimeTypeList())->GetIcon("TLeaf",kFALSE);
2587 if (!pic) pic = gClient->GetPicture("leaf_t.xpm");
2588 spic = gClient->GetMimeTypeList()->GetIcon("TLeaf",kTRUE);
2589 if (!spic) spic = gClient->GetPicture("leaf_t.xpm");
2590 entry = new TTVLVEntry(fLVContainer,pic,spic,textEntry,0,kLVSmallIcons);
2592 entry->SetToolTipText("Double-click to draw. Drag to X, Y, Z or scan box.");
2593 fLVContainer->AddThisItem(entry);
2594 entry->MapWindow();
2595 entry->SetAlias(textEntry->GetString());
2596 }
2597 }
2598 }
2599 }
2600
2601 TObjArray *branches = branch->GetListOfBranches();
2602 TBranch *branchDaughter = 0;
2603
2604 // loop all sub-branches
2605 for (ind=0; ind<branches->GetEntries(); ind++) {
2606 branchDaughter = (TBranch *)branches->UncheckedAt(ind);
2607 // map also all sub-branches
2608 MapBranch(branchDaughter, "", branchItem, listIt);
2609 }
2610}
2611
2612////////////////////////////////////////////////////////////////////////////////
2613/// Create new expression
2614
2616{
2618 const TGPicture *pic = gClient->GetPicture("expression_t.xpm");
2619 const TGPicture *spic = gClient->GetPicture("expression_t.xpm");
2620
2621 TTVLVEntry *entry = new TTVLVEntry(fLVContainer,pic,spic,
2622 new TGString(),0,kLVSmallIcons);
2624 fLVContainer->AddThisItem(entry);
2625 entry->MapWindow();
2626 entry->Empty();
2629 fListView->Layout();
2630 fNexpressions++;
2631}
2632
2633////////////////////////////////////////////////////////////////////////////////
2634/// Find parent tree of a clicked item.
2635
2637{
2638 if (!item) return;
2639 ULong_t *itemType = (ULong_t *)item->GetUserData();
2640 if (!itemType) return;
2641 TGListTreeItem *parent = 0;
2642 Int_t index;
2643 if (!(*itemType & kLTTreeType)) {
2644 parent = item->GetParent();
2645 SetParentTree(parent);
2646 } else {
2647 index = (Int_t)(*itemType >> 8);
2648 SwitchTree(index);
2649 }
2650}
2651
2652////////////////////////////////////////////////////////////////////////////////
2653/// Send a message on the status bar.
2654
2655void TTreeViewer::Message(const char* msg)
2656{
2657 fStatusBar->SetText(msg);
2658}
2659
2660////////////////////////////////////////////////////////////////////////////////
2661/// Put error/warning into TMsgBox and also forward to console.
2662
2663void TTreeViewer::DoError(int level, const char *location, const char *fmt, va_list va) const
2664{
2665 TObject::DoError(level, location, fmt, va);
2666
2667 // in case level will abort we will not come here...
2668
2669 static const int buf_size = 2048;
2670 char buf[buf_size], *bp;
2671
2672 int n = vsnprintf(buf, buf_size, fmt, va);
2673 // old vsnprintf's return -1 if string is truncated new ones return
2674 // total number of characters that would have been written
2675 if (n == -1 || n >= buf_size) {
2676 TObject::Warning("DoError", "Error message string truncated...");
2677 }
2678 if (level >= kSysError && level < kFatal)
2679 bp = Form("%s (%s)", buf, gSystem->GetError());
2680 else
2681 bp = buf;
2682
2683 const char *title = "";
2684 if (level == kInfo)
2685 title = "Info";
2686 if (level == kWarning)
2687 title = "Warning";
2688 if (level == kError)
2689 title = "Error";
2690 if (level == kSysError)
2691 title = "System Error";
2692
2693 new TGMsgBox(fClient->GetRoot(), this, title, bp, kMBIconExclamation);
2694}
2695
2696////////////////////////////////////////////////////////////////////////////////
2697/// Print the number of selected entries on status-bar.
2698
2700{
2701 if (!fTree) return;
2702 char * msg = new char[100];
2703 snprintf(msg,100, "First entry : %lld Last entry : %lld",
2705 Message(msg);
2706 delete[] msg;
2707}
2708
2709////////////////////////////////////////////////////////////////////////////////
2710/// Save current session as a C++ macro file.
2711
2712void TTreeViewer::SaveSource(const char* filename, Option_t *)
2713{
2714 if (!fTree) return;
2715 char quote = '"';
2716 std::ofstream out;
2717 Int_t lenfile = strlen(filename);
2718 char * fname;
2719 if (!lenfile) {
2720 fname = (char*)fSourceFile;
2721 lenfile = strlen(fname);
2722 } else {
2723 fname = (char*)filename;
2724 fSourceFile = filename;
2725 }
2726 // if filename is given, open this file, otherwise create a file
2727 // with a name : treeviewer.C
2728 if (lenfile) {
2729 out.open(fname, std::ios::out);
2730 } else {
2731 fname = new char[13];
2732 strlcpy(fname, "treeviewer.C",13);
2733 out.open(fname, std::ios::out);
2734 }
2735 if (!out.good ()) {
2736 printf("SaveSource cannot open file : %s\n", fname);
2737 fSourceFile = "treeviewer.C";
2738 if (!lenfile) delete [] fname;
2739 return;
2740 }
2741 // Write macro header and date/time stamp
2742 TDatime t;
2743 TString sname(fname);
2744 sname = sname.ReplaceAll(".C", "");
2745 out <<"void open_session(void *p = 0);"<<std::endl<<std::endl;
2746 out <<"void "<<sname.Data()<<"() {"<<std::endl;
2747 out <<"//=========Macro generated by ROOT version"<<gROOT->GetVersion()<<std::endl;
2748 out <<"//=========for tree "<<quote<<fTree->GetName()<<quote<<" ("<<t.AsString()<<")"<<std::endl;
2749 out <<"//===This macro can be opened from a TreeViewer session after loading"<<std::endl;
2750 out <<"//===the corresponding tree, or by running root with the macro name argument"<<std::endl<<std::endl;
2751 out <<" open_session();"<<std::endl;
2752 out <<"}"<<std::endl<<std::endl;
2753 out <<"void open_session(void *p = 0) {"<<std::endl;
2754 out <<" gSystem->Load("<<quote<<"libTreeViewer"<<quote<<");"<<std::endl;
2755 out <<" TTreeViewer *treeview = (TTreeViewer *) p;"<<std::endl;
2756 out <<" if (!treeview) treeview = new TTreeViewer();"<<std::endl;
2757 out <<" TTree *tv_tree = (TTree*)gROOT->FindObject("<<quote<<fTree->GetName()<<quote<<");"<<std::endl;
2758 out <<" TFile *tv_file = (TFile*)gROOT->GetListOfFiles()->FindObject("<<quote<<fFilename<<quote<<");"<<std::endl;
2759 out <<" if (!tv_tree) {"<<std::endl;
2760 out <<" if (!tv_file) tv_file = new TFile("<<quote<<fFilename<<quote<<");"<<std::endl;
2761 out <<" if (tv_file) tv_tree = (TTree*)tv_file->Get("<<quote<<fTree->GetName()<<quote<<");"<<std::endl;
2762 out <<" if(!tv_tree) {"<<std::endl;
2763 out <<" printf(\"Tree %s not found\", "<<quote<<fTree->GetName()<<quote<<");"<<std::endl;
2764 out <<" return;"<<std::endl;
2765 out <<" }"<<std::endl;
2766 out <<" }"<<std::endl<<std::endl;
2767 out <<" treeview->SetTreeName("<<quote<<fTree->GetName()<<quote<<");"<<std::endl;
2768 out <<" treeview->SetNexpressions("<<fNexpressions<<");"<<std::endl;
2769 // get expressions
2770 TTVLVEntry *item;
2771 out <<"// Set expressions on axis and cut"<<std::endl;
2772 out <<" TTVLVEntry *item;"<<std::endl;
2773 for (Int_t i=0; i<4; i++) {
2774 switch (i) {
2775 case 0:
2776 out <<"// X expression"<<std::endl;
2777 break;
2778 case 1:
2779 out <<"// Y expression"<<std::endl;
2780 break;
2781 case 2:
2782 out <<"// Z expression"<<std::endl;
2783 break;
2784 case 3:
2785 out <<"// Cut expression"<<std::endl;
2786 break;
2787 default:
2788 break;
2789 }
2790 item = ExpressionItem(i);
2791 out <<" item = treeview->ExpressionItem("<<i<<");"<<std::endl;
2792 out <<" item->SetExpression("<<quote<<item->GetTrueName()<<quote
2793 <<", "<<quote<<item->GetAlias()<<quote<<");"<<std::endl;
2794 }
2795 out <<"// Scan list"<<std::endl;
2796 item = ExpressionItem(4);
2797 out <<" item = treeview->ExpressionItem(4);"<<std::endl;
2798 out <<" item->SetExpression("<<quote<<item->GetTrueName()<<quote
2799 <<", "<<quote<<"Scan box"<<quote<<");"<<std::endl;
2800 out <<"// User defined expressions"<<std::endl;
2801 TString itemType;
2802 for (Int_t crt=5; crt<fNexpressions+5; crt++) {
2803 item = ExpressionItem(crt);
2804 if (item->IsCut())
2805 itemType = "kTRUE";
2806 else
2807 itemType = "kFALSE";
2808 out <<" item = treeview->ExpressionItem("<<crt<<");"<<std::endl;
2809 out <<" item->SetExpression("<<quote<<item->GetTrueName()<<quote
2810 <<", "<<quote<<item->GetAlias()<<quote<<", "<<itemType.Data()<<");"<<std::endl;
2811 }
2812 fSession->SaveSource(out);
2813 out <<"}"<<std::endl;
2814 out.close();
2815 printf("C++ Macro file: %s has been generated\n", fname);
2816 if (!lenfile) delete [] fname;
2817}
2818
2819////////////////////////////////////////////////////////////////////////////////
2820/// Makes current the tree at a given index in the list.
2821
2823{
2824 TTree *tree = (TTree *) fTreeList->At(index);
2825 if (!tree) {
2826 Warning("SwitchTree", "No tree found.");
2827 return kFALSE;
2828 }
2829 if ((tree == fTree) && (tree == fMappedTree)) return kFALSE; // nothing to switch
2830 std::string command;
2831 if (tree != fTree) {
2832 command = "tv__tree = (TTree *) tv__tree_list->At";
2833 command += Form("(%i)",index);
2834 ExecuteCommand(command.c_str());
2835 }
2836
2837 fTree = tree;
2840 command = "Current Tree : ";
2841 command += fTree->GetName();
2842 fLbl2->SetText(new TGString(command.c_str()));
2843 fTreeHdr->Layout();
2844 MapSubwindows();
2846 MapWindow();
2847 ///Resize(); //ia
2848 PrintEntries();
2849 return kTRUE;
2850}
2851
2852////////////////////////////////////////////////////////////////////////////////
2853/// Set record name
2854
2855void TTreeViewer::SetRecordName(const char *name)
2856{
2858}
2859
2860////////////////////////////////////////////////////////////////////////////////
2861/// Set current record
2862
2864{
2865 fCombo->Select(entry);
2866}
2867
2868////////////////////////////////////////////////////////////////////////////////
2869/// Set title of Histogram
2870
2871void TTreeViewer::SetHistogramTitle(const char *title)
2872{
2873 if (!gPad) return;
2874 TH1 *hist = (TH1*)gPad->GetListOfPrimitives()->FindObject(fBarHist->GetText());
2875 if (hist) {
2876 hist->SetTitle(title);
2877 gPad->Update();
2878 }
2879}
2880
2881////////////////////////////////////////////////////////////////////////////////
2882/// user defined command for current record
2883
2884void TTreeViewer::SetUserCode(const char *code, Bool_t autoexec)
2885{
2886 TTVRecord *rec = fSession->GetCurrent();
2887 if (rec) rec->SetUserCode(code, autoexec);
2888}
2889
2890////////////////////////////////////////////////////////////////////////////////
2891/// Updates combo box to current session entries.
2892
2894{
2895 TTVRecord *record;
2896 fCombo->RemoveEntries(0, 1000);
2897 for (Long64_t entry=0; entry<fSession->GetEntries(); entry++) {
2898 if ((record = fSession->GetRecord(entry)))
2899 fCombo->AddEntry(record->GetName(), entry);
2900 }
2901}
2902
2903////////////////////////////////////////////////////////////////////////////////
2904/// Updates current record to new X, Y, Z items.
2905
2906void TTreeViewer::UpdateRecord(const char *name)
2907{
2909}
2910
2911////////////////////////////////////////////////////////////////////////////////
2912/// This slot is called when button REFR is clicked
2913
2915{
2916 fTree->Refresh();
2918 Float_t max = (Float_t)fTree->GetEntries()-1;
2919 fSlider->SetRange(min,max);
2920 fSlider->SetPosition(min,max);
2921 ExecuteDraw();
2922}
@ kButtonPress
Definition: GuiTypes.h:59
EMouseButton
Definition: GuiTypes.h:213
@ kButton2
Definition: GuiTypes.h:213
@ kButton3
Definition: GuiTypes.h:213
@ kButton1
Definition: GuiTypes.h:213
R__EXTERN const char gTVHelpLayout[]
Definition: HelpTextTV.h:19
R__EXTERN const char gTVHelpStart[]
Definition: HelpTextTV.h:18
R__EXTERN const char gTVHelpDrawing[]
Definition: HelpTextTV.h:26
R__EXTERN const char gTVHelpOpenSave[]
Definition: HelpTextTV.h:20
R__EXTERN const char gTVHelpMacros[]
Definition: HelpTextTV.h:27
R__EXTERN const char gTVHelpEditExpressions[]
Definition: HelpTextTV.h:22
R__EXTERN const char gTVHelpAbout[]
Definition: HelpTextTV.h:17
R__EXTERN const char gTVHelpDraggingItems[]
Definition: HelpTextTV.h:21
R__EXTERN const char gTVHelpUserCommands[]
Definition: HelpTextTV.h:24
R__EXTERN const char gTVHelpSession[]
Definition: HelpTextTV.h:23
R__EXTERN const char gTVHelpContext[]
Definition: HelpTextTV.h:25
R__EXTERN const char gHelpAbout[]
Definition: HelpText.h:14
#define e(i)
Definition: RSha256.hxx:103
const Ssiz_t kNPOS
Definition: RtypesCore.h:111
int Int_t
Definition: RtypesCore.h:41
char Char_t
Definition: RtypesCore.h:29
unsigned int UInt_t
Definition: RtypesCore.h:42
const Bool_t kFALSE
Definition: RtypesCore.h:88
unsigned long ULong_t
Definition: RtypesCore.h:51
long Long_t
Definition: RtypesCore.h:50
bool Bool_t
Definition: RtypesCore.h:59
long long Long64_t
Definition: RtypesCore.h:69
float Float_t
Definition: RtypesCore.h:53
const Bool_t kTRUE
Definition: RtypesCore.h:87
const char Option_t
Definition: RtypesCore.h:62
#define ClassImp(name)
Definition: Rtypes.h:365
R__EXTERN TApplication * gApplication
Definition: TApplication.h:165
const Int_t kDoNotProcess
Definition: TBranch.h:53
#define gDirectory
Definition: TDirectory.h:218
const Int_t kError
Definition: TError.h:39
const Int_t kSysError
Definition: TError.h:41
const Int_t kFatal
Definition: TError.h:42
const Int_t kWarning
Definition: TError.h:38
const Int_t kInfo
Definition: TError.h:37
@ kButtonDown
Definition: TGButton.h:54
@ kButtonDisabled
Definition: TGButton.h:56
@ kButtonUp
Definition: TGButton.h:53
#define gClient
Definition: TGClient.h:166
@ kDoubleScaleBoth
@ kFDOpen
Definition: TGFileDialog.h:38
@ kFitHeight
Definition: TGFrame.h:66
@ kSunkenFrame
Definition: TGFrame.h:61
@ kVerticalFrame
Definition: TGFrame.h:59
@ kDoubleBorder
Definition: TGFrame.h:63
@ kFixedWidth
Definition: TGFrame.h:65
@ kHorizontalFrame
Definition: TGFrame.h:60
@ kLHintsRight
Definition: TGLayout.h:33
@ kLHintsExpandY
Definition: TGLayout.h:38
@ kLHintsLeft
Definition: TGLayout.h:31
@ kLHintsCenterY
Definition: TGLayout.h:35
@ kLHintsCenterX
Definition: TGLayout.h:32
@ kLHintsBottom
Definition: TGLayout.h:36
@ kLHintsTop
Definition: TGLayout.h:34
@ kLHintsExpandX
Definition: TGLayout.h:37
@ kLVSmallIcons
Definition: TGListView.h:41
@ kLVList
Definition: TGListView.h:42
@ kMBIconExclamation
Definition: TGMsgBox.h:35
XFontStruct * id
Definition: TGX11.cxx:108
char name[80]
Definition: TGX11.cxx:109
int nentries
Definition: THbookFile.cxx:89
#define gInterpreter
Definition: TInterpreter.h:553
#define gROOT
Definition: TROOT.h:414
char * Form(const char *fmt,...)
char * StrDup(const char *str)
Duplicate the string str.
Definition: TString.cxx:2490
R__EXTERN TSystem * gSystem
Definition: TSystem.h:560
static const char * gOpt2D[14]
EButtonIdentifiers
@ kBGFirst
@ kSLIDER
@ kSTOP
@ kBGLast
@ kRESET
@ kBGRecord
@ kCLOSE
@ kDRAW
@ kBGNext
@ kBGPrevious
static const char * gOptgen[16]
static const char * gOpenTypes[]
static const char * gOpt1D[12]
static const char * gMacroTypes[]
ERootTreeViewerCommands
@ kBarCommand
@ kHelpDragging
@ kFilePrint
@ kOptionsGeneral
@ kHelpContext
@ kFileClose
@ kHelpStart
@ kHelpCommands
@ kHelpMacros
@ kOptions2D
@ kFileQuit
@ kEditMacro
@ kHelpDrawing
@ kFileCanvas
@ kFileLoadLibrary
@ kFileSaveMacro
@ kFileBrowse
@ kFileOpenSession
@ kOptionsReset
@ kBarOption
@ kRunMacro
@ kAxis
@ kHelpEditing
@ kBarCut
@ kRunCommand
@ kHelpOpenSave
@ kHelpSession
@ kEditCut
@ kHelpAboutTV
@ kOptions1D
@ kEditEvent
@ kHelpLayout
@ kHelpAbout
@ kEditExpression
#define gPad
Definition: TVirtualPad.h:286
#define gVirtualX
Definition: TVirtualX.h:345
@ kWatch
Definition: TVirtualX.h:47
@ kPointer
Definition: TVirtualX.h:47
Int_t MK_MSG(EWidgetMessageTypes msg, EWidgetMessageTypes submsg)
Int_t GET_MSG(Long_t val)
@ kCT_SELCHANGED
@ kCM_COMBOBOX
@ kCM_MENU
@ kTE_ENTER
@ kCT_ITEMCLICK
@ kC_COMMAND
@ kCM_BUTTON
@ kC_TEXTENTRY
@ kC_LISTTREE
@ kC_VSLIDER
@ kCT_ITEMDBLCLICK
@ kC_CONTAINER
Int_t GET_SUBMSG(Long_t val)
#define snprintf
Definition: civetweb.c:1540
virtual void Terminate(Int_t status=0)
Terminate the application by call TSystem::Exit() unless application has been told to return from Run...
Class to manage histogram axis.
Definition: TAxis.h:30
A TTree is a list of TBranches.
Definition: TBranch.h:65
TObjArray * GetListOfBranches()
Definition: TBranch.h:214
Int_t GetNleaves() const
Definition: TBranch.h:217
TObjArray * GetListOfLeaves()
Definition: TBranch.h:215
This class provides an interface to context sensitive popup menus.
Definition: TContextMenu.h:40
virtual void Action(TObject *object, TMethod *method)
Action to be performed when this menu item is selected.
virtual void Popup(Int_t x, Int_t y, TObject *obj, TVirtualPad *c=0, TVirtualPad *p=0)
Popup context menu at given location in canvas c and pad p for selected object.
This class stores the date and time with a precision of one second in an unsigned 32 bit word (950130...
Definition: TDatime.h:37
const char * AsString() const
Return the date & time as a string (ctime() format).
Definition: TDatime.cxx:101
Describe directory structure in memory.
Definition: TDirectory.h:34
virtual Bool_t cd(const char *path=0)
Change current directory to "this" directory.
Definition: TDirectory.cxx:497
virtual TFile * GetFile() const
Definition: TDirectory.h:152
A TEventList object is a list of selected events (entries) in a TTree.
Definition: TEventList.h:31
A ROOT file is a suite of consecutive data records (TKey instances) with a well defined format.
Definition: TFile.h:48
A TFriendElement TF describes a TTree object TF in a file.
virtual TTree * GetTree()
Return pointer to friend TTree.
virtual void SetToolTipText(const char *text, Long_t delayms=400)
Set tool tip text associated with this button.
Definition: TGButton.cxx:395
virtual EButtonState GetState() const
Definition: TGButton.h:112
virtual void SetState(EButtonState state, Bool_t emit=kFALSE)
Set button state.
Definition: TGButton.cxx:185
virtual void SetContainer(TGFrame *f)
Definition: TGCanvas.h:232
TGViewPort * GetViewPort() const
Definition: TGCanvas.h:227
virtual void SetState(EButtonState state, Bool_t emit=kFALSE)
Set check button state.
Definition: TGButton.cxx:1200
const TGWindow * GetRoot() const
Returns current root (i.e.
Definition: TGClient.cxx:224
const TGPicture * GetPicture(const char *name)
Get picture from the picture pool.
Definition: TGClient.cxx:289
void NeedRedraw(TGWindow *w, Bool_t force=kFALSE)
Set redraw flags.
Definition: TGClient.cxx:372
virtual void AddEntry(TGString *s, Int_t id)
Definition: TGComboBox.h:106
virtual void RemoveEntries(Int_t from_ID, Int_t to_ID)
Definition: TGComboBox.h:125
virtual void Select(Int_t id, Bool_t emit=kTRUE)
Make the selected item visible in the combo box window and emit signals according to the second param...
Definition: TGComboBox.cxx:443
TGCompositeFrame(const TGCompositeFrame &)
virtual void AddFrame(TGFrame *f, TGLayoutHints *l=0)
Add frame to the composite frame using the specified layout hints.
Definition: TGFrame.cxx:1099
virtual UInt_t GetDefaultWidth() const
Definition: TGFrame.h:371
virtual void Layout()
Layout the elements of the composite frame.
Definition: TGFrame.cxx:1239
virtual TGDimension GetDefaultSize() const
std::cout << fWidth << "x" << fHeight << std::endl;
Definition: TGFrame.h:375
virtual void MapSubwindows()
Map all sub windows that are part of the composite frame.
Definition: TGFrame.cxx:1146
virtual UInt_t GetDefaultHeight() const
Definition: TGFrame.h:373
virtual void Associate(const TGWindow *w)
Definition: TGCanvas.h:99
virtual void RemoveItem(TGFrame *item)
Remove item from container.
Definition: TGCanvas.cxx:655
virtual const TGFrame * GetNextSelected(void **current)
Return the next selected item.
Definition: TGCanvas.cxx:676
virtual Int_t NumSelected() const
Definition: TGCanvas.h:114
virtual void RemoveAll()
Remove all items from the container.
Definition: TGCanvas.cxx:636
virtual Float_t GetMaxPosition() const
virtual Float_t GetMinPosition() const
virtual void SetRange(Float_t min, Float_t max)
virtual void SetPosition(Float_t min, Float_t max)
char * fFilename
Definition: TGFileDialog.h:61
const char ** fFileTypes
Definition: TGFileDialog.h:63
char * fIniDir
Definition: TGFileDialog.h:62
static Pixel_t GetWhitePixel()
Get white pixel value.
Definition: TGFrame.cxx:691
virtual UInt_t GetDefaultHeight() const
Definition: TGFrame.h:238
virtual void SetBackgroundColor(Pixel_t back)
Set background color (override from TGWindow base class).
Definition: TGFrame.cxx:294
virtual void DeleteWindow()
Delete window.
Definition: TGFrame.cxx:258
virtual void Resize(UInt_t w=0, UInt_t h=0)
Resize the frame.
Definition: TGFrame.cxx:587
virtual void SetWidth(UInt_t w)
Definition: TGFrame.h:293
virtual void MapWindow()
Definition: TGFrame.h:251
UInt_t GetWidth() const
Definition: TGFrame.h:271
virtual void SetHeight(UInt_t h)
Definition: TGFrame.h:294
void ShowPosition(Bool_t set=kTRUE, Bool_t percent=kTRUE, const char *format="%.2f")
Show postion text, either in percent or formatted according format.
void * GetUserData() const
Definition: TGListView.h:113
void SetUserData(void *userData)
Definition: TGListView.h:112
virtual void SetText(TGString *newText)
Set new text in label.
Definition: TGLabel.cxx:177
virtual const char * GetText() const =0
TGListTreeItem * GetParent() const
Definition: TGListTree.h:73
virtual void * GetUserData() const =0
void ClearHighlighted()
Un highlight items.
void AddItem(TGListTreeItem *parent, TGListTreeItem *item)
Add given item to list tree.
void OpenItem(TGListTreeItem *item)
Open item in list tree (i.e. show child items).
TGListTreeItem * GetSelected() const
Definition: TGListTree.h:397
TGListTreeItem * FindChildByName(TGListTreeItem *item, const char *name)
Find child of item by name.
void HighlightItem(TGListTreeItem *item)
Highlight item.
virtual void Layout()
Layout list view components (container and contents of container).
virtual void SetViewMode(EListViewMode viewMode)
Set list view mode.
virtual void SetContainer(TGFrame *f)
Set list view container.
virtual void SendCloseMessage()
Send close message to self.
Definition: TGFrame.cxx:1702
void SetWindowName(const char *name=0)
Set window name. This is typically done via the window manager.
Definition: TGFrame.cxx:1746
virtual void AddPopup(TGHotString *s, TGPopupMenu *menu, TGLayoutHints *l, TGPopupMenu *before=0)
Add popup menu to menu bar.
Definition: TGMenu.cxx:415
TGClient * fClient
Definition: TGObject.h:37
virtual void AddPopup(TGHotString *s, TGPopupMenu *popup, TGMenuEntry *before=0, const TGPicture *p=0)
Add a (cascading) popup menu to a popup menu.
Definition: TGMenu.cxx:1149
virtual Bool_t IsEntryChecked(Int_t id)
Return true if menu item is checked.
Definition: TGMenu.cxx:1835
virtual void AddEntry(TGHotString *s, Int_t id, void *ud=0, const TGPicture *p=0, TGMenuEntry *before=0)
Add a menu entry.
Definition: TGMenu.cxx:987
virtual void AddSeparator(TGMenuEntry *before=0)
Add a menu separator to the menu.
Definition: TGMenu.cxx:1057
virtual void CheckEntry(Int_t id)
Check a menu entry (i.e. add a check mark in front of it).
Definition: TGMenu.cxx:1772
virtual void DisableEntry(Int_t id)
Disable entry (disabled entries appear in a sunken relieve).
Definition: TGMenu.cxx:1714
virtual void UnCheckEntry(Int_t id)
Uncheck menu entry (i.e. remove check mark).
Definition: TGMenu.cxx:1797
virtual void Associate(const TGWindow *w)
Definition: TGMenu.h:219
void SetPosition(Float_t pos)
Set progress position between [min,max].
void SetFillType(EFillType type)
Set fill type.
virtual void SetBarColor(Pixel_t color)
Set progress bar color.
This class represent a specialized expression editor for TTVLVEntry 'true name' and 'alias' data memb...
TTVLVEntry * EditedEntry()
void InsertText(const char *text)
Insert text in text entry.
void GrabPointer()
Just focus the cursor inside.
void SetLabel(const char *title)
Set label of selection box.
void SetEntry(TTVLVEntry *entry)
Connect one entry.
static TGSelectBox * GetInstance()
Return the pointer to the instantiated singleton.
void Draw3DCorner(Bool_t corner)
Definition: TGStatusBar.h:67
virtual void SetText(TGString *text, Int_t partidx=0)
Set text in partition partidx in status bar.
const char * GetString() const
Definition: TGString.h:40
virtual void SetDefaultSize(UInt_t w, UInt_t h)
Set the default / minimal size of the widget.
const char * GetText() const
Definition: TGTextEntry.h:134
virtual void SetCursorPosition(Int_t pos)
Set the cursor position to newPos.
virtual void AppendText(const char *text)
Appends text to the end of text entry, clears the selection and moves the cursor to the end of the li...
virtual void SetToolTipText(const char *text, Long_t delayms=500)
Set tool tip text associated with this text entry.
virtual void SetText(const char *text, Bool_t emit=kTRUE)
Sets text entry to text, clears the selection and moves the cursor to the end of the line.
void Clear(Option_t *option="")
Clears up the text entry.
virtual Bool_t HandleButton(Event_t *event)
Handle mouse button event in text entry widget.
virtual void Associate(const TGWindow *w)
Definition: TGWidget.h:84
The TH1 histogram class.
Definition: TH1.h:56
virtual void SetTitle(const char *title)
See GetStatOverflows for more information.
Definition: TH1.cxx:6309
TAxis * GetZaxis()
Definition: TH1.h:318
TAxis * GetXaxis()
Get the behaviour adopted by the object about the statoverflows. See EStatOverflows for more informat...
Definition: TH1.h:316
TAxis * GetYaxis()
Definition: TH1.h:317
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2981
A TLeaf describes individual elements of a TBranch See TBranch structure in TTree.
Definition: TLeaf.h:49
A doubly linked list.
Definition: TList.h:44
virtual void Add(TObject *obj)
Definition: TList.h:87
virtual TObject * FindObject(const char *name) const
Find an object in this list using its name.
Definition: TList.cxx:575
virtual TObject * At(Int_t idx) const
Returns the object at position idx. Returns 0 if idx is out of range.
Definition: TList.cxx:354
virtual void Delete(Option_t *option="")
Remove all objects from the list AND delete all heap based objects.
Definition: TList.cxx:467
Each ROOT class (see TClass) has a linked list of methods.
Definition: TMethod.h:38
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
An array of TObjects.
Definition: TObjArray.h:37
Int_t GetEntries() const
Return the number of objects in array (i.e.
Definition: TObjArray.cxx:522
TObject * UncheckedAt(Int_t i) const
Definition: TObjArray.h:90
TObject * At(Int_t idx) const
Definition: TObjArray.h:166
Collectable string class.
Definition: TObjString.h:28
TString & String()
Definition: TObjString.h:48
Mother of all ROOT objects.
Definition: TObject.h:37
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition: TObject.h:172
virtual void DoError(int level, const char *location, const char *fmt, va_list va) const
Interface to ErrorHandler (protected).
Definition: TObject.cxx:841
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:866
virtual TObject * FindObject(const char *name) const
Must be redefined in derived classes.
Definition: TObject.cxx:321
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:443
The most important graphics class in the ROOT system.
Definition: TPad.h:29
TVirtualPad * cd(Int_t subpadnumber=0)
Set Current pad.
Definition: TPad.cxx:594
Bool_t Connect(const char *signal, const char *receiver_class, void *receiver, const char *slot)
Non-static method is used to connect from the signal of this object to the receiver slot.
Definition: TQObject.cxx:867
static const TString & GetBinDir()
Get the binary directory in the installation. Static utility function.
Definition: TROOT.cxx:2965
void SetText(const char *helpText)
Set help text from helpText buffer in TGTextView.
void Popup()
Show help dialog.
Sequenceable collection abstract base class.
Spider class.
Definition: TSpider.h:44
virtual void Draw(Option_t *options="")
Draw the spider.
Definition: TSpider.cxx:453
Basic string class.
Definition: TString.h:131
Bool_t EndsWith(const char *pat, ECaseCompare cmp=kExact) const
Return true if string ends with the specified string.
Definition: TString.cxx:2177
const char * Data() const
Definition: TString.h:364
TString & Remove(Ssiz_t pos)
Definition: TString.h:668
TString & Append(const char *cs)
Definition: TString.h:559
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:2311
virtual Int_t Exec(const char *shellcmd)
Execute a command.
Definition: TSystem.cxx:662
virtual const char * UnixPathName(const char *unixpathname)
Convert from a Unix pathname to a local pathname.
Definition: TSystem.cxx:1053
virtual const char * GetError()
Return system error string.
Definition: TSystem.cxx:259
This class represent the list view container for the.
void SelectItem(const char *name)
Select an item.
void SetViewer(TTreeViewer *viewer)
TList * ExpressionList()
Return the list of user-defined expressions.
const char * Cut()
Return the cut entry.
const char * Ez()
Return the expression on Z.
const char * Ey()
Return the expression on Y.
const char * ScanList()
Return the cut entry.
virtual void AddThisItem(TTVLVEntry *item)
TTVLVEntry * ExpressionItem(Int_t index)
Return the expression item at specific position.
void SetListView(TGListView *lv)
const char * Ex()
Return the expression on X.
void EmptyAll()
Clear all names and aliases for expression type items.
void RemoveNonStatic()
Remove all non-static items from the list view, except expressions.
This class represent entries that goes into the TreeViewer listview container.
TGItemContext * GetContext()
void SetSmallPic(const TGPicture *spic)
Set small picture.
void SetAlias(const char *alias)
void SetToolTipText(const char *text, Long_t delayms=1000)
Set tool tip text associated with this item.
void SetTrueName(const char *name)
const char * ConvertAliases()
Convert all aliases into true names.
Bool_t IsCut()
const char * GetAlias()
void Empty()
Clear all names and alias.
const char * GetTrueName()
I/O classes for TreeViewer session handling.
Definition: TTVSession.h:28
void SetUserCode(const char *code, Bool_t autoexec=kTRUE)
Definition: TTVSession.h:63
virtual const char * GetName() const
Returns name of object.
Definition: TTVSession.h:51
I/O classes for TreeViewer session handling.
Definition: TTVSession.h:69
TTVRecord * Last()
Definition: TTVSession.h:89
void SaveSource(std::ofstream &out)
Save the TTVSession in a C++ macro file.
Definition: TTVSession.cxx:253
TTVRecord * GetRecord(Int_t i)
Return record at index i.
Definition: TTVSession.cxx:186
void SetRecordName(const char *name)
Set record name.
Definition: TTVSession.cxx:209
void Show(TTVRecord *rec)
Display record rec.
Definition: TTVSession.cxx:242
TTVRecord * GetCurrent()
Definition: TTVSession.h:86
TTVRecord * Previous()
Definition: TTVSession.h:91
Int_t GetEntries()
Definition: TTVSession.h:85
TTVRecord * First()
Definition: TTVSession.h:88
void RemoveLastRecord()
Remove current record from list.
Definition: TTVSession.cxx:222
TTVRecord * Next()
Definition: TTVSession.h:90
TTVRecord * AddRecord(Bool_t fromFile=kFALSE)
Add a record.
Definition: TTVSession.cxx:159
void UpdateRecord(const char *name)
Updates current record according to new X, Y, Z settings.
Definition: TTVSession.cxx:270
Handles synchronous and a-synchronous timer events.
Definition: TTimer.h:51
virtual void TurnOff()
Remove timer from system timer list.
Definition: TTimer.cxx:229
virtual void TurnOn()
Add the timer to the system timer list.
Definition: TTimer.cxx:241
void Reset()
Reset the timer.
Definition: TTimer.cxx:157
Implement some of the functionality of the class TTree requiring access to extra libraries (Histogram...
Definition: TTreePlayer.h:37
A graphic user interface designed to handle ROOT trees and to take advantage of TTree class features.
Definition: TTreeViewer.h:56
TGLabel * fBarLbl2
Definition: TTreeViewer.h:118
void SetUserCode(const char *code, Bool_t autoexec=kTRUE)
TGSelectBox * fDialogBox
Definition: TTreeViewer.h:87
virtual ~TTreeViewer()
void SetCurrentRecord(Long64_t entry)
TGCheckButton * fBarRec
Definition: TTreeViewer.h:122
const TGPicture * fPicZ
Definition: TTreeViewer.h:90
TGHProgressBar * fProgressBar
Definition: TTreeViewer.h:137
void SetHistogramTitle(const char *title)
TTreeViewer(const char *treeName=0)
TTreeViewer default constructor.
const char * Ez()
Bool_t HandleTimer(TTimer *timer)
Execute action in response of a timer timing out.
void UpdateRecord(const char *name="new name")
TGTextEntry * fBarOption
Definition: TTreeViewer.h:124
TGLabel * fBLbl4
Definition: TTreeViewer.h:138
TGTextEntry * fBarCommand
Definition: TTreeViewer.h:123
TGPopupMenu * fOptions2D
Definition: TTreeViewer.h:111
TGListView * fListView
Definition: TTreeViewer.h:158
TGComboBox * fCombo
Definition: TTreeViewer.h:147
void EditExpression()
TGListTree * fLt
Definition: TTreeViewer.h:156
TList * fWidgets
Definition: TTreeViewer.h:161
void SetNexpressions(Int_t expr)
Bool_t fEnableCut
Definition: TTreeViewer.h:98
TGCheckButton * fBarH
Definition: TTreeViewer.h:120
TGPictureButton * fBGFirst
Definition: TTreeViewer.h:148
Cursor_t fDefaultCursor
Definition: TTreeViewer.h:93
void MapOptions(Long_t parm1)
TGLabel * fBLbl5
Definition: TTreeViewer.h:139
Bool_t IsScanRedirected()
virtual void CloseWindow()
Close and delete main frame.
TGCanvas * fTreeView
Definition: TTreeViewer.h:155
Int_t fTreeIndex
Definition: TTreeViewer.h:89
void NewExpression()
TGPictureButton * fSTOP
Definition: TTreeViewer.h:144
void SaveSource(const char *filename="", Option_t *option="")
Save the GUI main frame widget in a C++ macro file.
void SetSession(TTVSession *session)
void SetParentTree(TGListTreeItem *item)
const char * EmptyBrackets(const char *name)
TGMenuBar * fMenuBar
Definition: TTreeViewer.h:104
TGLabel * fLbl2
Definition: TTreeViewer.h:134
TTVLVEntry * ExpressionItem(Int_t index)
Cursor_t fWatchCursor
Definition: TTreeViewer.h:94
void AppendTree(TTree *tree)
Bool_t fScanMode
Definition: TTreeViewer.h:85
TGPopupMenu * fRunMenu
Definition: TTreeViewer.h:107
TGTextEntry * fBarListIn
Definition: TTreeViewer.h:140
void SetScanRedirect(Bool_t mode)
TTimer * fTimer
Definition: TTreeViewer.h:95
Long64_t Process(const char *filename, Option_t *option="", Long64_t nentries=TTree::kMaxEntries, Long64_t firstentry=0)
Bool_t SwitchTree(Int_t index)
TGCheckButton * fBarScan
Definition: TTreeViewer.h:121
void ExecuteDraw()
const TGPicture * fPicRefr
Definition: TTreeViewer.h:92
const char * fSourceFile
Definition: TTreeViewer.h:79
TGPopupMenu * fOptions1D
Definition: TTreeViewer.h:110
const TGPicture * fPicX
Definition: TTreeViewer.h:90
void SetRecordName(const char *name)
TGPictureButton * fDRAW
Definition: TTreeViewer.h:142
const char * Cut()
TString fLastOption
Definition: TTreeViewer.h:80
TList * fTreeList
Definition: TTreeViewer.h:88
const char * fFilename
Definition: TTreeViewer.h:78
@ kLTExpressionType
Definition: TTreeViewer.h:71
void RemoveItem()
TGTextEntry * fBarHist
Definition: TTreeViewer.h:125
TGLayoutHints * fBarLayout
Definition: TTreeViewer.h:115
TContextMenu * fContextMenu
Definition: TTreeViewer.h:86
void EmptyAll()
TGHorizontalFrame * fHpb
Definition: TTreeViewer.h:136
TGPopupMenu * fHelpMenu
Definition: TTreeViewer.h:112
TGPictureButton * fBGNext
Definition: TTreeViewer.h:151
Int_t fNexpressions
Definition: TTreeViewer.h:99
void SetTreeName(const char *treeName)
TGVerticalFrame * fV2
Definition: TTreeViewer.h:130
void Empty()
const char * En(Int_t n)
TGLayoutHints * fMenuBarItemLayout
Definition: TTreeViewer.h:102
Bool_t fStopMapping
Definition: TTreeViewer.h:97
void Message(const char *msg)
void ExecuteCommand(const char *command, Bool_t fast=kFALSE)
void BuildInterface()
void UpdateCombo()
TTree * fTree
Definition: TTreeViewer.h:76
const char * GetGrOpt()
TGCompositeFrame * fTreeHdr
Definition: TTreeViewer.h:131
TGLabel * fBarLbl3
Definition: TTreeViewer.h:119
Int_t Dimension()
Int_t MakeSelector(const char *selector=0)
void SetTree(TTree *tree)
TGDoubleVSlider * fSlider
Definition: TTreeViewer.h:128
TTree * fMappedTree
Definition: TTreeViewer.h:81
TGVerticalFrame * fV1
Definition: TTreeViewer.h:129
const TGPicture * fPicY
Definition: TTreeViewer.h:90
void DoRefresh()
const TGPicture * fPicDraw
Definition: TTreeViewer.h:91
TGPictureButton * fBGLast
Definition: TTreeViewer.h:152
const char * Ex()
TList * ExpressionList()
const TGPicture * fPicStop
Definition: TTreeViewer.h:91
void ActivateButtons(Bool_t first, Bool_t previous, Bool_t next, Bool_t last)
void MapTree(TTree *tree, TGListTreeItem *parent=0, Bool_t listIt=kTRUE)
void SetFile()
TGPopupMenu * fOptionsMenu
Definition: TTreeViewer.h:108
void RemoveLastRecord()
TGLayoutHints * fMenuBarHelpLayout
Definition: TTreeViewer.h:103
Bool_t fVarDraw
Definition: TTreeViewer.h:84
TGPictureButton * fBGPrevious
Definition: TTreeViewer.h:149
void PrintEntries()
TTVLVContainer * fLVContainer
Definition: TTreeViewer.h:159
TGLabel * fLbl1
Definition: TTreeViewer.h:133
void SetGrOpt(const char *option)
TGPopupMenu * fOptionsGen
Definition: TTreeViewer.h:109
TGLayoutHints * fMenuBarLayout
Definition: TTreeViewer.h:101
TGLabel * fBarLbl1
Definition: TTreeViewer.h:117
TGPictureButton * fREFR
Definition: TTreeViewer.h:145
TGHorizontalFrame * fHf
Definition: TTreeViewer.h:127
Bool_t fCounting
Definition: TTreeViewer.h:96
const char * Ey()
TGToolBar * fToolBar
Definition: TTreeViewer.h:114
TGHorizontalFrame * fBFrame
Definition: TTreeViewer.h:135
const char * ScanList()
TGPopupMenu * fEditMenu
Definition: TTreeViewer.h:106
TBranch * fMappedBranch
Definition: TTreeViewer.h:82
Bool_t ProcessMessage(Long_t msg, Long_t parm1, Long_t parm2)
TTVSession * fSession
Definition: TTreeViewer.h:77
Int_t fDimension
Definition: TTreeViewer.h:83
TGCompositeFrame * fListHdr
Definition: TTreeViewer.h:132
TGPopupMenu * fFileMenu
Definition: TTreeViewer.h:105
TGTextEntry * fBarListOut
Definition: TTreeViewer.h:141
TGTextButton * fReset
Definition: TTreeViewer.h:153
void ExecuteSpider()
TGPictureButton * fBGRecord
Definition: TTreeViewer.h:150
void MapBranch(TBranch *branch, const char *prefix="", TGListTreeItem *parent=0, Bool_t listIt=kTRUE)
void SetScanFileName(const char *name="")
TGStatusBar * fStatusBar
Definition: TTreeViewer.h:146
TGTextButton * fSPIDER
Definition: TTreeViewer.h:143
void DoError(int level, const char *location, const char *fmt, va_list va) const
Interface to ErrorHandler (protected).
A TTree represents a columnar dataset.
Definition: TTree.h:71
virtual void SetTimerInterval(Int_t msec=333)
Definition: TTree.h:581
virtual TBranch * GetBranch(const char *name)
Return pointer to the branch with the given name in this tree or its friends.
Definition: TTree.cxx:5075
virtual TObjArray * GetListOfLeaves()
Definition: TTree.h:428
TVirtualTreePlayer * GetPlayer()
Load the TTreePlayer (if not already done).
Definition: TTree.cxx:6083
virtual Long64_t Process(const char *filename, Option_t *option="", Long64_t nentries=kMaxEntries, Long64_t firstentry=0)
Process this tree executing the TSelector code in the specified filename.
Definition: TTree.cxx:7187
virtual void SetEventList(TEventList *list)
This function transfroms the given TEventList into a TEntryList The new TEntryList is owned by the TT...
Definition: TTree.cxx:8725
virtual Long64_t GetEntries() const
Definition: TTree.h:402
virtual Long64_t GetReadEntry() const
Definition: TTree.h:448
virtual TObjArray * GetListOfBranches()
Definition: TTree.h:427
TH1 * GetHistogram()
Definition: TTree.h:418
virtual Int_t MakeSelector(const char *selector=0, Option_t *option="")
Generate skeleton selector class for this tree.
Definition: TTree.cxx:6612
virtual void Refresh()
Refresh contents of this tree and its branches from the current status on disk.
Definition: TTree.cxx:7644
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
const Int_t n
Definition: legend1.C:16
UInt_t GetListOfMethods(TList &methods, TDirectory *dir=0)
Definition: tmvaglob.cxx:583
Short_t Abs(Short_t d)
Definition: TMathBase.h:120
Definition: file.py:1
Definition: first.py:1
Definition: tree.py:1
EGEventType fType
Definition: GuiTypes.h:174
REAL splitter
Definition: triangle.c:616