source: trunk/MagicSoft/Mars/mhist/MH.cc@ 1863

Last change on this file since 1863 was 1780, checked in by tbretz, 22 years ago
*** empty log message ***
File size: 21.8 KB
Line 
1/* ======================================================================== *\
2!
3! *
4! * This file is part of MARS, the MAGIC Analysis and Reconstruction
5! * Software. It is distributed to you in the hope that it can be a useful
6! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
7! * It is distributed WITHOUT ANY WARRANTY.
8! *
9! * Permission to use, copy, modify and distribute this software and its
10! * documentation for any purpose is hereby granted without fee,
11! * provided that the above copyright notice appear in all copies and
12! * that both that copyright notice and this permission notice appear
13! * in supporting documentation. It is provided "as is" without express
14! * or implied warranty.
15! *
16!
17!
18! Author(s): Thomas Bretz 07/2001 <mailto:tbretz@astro.uni-wuerzburg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2002
21!
22!
23\* ======================================================================== */
24
25//////////////////////////////////////////////////////////////////////////////
26// //
27// MH //
28// //
29// This is a base tasks for mars histograms. It defines a common interface //
30// for filling the histograms with events (MH::Fill) which is used by a //
31// common 'filler' And a SetupFill member function which may be used //
32// by MFillH. The idea is: //
33// 1) If your Histogram can become filled by one single container //
34// (like MHHillas) you overload MH::Fill and it gets called with //
35// a pointer to the container with which it should be filled. //
36// //
37// 2) You histogram needs several containers to get filled. Than you //
38// have to overload MH::SetupFill and get the necessary objects from //
39// the parameter list. Use this objects in Fill to fill your //
40// histogram. //
41// //
42// If you want to create your own histogram class the new class must be //
43// derived from MH (instead of the base MParContainer) and you must //
44// the fill function of MH. This is the function which is called to fill //
45// the histogram(s) by the data of a corresponding parameter container. //
46// //
47// Remark: the static member function (eg MakeDefCanvas) can be called //
48// from everywhere using: MH::MakeDefCanvas(...) //
49// //
50//////////////////////////////////////////////////////////////////////////////
51
52#include "MH.h"
53
54#include <TH1.h>
55#include <TH2.h>
56#include <TH3.h>
57#include <TGaxis.h>
58#include <TCanvas.h>
59#include <TLegend.h>
60#include <TPaveStats.h>
61
62#include "MLog.h"
63#include "MLogManip.h"
64
65#include "MParList.h"
66#include "MParContainer.h"
67
68#include "MBinning.h"
69
70ClassImp(MH);
71
72// --------------------------------------------------------------------------
73//
74// Default Constructor. It sets name and title only. Typically you won't
75// need to change this.
76//
77MH::MH(const char *name, const char *title)
78{
79 //
80 // set the name and title of this object
81 //
82 fName = name ? name : "MH";
83 fTitle = title ? title : "Base class for Mars histograms";
84}
85
86// --------------------------------------------------------------------------
87//
88// If you want to use the automatic filling of your derived class you
89// must overload this function. If it is not overloaded you cannot use
90// FillH with this class. The argument is a pointer to a container
91// in your paremeter list which is specified in the MFillH constructor.
92// If you are not going to use it you should at least add
93// Bool_t MH::Fill(const MParContainer *) { return kTRUE; }
94// to your class definition.
95//
96Bool_t MH::Fill(const MParContainer *par)
97{
98 *fLog << warn << GetDescriptor() << ": Fill not overloaded! Can't be used!" << endl;
99 return kFALSE;
100}
101
102// --------------------------------------------------------------------------
103//
104// This virtual function is ment as a generalized interface to retrieve
105// a pointer to a root histogram from the MH-derived class.
106//
107TH1 *MH::GetHistByName(const TString name)
108{
109 *fLog << warn << GetDescriptor() << ": GetHistByName not overloaded! Can't be used!" << endl;
110 return NULL;
111}
112
113// --------------------------------------------------------------------------
114//
115// This is a function which should replace the creation of default
116// canvases like root does. Because this is inconvinient in some aspects.
117// need to change this.
118// You can specify a name for the default canvas and a title. Also
119// width and height can be given.
120// MakeDefCanvas looks for a canvas with the given name. If now name is
121// given the DefCanvasName of root is used. If no such canvas is existing
122// it is created and returned. If such a canvas already exists a new canvas
123// with a name plus anumber is created (the number is calculated by the
124// number of all existing canvases plus one)
125//
126TCanvas *MH::MakeDefCanvas(TString name, const char *title,
127 const UInt_t w, const UInt_t h)
128{
129 const TList *list = (TList*)gROOT->GetListOfCanvases();
130
131 if (name.IsNull())
132 name = gROOT->GetDefCanvasName();
133
134 if (list->FindObject(name))
135 name += Form(" <%d>", list->GetSize()+1);
136
137 return new TCanvas(name, title, w, h);
138}
139
140// --------------------------------------------------------------------------
141//
142// This function works like MakeDefCanvas(name, title, w, h) but name
143// and title are retrieved from the given TObject.
144//
145TCanvas *MH::MakeDefCanvas(const TObject *obj,
146 const UInt_t w, const UInt_t h)
147{
148 return MakeDefCanvas(obj->GetName(), obj->GetTitle(), w, h);
149}
150
151// --------------------------------------------------------------------------
152//
153// Applies a given binning to a 1D-histogram
154//
155void MH::SetBinning(TH1 *h, const MBinning *binsx)
156{
157 //
158 // Another strange behaviour: TAxis::Set deletes the axis title!
159 //
160 TAxis &x = *h->GetXaxis();
161
162#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
163 TString xtitle = x.GetTitle();
164#endif
165
166 //
167 // This is a necessary workaround if one wants to set
168 // non-equidistant bins after the initialization
169 // TH1D::fNcells must be set correctly.
170 //
171 h->SetBins(binsx->GetNumBins(), 0, 1);
172
173 //
174 // Set the binning of the current histogram to the binning
175 // in one of the two given histograms
176 //
177 x.Set(binsx->GetNumBins(), binsx->GetEdges());
178#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
179 x.SetTitle(xtitle);
180#endif
181}
182
183// --------------------------------------------------------------------------
184//
185// Applies given binnings to the two axis of a 2D-histogram
186//
187void MH::SetBinning(TH2 *h, const MBinning *binsx, const MBinning *binsy)
188{
189 TAxis &x = *h->GetXaxis();
190 TAxis &y = *h->GetYaxis();
191
192 //
193 // Another strange behaviour: TAxis::Set deletes the axis title!
194 //
195#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
196 TString xtitle = x.GetTitle();
197 TString ytitle = y.GetTitle();
198#endif
199
200 //
201 // This is a necessary workaround if one wants to set
202 // non-equidistant bins after the initialization
203 // TH1D::fNcells must be set correctly.
204 //
205 h->SetBins(binsx->GetNumBins(), 0, 1,
206 binsy->GetNumBins(), 0, 1);
207
208 //
209 // Set the binning of the current histogram to the binning
210 // in one of the two given histograms
211 //
212 x.Set(binsx->GetNumBins(), binsx->GetEdges());
213 y.Set(binsy->GetNumBins(), binsy->GetEdges());
214#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
215 x.SetTitle(xtitle);
216 y.SetTitle(ytitle);
217#endif
218}
219
220// --------------------------------------------------------------------------
221//
222// Applies given binnings to the three axis of a 3D-histogram
223//
224void MH::SetBinning(TH3 *h, const MBinning *binsx, const MBinning *binsy, const MBinning *binsz)
225{
226 //
227 // Another strange behaviour: TAxis::Set deletes the axis title!
228 //
229 TAxis &x = *h->GetXaxis();
230 TAxis &y = *h->GetYaxis();
231 TAxis &z = *h->GetZaxis();
232
233#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
234 TString xtitle = x.GetTitle();
235 TString ytitle = y.GetTitle();
236 TString ztitle = z.GetTitle();
237#endif
238
239 //
240 // This is a necessary workaround if one wants to set
241 // non-equidistant bins after the initialization
242 // TH1D::fNcells must be set correctly.
243 //
244 h->SetBins(binsx->GetNumBins(), 0, 1,
245 binsy->GetNumBins(), 0, 1,
246 binsz->GetNumBins(), 0, 1);
247
248 //
249 // Set the binning of the current histogram to the binning
250 // in one of the two given histograms
251 //
252 x.Set(binsx->GetNumBins(), binsx->GetEdges());
253 y.Set(binsy->GetNumBins(), binsy->GetEdges());
254 z.Set(binsz->GetNumBins(), binsz->GetEdges());
255#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
256 x.SetTitle(xtitle);
257 y.SetTitle(ytitle);
258 z.SetTitle(ztitle);
259#endif
260}
261
262// --------------------------------------------------------------------------
263//
264// Applies given binning (the n+1 edges) to the axis of a 1D-histogram
265//
266void MH::SetBinning(TH1 *h, const TArrayD &binsx)
267{
268 MBinning bx;
269 bx.SetEdges(binsx);
270 SetBinning(h, &bx);
271}
272
273// --------------------------------------------------------------------------
274//
275// Applies given binning (the n+1 edges) to the two axis of a
276// 2D-histogram
277//
278void MH::SetBinning(TH2 *h, const TArrayD &binsx, const TArrayD &binsy)
279{
280 MBinning bx;
281 MBinning by;
282 bx.SetEdges(binsx);
283 by.SetEdges(binsy);
284 SetBinning(h, &bx, &by);
285}
286
287// --------------------------------------------------------------------------
288//
289// Applies given binning (the n+1 edges) to the three axis of a
290// 3D-histogram
291//
292void MH::SetBinning(TH3 *h, const TArrayD &binsx, const TArrayD &binsy, const TArrayD &binsz)
293{
294 MBinning bx;
295 MBinning by;
296 MBinning bz;
297 bx.SetEdges(binsx);
298 by.SetEdges(binsy);
299 bz.SetEdges(binsz);
300 SetBinning(h, &bx, &by, &bz);
301}
302
303// --------------------------------------------------------------------------
304//
305// Applies the binning of a TAxis (eg from a root histogram) to the axis
306// of a 1D-histogram
307//
308void MH::SetBinning(TH1 *h, const TAxis *binsx)
309{
310 const Int_t nx = binsx->GetNbins();
311
312 TArrayD bx(nx+1);
313 for (int i=0; i<nx; i++) bx[i] = binsx->GetBinLowEdge(i+1);
314 bx[nx] = binsx->GetXmax();
315
316 SetBinning(h, bx);
317}
318
319// --------------------------------------------------------------------------
320//
321// Applies the binnings of the TAxis' (eg from a root histogram) to the
322// two axis' of a 2D-histogram
323//
324void MH::SetBinning(TH2 *h, const TAxis *binsx, const TAxis *binsy)
325{
326 const Int_t nx = binsx->GetNbins();
327 const Int_t ny = binsy->GetNbins();
328
329 TArrayD bx(nx+1);
330 TArrayD by(ny+1);
331 for (int i=0; i<nx; i++) bx[i] = binsx->GetBinLowEdge(i+1);
332 for (int i=0; i<ny; i++) by[i] = binsy->GetBinLowEdge(i+1);
333 bx[nx] = binsx->GetXmax();
334 by[ny] = binsy->GetXmax();
335
336 SetBinning(h, bx, by);
337}
338
339// --------------------------------------------------------------------------
340//
341// Applies the binnings of the TAxis' (eg from a root histogram) to the
342// three axis' of a 3D-histogram
343//
344void MH::SetBinning(TH3 *h, const TAxis *binsx, const TAxis *binsy, const TAxis *binsz)
345{
346 const Int_t nx = binsx->GetNbins();
347 const Int_t ny = binsy->GetNbins();
348 const Int_t nz = binsz->GetNbins();
349
350 TArrayD bx(nx+1);
351 TArrayD by(ny+1);
352 TArrayD bz(nz+1);
353 for (int i=0; i<nx; i++) bx[i] = binsx->GetBinLowEdge(i+1);
354 for (int i=0; i<ny; i++) by[i] = binsy->GetBinLowEdge(i+1);
355 for (int i=0; i<nz; i++) bz[i] = binsz->GetBinLowEdge(i+1);
356 bx[nx] = binsx->GetXmax();
357 by[ny] = binsy->GetXmax();
358 bz[nz] = binsz->GetXmax();
359
360 SetBinning(h, bx, by, bz);
361}
362
363// --------------------------------------------------------------------------
364//
365// Applies the binnings of one root-histogram x to another one h
366// Both histograms must be of the same type: TH1, TH2 or TH3
367//
368void MH::SetBinning(TH1 *h, TH1 *x)
369{
370 if (h->InheritsFrom(TH3::Class()) && x->InheritsFrom(TH3::Class()))
371 {
372 SetBinning((TH3*)h, x->GetXaxis(), x->GetYaxis(), x->GetZaxis());
373 return;
374 }
375 if (h->InheritsFrom(TH3::Class()) || x->InheritsFrom(TH3::Class()))
376 return;
377 if (h->InheritsFrom(TH2::Class()) && x->InheritsFrom(TH2::Class()))
378 {
379 SetBinning((TH2*)h, x->GetXaxis(), x->GetYaxis());
380 return;
381 }
382 if (h->InheritsFrom(TH2::Class()) || x->InheritsFrom(TH2::Class()))
383 return;
384 if (h->InheritsFrom(TH1::Class()) && x->InheritsFrom(TH1::Class()))
385 {
386 SetBinning(h, x->GetXaxis());
387 return;
388 }
389}
390
391// --------------------------------------------------------------------------
392//
393// Multiplies all entries in a TArrayD by a float f
394//
395void MH::ScaleArray(TArrayD &bins, Double_t f)
396{
397 for (int i=0; i<bins.GetSize(); i++)
398 bins[i] *= f;
399}
400
401// --------------------------------------------------------------------------
402//
403// Scales the binning of a TAxis by a float f
404//
405TArrayD MH::ScaleAxis(TAxis &axe, Double_t f)
406{
407 TArrayD arr(axe.GetNbins()+1);
408
409 for (int i=1; i<=axe.GetNbins()+1; i++)
410 arr[i-1] = axe.GetBinLowEdge(i);
411
412 ScaleArray(arr, f);
413
414 return arr;
415}
416
417// --------------------------------------------------------------------------
418//
419// Scales the binning of one, two or three axis of a histogram by a float f
420//
421void MH::ScaleAxis(TH1 *h, Double_t fx, Double_t fy, Double_t fz)
422{
423 if (h->InheritsFrom(TH3::Class()))
424 {
425 SetBinning((TH3*)h,
426 ScaleAxis(*h->GetXaxis(), fx),
427 ScaleAxis(*h->GetYaxis(), fy),
428 ScaleAxis(*h->GetZaxis(), fz));
429 return;
430 }
431
432 if (h->InheritsFrom(TH2::Class()))
433 {
434 SetBinning((TH2*)h,
435 ScaleAxis(*h->GetXaxis(), fx),
436 ScaleAxis(*h->GetYaxis(), fy));
437 return;
438 }
439
440 if (h->InheritsFrom(TH1::Class()))
441 SetBinning(h, ScaleAxis(*h->GetXaxis(), fx));
442}
443
444// --------------------------------------------------------------------------
445//
446// Tries to find a MBinning container with the name "Binning"+name
447// in the given parameter list. If it was found it is applied to the
448// given histogram. This is only valid for 1D-histograms
449//
450Bool_t MH::ApplyBinning(const MParList &plist, TString name, TH1 *h)
451{
452 if (h->InheritsFrom(TH2::Class()) || h->InheritsFrom(TH3::Class()))
453 {
454 gLog << warn << "MH::ApplyBinning: '" << h->GetName() << "' is not a basic TH1 object... no binning applied." << endl;
455 return kFALSE;
456 }
457
458 const MBinning *bins = (MBinning*)plist.FindObject("Binning"+name);
459 if (!bins)
460 {
461 gLog << warn << "Object 'Binning" << name << "' [MBinning] not found... no binning applied." << endl;
462 return kFALSE;
463 }
464
465 SetBinning(h, bins);
466 return kTRUE;
467}
468
469void MH::FindGoodLimits(Int_t nbins, Int_t &newbins, Double_t &xmin, Double_t &xmax, Bool_t isInteger)
470{
471//*-*-*-*-*-*-*-*-*Find reasonable bin values*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
472//*-* ==========================
473
474 Double_t dx = 0.1*(xmax-xmin);
475 Double_t umin = xmin - dx;
476 Double_t umax = xmax + dx;
477
478 if (umin < 0 && xmin >= 0)
479 umin = 0;
480
481 if (umax > 0 && xmax <= 0)
482 umax = 0;
483
484 Double_t binlow =0;
485 Double_t binhigh =0;
486 Double_t binwidth=0;
487
488#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,01)
489 TGaxis::Optimize(umin, umax, nbins, binlow, binhigh, nbins, binwidth, "");
490#else
491 gLog << all << "*********************************************************" << endl;
492 gLog << all << "ERROR, because Abelardo removed TGaxis::Optimize from" << endl;
493 gLog << all << " MH::FindGoodLimits to be able to compile Mars with newer" << endl;
494 gLog << all << " root versions, whatever you are trying to do will fail!!!" << endl;
495 gLog << all << "*********************************************************" << endl;
496#endif
497
498 if (binwidth <= 0 || binwidth > 1.e+39)
499 {
500 xmin = -1;
501 xmax = 1;
502 }
503 else
504 {
505 xmin = binlow;
506 xmax = binhigh;
507 }
508
509 if (isInteger)
510 {
511 Int_t ixmin = (Int_t)xmin;
512 Int_t ixmax = (Int_t)xmax;
513 Double_t dxmin = (Double_t)ixmin;
514 Double_t dxmax = (Double_t)ixmax;
515
516 xmin = xmin<0 && xmin!=dxmin ? dxmin - 1 : dxmin;
517 xmax = xmax>0 && xmax!=dxmax ? dxmax + 1 : dxmax;
518
519 if (xmin>=xmax)
520 xmax = xmin+1;
521
522 Int_t bw = 1 + (Int_t)((xmax-xmin)/nbins);
523
524 nbins = (Int_t)((xmax-xmin)/bw);
525
526 if (xmin+nbins*bw < xmax)
527 {
528 nbins++;
529 xmax = xmin +nbins*bw;
530 }
531 }
532
533 newbins = nbins;
534}
535
536// --------------------------------------------------------------------------
537//
538// Returns the lowest entry in a histogram which is greater than gt (eg >0)
539//
540Double_t MH::GetMinimumGT(const TH1 &h, Double_t gt)
541{
542 Double_t min = FLT_MAX;
543
544 const TAxis &axe = *((TH1&)h).GetXaxis();
545
546 for (int i=1; i<=axe.GetNbins(); i++)
547 {
548 Double_t x = h.GetBinContent(i);
549 if (gt<x && x<min)
550 min = x;
551 }
552 return min;
553}
554
555// --------------------------------------------------------------------------
556//
557// Returns the bin center in a logarithmic scale. If the given bin
558// number is <1 it is set to 1. If it is =GetNbins() it is set to
559// GetNbins()
560//
561Double_t MH::GetBinCenterLog(const TAxis &axe, Int_t nbin)
562{
563 if (nbin>axe.GetNbins())
564 nbin = axe.GetNbins();
565
566 if (nbin<1)
567 nbin = 1;
568
569 const Double_t lo = axe.GetBinLowEdge(nbin);
570 const Double_t hi = axe.GetBinUpEdge(nbin);
571
572 const Double_t val = log10(lo) + log10(hi);
573
574 return pow(10, val/2);
575}
576
577// --------------------------------------------------------------------------
578//
579// Draws a copy of the two given histograms. Uses title as the pad title.
580// Also layout the two statistic boxes and a legend.
581//
582void MH::DrawCopy(const TH1 &hist1, const TH1 &hist2, const TString title)
583{
584 //
585 // Draw first histogram
586 //
587 TH1 *h1 = (TH1*)((TH1&)hist1).DrawCopy();
588 gPad->Update();
589
590 TPaveText *t = (TPaveText*)gPad->FindObject("title");
591 if (t)
592 {
593 t->SetName((TString)"MHTitle"); // rename object
594 t->Clear(); // clear old lines
595 t->AddText((TString)" "+title+" "); // add the new title
596 t->SetBit(kCanDelete); // make sure object is deleted
597
598 //
599 // FIXME: This is a stupid workaround to hide the redrawn
600 // (see THistPainter::PaintTitle) title
601 //
602 gPad->Modified(); // indicates a change
603 gPad->Update(); // recreates the original title
604 t->Pop(); // bring our title on top
605 }
606
607 //
608 // Rename first statistics box
609 //
610 TPaveStats &s1 = *(TPaveStats*)gPad->FindObject("stats");
611 s1.SetX1NDC(s1.GetX1NDC()-0.01);
612 s1.SetName("MHStat");
613
614 //
615 // Draw second histogram
616 //
617 TH1 *h2 = (TH1*)((TH1&)hist2).DrawCopy("sames");
618 gPad->Update();
619
620 //
621 // Set new position of second statistics box
622 //
623 TPaveStats &s2 = *(TPaveStats*)gPad->FindObject("stats");
624 s2.SetX1NDC(s1.GetX1NDC()-(s2.GetX2NDC()-s2.GetX1NDC())-0.01);
625 s2.SetX2NDC(s1.GetX1NDC()-0.01);
626
627 //
628 // Draw Legend
629 //
630 const Int_t n = s1.GetListOfLines()->GetSize();
631 const Double_t h = s1.GetY2NDC()-s1.GetY1NDC();
632 TLegend &l = *new TLegend(s1.GetX1NDC(), s1.GetY1NDC()-0.015-h*2/n,
633 s1.GetX2NDC(), s1.GetY1NDC()-0.01
634 );
635 l.AddEntry(h1, h1->GetTitle());
636 l.AddEntry(h2, h2->GetTitle());
637 l.SetTextSize(s2.GetTextSize());
638 l.SetTextFont(s2.GetTextFont());
639 l.SetBorderSize(s2.GetBorderSize());
640 l.SetBit(kCanDelete);
641
642 l.Draw();
643}
644
645// --------------------------------------------------------------------------
646//
647// Draws the two given histograms. Uses title as the pad title.
648// Also layout the two statistic boxes and a legend.
649//
650void MH::Draw(TH1 &hist1, TH1 &hist2, const TString title)
651{
652 //
653 // Draw first histogram
654 //
655 hist1.Draw();
656 gPad->Update();
657
658 TPaveText *t = (TPaveText*)gPad->FindObject("title");
659 if (t)
660 {
661 t->SetName((TString)"MHTitle"); // rename object
662 t->Clear(); // clear old lines
663 t->AddText((TString)" "+title+" "); // add the new title
664 t->SetBit(kCanDelete); // make sure object is deleted
665
666 //
667 // FIXME: This is a stupid workaround to hide the redrawn
668 // (see THistPainter::PaintTitle) title
669 //
670 gPad->Modified(); // indicates a change
671 gPad->Update(); // recreates the original title
672 t->Pop(); // bring our title on top
673 }
674
675 //
676 // Rename first statistics box
677 //
678 TPaveStats &s1 = *(TPaveStats*)gPad->FindObject("stats");
679 s1.SetX1NDC(s1.GetX1NDC()-0.01);
680 s1.SetName((TString)"Stat"+hist1.GetTitle());
681
682 //
683 // Draw second histogram
684 //
685 hist2.Draw("sames");
686
687 gPad->Update();
688
689 //
690 // Set new position of second statistics box
691 //
692 TPaveStats &s2 = *(TPaveStats*)gPad->FindObject("stats");
693 s2.SetX1NDC(s1.GetX1NDC()-(s2.GetX2NDC()-s2.GetX1NDC())-0.01);
694 s2.SetX2NDC(s1.GetX1NDC()-0.01);
695
696 //
697 // Draw Legend
698 //
699 TLegend &l = *new TLegend(s1.GetX1NDC(),
700 s1.GetY1NDC()-0.015-(s1.GetY2NDC()-s1.GetY1NDC())/2,
701 s1.GetX2NDC(),
702 s1.GetY1NDC()-0.01
703 );
704 l.AddEntry(&hist1, hist1.GetTitle());
705 l.AddEntry(&hist2, hist2.GetTitle());
706 l.SetTextSize(s2.GetTextSize());
707 l.SetTextFont(s2.GetTextFont());
708 l.SetBorderSize(s2.GetBorderSize());
709 l.SetBit(kCanDelete);
710 l.Draw();
711}
712
Note: See TracBrowser for help on using the repository browser.