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