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