Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
TParallelCoord.cxx
Go to the documentation of this file.
1 // @(#)root/treeviewer:$Id$
2 // Author: Bastien Dalla Piazza 02/08/2007
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2007, 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 #include "TParallelCoord.h"
13 #include "TParallelCoordVar.h"
14 #include "TParallelCoordRange.h"
15 
16 #include "Riostream.h"
17 #include "TROOT.h"
18 #include "TVirtualX.h"
19 #include "TPad.h"
20 #include "TPolyLine.h"
21 #include "TGraph.h"
22 #include "TPaveText.h"
23 #include "float.h"
24 #include "TMath.h"
25 #include "TBox.h"
26 #include "TH1.h"
27 #include "TStyle.h"
28 #include "TEntryList.h"
29 #include "TFrame.h"
30 #include "TTree.h"
31 #include "TTreePlayer.h"
32 #include "TSelectorDraw.h"
33 #include "TTreeFormula.h"
34 #include "TView.h"
35 #include "TRandom.h"
36 #include "TEnv.h"
37 #include "TCanvas.h"
38 #include "TGaxis.h"
39 #include "TFile.h"
40 
41 ClassImp(TParallelCoord);
42 
43 /** \class TParallelCoord
44 Parallel Coordinates class.
45 
46 The multidimensional system of Parallel coordinates is a common way of studying
47 high-dimensional geometry and visualizing multivariate problems. It has first
48 been proposed by A. Inselberg in 1981.
49 
50 To show a set of points in an n-dimensional space, a backdrop is drawn
51 consisting of n parallel lines. A point in n-dimensional space is represented as
52 a polyline with vertices on the parallel axes; the position of the vertex on the
53 i-th axis corresponds to the i-th coordinate of the point.
54 
55 This tool comes with a rather large gui in the editor. It is necessary to use
56 this editor in order to explore a data set, as explained below.
57 
58 ### Reduce cluttering:
59 
60 The main issue for parallel coordinates is the very high cluttering of the
61 output when dealing with large data set. Two techniques have been implemented to
62 bypass that so far:
63 
64  - Draw doted lines instead of plain lines with an adjustable dots spacing. A
65  slider to adjust the dots spacing is available in the editor.
66  - Sort the entries to display with a "weight cut". On each axis is drawn a
67  histogram describing the distribution of the data on the corresponding
68  variable. The "weight" of an entry is the sum of the bin content of each bin
69  the entry is going through. An entry going through the histograms peaks will
70  have a big weight wether an entry going randomly through the histograms will
71  have a rather small weight. Setting a cut on this weight allows to draw only
72  the most representative entries. A slider set the cut is also available in
73  the gui.
74 
75 ## Selections:
76 
77 Selections of specific entries can be defined over the data se using parallel
78 coordinates. With that representation, a selection is an ensemble of ranges
79 defined on the axes. Ranges defined on the same axis are conjugated with OR
80 (an entry must be in one or the other ranges to be selected). Ranges on
81 different axes are are conjugated with AND (an entry must be in all the ranges
82 to be selected). Several selections can be defined with different colors. It is
83 possible to generate an entry list from a given selection and apply it to the
84 tree using the editor ("Apply to tree" button).
85 
86 ## Axes:
87 
88 Options can be defined each axis separately using the right mouse click. These
89 options can be applied to every axes using the editor.
90 
91  - Axis width: If set to 0, the axis is simply a line. If higher, a color
92  histogram is drawn on the axis.
93  - Axis histogram height: If not 0, a usual bar histogram is drawn on the plot.
94 
95 The order in which the variables are drawn is essential to see the clusters. The
96 axes can be dragged to change their position. A zoom is also available. The
97 logarithm scale is also available by right clicking on the axis.
98 
99 ## Candle chart:
100 
101 TParallelCoord can also be used to display a candle chart. In that mode, every
102 variable is drawn in the same scale. The candle chart can be combined with the
103 parallel coordinates mode, drawing the candle sticks over the axes.
104 
105 ~~~ {.cpp}
106 {
107  TCanvas *c1 = new TCanvas("c1");
108  TFile *f = TFile::Open("cernstaff.root");
109  TTree *T = (TTree*)f->Get("T");
110  T->Draw("Age:Grade:Step:Cost:Division:Nation","","para");
111  TParallelCoord* para = (TParallelCoord*)gPad->GetListOfPrimitives()->FindObject("ParaCoord");
112  TParallelCoordVar* grade = (TParallelCoordVar*)para->GetVarList()->FindObject("Grade");
113  grade->AddRange(new TParallelCoordRange(grade,11.5,14));
114  para->AddSelection("less30");
115  para->GetCurrentSelection()->SetLineColor(kViolet);
116  TParallelCoordVar* age = (TParallelCoordVar*)para->GetVarList()->FindObject("Age");
117  age->AddRange(new TParallelCoordRange(age,21,30));
118 }
119 ~~~
120 
121 ### Some references:
122 
123  - Alfred Inselberg's Homepage <http://www.math.tau.ac.il/~aiisreal>, with
124  Visual Tutorial, History, Selected Publications and Applications.
125  - Almir Olivette Artero, Maria Cristina Ferreira de Oliveira, Haim Levkowitz,
126  "Uncovering Clusters in Crowded Parallel Coordinates Visualizations,"
127  infovis, pp. 81-88, IEEE Symposium on Information Visualization
128  (INFOVIS'04), 2004.
129 */
130 
131 ////////////////////////////////////////////////////////////////////////////////
132 /// Default constructor.
133 
134 TParallelCoord::TParallelCoord()
135  :TNamed()
136 {
137  Init();
138 }
139 
140 ////////////////////////////////////////////////////////////////////////////////
141 /// Constructor without a reference to a tree,
142 /// the datas must be added afterwards with
143 /// TParallelCoord::AddVariable(Double_t*,const char*).
144 
145 TParallelCoord::TParallelCoord(Long64_t nentries)
146 {
147  Init();
148  fNentries = nentries;
149  fCurrentN = fNentries;
150  fVarList = new TList();
151  fSelectList = new TList();
152  fCurrentSelection = new TParallelCoordSelect();
153  fSelectList->Add(fCurrentSelection);
154 }
155 
156 ////////////////////////////////////////////////////////////////////////////////
157 /// Normal constructor, the datas must be added afterwards
158 /// with TParallelCoord::AddVariable(Double_t*,const char*).
159 
160 TParallelCoord::TParallelCoord(TTree* tree, Long64_t nentries)
161  :TNamed("ParaCoord","ParaCoord")
162 {
163  Init();
164  Int_t estimate = tree->GetEstimate();
165  if (nentries>estimate) {
166  Warning("TParallelCoord","Call tree->SetEstimate(tree->GetEntries()) to display all the tree variables");
167  fNentries = estimate;
168  } else {
169  fNentries = nentries;
170  }
171  fCurrentN = fNentries;
172  fTree = tree;
173  fTreeName = fTree->GetName();
174  if (fTree->GetCurrentFile()) fTreeFileName = fTree->GetCurrentFile()->GetName();
175  else fTreeFileName = "";
176  fVarList = new TList();
177  fSelectList = new TList();
178  fCurrentSelection = new TParallelCoordSelect();
179  fSelectList->Add(fCurrentSelection);
180 }
181 
182 ////////////////////////////////////////////////////////////////////////////////
183 /// Destructor.
184 
185 TParallelCoord::~TParallelCoord()
186 {
187  if (fInitEntries != fCurrentEntries && fCurrentEntries != 0) delete fCurrentEntries;
188  if (fVarList) {
189  fVarList->Delete();
190  delete fVarList;
191  }
192  if (fSelectList) {
193  fSelectList->Delete();
194  delete fSelectList;
195  }
196  if (fCandleAxis) delete fCandleAxis;
197  SetDotsSpacing(0);
198 }
199 
200 ////////////////////////////////////////////////////////////////////////////////
201 /// Add a variable.
202 
203 void TParallelCoord::AddVariable(Double_t* val, const char* title)
204 {
205  ++fNvar;
206  fVarList->Add(new TParallelCoordVar(val,title,fVarList->GetSize(),this));
207  SetAxesPosition();
208 }
209 
210 ////////////////////////////////////////////////////////////////////////////////
211 /// Add a variable from an expression.
212 
213 void TParallelCoord::AddVariable(const char* varexp)
214 {
215  if(!fTree) return; // The tree from which one will get the data must be defined.
216 
217  // Select in the only the entries of this TParallelCoord.
218  TEntryList *list = GetEntryList(kFALSE);
219  fTree->SetEntryList(list);
220 
221  // ensure that there is only one variable given:
222 
223  TString exp = varexp;
224 
225  if (exp.Contains(':') || exp.Contains(">>") || exp.Contains("<<")) {
226  Warning("AddVariable","Only a single variable can be added at a time.");
227  return;
228  }
229  if (exp == ""){
230  Warning("AddVariable","Nothing to add");
231  return;
232  }
233 
234  Long64_t en = fTree->Draw(varexp,"","goff");
235  if (en<0) {
236  Warning("AddVariable","%s could not be evaluated",varexp);
237  return;
238  }
239 
240  AddVariable(fTree->GetV1(),varexp);
241 }
242 
243 ////////////////////////////////////////////////////////////////////////////////
244 /// Add a selection.
245 
246 void TParallelCoord::AddSelection(const char* title)
247 {
248  TParallelCoordSelect *sel = new TParallelCoordSelect(title);
249  fSelectList->Add(sel);
250  fCurrentSelection = sel;
251 }
252 
253 ////////////////////////////////////////////////////////////////////////////////
254 /// Apply the current selection to the tree.
255 
256 void TParallelCoord::ApplySelectionToTree()
257 {
258  if(!fTree) return;
259  if(fSelectList) {
260  if(fSelectList->GetSize() == 0) return;
261  if(fCurrentSelection == 0) fCurrentSelection = (TParallelCoordSelect*)fSelectList->First();
262  }
263  fCurrentEntries = GetEntryList();
264  fNentries = fCurrentEntries->GetN();
265  fCurrentFirst = 0;
266  fCurrentN = fNentries;
267  fTree->SetEntryList(fCurrentEntries);
268  TString varexp = "";
269  TIter next(fVarList);
270  TParallelCoordVar* var;
271  while ((var = (TParallelCoordVar*)next())) varexp.Append(Form(":%s",var->GetTitle()));
272  varexp.Remove(TString::kLeading,':');
273  TSelectorDraw* selector = (TSelectorDraw*)((TTreePlayer*)fTree->GetPlayer())->GetSelector();
274  fTree->Draw(varexp.Data(),"","goff");
275  next.Reset();
276  Int_t i = 0;
277  while ((var = (TParallelCoordVar*)next())) {
278  var->SetValues(fNentries, selector->GetVal(i));
279  ++i;
280  }
281  if (fSelectList) { // FIXME It would be better to update the selections by deleting
282  fSelectList->Delete(); // the meaningless ranges (selecting everything or nothing for example)
283  fCurrentSelection = 0; // after applying a new entrylist to the tree.
284  }
285  gPad->Modified();
286  gPad->Update();
287 }
288 
289 ////////////////////////////////////////////////////////////////////////////////
290 /// Call constructor and add the variables.
291 
292 void TParallelCoord::BuildParallelCoord(TSelectorDraw* selector, Bool_t candle)
293 {
294  TParallelCoord* pc = new TParallelCoord(selector->GetTree(),selector->GetNfill());
295  pc->SetBit(kCanDelete);
296  selector->SetObject(pc);
297  TString varexp = "";
298  for(Int_t i=0;i<selector->GetDimension();++i) {
299  if (selector->GetVal(i)) {
300  if (selector->GetVar(i)) {
301  pc->AddVariable(selector->GetVal(i),selector->GetVar(i)->GetTitle());
302  varexp.Append(Form(":%s",selector->GetVar(i)->GetTitle()));
303  }
304  }
305  }
306  varexp.Remove(TString::kLeading,':');
307  if (selector->GetSelect()) varexp.Append(Form("{%s}",selector->GetSelect()->GetTitle()));
308  pc->SetTitle(varexp.Data());
309  if (!candle) pc->Draw();
310  else pc->Draw("candle");
311 }
312 
313 ////////////////////////////////////////////////////////////////////////////////
314 /// Clean up the selections from the ranges which could have been deleted
315 /// when a variable has been deleted.
316 
317 void TParallelCoord::CleanUpSelections(TParallelCoordRange* range)
318 {
319  TIter next(fSelectList);
320  TParallelCoordSelect* select;
321  while ((select = (TParallelCoordSelect*)next())){
322  if(select->Contains(range)) select->Remove(range);
323  }
324 }
325 
326 ////////////////////////////////////////////////////////////////////////////////
327 /// Delete a selection.
328 
329 void TParallelCoord::DeleteSelection(TParallelCoordSelect * sel)
330 {
331  fSelectList->Remove(sel);
332  delete sel;
333  if(fSelectList->GetSize() == 0) fCurrentSelection = 0;
334  else fCurrentSelection = (TParallelCoordSelect*)fSelectList->At(0);
335 }
336 
337 ////////////////////////////////////////////////////////////////////////////////
338 /// Compute the distance from the TParallelCoord.
339 
340 Int_t TParallelCoord::DistancetoPrimitive(Int_t px, Int_t py)
341 {
342  if(!gPad) return 9999;
343 
344  TFrame *frame = gPad->GetFrame();
345 
346  Double_t x1,x2,y1,y2,xx,yy;
347 
348  x1 = frame->GetX1()+0.01;
349  x2 = frame->GetX2()-0.01;
350  y2 = frame->GetY2()-0.01;
351  y1 = frame->GetY1()+0.01;
352 
353  xx = gPad->AbsPixeltoX(px);
354  yy = gPad->AbsPixeltoY(py);
355 
356  if(xx>x1 && xx<x2 && yy>y1 && yy<y2) return 0;
357  else return 9999;
358 }
359 
360 ////////////////////////////////////////////////////////////////////////////////
361 /// Draw the parallel coordinates graph.
362 
363 void TParallelCoord::Draw(Option_t* option)
364 {
365  if (!GetTree()) return;
366  if (!fCurrentEntries) fCurrentEntries = fInitEntries;
367  Bool_t optcandle = kFALSE;
368  TString opt = option;
369  opt.ToLower();
370  if(opt.Contains("candle")) {
371  optcandle = kTRUE;
372  opt.ReplaceAll("candle","");
373  }
374  if(optcandle) {
375  SetBit(kPaintEntries,kFALSE);
376  SetBit(kCandleChart,kTRUE);
377  SetGlobalScale(kTRUE);
378  }
379 
380  if (gPad) {
381  if (!gPad->IsEditable()) gROOT->MakeDefCanvas();
382  } else gROOT->MakeDefCanvas();
383  TView *view = gPad->GetView();
384  if(view){
385  delete view;
386  gPad->SetView(0);
387  }
388  gPad->Clear();
389  if (!optcandle) {
390  if (gPad && gPad->IsA() == TCanvas::Class()
391  && !((TCanvas*)gPad)->GetShowEditor()) {
392  ((TCanvas*)gPad)->ToggleEditor();
393  ((TCanvas*)gPad)->ToggleEventStatus();
394  }
395  }
396 
397  gPad->SetBit(TGraph::kClipFrame,kTRUE);
398 
399  TFrame *frame = new TFrame(0.1,0.1,0.9,0.9);
400  frame->SetBorderSize(0);
401  frame->SetBorderMode(0);
402  frame->SetFillStyle(0);
403  frame->SetLineColor(gPad->GetFillColor());
404  frame->Draw();
405  AppendPad(option);
406  TPaveText *title = new TPaveText(0.05,0.95,0.35,1);
407  title->AddText(GetTitle());
408  title->Draw();
409  SetAxesPosition();
410  TIter next(fVarList);
411  TParallelCoordVar* var;
412  while ((var = (TParallelCoordVar*)next())) {
413  if(optcandle) {
414  var->SetBoxPlot(kTRUE);
415  var->SetHistogramHeight(0.5);
416  var->SetHistogramLineWidth(0);
417  }
418  }
419 
420  if (optcandle) {
421  if (TestBit(kVertDisplay)) fCandleAxis = new TGaxis(0.05,0.1,0.05,0.9,GetGlobalMin(),GetGlobalMax());
422  else fCandleAxis = new TGaxis(0.1,0.05,0.9,0.05,GetGlobalMin(),GetGlobalMax());
423  fCandleAxis->Draw();
424  }
425 
426  if (gPad && gPad->IsA() == TCanvas::Class())
427  ((TCanvas*)gPad)->Selected(gPad,this,1);
428 }
429 
430 ////////////////////////////////////////////////////////////////////////////////
431 /// Execute the corresponding entry.
432 
433 void TParallelCoord::ExecuteEvent(Int_t /*entry*/, Int_t /*px*/, Int_t /*py*/)
434 {
435  if (!gPad) return;
436  gPad->SetCursor(kHand);
437 }
438 
439 ////////////////////////////////////////////////////////////////////////////////
440 /// Return the selection currently being edited.
441 
442 TParallelCoordSelect* TParallelCoord::GetCurrentSelection()
443 {
444  if (!fSelectList) return 0;
445  if (!fCurrentSelection) {
446  fCurrentSelection = (TParallelCoordSelect*)fSelectList->First();
447  }
448  return fCurrentSelection;
449 }
450 
451 ////////////////////////////////////////////////////////////////////////////////
452 /// Get the whole entry list or one for a selection.
453 
454 TEntryList* TParallelCoord::GetEntryList(Bool_t sel)
455 {
456  if(!sel || fCurrentSelection->GetSize() == 0){ // If no selection is specified, return the entry list of all the entries.
457  return fInitEntries;
458  } else { // return the entry list corresponding to the current selection.
459  TEntryList *enlist = new TEntryList(fTree);
460  TIter next(fVarList);
461  for (Long64_t li=0;li<fNentries;++li) {
462  next.Reset();
463  Bool_t inrange=kTRUE;
464  TParallelCoordVar* var;
465  while((var = (TParallelCoordVar*)next())){
466  if(!var->Eval(li,fCurrentSelection)) inrange = kFALSE;
467  }
468  if(!inrange) continue;
469  enlist->Enter(fCurrentEntries->GetEntry(li));
470  }
471  return enlist;
472  }
473 }
474 
475 ////////////////////////////////////////////////////////////////////////////////
476 /// return the global maximum.
477 
478 Double_t TParallelCoord::GetGlobalMax()
479 {
480  Double_t gmax=-DBL_MAX;
481  TIter next(fVarList);
482  TParallelCoordVar* var;
483  while ((var = (TParallelCoordVar*)next())) {
484  if (gmax < var->GetCurrentMax()) gmax = var->GetCurrentMax();
485  }
486  return gmax;
487 }
488 
489 ////////////////////////////////////////////////////////////////////////////////
490 /// return the global minimum.
491 
492 Double_t TParallelCoord::GetGlobalMin()
493 {
494  Double_t gmin=DBL_MAX;
495  TIter next(fVarList);
496  TParallelCoordVar* var;
497  while ((var = (TParallelCoordVar*)next())) {
498  if (gmin > var->GetCurrentMin()) gmin = var->GetCurrentMin();
499  }
500  return gmin;
501 }
502 
503 ////////////////////////////////////////////////////////////////////////////////
504 /// get the binning of the histograms.
505 
506 Int_t TParallelCoord::GetNbins()
507 {
508  return ((TParallelCoordVar*)fVarList->First())->GetNbins();
509 }
510 
511 ////////////////////////////////////////////////////////////////////////////////
512 /// Get a selection from its title.
513 
514 TParallelCoordSelect* TParallelCoord::GetSelection(const char* title)
515 {
516  TIter next(fSelectList);
517  TParallelCoordSelect* sel;
518  while ((sel = (TParallelCoordSelect*)next()) && strcmp(title,sel->GetTitle())) { }
519  return sel;
520 }
521 
522 ////////////////////////////////////////////////////////////////////////////////
523 /// return the tree if fTree is defined. If not, the method try to load the tree
524 /// from fTreeFileName.
525 
526 TTree* TParallelCoord::GetTree()
527 {
528  if (fTree) return fTree;
529  if (fTreeFileName=="" || fTreeName=="") {
530  Error("GetTree","Cannot load the tree: no tree defined!");
531  return 0;
532  }
533  TFile *f = TFile::Open(fTreeFileName.Data());
534  if (!f) {
535  Error("GetTree","Tree file name : \"%s\" does not exist (Are you in the correct directory?).",fTreeFileName.Data());
536  return 0;
537  } else if (f->IsZombie()) {
538  Error("GetTree","while opening \"%s\".",fTreeFileName.Data());
539  return 0;
540  } else {
541  fTree = (TTree*)f->Get(fTreeName.Data());
542  if (!fTree) {
543  Error("GetTree","\"%s\" not found in \"%s\".", fTreeName.Data(), fTreeFileName.Data());
544  return 0;
545  } else {
546  fTree->SetEntryList(fCurrentEntries);
547  TString varexp = "";
548  TIter next(fVarList);
549  TParallelCoordVar* var;
550  while ((var = (TParallelCoordVar*)next())) varexp.Append(Form(":%s",var->GetTitle()));
551  varexp.Remove(TString::kLeading,':');
552  fTree->Draw(varexp.Data(),"","goff");
553  TSelectorDraw* selector = (TSelectorDraw*)((TTreePlayer*)fTree->GetPlayer())->GetSelector();
554  next.Reset();
555  Int_t i = 0;
556  while ((var = (TParallelCoordVar*)next())) {
557  var->SetValues(fNentries, selector->GetVal(i));
558  ++i;
559  }
560  return fTree;
561  }
562  }
563 }
564 
565 ////////////////////////////////////////////////////////////////////////////////
566 /// Get the variables values from its title.
567 
568 Double_t* TParallelCoord::GetVariable(const char* vartitle)
569 {
570  TIter next(fVarList);
571  TParallelCoordVar* var = 0;
572  while(((var = (TParallelCoordVar*)next()) != 0) && (var->GetTitle() != vartitle)) { }
573  if(!var) return 0;
574  else return var->GetValues();
575 }
576 
577 ////////////////////////////////////////////////////////////////////////////////
578 /// Get the variables values from its index.
579 
580 Double_t* TParallelCoord::GetVariable(Int_t i)
581 {
582  if(i<0 || (UInt_t)i>fNvar) return 0;
583  else return ((TParallelCoordVar*)fVarList->At(i))->GetValues();
584 }
585 
586 ////////////////////////////////////////////////////////////////////////////////
587 /// Initialise the data members of TParallelCoord.
588 
589 void TParallelCoord::Init()
590 {
591  fNentries = 0;
592  fVarList = 0;
593  fSelectList = 0;
594  SetBit(kVertDisplay,kTRUE);
595  SetBit(kCurveDisplay,kFALSE);
596  SetBit(kPaintEntries,kTRUE);
597  SetBit(kLiveUpdate,kFALSE);
598  SetBit(kGlobalScale,kFALSE);
599  SetBit(kCandleChart,kFALSE);
600  SetBit(kGlobalLogScale,kFALSE);
601  fTree = 0;
602  fCurrentEntries = 0;
603  fInitEntries = 0;
604  fCurrentSelection = 0;
605  fNvar = 0;
606  fDotsSpacing = 0;
607  fCurrentFirst = 0;
608  fCurrentN = 0;
609  fCandleAxis = 0;
610  fWeightCut = 0;
611  fLineWidth = 1;
612  fLineColor = kGreen-8;
613  fTreeName = "";
614  fTreeFileName = "";
615 }
616 
617 ////////////////////////////////////////////////////////////////////////////////
618 /// Paint the parallel coordinates graph.
619 
620 void TParallelCoord::Paint(Option_t* /*option*/)
621 {
622  if (!GetTree()) return;
623  gPad->Range(0,0,1,1);
624  TFrame *frame = gPad->GetFrame();
625  frame->SetLineColor(gPad->GetFillColor());
626  SetAxesPosition();
627  if(TestBit(kPaintEntries)){
628  PaintEntries(0);
629  TIter next(fSelectList);
630  TParallelCoordSelect* sel;
631  while((sel = (TParallelCoordSelect*)next())) {
632  if(sel->GetSize()>0 && sel->TestBit(TParallelCoordSelect::kActivated)) {
633  PaintEntries(sel);
634  }
635  }
636  }
637  gPad->RangeAxis(0,0,1,1);
638 
639  TIter nextVar(fVarList);
640  TParallelCoordVar* var=0;
641  while((var = (TParallelCoordVar*)nextVar())) {
642  var->Paint();
643  }
644 }
645 
646 ////////////////////////////////////////////////////////////////////////////////
647 /// Loop over the entries and paint them.
648 
649 void TParallelCoord::PaintEntries(TParallelCoordSelect* sel)
650 {
651  if (fVarList->GetSize() < 2) return;
652  Int_t i=0;
653  Long64_t n=0;
654 
655  Double_t *x = new Double_t[fNvar];
656  Double_t *y = new Double_t[fNvar];
657 
658  TGraph *gr = 0;
659  TPolyLine *pl = 0;
660  TAttLine *evline = 0;
661 
662  if (TestBit (kCurveDisplay)) {gr = new TGraph(fNvar); evline = (TAttLine*)gr;}
663  else {pl = new TPolyLine(fNvar); evline = (TAttLine*)pl;}
664 
665  if (fDotsSpacing == 0) evline->SetLineStyle(1);
666  else evline->SetLineStyle(11);
667  if (!sel){
668  evline->SetLineWidth(GetLineWidth());
669  evline->SetLineColor(GetLineColor());
670  } else {
671  evline->SetLineWidth(sel->GetLineWidth());
672  evline->SetLineColor(sel->GetLineColor());
673  }
674  TParallelCoordVar *var;
675 
676  TFrame *frame = gPad->GetFrame();
677  Double_t lx = ((frame->GetX2() - frame->GetX1())/(fNvar-1));
678  Double_t ly = ((frame->GetY2() - frame->GetY1())/(fNvar-1));
679  Double_t a,b;
680  TRandom r;
681 
682  for (n=fCurrentFirst; n<fCurrentFirst+fCurrentN; ++n) {
683  TListIter next(fVarList);
684  Bool_t inrange = kTRUE;
685  // Loop to check whenever the entry must be painted.
686  if (sel) {
687  while ((var = (TParallelCoordVar*)next())){
688  if (!var->Eval(n,sel)) inrange = kFALSE;
689  }
690  }
691  if (fWeightCut > 0) {
692  next.Reset();
693  Int_t entryweight = 0;
694  while ((var = (TParallelCoordVar*)next())) entryweight+=var->GetEntryWeight(n);
695  if (entryweight/(Int_t)fNvar < fWeightCut) inrange = kFALSE;
696  }
697  if(!inrange) continue;
698  i = 0;
699  next.Reset();
700  // Loop to set the polyline points.
701  while ((var = (TParallelCoordVar*)next())) {
702  var->GetEntryXY(n,x[i],y[i]);
703  ++i;
704  }
705  // beginning to paint the first point at a random distance
706  // to avoid artefacts when increasing the dots spacing.
707  if (fDotsSpacing != 0) {
708  if (TestBit(kVertDisplay)) {
709  a = (y[1]-y[0])/(x[1]-x[0]);
710  b = y[0]-a*x[0];
711  x[0] = x[0]+lx*r.Rndm();
712  y[0] = a*x[0]+b;
713  } else {
714  a = (x[1]-x[0])/(y[1]-y[0]);
715  b = x[0]-a*y[0];
716  y[0] = y[0]+ly*r.Rndm();
717  x[0] = a*y[0]+b;
718  }
719  }
720  if (pl) pl->PaintPolyLine(fNvar,x,y);
721  else gr->PaintGraph(fNvar,x,y,"C");
722  }
723 
724  if (pl) delete pl;
725  if (gr) delete gr;
726  delete [] x;
727  delete [] y;
728 }
729 
730 ////////////////////////////////////////////////////////////////////////////////
731 /// Delete a variable from the graph.
732 
733 void TParallelCoord::RemoveVariable(TParallelCoordVar *var)
734 {
735  fVarList->Remove(var);
736  fNvar = fVarList->GetSize();
737  SetAxesPosition();
738 }
739 
740 ////////////////////////////////////////////////////////////////////////////////
741 /// Delete the variable "vartitle" from the graph.
742 
743 Bool_t TParallelCoord::RemoveVariable(const char* vartitle)
744 {
745  TIter next(fVarList);
746  TParallelCoordVar* var=0;
747  while((var = (TParallelCoordVar*)next())) {
748  if (!strcmp(var->GetTitle(),vartitle)) break;
749  }
750  if(!var) {
751  Error("RemoveVariable","\"%s\" not a variable",vartitle);
752  return kFALSE;
753  } else {
754  RemoveVariable(var);
755  delete var;
756  return kTRUE;
757  }
758 }
759 
760 ////////////////////////////////////////////////////////////////////////////////
761 /// Reset the tree entry list to the initial one..
762 
763 void TParallelCoord::ResetTree()
764 {
765  if(!fTree) return;
766  fTree->SetEntryList(fInitEntries);
767  fCurrentEntries = fInitEntries;
768  fNentries = fCurrentEntries->GetN();
769  fCurrentFirst = 0;
770  fCurrentN = fNentries;
771  TString varexp = "";
772  TIter next(fVarList);
773  TParallelCoordVar* var;
774  while ((var = (TParallelCoordVar*)next())) varexp.Append(Form(":%s",var->GetTitle()));
775  varexp.Remove(TString::kLeading,':');
776  fTree->Draw(varexp.Data(),"","goff");
777  next.Reset();
778  TSelectorDraw* selector = (TSelectorDraw*)((TTreePlayer*)fTree->GetPlayer())->GetSelector();
779  Int_t i = 0;
780  while ((var = (TParallelCoordVar*)next())) {
781  var->SetValues(fNentries, selector->GetVal(i));
782  ++i;
783  }
784  if (fSelectList) { // FIXME It would be better to update the selections by deleting
785  fSelectList->Delete(); // the meaningless ranges (selecting everything or nothing for example)
786  fCurrentSelection = 0; // after applying a new entrylist to the tree.
787  }
788  gPad->Modified();
789  gPad->Update();
790 }
791 
792 ////////////////////////////////////////////////////////////////////////////////
793 /// Save the entry lists in a root file "filename.root".
794 
795 void TParallelCoord::SaveEntryLists(const char* filename, Bool_t overwrite)
796 {
797  TString sfile = filename;
798  if (sfile == "") sfile = Form("%s_parallelcoord_entries.root",fTree->GetName());
799 
800  TFile* f = TFile::Open(sfile.Data());
801  if (f) {
802  Warning("SaveEntryLists","%s already exists.", sfile.Data());
803  if (!overwrite) return;
804  else Warning("SaveEntryLists","Overwriting.");
805  f = new TFile(sfile.Data(),"RECREATE");
806  } else {
807  f = new TFile(sfile.Data(),"CREATE");
808  }
809  gDirectory = f;
810  fInitEntries->Write("initentries");
811  fCurrentEntries->Write("currententries");
812  Info("SaveEntryLists","File \"%s\" written.",sfile.Data());
813 }
814 
815 ////////////////////////////////////////////////////////////////////////////////
816 /// Save the TParallelCoord in a macro.
817 
818 void TParallelCoord::SavePrimitive(std::ostream & out, Option_t* options)
819 {
820  TString opt = options;
821  opt.ToLower();
822  //Bool_t overwrite = opt.Contains("overwrite"); // Is there a way to specify "options" when saving ?
823  // Save the entrylists.
824  const char* filename = Form("%s_parallelcoord_entries.root",fTree->GetName());
825  SaveEntryLists(filename,kTRUE); // FIXME overwriting by default.
826  SaveTree(fTreeFileName,kTRUE); // FIXME overwriting by default.
827  out<<" // Create a TParallelCoord."<<std::endl;
828  out<<" TFile *f = TFile::Open(\""<<fTreeFileName.Data()<<"\");"<<std::endl;
829  out<<" TTree* tree = (TTree*)f->Get(\""<<fTreeName.Data()<<"\");"<<std::endl;
830  out<<" TParallelCoord* para = new TParallelCoord(tree,"<<fNentries<<");"<<std::endl;
831  out<<" // Load the entrylists."<<std::endl;
832  out<<" TFile *entries = TFile::Open(\""<<filename<<"\");"<<std::endl;
833  out<<" TEntryList *currententries = (TEntryList*)entries->Get(\"currententries\");"<<std::endl;
834  out<<" tree->SetEntryList(currententries);"<<std::endl;
835  out<<" para->SetInitEntries((TEntryList*)entries->Get(\"initentries\"));"<<std::endl;
836  out<<" para->SetCurrentEntries(currententries);"<<std::endl;
837  TIter next(fSelectList);
838  TParallelCoordSelect* sel;
839  out<<" TParallelCoordSelect* sel;"<<std::endl;
840  out<<" para->GetSelectList()->Delete();"<<std::endl;
841  while ((sel = (TParallelCoordSelect*)next())) {
842  out<<" para->AddSelection(\""<<sel->GetTitle()<<"\");"<<std::endl;
843  out<<" sel = (TParallelCoordSelect*)para->GetSelectList()->Last();"<<std::endl;
844  out<<" sel->SetLineColor("<<sel->GetLineColor()<<");"<<std::endl;
845  out<<" sel->SetLineWidth("<<sel->GetLineWidth()<<");"<<std::endl;
846  }
847  TIter nextbis(fVarList);
848  TParallelCoordVar* var;
849  TString varexp = "";
850  while ((var = (TParallelCoordVar*)nextbis())) varexp.Append(Form(":%s",var->GetTitle()));
851  varexp.Remove(TString::kLeading,':');
852  out<<" tree->Draw(\""<<varexp.Data()<<"\",\"\",\"goff\");"<<std::endl;
853  out<<" TSelectorDraw* selector = (TSelectorDraw*)((TTreePlayer*)tree->GetPlayer())->GetSelector();"<<std::endl;
854  nextbis.Reset();
855  Int_t i=0;
856  out<<" TParallelCoordVar* var;"<<std::endl;
857  while ((var = (TParallelCoordVar*)nextbis())) {
858  out<<" //***************************************"<<std::endl;
859  out<<" // Create the axis \""<<var->GetTitle()<<"\"."<<std::endl;
860  out<<" para->AddVariable(selector->GetVal("<<i<<"),\""<<var->GetTitle()<<"\");"<<std::endl;
861  out<<" var = (TParallelCoordVar*)para->GetVarList()->Last();"<<std::endl;
862  var->SavePrimitive(out,"pcalled");
863  ++i;
864  }
865  out<<" //***************************************"<<std::endl;
866  out<<" // Set the TParallelCoord parameters."<<std::endl;
867  out<<" para->SetCurrentFirst("<<fCurrentFirst<<");"<<std::endl;
868  out<<" para->SetCurrentN("<<fCurrentN<<");"<<std::endl;
869  out<<" para->SetWeightCut("<<fWeightCut<<");"<<std::endl;
870  out<<" para->SetDotsSpacing("<<fDotsSpacing<<");"<<std::endl;
871  out<<" para->SetLineColor("<<GetLineColor()<<");"<<std::endl;
872  out<<" para->SetLineWidth("<<GetLineWidth()<<");"<<std::endl;
873  out<<" para->SetBit(TParallelCoord::kVertDisplay,"<<TestBit(kVertDisplay)<<");"<<std::endl;
874  out<<" para->SetBit(TParallelCoord::kCurveDisplay,"<<TestBit(kCurveDisplay)<<");"<<std::endl;
875  out<<" para->SetBit(TParallelCoord::kPaintEntries,"<<TestBit(kPaintEntries)<<");"<<std::endl;
876  out<<" para->SetBit(TParallelCoord::kLiveUpdate,"<<TestBit(kLiveUpdate)<<");"<<std::endl;
877  out<<" para->SetBit(TParallelCoord::kGlobalLogScale,"<<TestBit(kGlobalLogScale)<<");"<<std::endl;
878  if (TestBit(kGlobalScale)) out<<" para->SetGlobalScale(kTRUE);"<<std::endl;
879  if (TestBit(kCandleChart)) out<<" para->SetCandleChart(kTRUE);"<<std::endl;
880  if (TestBit(kGlobalLogScale)) out<<" para->SetGlobalLogScale(kTRUE);"<<std::endl;
881  out<<std::endl<<" para->Draw();"<<std::endl;
882 }
883 
884 ////////////////////////////////////////////////////////////////////////////////
885 /// Save the tree in a file if fTreeFileName == "".
886 
887 void TParallelCoord::SaveTree(const char* filename, Bool_t overwrite)
888 {
889  if (!(fTreeFileName=="")) return;
890  TString sfile = filename;
891  if (sfile == "") sfile = Form("%s.root",fTree->GetName());
892 
893  TFile* f = TFile::Open(sfile.Data());
894  if (f) {
895  Warning("SaveTree","%s already exists.", sfile.Data());
896  if (!overwrite) return;
897  else Warning("SaveTree","Overwriting.");
898  f = new TFile(sfile.Data(),"RECREATE");
899  } else {
900  f = new TFile(sfile.Data(),"CREATE");
901  }
902  gDirectory = f;
903  fTree->Write(fTreeName.Data());
904  fTreeFileName = sfile;
905  Info("SaveTree","File \"%s\" written.",sfile.Data());
906 }
907 
908 ////////////////////////////////////////////////////////////////////////////////
909 /// Update the position of the axes.
910 
911 void TParallelCoord::SetAxesPosition()
912 {
913  if(!gPad) return;
914  Bool_t vert = TestBit (kVertDisplay);
915  TFrame *frame = gPad->GetFrame();
916  if (fVarList->GetSize() > 1) {
917  if (vert) {
918  frame->SetX1(1.0/((Double_t)fVarList->GetSize()+1));
919  frame->SetX2(1-frame->GetX1());
920  frame->SetY1(0.1);
921  frame->SetY2(0.9);
922  gPad->RangeAxis(1.0/((Double_t)fVarList->GetSize()+1),0.1,1-frame->GetX1(),0.9);
923  } else {
924  frame->SetX1(0.1);
925  frame->SetX2(0.9);
926  frame->SetY1(1.0/((Double_t)fVarList->GetSize()+1));
927  frame->SetY2(1-frame->GetY1());
928  gPad->RangeAxis(0.1,1.0/((Double_t)fVarList->GetSize()+1),0.9,1-frame->GetY1());
929  }
930 
931  Double_t horSpace = (frame->GetX2() - frame->GetX1())/(fNvar-1);
932  Double_t verSpace = (frame->GetY2() - frame->GetY1())/(fNvar-1);
933  Int_t i=0;
934  TIter next(fVarList);
935 
936  TParallelCoordVar* var;
937  while((var = (TParallelCoordVar*)next())){
938  if (vert) var->SetX(gPad->GetFrame()->GetX1() + i*horSpace,TestBit(kGlobalScale));
939  else var->SetY(gPad->GetFrame()->GetY1() + i*verSpace,TestBit(kGlobalScale));
940  ++i;
941  }
942  } else if (fVarList->GetSize()==1) {
943  frame->SetX1(0.1);
944  frame->SetX2(0.9);
945  frame->SetY1(0.1);
946  frame->SetY2(0.9);
947  if (vert) ((TParallelCoordVar*)fVarList->First())->SetX(0.5,TestBit(kGlobalScale));
948  else ((TParallelCoordVar*)fVarList->First())->SetY(0.5,TestBit(kGlobalScale));
949  }
950 }
951 
952 ////////////////////////////////////////////////////////////////////////////////
953 /// Set the same histogram axis binning for all axis.
954 
955 void TParallelCoord::SetAxisHistogramBinning(Int_t n)
956 {
957  TIter next(fVarList);
958  TParallelCoordVar *var;
959  while((var = (TParallelCoordVar*)next())) var->SetHistogramBinning(n);
960 }
961 
962 ////////////////////////////////////////////////////////////////////////////////
963 /// Set the same histogram axis height for all axis.
964 
965 void TParallelCoord::SetAxisHistogramHeight(Double_t h)
966 {
967  TIter next(fVarList);
968  TParallelCoordVar *var;
969  while((var = (TParallelCoordVar*)next())) var->SetHistogramHeight(h);
970 }
971 
972 ////////////////////////////////////////////////////////////////////////////////
973 /// All axes in log scale.
974 
975 void TParallelCoord::SetGlobalLogScale(Bool_t lt)
976 {
977  if (lt == TestBit(kGlobalLogScale)) return;
978  SetBit(kGlobalLogScale,lt);
979  TIter next(fVarList);
980  TParallelCoordVar* var;
981  while ((var = (TParallelCoordVar*)next())) var->SetLogScale(lt);
982  if (TestBit(kGlobalScale)) SetGlobalScale(kTRUE);
983 }
984 
985 ////////////////////////////////////////////////////////////////////////////////
986 /// Constraint all axes to the same scale.
987 
988 void TParallelCoord::SetGlobalScale(Bool_t gl)
989 {
990  SetBit(kGlobalScale,gl);
991  if (fCandleAxis) {
992  delete fCandleAxis;
993  fCandleAxis = 0;
994  }
995  if (gl) {
996  Double_t min,max;
997  min = GetGlobalMin();
998  max = GetGlobalMax();
999  if (TestBit(kGlobalLogScale) && min<=0) min = 0.00001*max;
1000  if (TestBit(kVertDisplay)) {
1001  if (!TestBit(kGlobalLogScale)) fCandleAxis = new TGaxis(0.05,0.1,0.05,0.9,min,max);
1002  else fCandleAxis = new TGaxis(0.05,0.1,0.05,0.9,min,max,510,"G");
1003  } else {
1004  if (!TestBit(kGlobalLogScale)) fCandleAxis = new TGaxis(0.1,0.05,0.9,0.05,min,max);
1005  else fCandleAxis = new TGaxis(0.1,0.05,0.9,0.05,min,max,510,"G");
1006  }
1007  fCandleAxis->Draw();
1008  SetGlobalMin(min);
1009  SetGlobalMax(max);
1010  TIter next(fVarList);
1011  TParallelCoordVar* var;
1012  while ((var = (TParallelCoordVar*)next())) var->GetHistogram();
1013  }
1014  gPad->Modified();
1015  gPad->Update();
1016 }
1017 
1018 ////////////////////////////////////////////////////////////////////////////////
1019 /// Set the same histogram axis line width for all axis.
1020 
1021 void TParallelCoord::SetAxisHistogramLineWidth(Int_t lw)
1022 {
1023  TIter next(fVarList);
1024  TParallelCoordVar *var;
1025  while((var = (TParallelCoordVar*)next())) var->SetHistogramLineWidth(lw);
1026 }
1027 
1028 ////////////////////////////////////////////////////////////////////////////////
1029 /// Set a candle chart display.
1030 
1031 void TParallelCoord::SetCandleChart(Bool_t can)
1032 {
1033  SetBit(kCandleChart,can);
1034  SetGlobalScale(can);
1035  TIter next(fVarList);
1036  TParallelCoordVar* var;
1037  while ((var = (TParallelCoordVar*)next())) {
1038  var->SetBoxPlot(can);
1039  var->SetHistogramLineWidth(0);
1040  }
1041  if (fCandleAxis) delete fCandleAxis;
1042  fCandleAxis = 0;
1043  SetBit(kPaintEntries,!can);
1044  if (can) {
1045  if (TestBit(kVertDisplay)) fCandleAxis = new TGaxis(0.05,0.1,0.05,0.9,GetGlobalMin(),GetGlobalMax());
1046  else fCandleAxis = new TGaxis(0.1,0.05,0.9,0.05,GetGlobalMin(),GetGlobalMax());
1047  fCandleAxis->Draw();
1048  } else {
1049  if (fCandleAxis) {
1050  delete fCandleAxis;
1051  fCandleAxis = 0;
1052  }
1053  }
1054  gPad->Modified();
1055  gPad->Update();
1056 }
1057 
1058 ////////////////////////////////////////////////////////////////////////////////
1059 /// Set the first entry to be displayed.
1060 
1061 void TParallelCoord::SetCurrentFirst(Long64_t f)
1062 {
1063  if(f<0 || f>fNentries) return;
1064  fCurrentFirst = f;
1065  if(fCurrentFirst + fCurrentN > fNentries) fCurrentN = fNentries-fCurrentFirst;
1066  TIter next(fVarList);
1067  TParallelCoordVar* var;
1068  while ((var = (TParallelCoordVar*)next())) {
1069  var->GetMinMaxMean();
1070  var->GetHistogram();
1071  if (var->TestBit(TParallelCoordVar::kShowBox)) var->GetQuantiles();
1072  }
1073 }
1074 
1075 ////////////////////////////////////////////////////////////////////////////////
1076 /// Set the number of entry to be displayed.
1077 
1078 void TParallelCoord::SetCurrentN(Long64_t n)
1079 {
1080  if(n<=0) return;
1081  if(fCurrentFirst+n>fNentries) fCurrentN = fNentries-fCurrentFirst;
1082  else fCurrentN = n;
1083  TIter next(fVarList);
1084  TParallelCoordVar* var;
1085  while ((var = (TParallelCoordVar*)next())) {
1086  var->GetMinMaxMean();
1087  var->GetHistogram();
1088  if (var->TestBit(TParallelCoordVar::kShowBox)) var->GetQuantiles();
1089  }
1090 }
1091 
1092 ////////////////////////////////////////////////////////////////////////////////
1093 /// Set the selection being edited.
1094 
1095 TParallelCoordSelect* TParallelCoord::SetCurrentSelection(const char* title)
1096 {
1097  if (fCurrentSelection && fCurrentSelection->GetTitle() == title) return fCurrentSelection;
1098  TIter next(fSelectList);
1099  TParallelCoordSelect* sel;
1100  while((sel = (TParallelCoordSelect*)next()) && strcmp(sel->GetTitle(),title))
1101  if (sel) fCurrentSelection = sel;
1102  return sel;
1103 }
1104 
1105 ////////////////////////////////////////////////////////////////////////////////
1106 /// Set the selection being edited.
1107 
1108 void TParallelCoord::SetCurrentSelection(TParallelCoordSelect* sel)
1109 {
1110  if (fCurrentSelection == sel) return;
1111  fCurrentSelection = sel;
1112 }
1113 
1114 ////////////////////////////////////////////////////////////////////////////////
1115 /// Set dots spacing. Modify the line style 11.
1116 /// If the canvas support transparency dot spacing is ignored.
1117 
1118 void TParallelCoord::SetDotsSpacing(Int_t s)
1119 {
1120  if (gPad->GetCanvas()->SupportAlpha()) return;
1121  if (s == fDotsSpacing) return;
1122  fDotsSpacing = s;
1123  gStyle->SetLineStyleString(11,Form("%d %d",4,s*8));
1124 }
1125 
1126 ////////////////////////////////////////////////////////////////////////////////
1127 /// Set the entry lists of "para".
1128 
1129 void TParallelCoord::SetEntryList(TParallelCoord* para, TEntryList* enlist)
1130 {
1131  para->SetCurrentEntries(enlist);
1132  para->SetInitEntries(enlist);
1133 }
1134 
1135 ////////////////////////////////////////////////////////////////////////////////
1136 /// Force all variables to adopt the same max.
1137 
1138 void TParallelCoord::SetGlobalMax(Double_t max)
1139 {
1140  TIter next(fVarList);
1141  TParallelCoordVar* var;
1142  while ((var = (TParallelCoordVar*)next())) {
1143  var->SetCurrentMax(max);
1144  }
1145 }
1146 
1147 ////////////////////////////////////////////////////////////////////////////////
1148 /// Force all variables to adopt the same min.
1149 
1150 void TParallelCoord::SetGlobalMin(Double_t min)
1151 {
1152  TIter next(fVarList);
1153  TParallelCoordVar* var;
1154  while ((var = (TParallelCoordVar*)next())) {
1155  var->SetCurrentMin(min);
1156  }
1157 }
1158 
1159 ////////////////////////////////////////////////////////////////////////////////
1160 /// If true, the pad is updated while the motion of a dragged range.
1161 
1162 void TParallelCoord::SetLiveRangesUpdate(Bool_t on)
1163 {
1164  SetBit(kLiveUpdate,on);
1165  TIter next(fVarList);
1166  TParallelCoordVar* var;
1167  while((var = (TParallelCoordVar*)next())) var->SetLiveRangesUpdate(on);
1168 }
1169 
1170 ////////////////////////////////////////////////////////////////////////////////
1171 /// Set the vertical or horizontal display.
1172 
1173 void TParallelCoord::SetVertDisplay(Bool_t vert)
1174 {
1175  if (vert == TestBit (kVertDisplay)) return;
1176  SetBit(kVertDisplay,vert);
1177  if (!gPad) return;
1178  TFrame* frame = gPad->GetFrame();
1179  if (!frame) return;
1180  UInt_t ui = 0;
1181  Double_t horaxisspace = (frame->GetX2() - frame->GetX1())/(fNvar-1);
1182  Double_t veraxisspace = (frame->GetY2() - frame->GetY1())/(fNvar-1);
1183  TIter next(fVarList);
1184  TParallelCoordVar* var;
1185  while ((var = (TParallelCoordVar*)next())) {
1186  if (vert) var->SetX(frame->GetX1() + ui*horaxisspace,TestBit(kGlobalScale));
1187  else var->SetY(frame->GetY1() + ui*veraxisspace,TestBit(kGlobalScale));
1188  ++ui;
1189  }
1190  if (TestBit(kCandleChart)) {
1191  if (fCandleAxis) delete fCandleAxis;
1192  if (TestBit(kVertDisplay)) fCandleAxis = new TGaxis(0.05,0.1,0.05,0.9,GetGlobalMin(),GetGlobalMax());
1193  else fCandleAxis = new TGaxis(0.1,0.05,0.9,0.05,GetGlobalMin(),GetGlobalMax());
1194  fCandleAxis->Draw();
1195  }
1196  gPad->Modified();
1197  gPad->Update();
1198 }
1199 
1200 ////////////////////////////////////////////////////////////////////////////////
1201 /// Unzoom all variables.
1202 
1203 void TParallelCoord::UnzoomAll()
1204 {
1205  TIter next(fVarList);
1206  TParallelCoordVar* var;
1207  while((var = (TParallelCoordVar*)next())) var->Unzoom();
1208 }