Logo ROOT   6.08/07
Reference Guide
TLegend.cxx
Go to the documentation of this file.
1 // @(#)root/graf:$Id$
2 // Author: Matthew.Adam.Dobbs 06/09/99
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 #include <stdio.h>
13 
14 #include "TStyle.h"
15 #include "TLatex.h"
16 #include "TLine.h"
17 #include "TBox.h"
18 #include "TMarker.h"
19 #include "TLegend.h"
20 #include "TList.h"
21 #include "TVirtualPad.h"
22 #include "TMath.h"
23 #include "TROOT.h"
24 #include "TLegendEntry.h"
25 #include "Riostream.h"
26 #include "TMultiGraph.h"
27 #include "THStack.h"
28 
30 
31 /** \class TLegend
32 \ingroup BasicGraphics
33 
34 This class displays a legend box (TPaveText) containing several legend entries.
35 
36 Each legend entry is made of a reference to a ROOT object, a text label and an
37 option specifying which graphical attributes (marker/line/fill) should be
38 displayed.
39 
40 The following example shows how to create a legend. In this example the legend
41 contains a histogram, a function and a graph. The histogram is put in the legend
42 using its reference pointer whereas the graph and the function are added
43 using their names. Note that, because `TGraph` constructors do not have the
44 `TGraph` name as parameter, the graph name should be specified using the
45 `SetName` method.
46 
47 When an object is added by name, a scan is performed on the list of objects
48 contained in the current pad (`gPad`) and also in the possible
49 `TMultiGraph` and `THStack` present in the pad. If a matching
50 name is found, the corresponding object is added in the legend using its pointer.
51 
52 Begin_Macro(source)
53 {
54  TCanvas *c1 = new TCanvas("c1","c1",600,500);
55  gStyle->SetOptStat(0);
56 
57  TH1F *h1 = new TH1F("h1","TLegend Example",200,-10,10);
58  h1->FillRandom("gaus",30000);
59  h1->SetFillColor(kGreen);
60  h1->SetFillStyle(3003);
61  h1->Draw();
62 
63  TF1 *f1=new TF1("f1","1000*TMath::Abs(sin(x)/x)",-10,10);
64  f1->SetLineColor(kBlue);
65  f1->SetLineWidth(4);
66  f1->Draw("same");
67 
68  const Int_t n = 20;
69  Double_t x[n], y[n], ex[n], ey[n];
70  for (Int_t i=0;i<n;i++) {
71  x[i] = i*0.1;
72  y[i] = 1000*sin(x[i]+0.2);
73  x[i] = 17.8*x[i]-8.9;
74  ex[i] = 1.0;
75  ey[i] = 10.*i;
76  }
77  TGraphErrors *gr = new TGraphErrors(n,x,y,ex,ey);
78  gr->SetName("gr");
79  gr->SetLineColor(kRed);
80  gr->SetLineWidth(2);
81  gr->SetMarkerStyle(21);
82  gr->SetMarkerSize(1.3);
83  gr->SetMarkerColor(7);
84  gr->Draw("P");
85 
86  leg = new TLegend(0.1,0.7,0.48,0.9);
87  leg->SetHeader("The Legend Title","C"); // option "C" allows to center the header
88  leg->AddEntry(h1,"Histogram filled with random numbers","f");
89  leg->AddEntry("f1","Function abs(#frac{sin(x)}{x})","l");
90  leg->AddEntry("gr","Graph with error bars","lep");
91  leg->Draw();
92 
93  return c1;
94 }
95 End_Macro
96 
97 
98 `TLegend` inherits from `TAttText` therefore changing any
99 text attributes (text alignment, font, color...) on a legend will changed the
100 text attributes on each line.
101 
102 In particular it can be interesting to change the text alignement that way. In
103 order to have a base-line vertical alignment instead of a centered one simply do:
104 ~~~ {.cpp}
105  leg->SetTextAlign(13);
106 ~~~
107 or
108 ~~~ {.cpp}
109  leg->SetTextAlign(11);
110 ~~~
111 The default value of some `TLegend` attributes can be changed using
112 `gStyle`. The default settings are:
113 ~~~ {.cpp}
114  SetLegendBorderSize(1);
115  SetLegendFillColor(0);
116  SetLegendFont(42);
117  SetLegendTextSize(0.);
118 ~~~
119 The global attributes change the default values for the next created legends.
120 
121 Text attributes can be also changed individually on each legend entry:
122 ~~~ {.cpp}
123  TLegendEntry *le = leg->AddEntry(h1,"Histogram filled with random numbers","f");
124  le->SetTextColor(kBlue);;
125 ~~~
126 
127 Note that the `TPad` class has a method to build automatically a legend
128 for all objects in the pad. It is called `TPad::BuildLegend()`.
129 
130 Each item in the legend is added using the `AddEntry` method. This
131 method defines the object to be added (by reference or name), the label
132 associated to this object and an option which a combination of:
133 
134  - L: draw line associated with TAttLine if obj inherits from TAttLine
135  - P: draw polymarker associated with TAttMarker if obj inherits from TAttMarker
136  - F: draw a box with fill associated wit TAttFill if obj inherits TAttFill
137  - E: draw vertical error bar
138 
139 As shown in the following example, passing a NULL pointer as first parameter in
140 `AddEntry` is also valid. This allows to add text or blank lines in a
141 legend.
142 
143 Begin_Macro(source)
144 {
145  TCanvas *c2 = new TCanvas("c2","c2",500,300);
146 
147  TLegend* leg = new TLegend(0.2, 0.2, .8, .8);
148  TH1* h = new TH1F("", "", 1, 0, 1);
149 
150  leg->AddEntry(h, "Histogram \"h\"", "l");
151  leg->AddEntry((TObject*)0, "", "");
152  leg->AddEntry((TObject*)0, "Some text", "");
153  leg->AddEntry((TObject*)0, "", "");
154  leg->AddEntry(h, "Histogram \"h\" again", "l");
155 
156  leg->Draw();
157  return c2;
158 }
159 End_Macro
160 
161 It is possible to draw the legend entries over several columns using
162 the method `SetNColumns()` like in the following example.
163 
164 Begin_Macro(source)
165 {
166  TCanvas *c3 = new TCanvas("c2","c2",500,300);
167 
168  TLegend* leg = new TLegend(0.2, 0.2, .8, .8);
169  TH1* h = new TH1F("", "", 1, 0, 1);
170 
171  leg-> SetNColumns(2);
172 
173  leg->AddEntry(h, "Column 1 line 1", "l");
174  leg->AddEntry(h, "Column 2 line 1", "l");
175  leg->AddEntry(h, "Column 1 line 2", "l");
176  leg->AddEntry(h, "Column 2 line 2", "l");
177 
178  leg->Draw();
179  return c3;
180 }
181 End_Macro
182 */
183 
184 ////////////////////////////////////////////////////////////////////////////////
185 /// Default constructor.
186 
187 TLegend::TLegend(): TPave(), TAttText(12,0,1,gStyle->GetLegendFont(),0)
188 {
189  fPrimitives = 0;
190  SetDefaults();
193 }
194 
195 ////////////////////////////////////////////////////////////////////////////////
196 /// Normal constructor.
197 ///
198 /// A TLegend is a Pave with several TLegendEntry(s).
199 /// x1,y1,x2,y2 are the coordinates of the Legend in the current pad
200 /// (in normalised coordinates by default)
201 /// "header" is the title that will be displayed at the top of the legend
202 /// it is treated like a regular entry and supports TLatex. The default
203 /// is no header (header = 0).
204 /// The options are the same as for TPave Default = "brNDC"
205 
207  const char *header, Option_t *option)
208  :TPave(x1,y1,x2,y2,4,option), TAttText(12,0,1,gStyle->GetLegendFont(),0)
209 {
210  fPrimitives = new TList;
211  if ( header && strlen(header) > 0) {
212  TLegendEntry *headerEntry = new TLegendEntry( 0, header, "h" );
213  headerEntry->SetTextAlign(0);
214  headerEntry->SetTextAngle(0);
215  headerEntry->SetTextColor(0);
216  headerEntry->SetTextFont(gStyle->GetLegendFont());
217  headerEntry->SetTextSize(0);
218  fPrimitives->AddFirst(headerEntry);
219  }
220  SetDefaults();
223 }
224 
225 ////////////////////////////////////////////////////////////////////////////////
226 /// Copy constructor.
227 
228 TLegend::TLegend( const TLegend &legend ) : TPave(legend), TAttText(legend),
229  fPrimitives(0)
230 {
231  if (legend.fPrimitives) {
232  fPrimitives = new TList();
233  TListIter it(legend.fPrimitives);
234  while (TLegendEntry *e = (TLegendEntry *)it.Next()) {
235  TLegendEntry *newentry = new TLegendEntry(*e);
236  fPrimitives->Add(newentry);
237  }
238  }
239  ((TLegend&)legend).Copy(*this);
240 }
241 
242 ////////////////////////////////////////////////////////////////////////////////
243 /// Assignment operator.
244 
246 {
247  if(this!=&lg) {
248  TPave::operator=(lg);
249  TAttText::operator=(lg);
252  fMargin=lg.fMargin;
253  fNColumns=lg.fNColumns;
254  }
255  return *this;
256 }
257 
258 ////////////////////////////////////////////////////////////////////////////////
259 /// Default destructor.
260 
262 {
264  delete fPrimitives;
265  fPrimitives = 0;
266 }
267 
268 ////////////////////////////////////////////////////////////////////////////////
269 /// Add a new entry to this legend. "obj" is the object to be represented.
270 /// "label" is the text you wish to associate with obj in the legend.
271 /// If "label" is null or empty, the title of the object will be used.
272 ///
273 /// Options are:
274 ///
275 /// - L: draw line associated with TAttLine if obj inherits from TAttLine
276 /// - P: draw polymarker associated with TAttMarker if obj inherits from TAttMarker
277 /// - F: draw a box with fill associated wit TAttFill if obj inherits TAttFill
278 /// - E: draw vertical error bar if option "L" is also specified
279 
280 TLegendEntry *TLegend::AddEntry(const TObject *obj, const char *label, Option_t *option)
281 {
282  const char *lab = label;
283 
284  if (obj && (!label || strlen(label)==0)) lab = obj->GetTitle();
285  TLegendEntry *newentry = new TLegendEntry( obj, lab, option );
286  if ( !fPrimitives ) fPrimitives = new TList;
287  fPrimitives->Add(newentry);
288  return newentry;
289 }
290 
291 ////////////////////////////////////////////////////////////////////////////////
292 /// Add a new entry to this legend. "name" is the name of an object in the pad to
293 /// be represented label is the text you wish to associate with obj in the legend
294 /// if label is null or empty, the title of the object will be used.
295 ///
296 /// Options are:
297 ///
298 /// - L: draw line associated with TAttLine if obj inherits from TAttLine
299 /// - P: draw polymarker associated with TAttMarker if obj inherits from TAttMarker
300 /// - F: draw a box with fill associated wit TAttFill if obj inherits TAttFill
301 /// - E: draw vertical error bar if option "L" is also specified
302 
303 TLegendEntry *TLegend::AddEntry(const char *name, const char *label, Option_t *option)
304 {
305  if (!gPad) {
306  Error("AddEntry", "need to create a canvas first");
307  return 0;
308  }
309 
310  TObject *obj = gPad->FindObject(name);
311 
312  // If the object "name" has not been found, the following code tries to
313  // find it in TMultiGraph or THStack possibly present in the current pad.
314  if (!obj) {
315  TList *lop = gPad->GetListOfPrimitives();
316  if (lop) {
317  TObject *o=0;
318  TIter next(lop);
319  while( (o=next()) ) {
320  if ( o->InheritsFrom(TMultiGraph::Class() ) ) {
321  TList * grlist = ((TMultiGraph *)o)->GetListOfGraphs();
322  obj = grlist->FindObject(name);
323  if (obj) break;
324  }
325  if ( o->InheritsFrom(THStack::Class() ) ) {
326  TList * hlist = ((THStack *)o)->GetHists();
327  obj = hlist->FindObject(name);
328  if (obj) break;
329  }
330  }
331  }
332  }
333 
334  return AddEntry( obj, label, option );
335 }
336 
337 ////////////////////////////////////////////////////////////////////////////////
338 /// Clear all entries in this legend, including the header.
339 
341 {
342  if (!fPrimitives) return;
343  fPrimitives->Delete();
344 }
345 
346 ////////////////////////////////////////////////////////////////////////////////
347 /// Copy this legend into "obj".
348 
349 void TLegend::Copy( TObject &obj ) const
350 {
351  TPave::Copy(obj);
352  TAttText::Copy((TLegend&)obj);
353  ((TLegend&)obj).fEntrySeparation = fEntrySeparation;
354  ((TLegend&)obj).fMargin = fMargin;
355  ((TLegend&)obj).fNColumns = fNColumns;
356 }
357 
358 ////////////////////////////////////////////////////////////////////////////////
359 /// Delete entry at the mouse position.
360 
362 {
363  if ( !fPrimitives ) return;
364  TLegendEntry* entry = GetEntry(); // get entry pointed by the mouse
365  if ( !entry ) return;
366  fPrimitives->Remove(entry);
367  delete entry;
368 }
369 
370 ////////////////////////////////////////////////////////////////////////////////
371 /// Draw this legend with its current attributes.
372 
373 void TLegend::Draw( Option_t *option )
374 {
375  AppendPad(option);
376 }
377 
378 ////////////////////////////////////////////////////////////////////////////////
379 /// Edit the fill attributes for the entry pointed by the mouse.
380 
382 {
383  TLegendEntry* entry = GetEntry(); // get entry pointed by the mouse
384  if ( !entry ) return;
385  gROOT->SetSelectedPrimitive( entry );
386  entry->SetFillAttributes();
387 }
388 
389 ////////////////////////////////////////////////////////////////////////////////
390 /// Edit the line attributes for the entry pointed by the mouse.
391 
393 {
394  TLegendEntry* entry = GetEntry(); // get entry pointed by the mouse
395  if ( !entry ) return;
396  gROOT->SetSelectedPrimitive( entry );
397  entry->SetLineAttributes();
398 }
399 
400 ////////////////////////////////////////////////////////////////////////////////
401 /// Edit the marker attributes for the entry pointed by the mouse.
402 
404 {
405  TLegendEntry* entry = GetEntry(); // get entry pointed by the mouse
406  if ( !entry ) return;
407  gROOT->SetSelectedPrimitive( entry );
408  entry->SetMarkerAttributes();
409 }
410 
411 ////////////////////////////////////////////////////////////////////////////////
412 /// Edit the text attributes for the entry pointed by the mouse.
413 
415 {
416  TLegendEntry* entry = GetEntry(); // get entry pointed by the mouse
417  if ( !entry ) return;
418  gROOT->SetSelectedPrimitive( entry );
419  entry->SetTextAttributes();
420 }
421 
422 ////////////////////////////////////////////////////////////////////////////////
423 /// Get entry pointed to by the mouse.
424 /// This method is mostly a tool for other methods inside this class.
425 
427 {
428  if (!gPad) {
429  Error("GetEntry", "need to create a canvas first");
430  return 0;
431  }
432 
433  Int_t nRows = GetNRows();
434  if ( nRows == 0 ) return 0;
435 
436  Double_t ymouse = gPad->AbsPixeltoY(gPad->GetEventY())-fY1;
437  Double_t yspace = (fY2 - fY1)/nRows;
438 
439  Int_t nColumns = GetNColumns();
440  Double_t xmouse = gPad->AbsPixeltoX(gPad->GetEventX())-fX1;
441  Double_t xspace = 0.;
442  if (nColumns > 0) xspace = (fX2 - fX1)/nColumns;
443 
444  Int_t ix = 1;
445  if (xspace > 0.) ix = (Int_t)(xmouse/xspace)+1;
446  if (ix > nColumns) ix = nColumns;
447  if (ix < 1) ix = 1;
448 
449  Int_t iy = nRows-(Int_t)(ymouse/yspace);
450  if (iy > nRows) iy = nRows;
451  if (iy < 1) iy = 1;
452 
453  Int_t nloops = TMath::Min(ix+(nColumns*(iy-1)), fPrimitives->GetSize());
454 
455  TIter next(fPrimitives);
456  TLegendEntry *entry = 0;
457 
458  for (Int_t i=1; i<= nloops; i++) entry = (TLegendEntry *)next();
459 
460  return entry;
461 }
462 
463 ////////////////////////////////////////////////////////////////////////////////
464 /// Returns the header, which is the title that appears at the top
465 /// of the legend.
466 
467 const char *TLegend::GetHeader() const
468 {
469  if ( !fPrimitives ) return 0;
470  TIter next(fPrimitives);
471  TLegendEntry *first; // header is always the first entry
472  if (( first = (TLegendEntry*)next() )) {
473  TString opt = first->GetOption();
474  opt.ToLower();
475  if ( opt.Contains("h") ) return first->GetLabel();
476  }
477  return 0;
478 }
479 
480 ////////////////////////////////////////////////////////////////////////////////
481 /// Add a new entry before the entry at the mouse position.
482 
483 void TLegend::InsertEntry( const char* objectName, const char* label, Option_t* option)
484 {
485  if (!gPad) {
486  Error("InsertEntry", "need to create a canvas first");
487  return;
488  }
489 
490  TLegendEntry* beforeEntry = GetEntry(); // get entry pointed by the mouse
491  TObject *obj = gPad->FindObject( objectName );
492 
493  // note either obj OR beforeEntry may be zero at this point
494 
495  TLegendEntry *newentry = new TLegendEntry( obj, label, option );
496 
497  if ( !fPrimitives ) fPrimitives = new TList;
498  if ( beforeEntry ) {
499  fPrimitives->AddBefore( (TObject*)beforeEntry, (TObject*)newentry );
500  } else {
501  fPrimitives->Add((TObject*)newentry);
502  }
503 }
504 
505 ////////////////////////////////////////////////////////////////////////////////
506 /// Paint this legend with its current attributes.
507 
508 void TLegend::Paint( Option_t* option )
509 {
512  PaintPrimitives();
513 }
514 
515 ////////////////////////////////////////////////////////////////////////////////
516 /// Get the number of rows.
517 
519 {
520  Int_t nEntries = 0;
521  if ( fPrimitives ) nEntries = fPrimitives->GetSize();
522  if ( nEntries == 0 ) return 0;
523 
524  Int_t nRows;
525  if(GetHeader() != NULL) nRows = 1 + (Int_t) TMath::Ceil((Double_t) (nEntries-1)/fNColumns);
526  else nRows = (Int_t) TMath::Ceil((Double_t) nEntries/fNColumns);
527 
528  return nRows;
529 }
530 
531 ////////////////////////////////////////////////////////////////////////////////
532 /// Set the number of columns for the legend. The header, if set, is given
533 /// its own row. After that, every nColumns entries are inserted into the
534 /// same row. For example, if one calls legend.SetNColumns(2), and there
535 /// is no header, then the first two TObjects added to the legend will be
536 /// in the first row, the next two will appear in the second row, and so on.
537 
539 {
540  if(nColumns < 1) {
541  Warning("TLegend::SetNColumns", "illegal value nColumns = %d; keeping fNColumns = %d", nColumns, fNColumns);
542  return;
543  }
544  fNColumns = nColumns;
545 }
546 
547 ////////////////////////////////////////////////////////////////////////////////
548 /// Paint the entries (list of primitives) for this legend.
549 
551 {
552  Int_t nRows = GetNRows();
553  if ( nRows == 0 ) return;
554 
555  // Evaluate text size as a function of the number of entries
556  // taking into account their real size after drawing latex
557  // Note: in pixel coords y1 > y2=0, but x2 > x1=0
558  // in NDC y2 > y1, and x2 > x1
559 
560  Double_t x1 = fX1NDC;
561  Double_t y1 = fY1NDC;
562  Double_t x2 = fX2NDC;
563  Double_t y2 = fY2NDC;
564  Double_t margin = fMargin*( x2-x1 )/fNColumns;
565  Double_t boxwidth = margin;
566  Double_t boxw = boxwidth*0.35;
567  Double_t yspace = (y2-y1)/nRows;
568  Double_t yspace2 = yspace/2.;
569  Double_t textsize = GetTextSize();
570  Double_t save_textsize = textsize;
571  if (textsize==0.) {
573  textsize = GetTextSize();
574  }
575  Bool_t autosize = kFALSE;
576  Double_t* columnWidths = new Double_t[fNColumns];
577  memset(columnWidths, 0, fNColumns*sizeof(Double_t));
578 
579  if ( textsize == 0 ) {
580  autosize = kTRUE;
581  textsize = ( 1. - fEntrySeparation ) * yspace;
582 
583  // find the max width and height (in pad coords) of one latex entry label
584  Double_t maxentrywidth = 0, maxentryheight = 0;
585  TIter nextsize(fPrimitives);
586  TLegendEntry *entrysize;
587  Int_t iColumn = 0;
588  while (( entrysize = (TLegendEntry *)nextsize() )) {
589  TLatex entrytex( 0, 0, entrysize->GetLabel() );
590  entrytex.SetNDC();
591  Style_t tfont = entrysize->GetTextFont();
592  if (tfont == 0) tfont = GetTextFont();
593  if (tfont%10 == 3) --tfont;
594  entrytex.SetTextFont(tfont);
595  entrytex.SetTextSize(textsize);
596  if ( entrytex.GetYsize() > maxentryheight ) {
597  maxentryheight = entrytex.GetYsize();
598  }
599  TString opt = entrysize->GetOption();
600  opt.ToLower();
601  if ( opt.Contains("h") ) {
602  if ( entrytex.GetXsize() > maxentrywidth ) {
603  maxentrywidth = entrytex.GetXsize();
604  }
605  } else {
606  if ( entrytex.GetXsize() > columnWidths[iColumn] ) {
607  columnWidths[iColumn] = entrytex.GetXsize();
608  }
609  iColumn++;
610  iColumn %= fNColumns;
611  }
612  Double_t tmpMaxWidth = 0.0;
613  for(int i=0; i<fNColumns; i++) tmpMaxWidth += columnWidths[i];
614  if ( tmpMaxWidth > maxentrywidth) maxentrywidth = tmpMaxWidth;
615  }
616  // make sure all labels fit in the allotted space
617  Double_t tmpsize_h = maxentryheight /(gPad->GetY2() - gPad->GetY1());
618  textsize = TMath::Min( textsize, tmpsize_h );
619  Double_t tmpsize_w = textsize*(fX2-fX1)*(1.0-fMargin)/maxentrywidth;
620  if(fNColumns > 1) tmpsize_w = textsize*(fX2-fX1)*(1.0-fMargin-fColumnSeparation)/maxentrywidth;
621  textsize = TMath::Min( textsize, tmpsize_w );
622  SetTextSize( textsize );
623  }
624 
625  // Update column widths, put into NDC units
626  // block off this section of code to make sure all variables are local:
627  // don't want to ruin initialisation of these variables later on
628  {
629  TIter next(fPrimitives);
630  TLegendEntry *entry;
631  Int_t iColumn = 0;
632  memset(columnWidths, 0, fNColumns*sizeof(Double_t));
633  while (( entry = (TLegendEntry *)next() )) {
634  TLatex entrytex( 0, 0, entry->GetLabel() );
635  entrytex.SetNDC();
636  Style_t tfont = entry->GetTextFont();
637  if (tfont == 0) tfont = GetTextFont();
638  if (autosize && tfont%10 == 3) --tfont;
639  entrytex.SetTextFont(tfont);
640  if(entry->GetTextSize() == 0) entrytex.SetTextSize(textsize);
641  TString opt = entry->GetOption();
642  opt.ToLower();
643  if (!opt.Contains("h")) {
644  if ( entrytex.GetXsize() > columnWidths[iColumn] ) {
645  columnWidths[iColumn] = entrytex.GetXsize();
646  }
647  iColumn++;
648  iColumn %= fNColumns;
649  }
650  }
651  double totalWidth = 0.0;
652  for(int i=0; i<fNColumns; i++) totalWidth += columnWidths[i];
653  if(fNColumns > 1) totalWidth /= (1.0-fMargin-fColumnSeparation);
654  else totalWidth /= (1.0 - fMargin);
655  for(int i=0; i<fNColumns; i++) {
656  columnWidths[i] = columnWidths[i]/totalWidth*(x2-x1) + margin;
657  }
658  }
659 
660  Double_t ytext = y2 + 0.5*yspace; // y-location of 0th entry
661 
662  // iterate over and paint all the TLegendEntries
663  TIter next(fPrimitives);
664  TLegendEntry *entry;
665  Int_t iColumn = 0;
666  while (( entry = (TLegendEntry *)next() )) {
667  if(iColumn == 0) ytext -= yspace;
668 
669  // Draw Label in Latexmargin
670 
671  Short_t talign = entry->GetTextAlign();
672  Float_t tangle = entry->GetTextAngle();
673  Color_t tcolor = entry->GetTextColor();
674  Style_t tfont = entry->GetTextFont();
675  Size_t tsize = entry->GetTextSize();
676  // if the user hasn't set a parameter, then set it to the TLegend value
677  if (talign == 0) entry->SetTextAlign(GetTextAlign());
678  if (tangle == 0) entry->SetTextAngle(GetTextAngle());
679  if (tcolor == 0) entry->SetTextColor(GetTextColor());
680  if (tfont == 0) {
681  tfont = GetTextFont();
682  if (autosize && tfont%10 == 3) --tfont;
683  entry->SetTextFont(tfont);
684  }
685  if (tsize == 0) entry->SetTextSize(GetTextSize());
686  // set x,y according to the requested alignment
687  Double_t x=0,y=0;
688  Int_t halign = entry->GetTextAlign()/10;
689  Double_t entrymargin = margin;
690  // for the header the margin is near zero
691  TString opt = entry->GetOption();
692  opt.ToLower();
693  x1 = fX1NDC;
694  x2 = fX2NDC;
695  if ( opt.Contains("h") ) entrymargin = margin/10.;
696  else if (fNColumns > 1) {
697  for(int i=0; i<iColumn; i++) x1 += columnWidths[i] + fColumnSeparation*(fX2NDC-fX1NDC)/(fNColumns-1);
698  x2 = x1 + columnWidths[iColumn];
699  iColumn++;
700  iColumn %= fNColumns;
701  }
702  if (halign == 1) x = x1 + entrymargin;
703  if (halign == 2) x = 0.5*( (x1+entrymargin) + x2 );
704  if (halign == 3) x = x2 - entrymargin/10.;
705  Int_t valign = entry->GetTextAlign()%10;
706 
707  if (valign == 1) y = ytext - (1. - fEntrySeparation)* yspace2;
708  if (valign == 3) y = ytext + (1. - fEntrySeparation)* yspace2;
709 
710  // The vertical alignment "centered" is treated in a special way
711  // to ensure a better spacing between lines.
712  if (valign == 2) {
713  Float_t tsizepad = textsize;
714  if (tfont%10 == 3) tsizepad = (gPad->AbsPixeltoY(0) - gPad->AbsPixeltoY(textsize))/(gPad->GetY2() - gPad->GetY1());
715  if (yspace2 < tsizepad) {
716  entry->SetTextAlign(10*halign+1);
717  y = ytext - (1. - fEntrySeparation)* yspace2/2.;
718  } else {
719  y = ytext;
720  }
721  }
722 
723  TLatex entrytex( x, y, entry->GetLabel() );
724  entrytex.SetNDC();
725  entry->TAttText::Copy(entrytex);
726  entrytex.Paint();
727 
728  // reset attributes back to their original values
729  entry->SetTextAlign(talign);
730  entry->SetTextAngle(tangle);
731  entry->SetTextColor(tcolor);
732  entry->SetTextFont(tfont);
733  entry->SetTextSize(tsize);
734 
735  // define x,y as the center of the symbol for this entry
736  Double_t xsym = x1 + margin/2.;
737  Double_t ysym = ytext;
738 
739  TObject *eobj = entry->GetObject();
740 
741  // Draw fill pattern (in a box)
742 
743  if ( opt.Contains("f")) {
744  if (eobj && eobj->InheritsFrom(TAttFill::Class())) {
745  dynamic_cast<TAttFill*>(eobj)->Copy(*entry);
746  }
747 
748  // box total height is yspace*0.7
749  entry->TAttFill::Modify();
750  Double_t xf[4],yf[4];
751  xf[0] = xsym - boxw;
752  yf[0] = ysym - yspace*0.35;
753  xf[1] = xsym + boxw;
754  yf[1] = yf[0];
755  xf[2] = xf[1];
756  yf[2] = ysym + yspace*0.35;
757  xf[3] = xf[0];
758  yf[3] = yf[2];
759  for (Int_t i=0;i<4;i++) {
760  xf[i] = gPad->GetX1() + xf[i]*(gPad->GetX2()-gPad->GetX1());
761  yf[i] = gPad->GetY1() + yf[i]*(gPad->GetY2()-gPad->GetY1());
762  }
763  gPad->PaintFillArea(4,xf,yf);
764  }
765 
766  // Draw line
767 
768  if ( opt.Contains("l") || opt.Contains("f")) {
769 
770  if (eobj && eobj->InheritsFrom(TAttLine::Class())) {
771  dynamic_cast<TAttLine*>(eobj)->Copy(*entry);
772  }
773  // line total length (in x) is margin*0.8
774  TLine entryline( xsym - boxw, ysym, xsym + boxw, ysym );
775  entryline.SetBit(TLine::kLineNDC);
776  entry->TAttLine::Copy(entryline);
777  // if the entry is filled, then surround the box with the line instead
778  if ( opt.Contains("f") && !opt.Contains("l")) {
779  // box total height is yspace*0.7
780  boxwidth = yspace*
781  (gPad->GetX2()-gPad->GetX1())/(gPad->GetY2()-gPad->GetY1());
782  if ( boxwidth > margin ) boxwidth = margin;
783 
784  entryline.PaintLineNDC( xsym - boxw, ysym + yspace*0.35,
785  xsym + boxw, ysym + yspace*0.35);
786  entryline.PaintLineNDC( xsym - boxw, ysym - yspace*0.35,
787  xsym + boxw, ysym - yspace*0.35);
788  entryline.PaintLineNDC( xsym + boxw, ysym - yspace*0.35,
789  xsym + boxw, ysym + yspace*0.35);
790  entryline.PaintLineNDC( xsym - boxw, ysym - yspace*0.35,
791  xsym - boxw, ysym + yspace*0.35);
792  } else {
793  entryline.Paint();
794  if (opt.Contains("e")) {
795  entryline.PaintLineNDC( xsym, ysym - yspace*0.30,
796  xsym, ysym + yspace*0.30);
797  }
798  }
799  }
800 
801  // Draw error only
802 
803  if (opt.Contains("e") && !(opt.Contains("l") || opt.Contains("f"))) {
804  if (eobj && eobj->InheritsFrom(TAttLine::Class())) {
805  dynamic_cast<TAttLine*>(eobj)->Copy(*entry);
806  }
807  TLine entryline(xsym, ysym - yspace*0.30,
808  xsym, ysym + yspace*0.30);
809  entryline.SetBit(TLine::kLineNDC);
810  entry->TAttLine::Copy(entryline);
811  entryline.Paint();
812  }
813 
814  // Draw Polymarker
815 
816  if ( opt.Contains("p")) {
817 
818  if (eobj && eobj->InheritsFrom(TAttMarker::Class())) {
819  dynamic_cast<TAttMarker*>(eobj)->Copy(*entry);
820  }
821  TMarker entrymarker( xsym, ysym, 0 );
822  entrymarker.SetNDC();
823  entry->TAttMarker::Copy(entrymarker);
824  entrymarker.Paint();
825  }
826  }
827 
828  SetTextSize(save_textsize);
829  delete [] columnWidths;
830 }
831 
832 ////////////////////////////////////////////////////////////////////////////////
833 /// Dump this TLegend and its contents.
834 
835 void TLegend::Print( Option_t* option ) const
836 {
837  TPave::Print( option );
838  if (fPrimitives) fPrimitives->Print();
839 }
840 
841 ////////////////////////////////////////////////////////////////////////////////
842 /// Reset the legend entries pointing to "obj".
843 
845 {
846  TIter next(fPrimitives);
847  TLegendEntry *entry;
848  while (( entry = (TLegendEntry *)next() )) {
849  if (entry->GetObject() == obj) entry->SetObject((TObject*)0);
850  }
851 }
852 
853 ////////////////////////////////////////////////////////////////////////////////
854 /// Save this legend as C++ statements on output stream out
855 /// to be used with the SaveAs .C option.
856 
857 void TLegend::SavePrimitive(std::ostream &out, Option_t* )
858 {
859 
860  out << " " << std::endl;
861  char quote = '"';
862  if ( gROOT->ClassSaved( TLegend::Class() ) ) {
863  out << " ";
864  } else {
865  out << " TLegend *";
866  }
867  // note, we can always use NULL header, since its included in primitives
868  out << "leg = new TLegend("<<GetX1NDC()<<","<<GetY1NDC()<<","
869  <<GetX2NDC()<<","<<GetY2NDC()<<","
870  << "NULL" << "," <<quote<< fOption <<quote<<");" << std::endl;
871  if (fBorderSize != 4) {
872  out<<" leg->SetBorderSize("<<fBorderSize<<");"<<std::endl;
873  }
874  SaveTextAttributes(out,"leg",12,0,1,42,0);
875  SaveLineAttributes(out,"leg",-1,-1,-1);
876  SaveFillAttributes(out,"leg",-1,-1);
877  if ( fPrimitives ) {
878  TIter next(fPrimitives);
879  TLegendEntry *entry;
880  while (( entry = (TLegendEntry *)next() )) entry->SaveEntry(out,"leg");
881  }
882  out << " leg->Draw();"<<std::endl;
883 }
884 
885 ////////////////////////////////////////////////////////////////////////////////
886 /// Edit the label of the entry pointed to by the mouse.
887 
888 void TLegend::SetEntryLabel( const char* label )
889 {
890  TLegendEntry* entry = GetEntry(); // get entry pointed by the mouse
891  if ( entry ) entry->SetLabel( label );
892 }
893 
894 ////////////////////////////////////////////////////////////////////////////////
895 /// Edit the option of the entry pointed to by the mouse.
896 
898 {
899  TLegendEntry* entry = GetEntry(); // get entry pointed by the mouse
900  if ( entry ) entry->SetOption( option );
901 }
902 
903 ////////////////////////////////////////////////////////////////////////////////
904 /// Sets the header, which is the "title" that appears at the top of the legend.
905 /// If `option` contains `C`, the title will be centered.
906 
907 void TLegend::SetHeader( const char *header, Option_t* option )
908 {
909  TString opt;
910 
911  if ( !fPrimitives ) fPrimitives = new TList;
912  TIter next(fPrimitives);
913  TLegendEntry *first; // header is always the first entry
914  if (( first = (TLegendEntry*)next() )) {
915  opt = first->GetOption();
916  opt.ToLower();
917  if ( opt.Contains("h") ) {
918  first->SetLabel(header);
919  opt = option;
920  opt.ToLower();
921  if ( opt.Contains("c") ) first->SetTextAlign(22);
922  else first->SetTextAlign(0);
923  return;
924  }
925  }
926  first = new TLegendEntry( 0, header, "h" );
927  opt = option;
928  opt.ToLower();
929  if ( opt.Contains("c") ) first->SetTextAlign(22);
930  else first->SetTextAlign(0);
931  first->SetTextAngle(0);
932  first->SetTextColor(0);
933  first->SetTextFont(GetTextFont()); // default font is TLegend font for the header
934  first->SetTextSize(0);
935  fPrimitives->AddFirst((TObject*)first);
936 }
virtual void Clear(Option_t *option="")
Clear all entries in this legend, including the header.
Definition: TLegend.cxx:340
Style_t GetLegendFont() const
Definition: TStyle.h:206
TBox & operator=(const TBox &)
Assignment operator.
Definition: TBox.cxx:92
virtual void SetLineAttributes()
Invoke the DialogCanvas Line attributes.
Definition: TAttLine.cxx:280
virtual void Delete(Option_t *option="")
Remove all objects from the list AND delete all heap based objects.
Definition: TList.cxx:405
virtual void Print(Option_t *option="") const
Dump this pave with its attributes.
Definition: TPave.cxx:593
virtual void SetFillAttributes()
Invoke the DialogCanvas Fill attributes.
Definition: TAttFill.cxx:249
Int_t fNColumns
Number of columns in the legend.
Definition: TLegend.h:81
This class displays a legend box (TPaveText) containing several legend entries.
Definition: TLegend.h:27
The Histogram stack class.
Definition: THStack.h:35
virtual void SetNDC(Bool_t isNDC=kTRUE)
Set NDC mode on if isNDC = kTRUE, off otherwise.
Definition: TMarker.cxx:306
short Style_t
Definition: RtypesCore.h:76
virtual Color_t GetTextColor() const
Return the text color.
Definition: TAttText.h:40
float Float_t
Definition: RtypesCore.h:53
float Size_t
Definition: RtypesCore.h:83
virtual Short_t GetTextAlign() const
Return the text alignment.
Definition: TAttText.h:38
const char Option_t
Definition: RtypesCore.h:62
virtual void EditEntryAttFill()
Edit the fill attributes for the entry pointed by the mouse.
Definition: TLegend.cxx:381
R__EXTERN TStyle * gStyle
Definition: TStyle.h:418
virtual void Draw(Option_t *option="")
Draw this legend with its current attributes.
Definition: TLegend.cxx:373
Color_t GetLegendFillColor() const
Definition: TStyle.h:205
virtual Float_t GetTextAngle() const
Return the text angle.
Definition: TAttText.h:39
Double_t fY2
Y of 2nd point.
Definition: TBox.h:45
TLegend & operator=(const TLegend &)
Assignment operator.
Definition: TLegend.cxx:245
virtual void EditEntryAttText()
Edit the text attributes for the entry pointed by the mouse.
Definition: TLegend.cxx:414
A TMultiGraph is a collection of TGraph (or derived) objects.
Definition: TMultiGraph.h:37
virtual void AddFirst(TObject *obj)
Add object at the beginning of the list.
Definition: TList.cxx:93
Storage class for one entry of a TLegend.
Definition: TLegendEntry.h:36
virtual void SetTextAttributes()
Invoke the DialogCanvas Text attributes.
Definition: TAttText.cxx:375
#define gROOT
Definition: TROOT.h:364
virtual void EditEntryAttMarker()
Edit the marker attributes for the entry pointed by the mouse.
Definition: TLegend.cxx:403
void Copy(TObject &pave) const
Copy this pave to pave.
Definition: TPave.cxx:167
Basic string class.
Definition: TString.h:137
Manages Markers.
Definition: TMarker.h:31
virtual void SetMarkerAttributes()
Invoke the DialogCanvas Marker attributes.
Definition: TAttMarker.cxx:249
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:170
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1089
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual void SaveEntry(std::ostream &out, const char *name)
Save this TLegendEntry as C++ statements on output stream out to be used with the SaveAs ...
Width_t GetLegendBorderSize() const
Definition: TStyle.h:204
const char * Class
Definition: TXMLSetup.cxx:64
TObject * Next()
Return next object in the list. Returns 0 when no more objects in list.
Definition: TList.cxx:933
TLegend * legend
Definition: pirndm.C:35
virtual void DeleteEntry()
Delete entry at the mouse position.
Definition: TLegend.cxx:361
Iterator of linked list.
Definition: TList.h:187
Double_t GetLegendTextSize() const
Definition: TStyle.h:207
virtual void SetEntryLabel(const char *label)
Edit the label of the entry pointed to by the mouse.
Definition: TLegend.cxx:888
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:739
virtual TObject * FindObject(const char *name) const
Find an object in this list using its name.
Definition: TList.cxx:497
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition: TObject.cxx:165
A TBox with a bordersize and a shadow option.
Definition: TPave.h:23
virtual void ConvertNDCtoPad()
Convert pave coordinates from NDC to Pad coordinates.
Definition: TPave.cxx:121
Double_t fY1
Y of 1st point.
Definition: TBox.h:43
Marker Attributes class.
Definition: TAttMarker.h:24
virtual void SetTextFont(Font_t tfont=62)
Set the text font.
Definition: TAttText.h:51
virtual void Print(Option_t *option="") const
Dump this TLegend and its contents.
Definition: TLegend.cxx:835
static const double x2[5]
virtual void PaintLineNDC(Double_t u1, Double_t v1, Double_t u2, Double_t v2)
Draw this line with new coordinates in NDC.
Definition: TLine.cxx:389
virtual void SetObject(TObject *obj)
(re)set the obj pointed to by this entry
Fill Area Attributes class.
Definition: TAttFill.h:24
Double_t x[n]
Definition: legend1.C:17
Int_t GetNColumns() const
Definition: TLegend.h:55
Double_t fX1NDC
X1 point in NDC coordinates.
Definition: TPave.h:26
virtual void SaveLineAttributes(std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1, Int_t widdef=1)
Save line attributes as C++ statement(s) on output stream out.
Definition: TAttLine.cxx:260
virtual Float_t GetTextSize() const
Return the text size.
Definition: TAttText.h:42
virtual void InsertEntry(const char *objectName="", const char *label="", Option_t *option="lpf")
Add a new entry before the entry at the mouse position.
Definition: TLegend.cxx:483
To draw Mathematical Formula.
Definition: TLatex.h:25
virtual const char * GetHeader() const
Returns the header, which is the title that appears at the top of the legend.
Definition: TLegend.cxx:467
virtual Option_t * GetOption() const
Definition: TLegendEntry.h:46
virtual void Paint(Option_t *option="")
Paint this line with its current attributes.
Definition: TLine.cxx:371
virtual void SetNDC(Bool_t isNDC=kTRUE)
Set NDC mode on if isNDC = kTRUE, off otherwise.
Definition: TText.cxx:813
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save this legend as C++ statements on output stream out to be used with the SaveAs ...
Definition: TLegend.cxx:857
TLegendEntry * GetEntry() const
Get entry pointed to by the mouse.
Definition: TLegend.cxx:426
c SetBorderSize(2)
short Color_t
Definition: RtypesCore.h:79
virtual void SetTextAlign(Short_t align=11)
Set the text alignment.
Definition: TAttText.h:47
virtual void SetLabel(const char *label="")
Definition: TLegendEntry.h:49
A doubly linked list.
Definition: TList.h:47
void Copy(TAttText &atttext) const
Copy this text attributes to a new TAttText.
Definition: TAttText.cxx:294
virtual void Paint(Option_t *option="")
Paint this marker with its current attributes.
Definition: TMarker.cxx:247
virtual void EditEntryAttLine()
Edit the line attributes for the entry pointed by the mouse.
Definition: TLegend.cxx:392
Double_t GetX1NDC() const
Definition: TPave.h:60
void SetNColumns(Int_t nColumns)
Set the number of columns for the legend.
Definition: TLegend.cxx:538
virtual void SaveTextAttributes(std::ostream &out, const char *name, Int_t alidef=12, Float_t angdef=0, Int_t coldef=1, Int_t fondef=61, Float_t sizdef=1)
Save text attributes as C++ statement(s) on output stream out.
Definition: TAttText.cxx:347
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition: TAttFill.h:42
Float_t fColumnSeparation
Separation between columns, as a fraction of The space allowed to one column.
Definition: TLegend.h:82
virtual Font_t GetTextFont() const
Return the text font.
Definition: TAttText.h:41
virtual TObject * Remove(TObject *obj)
Remove object from the list.
Definition: TList.cxx:675
Int_t GetNRows() const
Get the number of rows.
Definition: TLegend.cxx:518
Double_t fX2
X of 2nd point.
Definition: TBox.h:44
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:488
virtual void SetTextAngle(Float_t tangle=0)
Set the text angle.
Definition: TAttText.h:48
Text Attributes class.
Definition: TAttText.h:24
h1 SetFillColor(kGreen)
virtual TObject * GetObject() const
Definition: TLegendEntry.h:45
virtual void SaveFillAttributes(std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1001)
Save fill attributes as C++ statement(s) on output stream out.
Definition: TAttFill.cxx:232
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:925
A simple line.
Definition: TLine.h:33
virtual void SetOption(Option_t *option="lpf")
Definition: TLegendEntry.h:52
Double_t fX1
X of 1st point.
Definition: TBox.h:42
Use NDC coordinates.
Definition: TLine.h:44
short Short_t
Definition: RtypesCore.h:35
TLegend()
Default constructor.
Definition: TLegend.cxx:187
Double_t GetY1NDC() const
Definition: TPave.h:62
virtual void AddBefore(const TObject *before, TObject *obj)
Insert object before object before in the list.
Definition: TList.cxx:173
Int_t fBorderSize
window box bordersize in pixels
Definition: TPave.h:30
Float_t fEntrySeparation
Separation between entries, as a fraction of The space allocated to one entry.
Definition: TLegend.h:77
Double_t GetY2NDC() const
Definition: TPave.h:63
virtual const char * GetLabel() const
Definition: TLegendEntry.h:44
virtual void Paint(Option_t *option="")
This method must be overridden if a class wants to paint itself.
Definition: TObject.cxx:564
Double_t fY2NDC
Y2 point in NDC coordinates.
Definition: TPave.h:29
static const double x1[5]
#define ClassImp(name)
Definition: Rtypes.h:279
double Double_t
Definition: RtypesCore.h:55
TLegendEntry * AddEntry(const TObject *obj, const char *label="", Option_t *option="lpf")
Add a new entry to this legend.
Definition: TLegend.cxx:280
Double_t GetX2NDC() const
Definition: TPave.h:61
Double_t y[n]
Definition: legend1.C:17
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:567
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
Definition: TRolke.cxx:630
virtual void RecursiveRemove(TObject *obj)
Reset the legend entries pointing to "obj".
Definition: TLegend.cxx:844
Double_t fX2NDC
X2 point in NDC coordinates.
Definition: TPave.h:28
Mother of all ROOT objects.
Definition: TObject.h:37
void SetDefaults()
Definition: TLegend.h:64
TList * fPrimitives
List of TLegendEntries.
Definition: TLegend.h:76
Float_t fMargin
Fraction of total width used for symbol.
Definition: TLegend.h:80
virtual const char * GetTitle() const
Returns title of object.
Definition: TObject.cxx:460
virtual ~TLegend()
Default destructor.
Definition: TLegend.cxx:261
virtual void Add(TObject *obj)
Definition: TList.h:81
Double_t Ceil(Double_t x)
Definition: TMath.h:467
#define NULL
Definition: Rtypes.h:82
#define gPad
Definition: TVirtualPad.h:289
virtual void SetTextColor(Color_t tcolor=1)
Set the text color.
Definition: TAttText.h:49
Double_t fY1NDC
Y1 point in NDC coordinates.
Definition: TPave.h:27
virtual void Paint(Option_t *option="")
Paint this legend with its current attributes.
Definition: TLegend.cxx:508
Definition: first.py:1
virtual void PaintPave(Double_t x1, Double_t y1, Double_t x2, Double_t y2, Int_t bordersize=4, Option_t *option="br")
Draw this pave with new coordinates.
Definition: TPave.cxx:294
virtual void SetTextSize(Float_t tsize=1)
Set the text size.
Definition: TAttText.h:52
virtual void Print(Option_t *option="") const
Default print for collections, calls Print(option, 1).
TString fOption
Pave style.
Definition: TPave.h:34
virtual Int_t GetSize() const
Definition: TCollection.h:95
const Bool_t kTRUE
Definition: Rtypes.h:91
virtual void Copy(TObject &obj) const
Copy this legend into "obj".
Definition: TLegend.cxx:349
Line Attributes class.
Definition: TAttLine.h:24
char name[80]
Definition: TGX11.cxx:109
virtual void SetBorderSize(Int_t bordersize=4)
Definition: TPave.h:74
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:911
virtual void PaintPrimitives()
Paint the entries (list of primitives) for this legend.
Definition: TLegend.cxx:550
virtual void SetEntryOption(Option_t *option)
Edit the option of the entry pointed to by the mouse.
Definition: TLegend.cxx:897
Int_t GetBorderSize() const
Definition: TPave.h:55
virtual void SetHeader(const char *header="", Option_t *option="")
Sets the header, which is the "title" that appears at the top of the legend.
Definition: TLegend.cxx:907