#include <stdio.h>
#include "TStyle.h"
#include "TLatex.h"
#include "TLine.h"
#include "TBox.h"
#include "TMarker.h"
#include "TLegend.h"
#include "TList.h"
#include "TVirtualPad.h"
#include "TMath.h"
#include "TROOT.h"
#include "TLegendEntry.h"
#include "Riostream.h"
#include "TMultiGraph.h"
#include "THStack.h"
ClassImp(TLegend)
TLegend::TLegend(): TPave(), TAttText()
{
fPrimitives = 0;
SetDefaults();
}
TLegend::TLegend( Double_t x1, Double_t y1,Double_t x2, Double_t y2,
const char *header, Option_t *option)
:TPave(x1,y1,x2,y2,4,option), TAttText(12,0,1,gStyle->GetTextFont(),0)
{
fPrimitives = new TList;
if ( header && strlen(header) > 0) {
TLegendEntry *headerEntry = new TLegendEntry( 0, header, "h" );
headerEntry->SetTextAlign(0);
headerEntry->SetTextAngle(0);
headerEntry->SetTextColor(0);
headerEntry->SetTextFont(62);
headerEntry->SetTextSize(0);
fPrimitives->AddFirst(headerEntry);
}
SetDefaults();
SetBorderSize(gStyle->GetLegendBorderSize());
}
TLegend::TLegend( const TLegend &legend ) : TPave(legend), TAttText(legend),
fPrimitives(0)
{
if (legend.fPrimitives) {
fPrimitives = new TList();
TListIter it(legend.fPrimitives);
while (TLegendEntry *e = (TLegendEntry *)it.Next()) {
TLegendEntry *newentry = new TLegendEntry(*e);
fPrimitives->Add(newentry);
}
}
((TLegend&)legend).Copy(*this);
}
TLegend& TLegend::operator=(const TLegend &lg)
{
if(this!=&lg) {
TPave::operator=(lg);
TAttText::operator=(lg);
fPrimitives=lg.fPrimitives;
fEntrySeparation=lg.fEntrySeparation;
fMargin=lg.fMargin;
fNColumns=lg.fNColumns;
}
return *this;
}
TLegend::~TLegend()
{
if (fPrimitives) fPrimitives->Delete();
delete fPrimitives;
fPrimitives = 0;
}
TLegendEntry *TLegend::AddEntry(const TObject *obj, const char *label, Option_t *option)
{
const char *lab = label;
if ((obj && !label) || strlen(label)==0) lab = obj->GetTitle();
TLegendEntry *newentry = new TLegendEntry( obj, lab, option );
if ( !fPrimitives ) fPrimitives = new TList;
fPrimitives->Add(newentry);
return newentry;
}
TLegendEntry *TLegend::AddEntry(const char *name, const char *label, Option_t *option)
{
TObject *obj = gPad->FindObject(name);
if (!obj) {
TList *lop = gPad->GetListOfPrimitives();
if (lop) {
TObject *o=0;
TIter next(lop);
while( (o=next()) ) {
if ( o->InheritsFrom("TMultiGraph" ) ) {
TList * grlist = ((TMultiGraph *)o)->GetListOfGraphs();
obj = grlist->FindObject(name);
if (obj) continue;
}
if ( o->InheritsFrom("THStack" ) ) {
TList * hlist = ((THStack *)o)->GetHists();
obj = hlist->FindObject(name);
if (obj) continue;
}
}
}
}
return AddEntry( obj, label, option );
}
void TLegend::Clear( Option_t *)
{
if (!fPrimitives) return;
fPrimitives->Delete();
}
void TLegend::Copy( TObject &obj ) const
{
TPave::Copy(obj);
TAttText::Copy((TLegend&)obj);
((TLegend&)obj).fEntrySeparation = fEntrySeparation;
((TLegend&)obj).fMargin = fMargin;
((TLegend&)obj).fNColumns = fNColumns;
}
void TLegend::DeleteEntry()
{
if ( !fPrimitives ) return;
TLegendEntry* entry = GetEntry();
if ( !entry ) return;
fPrimitives->Remove(entry);
delete entry;
}
void TLegend::Draw( Option_t *option )
{
AppendPad(option);
}
void TLegend::EditEntryAttFill()
{
TLegendEntry* entry = GetEntry();
if ( !entry ) return;
gROOT->SetSelectedPrimitive( entry );
entry->SetFillAttributes();
}
void TLegend::EditEntryAttLine()
{
TLegendEntry* entry = GetEntry();
if ( !entry ) return;
gROOT->SetSelectedPrimitive( entry );
entry->SetLineAttributes();
}
void TLegend::EditEntryAttMarker()
{
TLegendEntry* entry = GetEntry();
if ( !entry ) return;
gROOT->SetSelectedPrimitive( entry );
entry->SetMarkerAttributes();
}
void TLegend::EditEntryAttText()
{
TLegendEntry* entry = GetEntry();
if ( !entry ) return;
gROOT->SetSelectedPrimitive( entry );
entry->SetTextAttributes();
}
TLegendEntry *TLegend::GetEntry() const
{
Int_t nRows = GetNRows();
if ( nRows == 0 ) return 0;
Double_t ymouse = gPad->AbsPixeltoY(gPad->GetEventY());
Double_t yspace = (fY2 - fY1)/nRows;
Double_t ybottomOfEntry = fY2;
TIter next(fPrimitives);
TLegendEntry *entry;
while (( entry = (TLegendEntry *)next() )) {
ybottomOfEntry -= yspace;
if ( ybottomOfEntry < ymouse ) return entry;
}
return 0;
}
const char *TLegend::GetHeader() const
{
if ( !fPrimitives ) return 0;
TIter next(fPrimitives);
TLegendEntry *first;
if (( first = (TLegendEntry*)next() )) {
TString opt = first->GetOption();
opt.ToLower();
if ( opt.Contains("h") ) return first->GetLabel();
}
return 0;
}
void TLegend::InsertEntry( const char* objectName, const char* label, Option_t* option)
{
TLegendEntry* beforeEntry = GetEntry();
TObject *obj = gPad->FindObject( objectName );
TLegendEntry *newentry = new TLegendEntry( obj, label, option );
if ( !fPrimitives ) fPrimitives = new TList;
if ( beforeEntry ) {
fPrimitives->AddBefore( (TObject*)beforeEntry, (TObject*)newentry );
} else {
fPrimitives->Add((TObject*)newentry);
}
}
void TLegend::Paint( Option_t* option )
{
TPave::ConvertNDCtoPad();
TPave::PaintPave(fX1,fY1,fX2,fY2,GetBorderSize(),option);
PaintPrimitives();
}
Int_t TLegend::GetNRows() const
{
Int_t nEntries = 0;
if ( fPrimitives ) nEntries = fPrimitives->GetSize();
if ( nEntries == 0 ) return 0;
Int_t nRows;
if(GetHeader() != NULL) nRows = 1 + (Int_t) TMath::Ceil((Double_t) (nEntries-1)/fNColumns);
else nRows = (Int_t) TMath::Ceil((Double_t) nEntries/fNColumns);
return nRows;
}
void TLegend::SetNColumns(Int_t nColumns)
{
if(nColumns < 1) {
Warning("TLegend::SetNColumns", "illegal value nColumns = %d; keeping fNColumns = %d", nColumns, fNColumns);
return;
}
fNColumns = nColumns;
}
void TLegend::PaintPrimitives()
{
Int_t nRows = GetNRows();
if ( nRows == 0 ) return;
Double_t x1 = fX1NDC;
Double_t y1 = fY1NDC;
Double_t x2 = fX2NDC;
Double_t y2 = fY2NDC;
Double_t margin = fMargin*( x2-x1 )/fNColumns;
Double_t boxwidth = margin;
Double_t boxw = boxwidth*0.35;
Double_t yspace = (y2-y1)/nRows;
Double_t textsize = GetTextSize();
Double_t save_textsize = textsize;
Double_t* columnWidths = new Double_t[fNColumns];
memset(columnWidths, 0, fNColumns*sizeof(Double_t));
if ( textsize == 0 ) {
textsize = ( 1. - fEntrySeparation ) * yspace;
Double_t maxentrywidth = 0, maxentryheight = 0;
TIter nextsize(fPrimitives);
TLegendEntry *entrysize;
Int_t iColumn = 0;
while (( entrysize = (TLegendEntry *)nextsize() )) {
TLatex entrytex( 0, 0, entrysize->GetLabel() );
entrytex.SetNDC();
Style_t tfont = entrysize->GetTextFont();
if (tfont == 0) tfont = GetTextFont();
entrytex.SetTextFont(tfont);
entrytex.SetTextSize(textsize);
if ( entrytex.GetYsize() > maxentryheight ) {
maxentryheight = entrytex.GetYsize();
}
TString opt = entrysize->GetOption();
opt.ToLower();
if ( opt.Contains("h") ) {
if ( entrytex.GetXsize() > maxentrywidth ) {
maxentrywidth = entrytex.GetXsize();
}
} else {
if ( entrytex.GetXsize() > columnWidths[iColumn] ) {
columnWidths[iColumn] = entrytex.GetXsize();
}
iColumn++;
iColumn %= fNColumns;
}
Double_t tmpMaxWidth = 0.0;
for(int i=0; i<fNColumns; i++) tmpMaxWidth += columnWidths[i];
if ( tmpMaxWidth > maxentrywidth) maxentrywidth = tmpMaxWidth;
}
Double_t tmpsize_h = maxentryheight /(gPad->GetY2() - gPad->GetY1());
textsize = TMath::Min( textsize, tmpsize_h );
Double_t tmpsize_w = textsize*(fX2-fX1)*(1.0-fMargin)/maxentrywidth;
if(fNColumns > 1) tmpsize_w = textsize*(fX2-fX1)*(1.0-fMargin-fColumnSeparation)/maxentrywidth;
textsize = TMath::Min( textsize, tmpsize_w );
SetTextSize( textsize );
}
{
TIter next(fPrimitives);
TLegendEntry *entry;
Int_t iColumn = 0;
memset(columnWidths, 0, fNColumns*sizeof(Double_t));
while (( entry = (TLegendEntry *)next() )) {
TLatex entrytex( 0, 0, entry->GetLabel() );
entrytex.SetNDC();
Style_t tfont = entry->GetTextFont();
if (tfont == 0) tfont = GetTextFont();
entrytex.SetTextFont(tfont);
if(entry->GetTextSize() == 0) entrytex.SetTextSize(textsize);
TString opt = entry->GetOption();
opt.ToLower();
if (!opt.Contains("h")) {
if ( entrytex.GetXsize() > columnWidths[iColumn] ) {
columnWidths[iColumn] = entrytex.GetXsize();
}
iColumn++;
iColumn %= fNColumns;
}
}
double totalWidth = 0.0;
for(int i=0; i<fNColumns; i++) totalWidth += columnWidths[i];
if(fNColumns > 1) totalWidth /= (1.0-fMargin-fColumnSeparation);
else totalWidth /= (1.0 - fMargin);
for(int i=0; i<fNColumns; i++) {
columnWidths[i] = columnWidths[i]/totalWidth*(x2-x1) + margin;
}
}
Double_t ytext = y2 + 0.5*yspace;
TIter next(fPrimitives);
TLegendEntry *entry;
Int_t iColumn = 0;
while (( entry = (TLegendEntry *)next() )) {
if(iColumn == 0) ytext -= yspace;
Short_t talign = entry->GetTextAlign();
Float_t tangle = entry->GetTextAngle();
Color_t tcolor = entry->GetTextColor();
Style_t tfont = entry->GetTextFont();
Size_t tsize = entry->GetTextSize();
if (talign == 0) entry->SetTextAlign(GetTextAlign());
if (tangle == 0) entry->SetTextAngle(GetTextAngle());
if (tcolor == 0) entry->SetTextColor(GetTextColor());
if (tfont == 0) entry->SetTextFont(GetTextFont());
if (tsize == 0) entry->SetTextSize(GetTextSize());
Double_t x=0,y=0;
Int_t halign = entry->GetTextAlign()/10;
Double_t entrymargin = margin;
TString opt = entry->GetOption();
opt.ToLower();
x1 = fX1NDC;
x2 = fX2NDC;
if ( opt.Contains("h") ) entrymargin = margin/10.;
else if (fNColumns > 1) {
for(int i=0; i<iColumn; i++) x1 += columnWidths[i] + fColumnSeparation*(fX2NDC-fX1NDC)/(fNColumns-1);
x2 = x1 + columnWidths[iColumn];
iColumn++;
iColumn %= fNColumns;
}
if (halign == 1) x = x1 + entrymargin;
if (halign == 2) x = 0.5*( (x1+entrymargin) + x2 );
if (halign == 3) x = x2 - entrymargin/10.;
Int_t valign = entry->GetTextAlign()%10;
if (valign == 1) y = ytext - (1. - fEntrySeparation)* yspace/2.;
if (valign == 2) y = ytext;
if (valign == 3) y = ytext + (1. - fEntrySeparation)* yspace/2.;
TLatex entrytex( x, y, entry->GetLabel() );
entrytex.SetNDC();
entry->TAttText::Copy(entrytex);
entrytex.Paint();
entry->SetTextAlign(talign);
entry->SetTextAngle(tangle);
entry->SetTextColor(tcolor);
entry->SetTextFont(tfont);
entry->SetTextSize(tsize);
Double_t xsym = x1 + margin/2.;
Double_t ysym = ytext;
TObject *eobj = entry->GetObject();
if ( opt.Contains("f")) {
if (eobj && eobj->InheritsFrom(TAttFill::Class())) {
dynamic_cast<TAttFill*>(eobj)->Copy(*entry);
}
entry->TAttFill::Modify();
Double_t xf[4],yf[4];
xf[0] = xsym - boxw;
yf[0] = ysym - yspace*0.35;
xf[1] = xsym + boxw;
yf[1] = yf[0];
xf[2] = xf[1];
yf[2] = ysym + yspace*0.35;
xf[3] = xf[0];
yf[3] = yf[2];
for (Int_t i=0;i<4;i++) {
xf[i] = gPad->GetX1() + xf[i]*(gPad->GetX2()-gPad->GetX1());
yf[i] = gPad->GetY1() + yf[i]*(gPad->GetY2()-gPad->GetY1());
}
gPad->PaintFillArea(4,xf,yf);
}
if ( opt.Contains("l") || opt.Contains("f")) {
if (eobj && eobj->InheritsFrom(TAttLine::Class())) {
dynamic_cast<TAttLine*>(eobj)->Copy(*entry);
}
TLine entryline( xsym - boxw, ysym, xsym + boxw, ysym );
entryline.SetBit(TLine::kLineNDC);
entry->TAttLine::Copy(entryline);
if ( opt.Contains("f") && !opt.Contains("l")) {
boxwidth = yspace*
(gPad->GetX2()-gPad->GetX1())/(gPad->GetY2()-gPad->GetY1());
if ( boxwidth > margin ) boxwidth = margin;
entryline.PaintLineNDC( xsym - boxw, ysym + yspace*0.35,
xsym + boxw, ysym + yspace*0.35);
entryline.PaintLineNDC( xsym - boxw, ysym - yspace*0.35,
xsym + boxw, ysym - yspace*0.35);
entryline.PaintLineNDC( xsym + boxw, ysym - yspace*0.35,
xsym + boxw, ysym + yspace*0.35);
entryline.PaintLineNDC( xsym - boxw, ysym - yspace*0.35,
xsym - boxw, ysym + yspace*0.35);
} else {
entryline.Paint();
if (opt.Contains("e")) {
entryline.PaintLineNDC( xsym, ysym - yspace*0.30,
xsym, ysym + yspace*0.30);
}
}
}
if ( opt.Contains("p")) {
if (eobj && eobj->InheritsFrom(TAttMarker::Class())) {
dynamic_cast<TAttMarker*>(eobj)->Copy(*entry);
}
TMarker entrymarker( xsym, ysym, 0 );
entrymarker.SetNDC();
entry->TAttMarker::Copy(entrymarker);
entrymarker.Paint();
}
}
SetTextSize(save_textsize);
}
void TLegend::Print( Option_t* option ) const
{
TPave::Print( option );
if (fPrimitives) fPrimitives->Print();
}
void TLegend::RecursiveRemove(TObject *obj)
{
TIter next(fPrimitives);
TLegendEntry *entry;
while (( entry = (TLegendEntry *)next() )) {
if (entry->GetObject() == obj) entry->SetObject((TObject*)0);
}
}
void TLegend::SavePrimitive(ostream &out, Option_t* )
{
out << " " << endl;
char quote = '"';
if ( gROOT->ClassSaved( TLegend::Class() ) ) {
out << " ";
} else {
out << " TLegend *";
}
out << "leg = new TLegend("<<GetX1NDC()<<","<<GetY1NDC()<<","
<<GetX2NDC()<<","<<GetY2NDC()<<","
<< "NULL" << "," <<quote<< fOption <<quote<<");" << endl;
if (fBorderSize != 4) {
out<<" leg->SetBorderSize("<<fBorderSize<<");"<<endl;
}
SaveTextAttributes(out,"leg",12,0,1,42,0);
SaveLineAttributes(out,"leg",-1,-1,-1);
SaveFillAttributes(out,"leg",-1,-1);
if ( fPrimitives ) {
TIter next(fPrimitives);
TLegendEntry *entry;
while (( entry = (TLegendEntry *)next() )) entry->SaveEntry(out,"leg");
}
out << " leg->Draw();"<<endl;
}
void TLegend::SetEntryLabel( const char* label )
{
TLegendEntry* entry = GetEntry();
if ( entry ) entry->SetLabel( label );
}
void TLegend::SetEntryOption( Option_t* option )
{
TLegendEntry* entry = GetEntry();
if ( entry ) entry->SetOption( option );
}
void TLegend::SetHeader( const char *header )
{
if ( !fPrimitives ) fPrimitives = new TList;
TIter next(fPrimitives);
TLegendEntry *first;
if (( first = (TLegendEntry*)next() )) {
TString opt = first->GetOption();
opt.ToLower();
if ( opt.Contains("h") ) {
first->SetLabel(header);
return;
}
}
first = new TLegendEntry( 0, header, "h" );
first->SetTextAlign(0);
first->SetTextAngle(0);
first->SetTextColor(0);
first->SetTextFont(GetTextFont());
first->SetTextSize(0);
fPrimitives->AddFirst((TObject*)first);
}