Logo ROOT   6.18/05
Reference Guide
TPie.cxx
Go to the documentation of this file.
1// @(#)root/graf:$Id$
2// Author: Guido Volpi, Olivier Couet 03/11/2006
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 "TPie.h"
13#include "TPieSlice.h"
14
15#include <Riostream.h>
16#include <TROOT.h>
17#include <TVirtualPad.h>
18#include <TArc.h>
19#include <TLegend.h>
20#include <TMath.h>
21#include <TStyle.h>
22#include <TLatex.h>
23#include <TPaveText.h>
24#include <TH1.h>
25#include <TColor.h>
26
28
29/** \class TPie
30\ingroup BasicGraphics
31
32Draw a Pie Chart,
33
34Example:
35
36Begin_Macro(source)
37../../../tutorials/graphics/piechart.C
38End_Macro
39*/
40
41Double_t gX = 0; // Temporary pie X position.
42Double_t gY = 0; // Temporary pie Y position.
43Double_t gRadius = 0; // Temporary pie Radius of the TPie.
44Double_t gRadiusOffset = 0; // Temporary slice's radial offset.
45Double_t gAngularOffset = 0; // Temporary slice's angular offset.
46Bool_t gIsUptSlice = kFALSE; // True if a slice in the TPie should
47 // be updated.
48Int_t gCurrent_slice = -1;// Current slice under mouse.
49Double_t gCurrent_phi1 = 0; // Phimin of the current slice.
50Double_t gCurrent_phi2 = 0; // Phimax of the current slice.
51Double_t gCurrent_rad = 0; // Current distance from the vertex of the
52 // current slice.
53Double_t gCurrent_x = 0; // Current x in the pad metric.
54Double_t gCurrent_y = 0; // Current y in the pad metric.
55Double_t gCurrent_ang = 0; // Current angular, within current_phi1
56 // and current_phi2.
57
58////////////////////////////////////////////////////////////////////////////////
59/// Default constructor.
60
62{
63 Init(1, 0, 0.5, 0.5, 0.4);
64}
65
66////////////////////////////////////////////////////////////////////////////////
67/// This constructor creates a pie chart when only the number of
68/// the slices is known. The number of slices is fixed.
69
70TPie::TPie(const char *name, const char *title, Int_t npoints) :
71 TNamed(name,title)
72{
73 Init(npoints, 0, 0.5, 0.5, 0.4);
74}
75
76////////////////////////////////////////////////////////////////////////////////
77/// Normal constructor. The 1st and 2nd parameters are the name of the object
78/// and its title.
79///
80/// The number of points passed at this point is used to allocate the memory.
81///
82/// Slices values are given as Double_t.
83///
84/// The 4th elements is an array containing, in double precision format,
85/// the value of each slice. It is also possible to specify the filled color
86/// of each slice. If the color array is not specified the slices are colored
87/// using a color sequence in the standard palette.
88
89TPie::TPie(const char *name, const char *title,
90 Int_t npoints, Double_t *vals,
91 Int_t *colors, const char *lbls[]) : TNamed(name,title)
92{
93 Init(npoints, 0, 0.5, 0.5, 0.4);
94 for (Int_t i=0; i<fNvals; ++i) fPieSlices[i]->SetValue(vals[i]);
95
97 SetLabels(lbls);
98}
99
100////////////////////////////////////////////////////////////////////////////////
101/// Normal constructor (Float_t).
102
103TPie::TPie(const char *name,
104 const char *title,
105 Int_t npoints, Float_t *vals,
106 Int_t *colors, const char *lbls[]) : TNamed(name,title)
107{
108 Init(npoints, 0, 0.5, 0.5, 0.4);
109 for (Int_t i=0; i<fNvals; ++i) fPieSlices[i]->SetValue(vals[i]);
110
112 SetLabels(lbls);
113}
114
115////////////////////////////////////////////////////////////////////////////////
116/// Constructor from a TH1
117
118TPie::TPie(const TH1 *h) : TNamed(h->GetName(),h->GetTitle())
119{
120 Int_t i;
121
122 const TAxis *axis = h->GetXaxis();
123 Int_t first = axis->GetFirst();
124 Int_t last = axis->GetLast();
125 Int_t np = last-first+1;
126 Init(np, 0, 0.5, 0.5, 0.4);
127
128 for (i=first; i<=last; ++i) fPieSlices[i-first]->SetValue(h->GetBinContent(i));
129 if (axis->GetLabels()) {
130 for (i=first; i<=last; ++i) fPieSlices[i-first]->SetTitle(axis->GetBinLabel(i));
131 } else {
132 SetLabelFormat("%val");
133 }
134 SetTextSize(axis->GetLabelSize());
136 SetTextFont(axis->GetLabelFont());
137}
138
139////////////////////////////////////////////////////////////////////////////////
140/// Copy constructor.
141
142TPie::TPie(const TPie &cpy) : TNamed(cpy), TAttText(cpy)
143{
144 Init(cpy.fNvals, cpy.fAngularOffset, cpy.fX, cpy.fY, cpy.fRadius);
145
146 for (Int_t i=0;i<fNvals;++i) {
147 fPieSlices[i] = cpy.fPieSlices[i];
148 }
149}
150
151////////////////////////////////////////////////////////////////////////////////
152/// Destructor.
153
155{
156 if (fNvals>0) {
157 for (int i=0; i<fNvals; ++i) delete fPieSlices[i];
158 delete [] fPieSlices;
159 }
160
161 if (fSlices) delete [] fSlices;
162 if (fLegend) delete fLegend;
163}
164
165////////////////////////////////////////////////////////////////////////////////
166/// Evaluate the distance to the chart in gPad.
167
169{
170 Int_t dist = 9999;
171
173 if ( gCurrent_slice>=0 ) {
174 if (gCurrent_rad<=fRadius) {
175 dist = 0;
176 }
177 }
178
179 return dist;
180}
181
182////////////////////////////////////////////////////////////////////////////////
183/// Returns the slice number at the pixel position (px,py).
184/// Returns -1 if no slice is picked.
185///
186/// Used by DistancetoPrimitive.
187
189{
190 MakeSlices();
191
192 Int_t result(-1);
193
194 // coordinates
195 Double_t xx = gPad->AbsPixeltoX(px); //gPad->PadtoX(gPad->AbsPixeltoX(px));
196 Double_t yy = gPad->AbsPixeltoY(py); //gPad->PadtoY(gPad->AbsPixeltoY(py));
197
198 // XY metric
199 Double_t radX = fRadius;
200 Double_t radY = fRadius;
201 Double_t radXY = 1.;
202 if (fIs3D==kTRUE) {
203 radXY = TMath::Sin(fAngle3D/180.*TMath::Pi());
204 radY = radXY*radX;
205 }
206
207 Double_t phimin;
208 Double_t cphi;
209 Double_t phimax;
210
211 Float_t dPxl = (gPad->PixeltoY(0)-gPad->PixeltoY(1))/radY;
212 for (Int_t i=0;i<fNvals;++i) {
214
215 if (gIsUptSlice && gCurrent_slice!=i) continue;
216
217 // Angles' values for this slice
218 phimin = fSlices[2*i ]*TMath::Pi()/180.;
219 cphi = fSlices[2*i+1]*TMath::Pi()/180.;
220 phimax = fSlices[2*i+2]*TMath::Pi()/180.;
221
222 Double_t radOffset = fPieSlices[i]->GetRadiusOffset();
223
224 Double_t dx = (xx-fX-radOffset*TMath::Cos(cphi))/radX;
225 Double_t dy = (yy-fY-radOffset*TMath::Sin(cphi)*radXY)/radY;
226
227 if (TMath::Abs(dy)<dPxl) dy = dPxl;
228
229 Double_t ang = TMath::ATan2(dy,dx);
230 if (ang<0) ang += TMath::TwoPi();
231
232 Double_t dist = TMath::Sqrt(dx*dx+dy*dy);
233
234 if ( ((ang>=phimin && ang <= phimax) || (phimax>TMath::TwoPi() &&
235 ang+TMath::TwoPi()>=phimin && ang+TMath::TwoPi()<phimax)) &&
236 dist<=1.) { // if true the pointer is in the slice region
237
238 gCurrent_x = dx;
239 gCurrent_y = dy;
240 gCurrent_ang = ang;
241 gCurrent_phi1 = phimin;
242 gCurrent_phi2 = phimax;
244
245 if (dist<.95 && dist>.65) {
246 Double_t range = phimax-phimin;
247 Double_t lang = ang-phimin;
248 Double_t rang = phimax-ang;
249 if (lang<0) lang += TMath::TwoPi();
250 else if (lang>=TMath::TwoPi()) lang -= TMath::TwoPi();
251 if (rang<0) rang += TMath::TwoPi();
252 else if (rang>=TMath::TwoPi()) rang -= TMath::TwoPi();
253
254 if (lang/range<.25 || rang/range<.25) {
256 result = -1;
257 }
258 else result = i;
259 } else {
260 result = i;
261 }
262
263 break;
264 }
265 }
266 return result;
267}
268
269////////////////////////////////////////////////////////////////////////////////
270/// Draw the pie chart.
271///
272/// The possible options are listed in the TPie::Paint() method.
273
274void TPie::Draw(Option_t *option)
275{
276 TString soption(option);
277 soption.ToLower();
278
279 if (soption.Length()==0) soption = "l";
280
281 if (gPad) {
282 if (!gPad->IsEditable()) gROOT->MakeDefCanvas();
283 if (!soption.Contains("same")) {
284 gPad->Clear();
285 gPad->Range(0.,0.,1.,1.);
286 }
287 }
288
289 for (Int_t i=0;i<fNvals;++i) fPieSlices[i]->AppendPad();
290 AppendPad(soption.Data());
291}
292
293////////////////////////////////////////////////////////////////////////////////
294/// This method is for internal use. It is used by Execute event to draw the
295/// outline of "this" TPie. Used when the opaque movements are not permitted.
296
298{
299 MakeSlices();
300
301 // XY metric
302 Double_t radXY = 1.;
303 if (fIs3D) {
304 radXY = TMath::Sin(fAngle3D/180.*TMath::Pi());
305 }
306
307 for (Int_t i=0;i<fNvals&&fIs3D==kTRUE;++i) {
308 Float_t minphi = (fSlices[i*2]+gAngularOffset+.5)*TMath::Pi()/180.;
309 Float_t avgphi = (fSlices[i*2+1]+gAngularOffset)*TMath::Pi()/180.;
310 Float_t maxphi = (fSlices[i*2+2]+gAngularOffset-.5)*TMath::Pi()/180.;
311
313 Double_t x0 = gX+radOffset*TMath::Cos(avgphi);
314 Double_t y0 = gY+radOffset*TMath::Sin(avgphi)*radXY-fHeight;
315
316 gVirtualX->DrawLine( gPad->XtoAbsPixel(x0), gPad->YtoAbsPixel(y0),
317 gPad->XtoAbsPixel(x0+gRadius*TMath::Cos(minphi)),
318 gPad->YtoAbsPixel(y0+gRadius*TMath::Sin(minphi)*radXY) );
319
320 Int_t ndiv = 10;
321 Double_t dphi = (maxphi-minphi)/ndiv;
322
323 if (dphi>.15) ndiv = (Int_t) ((maxphi-minphi)/.15);
324 dphi = (maxphi-minphi)/ndiv;
325
326 // Loop to draw the arc
327 for (Int_t j=0;j<ndiv;++j) {
328 Double_t phi = minphi+dphi*j;
329 gVirtualX->DrawLine( gPad->XtoAbsPixel(x0+gRadius*TMath::Cos(phi)),
330 gPad->YtoAbsPixel(y0+gRadius*TMath::Sin(phi)*radXY),
331 gPad->XtoAbsPixel(x0+gRadius*TMath::Cos(phi+dphi)),
332 gPad->YtoAbsPixel(y0+gRadius*TMath::Sin(phi+dphi)*radXY));
333 }
334
335 gVirtualX->DrawLine( gPad->XtoAbsPixel(x0+gRadius*TMath::Cos(maxphi)),
336 gPad->YtoAbsPixel(y0+gRadius*TMath::Sin(maxphi)*radXY),
337 gPad->XtoAbsPixel(x0), gPad->YtoAbsPixel(y0) );
338
339 gVirtualX->DrawLine(gPad->XtoAbsPixel(x0),
340 gPad->YtoAbsPixel(y0),
341 gPad->XtoAbsPixel(x0),
342 gPad->YtoAbsPixel(y0+fHeight));
343 gVirtualX->DrawLine(gPad->XtoAbsPixel(x0+gRadius*TMath::Cos(minphi)),
344 gPad->YtoAbsPixel(y0+gRadius*TMath::Sin(minphi)*radXY),
345 gPad->XtoAbsPixel(x0+gRadius*TMath::Cos(minphi)),
346 gPad->YtoAbsPixel(y0+gRadius*TMath::Sin(minphi)*radXY+fHeight));
347 gVirtualX->DrawLine(gPad->XtoAbsPixel(x0+gRadius*TMath::Cos(maxphi)),
348 gPad->YtoAbsPixel(y0+gRadius*TMath::Sin(maxphi)*radXY),
349 gPad->XtoAbsPixel(x0+gRadius*TMath::Cos(maxphi)),
350 gPad->YtoAbsPixel(y0+gRadius*TMath::Sin(maxphi)*radXY+fHeight));
351 }
352
353
354 // Loop over slices
355 for (Int_t i=0;i<fNvals;++i) {
356 Float_t minphi = (fSlices[i*2]+gAngularOffset+.5)*TMath::Pi()/180.;
357 Float_t avgphi = (fSlices[i*2+1]+gAngularOffset)*TMath::Pi()/180.;
358 Float_t maxphi = (fSlices[i*2+2]+gAngularOffset-.5)*TMath::Pi()/180.;
359
361 Double_t x0 = gX+radOffset*TMath::Cos(avgphi);
362 Double_t y0 = gY+radOffset*TMath::Sin(avgphi)*radXY;
363
364 gVirtualX->DrawLine( gPad->XtoAbsPixel(x0), gPad->YtoAbsPixel(y0),
365 gPad->XtoAbsPixel(x0+gRadius*TMath::Cos(minphi)),
366 gPad->YtoAbsPixel(y0+gRadius*TMath::Sin(minphi)*radXY) );
367
368
369 Int_t ndiv = 10;
370 Double_t dphi = (maxphi-minphi)/ndiv;
371
372 if (dphi>.15) ndiv = (Int_t) ((maxphi-minphi)/.15);
373 dphi = (maxphi-minphi)/ndiv;
374
375 // Loop to draw the arc
376 for (Int_t j=0;j<ndiv;++j) {
377 Double_t phi = minphi+dphi*j;
378 gVirtualX->DrawLine( gPad->XtoAbsPixel(x0+gRadius*TMath::Cos(phi)),
379 gPad->YtoAbsPixel(y0+gRadius*TMath::Sin(phi)*radXY),
380 gPad->XtoAbsPixel(x0+gRadius*TMath::Cos(phi+dphi)),
381 gPad->YtoAbsPixel(y0+gRadius*TMath::Sin(phi+dphi)*radXY));
382 }
383
384 gVirtualX->DrawLine( gPad->XtoAbsPixel(x0+gRadius*TMath::Cos(maxphi)),
385 gPad->YtoAbsPixel(y0+gRadius*TMath::Sin(maxphi)*radXY),
386 gPad->XtoAbsPixel(x0), gPad->YtoAbsPixel(y0) );
387 }
388}
389
390////////////////////////////////////////////////////////////////////////////////
391/// Execute the mouse events.
392
394{
395 if (!gPad) return;
396 if (!gPad->IsEditable() && event != kMouseEnter) return;
397
398 if (gCurrent_slice<=-10) {
399 gPad->SetCursor(kCross);
400 return;
401 }
402
403 MakeSlices();
404
405 static bool isMovingPie(kFALSE);
406 static bool isMovingSlice(kFALSE);
407 static bool isResizing(kFALSE);
408 static bool isRotating(kFALSE);
409 static bool onBorder(kFALSE);
410 bool isRedrawing(kFALSE);
411 static Int_t prev_event(-1);
412 static Int_t oldpx, oldpy;
413
414 // Portion of pie considered as "border"
415 const Double_t dr = gPad->PixeltoX(3);
416 const Double_t minRad = gPad->PixeltoX(10);
417
418 // Angular divisions in radial direction
419 const Double_t angstep1 = 0.5*TMath::PiOver4();
420 const Double_t angstep2 = 1.5*TMath::PiOver4();
421 const Double_t angstep3 = 2.5*TMath::PiOver4();
422 const Double_t angstep4 = 3.5*TMath::PiOver4();
423 const Double_t angstep5 = 4.5*TMath::PiOver4();
424 const Double_t angstep6 = 5.5*TMath::PiOver4();
425 const Double_t angstep7 = 6.5*TMath::PiOver4();
426 const Double_t angstep8 = 7.5*TMath::PiOver4();
427
428 // XY metric
429 Double_t radXY = 1.;
430 if (fIs3D==kTRUE) {
431 radXY = TMath::Sin(fAngle3D/180.*TMath::Pi());
432 }
433
434 Int_t dx, dy;
435 Double_t mdx, mdy;
436
437 switch(event) {
438 case kArrowKeyPress:
439 case kButton1Down:
440 // Change cursor to show pie's movement.
441 gVirtualX->SetLineColor(1);
442 gVirtualX->SetLineWidth(2);
443
444 // Current center and radius.
445 gX = fX;
446 gY = fY;
449 gAngularOffset = 0;
451
452 prev_event = kButton1Down;
453
454 case kMouseMotion:
455 if (gCurrent_rad>=fRadius-2.*dr && gCurrent_rad<=fRadius+dr
456 && !isMovingPie && !isMovingSlice && !isResizing) {
457 if (gCurrent_ang>=angstep8 || gCurrent_ang<angstep1)
458 gPad->SetCursor(kRightSide);
459 else if (gCurrent_ang>=angstep1 && gCurrent_ang<angstep2)
460 gPad->SetCursor(kTopRight);
461 else if (gCurrent_ang>=angstep2 && gCurrent_ang<angstep3)
462 gPad->SetCursor(kTopSide);
463 else if (gCurrent_ang>=angstep3 && gCurrent_ang<angstep4)
464 gPad->SetCursor(kTopLeft);
465 else if (gCurrent_ang>=angstep4 && gCurrent_ang<=angstep5)
466 gPad->SetCursor(kLeftSide);
467 else if (gCurrent_ang>=angstep5 && gCurrent_ang<angstep6)
468 gPad->SetCursor(kBottomLeft);
469 else if (gCurrent_ang>=angstep6 && gCurrent_ang<angstep7)
470 gPad->SetCursor(kBottomSide);
471 else if (gCurrent_ang>=angstep7 && gCurrent_ang<angstep8)
472 gPad->SetCursor(kBottomRight);
473 onBorder = kTRUE;
474 } else {
475 onBorder = kFALSE;
476 if (gCurrent_rad>fRadius*.6) {
477 gPad->SetCursor(kPointer);
478 } else if (gCurrent_rad<=fRadius*.3) {
479 gPad->SetCursor(kHand);
480 } else if (gCurrent_rad<=fRadius*.6 && gCurrent_rad>=fRadius*.3) {
481 gPad->SetCursor(kRotate);
482 }
483 }
484 oldpx = px;
485 oldpy = py;
486 if (isMovingPie || isMovingSlice) gPad->SetCursor(kMove);
487 break;
488
489 case kArrowKeyRelease:
490 case kButton1Motion:
491 if (!isMovingSlice || !isMovingPie || !isResizing || !isRotating) {
492 if (prev_event==kButton1Down) {
493 if (onBorder) {
494 isResizing = kTRUE;
495 } else if (gCurrent_rad>=fRadius*.6 && gCurrent_slice>=0) {
496 isMovingSlice = kTRUE;
497 } else if (gCurrent_rad<=fRadius*.3) {
498 isMovingPie = kTRUE;
499 } else if (gCurrent_rad<fRadius*.6 && gCurrent_rad>fRadius*.3) {
500 isRotating = kTRUE;
501 }
502 }
503 }
504
505 dx = px-oldpx;
506 dy = py-oldpy;
507
508 mdx = gPad->PixeltoX(dx);
509 mdy = gPad->PixeltoY(dy);
510
511 if (isMovingPie || isMovingSlice) {
512 gPad->SetCursor(kMove);
513 if (isMovingSlice) {
514 Float_t avgphi = fSlices[gCurrent_slice*2+1]*TMath::Pi()/180.;
515
516 if (!gPad->OpaqueMoving()) DrawGhost();
517
518 gRadiusOffset += TMath::Cos(avgphi)*mdx +TMath::Sin(avgphi)*mdy/radXY;
519 if (gRadiusOffset<0) gRadiusOffset = .0;
521
522 if (!gPad->OpaqueMoving()) DrawGhost();
523 } else {
524 if (!gPad->OpaqueMoving()) DrawGhost();
525
526 gX += mdx;
527 gY += mdy;
528
529 if (!gPad->OpaqueMoving()) DrawGhost();
530 }
531 } else if (isResizing) {
532 if (!gPad->OpaqueResizing()) DrawGhost();
533
535 if (gRadius+dr1>=minRad) {
536 gRadius += dr1;
537 } else {
538 gRadius = minRad;
539 }
540
541 if (!gPad->OpaqueResizing()) DrawGhost();
542 } else if (isRotating) {
543 if (!gPad->OpaqueMoving()) DrawGhost();
544
545 Double_t xx = gPad->AbsPixeltoX(px);
546 Double_t yy = gPad->AbsPixeltoY(py);
547
548 Double_t dx1 = xx-gX;
549 Double_t dy1 = yy-gY;
550
551 Double_t ang = TMath::ATan2(dy1,dx1);
552 if (ang<0) ang += TMath::TwoPi();
553
555
556 if (!gPad->OpaqueMoving()) DrawGhost();
557 }
558
559 oldpx = px;
560 oldpy = py;
561
562 if ( ((isMovingPie || isMovingSlice || isRotating) && gPad->OpaqueMoving()) ||
563 (isResizing && gPad->OpaqueResizing()) ) {
564 isRedrawing = kTRUE;
565 event = kButton1Up;
566 }
567 else break;
568
569 case kButton1Up:
570 if (!isRedrawing) {
571 prev_event = kButton1Up;
573 }
574
575 if (gROOT->IsEscaped()) {
576 gROOT->SetEscape(kFALSE);
578 isRedrawing = kFALSE;
579 break;
580 }
581
582 fX = gX;
583 fY = gY;
587
588 if (isRedrawing && (isMovingPie || isMovingSlice)) gPad->SetCursor(kMove);
589
590 if (isMovingPie) isMovingPie = kFALSE;
591 if (isMovingSlice) isMovingSlice = kFALSE;
592 if (isResizing) isResizing = kFALSE;
593 if (isRotating) {
594 isRotating = kFALSE;
595 // this is important mainly when OpaqueMoving==kTRUE
597 }
598
599 gPad->Modified(kTRUE);
600
601
602 isRedrawing = kFALSE;
604
605 gVirtualX->SetLineColor(-1);
606 gVirtualX->SetLineWidth(-1);
607
608 break;
609 case kButton1Locate:
610
611 ExecuteEvent(kButton1Down, px, py);
612
613 while (1) {
614 px = py = 0;
615 event = gVirtualX->RequestLocator(1, 1, px, py);
616
618
619 if (event != -1) { // button is released
620 ExecuteEvent(kButton1Up, px, py);
621 return;
622 }
623 }
624 break;
625
626 case kMouseEnter:
627 break;
628
629 default:
630 // unknown event
631 break;
632 }
633}
634
635////////////////////////////////////////////////////////////////////////////////
636/// Returns the label of the entry number "i".
637
639{
640 return GetSlice(i)->GetTitle();
641}
642
643////////////////////////////////////////////////////////////////////////////////
644/// Return the color of the slice number "i".
645
647{
648 return GetSlice(i)->GetFillColor();
649}
650
651////////////////////////////////////////////////////////////////////////////////
652/// Return the style use to fill the slice number "i".
653
655{
656 return GetSlice(i)->GetFillStyle();
657}
658
659////////////////////////////////////////////////////////////////////////////////
660/// Return the line color used to outline thi "i" slice
661
663{
664 return GetSlice(i)->GetLineColor();
665}
666
667////////////////////////////////////////////////////////////////////////////////
668/// Return the style used to outline thi "i" slice
669
671{
672 return GetSlice(i)->GetLineStyle();
673}
674
675////////////////////////////////////////////////////////////////////////////////
676/// Return the line width used to outline thi "i" slice
677
679{
680 return GetSlice(i)->GetLineWidth();
681}
682
683////////////////////////////////////////////////////////////////////////////////
684/// Return the radial offset's value for the slice number "i".
685
687{
688 return GetSlice(i)->GetRadiusOffset();
689}
690
691////////////////////////////////////////////////////////////////////////////////
692/// Return the value associated with the slice number "i".
693
695{
696 return GetSlice(i)->GetValue();
697}
698
699////////////////////////////////////////////////////////////////////////////////
700/// If created before by Paint option or by MakeLegend method return
701/// the pointer to the legend, otherwise return 0;
702
704{
705 return fLegend;
706}
707
708////////////////////////////////////////////////////////////////////////////////
709/// Return the reference to the slice of index 'id'. There are no controls
710/// of memory corruption, be carefull.
711
713{
714 return fPieSlices[id];
715}
716
717////////////////////////////////////////////////////////////////////////////////
718/// Common initialization for all constructors.
719/// This is a private function called to allocate the memory.
720
722{
724
725 fAngularOffset = ao;
726 fX = x;
727 fY = y;
728 fRadius = r;
729 fNvals = np;
730 fSum = 0.;
731 fSlices = 0;
732 fLegend = 0;
733 fHeight = 0.08;
734 fAngle3D = 30;
735 fIs3D = kFALSE;
736
738
740
741 for (Int_t i=0;i<fNvals;++i) {
742 TString tmplbl = "Slice";
743 tmplbl += i;
744 fPieSlices[i] = new TPieSlice(tmplbl.Data(), tmplbl.Data(), this);
750 fPieSlices[i]->SetFillStyle(1001);
751 }
752
753 fLabelFormat = "%txt";
754 fFractionFormat = "%3.2f";
755 fValueFormat = "%4.2f";
756 fPercentFormat = "%3.1f";
757}
758
759////////////////////////////////////////////////////////////////////////////////
760/// This method create a legend that explains the contents
761/// of the slice for this pie-chart.
762///
763/// The parameter passed reppresents the option passed to shown the slices,
764/// see TLegend::AddEntry() for further details.
765///
766/// The pointer of the TLegend is returned.
767
768TLegend* TPie::MakeLegend(Double_t x1, Double_t y1, Double_t x2, Double_t y2, const char *leg_header)
769{
770 if (!fLegend) fLegend = new TLegend(x1,y1,x2,y2,leg_header);
771 else fLegend->Clear();
772
773 for (Int_t i=0;i<fNvals;++i) {
775 }
776
777 if (gPad) fLegend->Draw();
778
779 return fLegend;
780}
781
782////////////////////////////////////////////////////////////////////////////////
783/// Paint a Pie chart in a canvas.
784/// The possible option are:
785///
786/// - "R" Print the labels along the central "R"adius of slices.
787/// - "T" Print the label in a direction "T"angent to circle that describes
788/// the TPie.
789/// - "SC" Paint the the labels with the "S"ame "C"olor as the slices.
790/// - "3D" Draw the pie-chart with a pseudo 3D effect.
791/// - "NOL" No OutLine: Don't draw the slices' outlines, any property over the
792/// slices' line is ignored.
793/// - ">" Sort the slices in increasing order.
794/// - "<" Sort the slices in decreasing order.
795///
796/// After the use of > or < options the internal order of the TPieSlices
797/// is changed.
798///
799/// Other options changing the labels' format are described in
800/// TPie::SetLabelFormat().
801
803{
804 MakeSlices();
805
806 TString soption(option);
807
808 bool optionSame(kFALSE);
809
810 // if true the lines around the slices are drawn, if false not
811 Bool_t optionLine(kTRUE);
812
813 // if true the labels' colors are the same as the slices' colors
814 Bool_t optionSameColor(kFALSE);
815
816 // For the label orientation there are 3 possibilities:
817 // 0: horizontal
818 // 1: radial
819 // 2: tangent
820 Int_t lblor(0);
821
822 // Parse the options
823 Int_t idx;
824 // Paint the TPie in an existing canvas
825 if ( (idx=soption.Index("same"))>=0 ) {
826 optionSame = kTRUE;
827 soption.Remove(idx,4);
828 }
829
830 if ( (idx=soption.Index("nol"))>=0 ) {
831 optionLine = kFALSE;
832 soption.Remove(idx,3);
833 }
834
835 if ( (idx=soption.Index("sc"))>=0 ) {
836 optionSameColor = kTRUE;
837 soption.Remove(idx,2);
838 }
839
840 // check if is active the pseudo-3d
841 if ( (idx=soption.Index("3d"))>=0 ) {
842 fIs3D = kTRUE;
843 soption.Remove(idx,2);
844 } else {
845 fIs3D = kFALSE;
846 }
847
848 // seek if have to draw the labels around the pie chart
849 if ( (idx=soption.Index("t"))>=0 ) {
850 lblor = 2;
851 soption.Remove(idx,1);
852 }
853
854 // Seek if have to paint the labels along the radii
855 if ( (idx=soption.Index("r"))>=0 ) {
856 lblor = 1;
857 soption.Remove(idx,1);
858 }
859
860 // Seeks if has to paint sort the slices in increasing mode
861 if ( (idx=soption.Index(">"))>=0 ) {
863 soption.Remove(idx,1);
864 }
865
866 // Seeks if has to paint sort the slices in decreasing mode
867 if ( (idx=soption.Index("<"))>=0 ) {
869 soption.Remove(idx,1);
870 }
871
872 if (fNvals<=0) {
873 Warning("Paint","No vals");
874 return;
875 }
876
877 if (!fPieSlices) {
878 Warning("Paint","No valid arrays of values");
879 return;
880 }
881
882 // Check if gPad exists and define the drawing range.
883 if (!gPad) return;
884
885 // Objects useful to draw labels and slices
886 TLatex *textlabel = new TLatex();
887 TArc *arc = new TArc();
888 TLine *line = new TLine();
889
890 // XY metric
891 Double_t radX = fRadius;
892 Double_t radY = fRadius;
893 Double_t radXY = 1.;
894
895 if (fIs3D) {
896 radXY = TMath::Sin(fAngle3D/180.*TMath::Pi());
897 radY = fRadius*radXY;
898 }
899
900 // Draw the slices.
901 Int_t pixelHeight = gPad->YtoPixel(0)-gPad->YtoPixel(fHeight);
902 for (Int_t pi=0;pi<pixelHeight&&fIs3D==kTRUE; ++pi) { // loop for pseudo-3d effect
903 for (Int_t i=0;i<fNvals;++i) {
904 // draw the arc
905 // set the color of the next slice
906 if (pi>0) {
907 arc->SetFillStyle(0);
908 arc->SetLineColor(TColor::GetColorDark((fPieSlices[i]->GetFillColor())));
909 } else {
910 arc->SetFillStyle(0);
911 if (optionLine==kTRUE) {
912 arc->SetLineColor(fPieSlices[i]->GetLineColor());
913 arc->SetLineStyle(fPieSlices[i]->GetLineStyle());
914 arc->SetLineWidth(fPieSlices[i]->GetLineWidth());
915 } else {
916 arc->SetLineWidth(1);
917 arc->SetLineColor(TColor::GetColorDark((fPieSlices[i]->GetFillColor())));
918 }
919 }
920 // Paint the slice
921 Float_t aphi = fSlices[2*i+1]*TMath::Pi()/180.;
922
924 Double_t ay = fY+TMath::Sin(aphi)*fPieSlices[i]->GetRadiusOffset()*radXY+gPad->PixeltoY(pixelHeight-pi);
925
926 arc->PaintEllipse(ax, ay, radX, radY, fSlices[2*i],
927 fSlices[2*i+2], 0.);
928
929 if (optionLine==kTRUE) {
930 line->SetLineColor(fPieSlices[i]->GetLineColor());
931 line->SetLineStyle(fPieSlices[i]->GetLineStyle());
932 line->SetLineWidth(fPieSlices[i]->GetLineWidth());
933 line->PaintLine(ax,ay,ax,ay);
934
935 Double_t x0, y0;
936 x0 = ax+radX*TMath::Cos(fSlices[2*i]/180.*TMath::Pi());
937 y0 = ay+radY*TMath::Sin(fSlices[2*i]/180.*TMath::Pi());
938 line->PaintLine(x0,y0,x0,y0);
939
940 x0 = ax+radX*TMath::Cos(fSlices[2*i+2]/180.*TMath::Pi());
941 y0 = ay+radY*TMath::Sin(fSlices[2*i+2]/180.*TMath::Pi());
942 line->PaintLine(x0,y0,x0,y0);
943 }
944 }
945 } // end loop for pseudo-3d effect
946
947 for (Int_t i=0;i<fNvals;++i) { // loop for the piechart
948 // Set the color of the next slice
949 arc->SetFillColor(fPieSlices[i]->GetFillColor());
950 arc->SetFillStyle(fPieSlices[i]->GetFillStyle());
951 if (optionLine==kTRUE) {
952 arc->SetLineColor(fPieSlices[i]->GetLineColor());
953 arc->SetLineStyle(fPieSlices[i]->GetLineStyle());
954 arc->SetLineWidth(fPieSlices[i]->GetLineWidth());
955 } else {
956 arc->SetLineWidth(1);
957 arc->SetLineColor(fPieSlices[i]->GetFillColor());
958 }
959
960 // Paint the slice
961 Float_t aphi = fSlices[2*i+1]*TMath::Pi()/180.;
962
964 Double_t ay = fY+TMath::Sin(aphi)*fPieSlices[i]->GetRadiusOffset()*radXY;
965 arc->PaintEllipse(ax, ay, radX, radY, fSlices[2*i],
966 fSlices[2*i+2], 0.);
967
968 } // end loop to draw the slices
969
970
971 // Set the font
972 textlabel->SetTextFont(GetTextFont());
973 textlabel->SetTextSize(GetTextSize());
974 textlabel->SetTextColor(GetTextColor());
975
976 // Loop to place the labels.
977 for (Int_t i=0;i<fNvals;++i) {
978 Float_t aphi = fSlices[2*i+1]*TMath::Pi()/180.;
979 //aphi = TMath::ATan2(TMath::Sin(aphi)*radXY,TMath::Cos(aphi));
980
981 Float_t label_off = fLabelsOffset;
982
983
984 // Paint the text in the pad
985 TString tmptxt = fLabelFormat;
986
987 tmptxt.ReplaceAll("%txt",fPieSlices[i]->GetTitle());
988 tmptxt.ReplaceAll("%val",Form(fValueFormat.Data(),fPieSlices[i]->GetValue()));
989 tmptxt.ReplaceAll("%frac",Form(fFractionFormat.Data(),fPieSlices[i]->GetValue()/fSum));
990 tmptxt.ReplaceAll("%perc",Form(Form("%s %s",fPercentFormat.Data(),"%s"),(fPieSlices[i]->GetValue()/fSum)*100,"%"));
991
992 textlabel->SetTitle(tmptxt.Data());
993 Double_t h = textlabel->GetYsize();
994 Double_t w = textlabel->GetXsize();
995
996 Float_t lx = fX+(fRadius+fPieSlices[i]->GetRadiusOffset()+label_off)*TMath::Cos(aphi);
997 Float_t ly = fY+(fRadius+fPieSlices[i]->GetRadiusOffset()+label_off)*TMath::Sin(aphi)*radXY;
998
999 Double_t lblang = 0;
1000
1001 if (lblor==1) { // radial direction for the label
1002 aphi = TMath::ATan2(TMath::Sin(aphi)*radXY,TMath::Cos(aphi));
1003 lblang += aphi;
1004 if (lblang<=0) lblang += TMath::TwoPi();
1005 if (lblang>TMath::TwoPi()) lblang-= TMath::TwoPi();
1006
1007 lx += h/2.*TMath::Sin(lblang);
1008 ly -= h/2.*TMath::Cos(lblang);
1009
1010 // This control prevent text up-side
1011 if (lblang>TMath::PiOver2() && lblang<=3.*TMath::PiOver2()) {
1012 lx += w*TMath::Cos(lblang)-h*TMath::Sin(lblang);
1013 ly += w*TMath::Sin(lblang)+h*TMath::Cos(lblang);
1014 lblang -= TMath::Pi();
1015 }
1016 } else if (lblor==2) { // tangential direction of the labels
1017 aphi -=TMath::PiOver2();
1018 aphi = TMath::ATan2(TMath::Sin(aphi)*radXY,TMath::Cos(aphi));
1019 lblang += aphi;//-TMath::PiOver2();
1020 if (lblang<0) lblang+=TMath::TwoPi();
1021
1022 lx -= w/2.*TMath::Cos(lblang);
1023 ly -= w/2.*TMath::Sin(lblang);
1024
1025 if (lblang>TMath::PiOver2() && lblang<3.*TMath::PiOver2()) {
1026 lx += w*TMath::Cos(lblang)-h*TMath::Sin(lblang);
1027 ly += w*TMath::Sin(lblang)+h*TMath::Cos(lblang);
1028 lblang -= TMath::Pi();
1029 }
1030 } else { // horizontal labels (default direction)
1031 aphi = TMath::ATan2(TMath::Sin(aphi)*radXY,TMath::Cos(aphi));
1032 if (aphi>TMath::PiOver2() || aphi<=-TMath::PiOver2()) lx -= w;
1033 if (aphi<0) ly -= h;
1034 }
1035
1036 Float_t rphi = TMath::ATan2((ly-fY)*radXY,lx-fX);
1037 if (rphi < 0 && fIs3D && label_off>=0.)
1038 ly -= fHeight;
1039
1040 if (optionSameColor) textlabel->SetTextColor((fPieSlices[i]->GetFillColor()));
1041 textlabel->PaintLatex(lx,ly,
1042 lblang*180/TMath::Pi()+GetTextAngle(),
1043 GetTextSize(), tmptxt.Data());
1044 }
1045
1046 delete arc;
1047 delete line;
1048 delete textlabel;
1049
1050 if (optionSame) return;
1051
1052 // Draw title
1053 TPaveText *title = 0;
1054 TObject *obj;
1055 if ((obj = gPad->GetListOfPrimitives()->FindObject("title"))) {
1056 title = dynamic_cast<TPaveText*>(obj);
1057 }
1058
1059 // Check the OptTitle option
1060 if (strlen(GetTitle()) == 0 || gStyle->GetOptTitle() <= 0) {
1061 if (title) delete title;
1062 return;
1063 }
1064
1065 // Height and width of the title
1066 Double_t ht = gStyle->GetTitleH();
1067 Double_t wt = gStyle->GetTitleW();
1068 if (ht<=0) ht = 1.1*gStyle->GetTitleFontSize();
1069 if (ht<=0) ht = 0.05; // minum height
1070 if (wt<=0) { // eval the width of the title
1071 TLatex l;
1072 l.SetTextSize(ht);
1073 l.SetTitle(GetTitle());
1074 // adjustment in case the title has several lines (#splitline)
1075 ht = TMath::Max(ht, 1.2*l.GetYsize()/(gPad->GetY2() - gPad->GetY1()));
1076 Double_t wndc = l.GetXsize()/(gPad->GetX2() - gPad->GetX1());
1077 wt = TMath::Min(0.7, 0.02+wndc);
1078 }
1079
1080 if (title) {
1081 TText *t0 = (TText*)title->GetLine(0);
1082 if (t0) {
1083 if (!strcmp(t0->GetTitle(),GetTitle())) return;
1084 t0->SetTitle(GetTitle());
1085 if (wt > 0) title->SetX2NDC(title->GetX1NDC()+wt);
1086 }
1087 return;
1088 }
1089
1090 Int_t talh = gStyle->GetTitleAlign()/10;
1091 if (talh < 1) talh = 1; else if (talh > 3) talh = 3;
1092 Int_t talv = gStyle->GetTitleAlign()%10;
1093 if (talv < 1) talv = 1; else if (talv > 3) talv = 3;
1094 Double_t xpos, ypos;
1095 xpos = gStyle->GetTitleX();
1096 ypos = gStyle->GetTitleY();
1097 if (talh == 2) xpos = xpos-wt/2.;
1098 if (talh == 3) xpos = xpos-wt;
1099 if (talv == 2) ypos = ypos+ht/2.;
1100 if (talv == 1) ypos = ypos+ht;
1101
1102 title = new TPaveText(xpos,ypos-ht,xpos+wt,ypos,"blNDC");
1105 title->SetName("title");
1106
1109 title->SetTextFont(gStyle->GetTitleFont(""));
1110 if (gStyle->GetTitleFont("")%10 > 2)
1112 title->AddText(GetTitle());
1113
1114 title->SetBit(kCanDelete);
1115
1116 title->Draw();
1117 title->Paint();
1118}
1119
1120////////////////////////////////////////////////////////////////////////////////
1121/// Save primitive as a C++ statement(s) on output stream out
1122
1123void TPie::SavePrimitive(std::ostream &out, Option_t *option)
1124{
1125 out << " " << std::endl;
1126 if (gROOT->ClassSaved(TPie::Class())) {
1127 out << " ";
1128 } else {
1129 out << " TPie *";
1130 }
1131 out << GetName() << " = new TPie(\"" << GetName() << "\", \"" << GetTitle()
1132 << "\", " << fNvals << ");" << std::endl;
1133 out << " " << GetName() << "->SetCircle(" << fX << ", " << fY << ", "
1134 << fRadius << ");" << std::endl;
1135 out << " " << GetName() << "->SetValueFormat(\"" << GetValueFormat()
1136 << "\");" << std::endl;
1137 out << " " << GetName() << "->SetLabelFormat(\"" << GetLabelFormat()
1138 << "\");" << std::endl;
1139 out << " " << GetName() << "->SetPercentFormat(\"" << GetPercentFormat()
1140 << "\");" << std::endl;
1141 out << " " << GetName() << "->SetLabelsOffset(" << GetLabelsOffset()
1142 << ");" << std::endl;
1143 out << " " << GetName() << "->SetAngularOffset(" << GetAngularOffset()
1144 << ");" << std::endl;
1145 out << " " << GetName() << "->SetTextAngle(" << GetTextAngle() << ");" << std::endl;
1146 out << " " << GetName() << "->SetTextColor(" << GetTextColor() << ");" << std::endl;
1147 out << " " << GetName() << "->SetTextFont(" << GetTextFont() << ");" << std::endl;
1148 out << " " << GetName() << "->SetTextSize(" << GetTextSize() << ");" << std::endl;
1149
1150
1151 // Save the values for the slices
1152 for (Int_t i=0;i<fNvals;++i) {
1153 out << " " << GetName() << "->GetSlice(" << i << ")->SetTitle(\""
1154 << fPieSlices[i]->GetTitle() << "\");" << std::endl;
1155 out << " " << GetName() << "->GetSlice(" << i << ")->SetValue("
1156 << fPieSlices[i]->GetValue() << ");" << std::endl;
1157 out << " " << GetName() << "->GetSlice(" << i << ")->SetRadiusOffset("
1158 << fPieSlices[i]->GetRadiusOffset() << ");" << std::endl;
1159 out << " " << GetName() << "->GetSlice(" << i << ")->SetFillColor("
1160 << fPieSlices[i]->GetFillColor() << ");" << std::endl;
1161 out << " " << GetName() << "->GetSlice(" << i << ")->SetFillStyle("
1162 << fPieSlices[i]->GetFillStyle() << ");" << std::endl;
1163 out << " " << GetName() << "->GetSlice(" << i << ")->SetLineColor("
1164 << fPieSlices[i]->GetLineColor() << ");" << std::endl;
1165 out << " " << GetName() << "->GetSlice(" << i << ")->SetLineStyle("
1166 << fPieSlices[i]->GetLineStyle() << ");" << std::endl;
1167 out << " " << GetName() << "->GetSlice(" << i << ")->SetLineWidth("
1168 << fPieSlices[i]->GetLineWidth() << ");" << std::endl;
1169 }
1170
1171 out << " " << GetName() << "->Draw(\"" << option << "\");" << std::endl;
1172}
1173
1174////////////////////////////////////////////////////////////////////////////////
1175/// Set the value of for the pseudo 3D view angle, in degree.
1176/// The range of the permitted values is: [0,90]
1177
1179 // check if val is in the permitted range
1180 while (val>360.) val -= 360.;
1181 while (val<0) val += 360.;
1182 if (val>=90 && val<180) val = 180-val;
1183 else if (val>=180 && val<=360) val = 360-val;
1184
1185 fAngle3D = val;
1186}
1187
1188////////////////////////////////////////////////////////////////////////////////
1189/// Set the global angular offset for slices in degree [0,360]
1190
1192{
1193 fAngularOffset = offset;
1194
1195 while (fAngularOffset>=360.) fAngularOffset -= 360.;
1196 while (fAngularOffset<0.) fAngularOffset += 360.;
1197
1199}
1200
1201////////////////////////////////////////////////////////////////////////////////
1202/// Set the coordinates of the circle that describe the pie:
1203/// - the 1st and the 2nd arguments are the x and y center's coordinates.
1204/// - the 3rd value is the pie-chart's radius.
1205///
1206/// All the coordinates are in NDC space.
1207
1209{
1210 fX = x;
1211 fY = y;
1212 fRadius = rad;
1213}
1214
1215////////////////////////////////////////////////////////////////////////////////
1216/// Set slice number "i" label. The first parameter is the index of the slice,
1217/// the other is the label text.
1218
1219void TPie::SetEntryLabel(Int_t i, const char *text)
1220{
1221 // Set the Label of a single slice
1222 if (i>=0 && i<fNvals) fPieSlices[i]->SetTitle(text);
1223}
1224
1225////////////////////////////////////////////////////////////////////////////////
1226/// Set the color for the slice's outline. "i" is the slice number.
1227
1229{
1230 if (i>=0 && i<fNvals) fPieSlices[i]->SetLineColor(color);
1231}
1232
1233////////////////////////////////////////////////////////////////////////////////
1234/// Set the style for the slice's outline. "i" is the slice number.
1235
1237{
1238 if (i>=0 && i<fNvals) fPieSlices[i]->SetLineStyle(style);
1239}
1240
1241////////////////////////////////////////////////////////////////////////////////
1242/// Set the width of the slice's outline. "i" is the slice number.
1243
1245{
1246 if (i>=0 && i<fNvals) fPieSlices[i]->SetLineWidth(width);
1247}
1248
1249////////////////////////////////////////////////////////////////////////////////
1250/// Set the color for the slice "i".
1251
1253{
1254 if (i>=0 && i<fNvals) fPieSlices[i]->SetFillColor(color);
1255}
1256
1257////////////////////////////////////////////////////////////////////////////////
1258/// Set the fill style for the "i" slice
1259
1261{
1262 if (i>=0 && i<fNvals) fPieSlices[i]->SetFillStyle(style);
1263}
1264
1265////////////////////////////////////////////////////////////////////////////////
1266/// Set the distance, in the direction of the radius of the slice.
1267
1269{
1270 if (i>=0 && i<fNvals) fPieSlices[i]->SetRadiusOffset(shift);
1271}
1272
1273////////////////////////////////////////////////////////////////////////////////
1274/// Set the value of a slice
1275
1277{
1278 if (i>=0 && i<fNvals) fPieSlices[i]->SetValue(val);
1279
1281}
1282
1283////////////////////////////////////////////////////////////////////////////////
1284/// Set the fill colors for all the TPie's slices.
1285
1287{
1288 if (!colors) return;
1289 for (Int_t i=0;i<fNvals;++i) fPieSlices[i]->SetFillColor(colors[i]);
1290}
1291
1292////////////////////////////////////////////////////////////////////////////////
1293/// Set the height, in pixel, for the piechart if is drawn using
1294/// the pseudo-3d mode.
1295///
1296/// The default value is 20 pixel.
1297
1299{
1300 fHeight = val;
1301}
1302
1303////////////////////////////////////////////////////////////////////////////////
1304/// This method is used to customize the label format. The format string
1305/// must contain one of these modifiers:
1306///
1307/// - %txt : to print the text label associated with the slice
1308/// - %val : to print the numeric value of the slice
1309/// - %frac : to print the relative fraction of this slice
1310/// - %perc : to print the % of this slice
1311///
1312/// ex. : mypie->SetLabelFormat("%txt (%frac)");
1313
1314void TPie::SetLabelFormat(const char *fmt)
1315{
1316 fLabelFormat = fmt;
1317}
1318
1319////////////////////////////////////////////////////////////////////////////////
1320/// Set numeric format in the label, is used if the label format
1321/// there is the modifier %frac, in this case the value is printed
1322/// using this format.
1323///
1324/// The numeric format use the standard C modifier used in stdio library:
1325/// %f, %2.1$, %e... for further documentation you can use the printf
1326/// man page ("man 3 printf" on linux)
1327///
1328/// Example:
1329/// ~~~ {.cpp}
1330/// mypie->SetLabelFormat("%txt (%frac)");
1331/// mypie->SetFractionFormat("2.1f");
1332/// ~~~
1333
1334void TPie::SetFractionFormat(const char *fmt)
1335{
1336 fFractionFormat = fmt;
1337}
1338
1339////////////////////////////////////////////////////////////////////////////////
1340/// Set the labels for all the slices.
1341
1342void TPie::SetLabels(const char *lbls[])
1343{
1344 if (!lbls) return;
1345 for (Int_t i=0;i<fNvals;++i) fPieSlices[i]->SetTitle(lbls[i]);
1346}
1347
1348////////////////////////////////////////////////////////////////////////////////
1349/// Set the distance between the label end the external line of the TPie.
1350
1352{
1353 fLabelsOffset = labelsoffset;
1354}
1355
1356////////////////////////////////////////////////////////////////////////////////
1357/// Set the numeric format for the percent value of a slice, default: %3.1f
1358
1359void TPie::SetPercentFormat(const char *fmt)
1360{
1361 fPercentFormat = fmt;
1362}
1363
1364////////////////////////////////////////////////////////////////////////////////
1365/// Set the pie chart's radius' value.
1366
1368{
1369 if (rad>0) {
1370 fRadius = rad;
1371 } else {
1372 Warning("SetRadius",
1373 "It's not possible set the radius to a negative value");
1374 }
1375}
1376
1377////////////////////////////////////////////////////////////////////////////////
1378/// Set the numeric format the slices' values.
1379/// Used by %val (see SetLabelFormat()).
1380
1381void TPie::SetValueFormat(const char *fmt)
1382{
1383 fValueFormat = fmt;
1384}
1385
1386////////////////////////////////////////////////////////////////////////////////
1387/// Set X value.
1388
1390{
1391 fX = x;
1392}
1393
1394////////////////////////////////////////////////////////////////////////////////
1395/// Set Y value.
1396
1398{
1399 fY = y;
1400}
1401
1402////////////////////////////////////////////////////////////////////////////////
1403/// Make the slices.
1404/// If they already exist it does nothing unless force=kTRUE.
1405
1407{
1408 if (fSlices && !force) return;
1409
1410 fSum = .0;
1411
1412 for (Int_t i=0;i<fNvals;++i) {
1413 if (fPieSlices[i]->GetValue()<0) {
1414 Warning("MakeSlices",
1415 "Negative values in TPie, absolute value will be used");
1416 fPieSlices[i]->SetValue(-1.*fPieSlices[i]->GetValue());
1417 }
1418 fSum += fPieSlices[i]->GetValue();
1419 }
1420
1421 if (fSum<=.0) return;
1422
1423 if (!fSlices) fSlices = new Float_t[2*fNvals+1];
1424
1425 // Compute the slices size and position (2 angles for each slice)
1427 for (Int_t i=0;i<fNvals;++i) {
1428 Float_t dphi = fPieSlices[i]->GetValue()/fSum*360.;
1429 fSlices[2*i+1] = fSlices[2*i]+dphi/2.;
1430 fSlices[2*i+2] = fSlices[2*i]+dphi;
1431 }
1432}
1433
1434////////////////////////////////////////////////////////////////////////////////
1435/// This method, mainly intended for internal use, ordered the slices according their values.
1436/// The default (amode=kTRUE) is increasing order, but is also possible in decreasing order (amode=kFALSE).
1437///
1438/// If the merge_threshold>0 the slice that contains a quantity smaller than merge_threshold are merged
1439/// together
1440
1441void TPie::SortSlices(Bool_t amode, Float_t merge_threshold)
1442{
1443
1444 // main loop to order, bubble sort, the array
1445 Bool_t isDone = kFALSE;
1446
1447 while (isDone==kFALSE) {
1448 isDone = kTRUE;
1449
1450 for (Int_t i=0;i<fNvals-1;++i) { // loop over the values
1451 if ( (amode && (fPieSlices[i]->GetValue()>fPieSlices[i+1]->GetValue())) ||
1452 (!amode && (fPieSlices[i]->GetValue()<fPieSlices[i+1]->GetValue()))
1453 )
1454 {
1455 // exchange the order
1456 TPieSlice *tmpcpy = fPieSlices[i];
1457 fPieSlices[i] = fPieSlices[i+1];
1458 fPieSlices[i+1] = tmpcpy;
1459
1460 isDone = kFALSE;
1461 }
1462 } // end loop the values
1463 } // end main ordering loop
1464
1465 if (merge_threshold>0) {
1466 // merge smallest slices
1467 TPieSlice *merged_slice = new TPieSlice("merged","other",this);
1468 merged_slice->SetRadiusOffset(0.);
1469 merged_slice->SetLineColor(1);
1470 merged_slice->SetLineStyle(1);
1471 merged_slice->SetLineWidth(1);
1472 merged_slice->SetFillColor(gStyle->GetColorPalette( (amode ? 0 : fNvals-1) ));
1473 merged_slice->SetFillStyle(1001);
1474
1475 if (amode) {
1476 // search slices under the threshold
1477 Int_t iMerged = 0;
1478 for (;iMerged<fNvals&&fPieSlices[iMerged]->GetValue()<merge_threshold;++iMerged) {
1479 merged_slice->SetValue( merged_slice->GetValue()+fPieSlices[iMerged]->GetValue() );
1480 }
1481
1482 // evaluate number of valid slices
1483 if (iMerged<=1) { // no slices to merge
1484 delete merged_slice;
1485 }
1486 else { // write a new array with the right dimension
1487 Int_t old_fNvals = fNvals;
1488 fNvals = fNvals-iMerged+1;
1489 TPieSlice **new_array = new TPieSlice*[fNvals];
1490 new_array[0] = merged_slice;
1491 for (Int_t i=0;i<old_fNvals;++i) {
1492 if (i<iMerged) delete fPieSlices[i];
1493 else new_array[i-iMerged+1] = fPieSlices[i];
1494 }
1495 delete [] fPieSlices;
1496 fPieSlices = new_array;
1497 }
1498 }
1499 else {
1500 Int_t iMerged = fNvals-1;
1501 for (;iMerged>=0&&fPieSlices[iMerged]->GetValue()<merge_threshold;--iMerged) {
1502 merged_slice->SetValue( merged_slice->GetValue()+fPieSlices[iMerged]->GetValue() );
1503 }
1504
1505 // evaluate number of valid slices
1506 Int_t nMerged = fNvals-1-iMerged;
1507 if (nMerged<=1) { // no slices to merge
1508 delete merged_slice;
1509 }
1510 else { // write a new array with the right dimension
1511 Int_t old_fNvals = fNvals;
1512 fNvals = fNvals-nMerged+1;
1513 TPieSlice **new_array = new TPieSlice*[fNvals];
1514 new_array[fNvals-1] = merged_slice;
1515 for (Int_t i=old_fNvals-1;i>=0;--i) {
1516 if (i>iMerged) delete fPieSlices[i];
1517 else new_array[i-nMerged-1] = fPieSlices[i];
1518 }
1519 delete [] fPieSlices;
1520 fPieSlices = new_array;
1521 }
1522
1523 }
1524 }
1525
1527}
1528
@ kMouseMotion
Definition: Buttons.h:23
@ kArrowKeyRelease
Definition: Buttons.h:21
@ kButton1Motion
Definition: Buttons.h:20
@ kButton1Up
Definition: Buttons.h:19
@ kArrowKeyPress
Definition: Buttons.h:21
@ kButton1Down
Definition: Buttons.h:17
@ kButton1Locate
Definition: Buttons.h:22
@ kMouseEnter
Definition: Buttons.h:23
void Class()
Definition: Class.C:29
ROOT::R::TRInterface & r
Definition: Object.C:4
#define h(i)
Definition: RSha256.hxx:106
static const double x2[5]
static const double x1[5]
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: RtypesCore.h:88
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
float Float_t
Definition: RtypesCore.h:53
const Bool_t kTRUE
Definition: RtypesCore.h:87
const char Option_t
Definition: RtypesCore.h:62
#define ClassImp(name)
Definition: Rtypes.h:365
include TDocParser_001 C image html pict1_TDocParser_001 png width
Definition: TDocParser.cxx:121
XFontStruct * id
Definition: TGX11.cxx:108
char name[80]
Definition: TGX11.cxx:109
Double_t gCurrent_ang
Definition: TPie.cxx:55
Double_t gX
Definition: TPie.cxx:41
Double_t gCurrent_phi1
Definition: TPie.cxx:49
Double_t gCurrent_phi2
Definition: TPie.cxx:50
Double_t gCurrent_y
Definition: TPie.cxx:54
Double_t gRadiusOffset
Definition: TPie.cxx:44
Double_t gAngularOffset
Definition: TPie.cxx:45
Double_t gCurrent_x
Definition: TPie.cxx:53
Double_t gCurrent_rad
Definition: TPie.cxx:51
Double_t gRadius
Definition: TPie.cxx:43
Int_t gCurrent_slice
Definition: TPie.cxx:48
Bool_t gIsUptSlice
Definition: TPie.cxx:46
Double_t gY
Definition: TPie.cxx:42
#define gROOT
Definition: TROOT.h:414
char * Form(const char *fmt,...)
R__EXTERN TStyle * gStyle
Definition: TStyle.h:406
#define gPad
Definition: TVirtualPad.h:286
#define gVirtualX
Definition: TVirtualX.h:345
@ kRightSide
Definition: TVirtualX.h:45
@ kBottomSide
Definition: TVirtualX.h:45
@ kTopLeft
Definition: TVirtualX.h:44
@ kBottomRight
Definition: TVirtualX.h:44
@ kTopSide
Definition: TVirtualX.h:45
@ kLeftSide
Definition: TVirtualX.h:45
@ kMove
Definition: TVirtualX.h:46
@ kTopRight
Definition: TVirtualX.h:44
@ kBottomLeft
Definition: TVirtualX.h:44
@ kHand
Definition: TVirtualX.h:46
@ kCross
Definition: TVirtualX.h:46
@ kRotate
Definition: TVirtualX.h:46
@ kPointer
Definition: TVirtualX.h:47
Color * colors
Definition: X3DBuffer.c:21
Create an Arc.
Definition: TArc.h:28
virtual Color_t GetLabelColor() const
Definition: TAttAxis.h:38
virtual Style_t GetLabelFont() const
Definition: TAttAxis.h:39
virtual Float_t GetLabelSize() const
Definition: TAttAxis.h:41
virtual Color_t GetFillColor() const
Return the fill area color.
Definition: TAttFill.h:30
virtual Style_t GetFillStyle() const
Return the fill area style.
Definition: TAttFill.h:31
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition: TAttFill.h:37
virtual void SetFillStyle(Style_t fstyle)
Set the fill area style.
Definition: TAttFill.h:39
virtual Color_t GetLineColor() const
Return the line color.
Definition: TAttLine.h:33
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition: TAttLine.h:42
virtual Width_t GetLineWidth() const
Return the line width.
Definition: TAttLine.h:35
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition: TAttLine.h:43
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:40
virtual Style_t GetLineStyle() const
Return the line style.
Definition: TAttLine.h:34
Text Attributes class.
Definition: TAttText.h:18
virtual Float_t GetTextSize() const
Return the text size.
Definition: TAttText.h:36
virtual Font_t GetTextFont() const
Return the text font.
Definition: TAttText.h:35
virtual Color_t GetTextColor() const
Return the text color.
Definition: TAttText.h:34
virtual Float_t GetTextAngle() const
Return the text angle.
Definition: TAttText.h:33
virtual void SetTextColor(Color_t tcolor=1)
Set the text color.
Definition: TAttText.h:43
virtual void SetTextFont(Font_t tfont=62)
Set the text font.
Definition: TAttText.h:45
virtual void SetTextSize(Float_t tsize=1)
Set the text size.
Definition: TAttText.h:46
Class to manage histogram axis.
Definition: TAxis.h:30
const char * GetBinLabel(Int_t bin) const
Return label for bin.
Definition: TAxis.cxx:426
Int_t GetLast() const
Return last bin on the axis i.e.
Definition: TAxis.cxx:455
Int_t GetFirst() const
Return first bin on the axis i.e.
Definition: TAxis.cxx:444
THashList * GetLabels() const
Definition: TAxis.h:117
static Int_t GetColorDark(Int_t color)
Static function: Returns the dark color number corresponding to n If the TColor object does not exist...
Definition: TColor.cxx:1935
virtual void PaintEllipse(Double_t x1, Double_t y1, Double_t r1, Double_t r2, Double_t phimin, Double_t phimax, Double_t theta, Option_t *option="")
Draw this ellipse with new coordinates.
Definition: TEllipse.cxx:526
The TH1 histogram class.
Definition: TH1.h:56
To draw Mathematical Formula.
Definition: TLatex.h:18
Double_t GetXsize()
Return size of the formula along X in pad coordinates.
Definition: TLatex.cxx:2510
Double_t GetYsize()
Return size of the formula along Y in pad coordinates.
Definition: TLatex.cxx:2597
virtual void PaintLatex(Double_t x, Double_t y, Double_t angle, Double_t size, const char *text)
Main drawing function.
Definition: TLatex.cxx:2038
This class displays a legend box (TPaveText) containing several legend entries.
Definition: TLegend.h:23
virtual void Clear(Option_t *option="")
Clear all entries in this legend, including the header.
Definition: TLegend.cxx:390
TLegendEntry * AddEntry(const TObject *obj, const char *label="", Option_t *option="lpf")
Add a new entry to this legend.
Definition: TLegend.cxx:330
virtual void Draw(Option_t *option="")
Draw this legend with its current attributes.
Definition: TLegend.cxx:423
A simple line.
Definition: TLine.h:23
virtual void PaintLine(Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Draw this line with new coordinates.
Definition: TLine.cxx:375
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:29
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition: TNamed.cxx:164
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:48
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
Mother of all ROOT objects.
Definition: TObject.h:37
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:866
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition: TObject.cxx:105
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:694
@ kCanDelete
if object in a list can be deleted
Definition: TObject.h:58
A Pave (see TPave) with text, lines or/and boxes inside.
Definition: TPaveText.h:21
virtual TText * AddText(Double_t x1, Double_t y1, const char *label)
Add a new Text line to this pavetext at given coordinates.
Definition: TPaveText.cxx:182
virtual void Draw(Option_t *option="")
Draw this pavetext with its current attributes.
Definition: TPaveText.cxx:233
virtual void Paint(Option_t *option="")
Paint this pavetext with its current attributes.
Definition: TPaveText.cxx:410
virtual TText * GetLine(Int_t number) const
Get Pointer to line number in this pavetext.
Definition: TPaveText.cxx:274
virtual void SetName(const char *name="")
Definition: TPave.h:75
virtual void SetBorderSize(Int_t bordersize=4)
Definition: TPave.h:73
Double_t GetX1NDC() const
Definition: TPave.h:59
virtual void SetX2NDC(Double_t x2)
Definition: TPave.h:79
A slice of a piechart, see the TPie class.
Definition: TPieSlice.h:18
void SetValue(Double_t)
Set the value for this slice.
Definition: TPieSlice.cxx:108
void SetIsActive(Bool_t is)
Definition: TPieSlice.h:37
void SetRadiusOffset(Double_t)
Set the radial offset of this slice.
Definition: TPieSlice.cxx:98
Double_t GetRadiusOffset()
return the value of the offset in radial direction for this slice.
Definition: TPieSlice.cxx:75
Double_t GetValue()
Return the value of this slice.
Definition: TPieSlice.cxx:83
Draw a Pie Chart,.
Definition: TPie.h:23
void SetEntryVal(Int_t, Double_t)
Set the value of a slice.
Definition: TPie.cxx:1276
TPieSlice * GetSlice(Int_t i)
Return the reference to the slice of index 'id'.
Definition: TPie.cxx:712
Float_t fLabelsOffset
LabelsOffset offset of label.
Definition: TPie.h:37
void SetAngularOffset(Double_t)
Set the global angular offset for slices in degree [0,360].
Definition: TPie.cxx:1191
Double_t GetEntryVal(Int_t)
Return the value associated with the slice number "i".
Definition: TPie.cxx:694
TString fLabelFormat
Format format of the slices' label.
Definition: TPie.h:38
Float_t * fSlices
!Subdivisions of the slices
Definition: TPie.h:29
void SetX(Double_t)
Set X value.
Definition: TPie.cxx:1389
Double_t fAngularOffset
Offset angular offset for the first slice.
Definition: TPie.h:36
void SortSlices(Bool_t amode=kTRUE, Float_t merge_thresold=.0)
This method, mainly intended for internal use, ordered the slices according their values.
Definition: TPie.cxx:1441
void SetEntryFillStyle(Int_t, Int_t)
Set the fill style for the "i" slice.
Definition: TPie.cxx:1260
TPie()
Default constructor.
Definition: TPie.cxx:61
void SetFractionFormat(const char *)
Set numeric format in the label, is used if the label format there is the modifier frac,...
Definition: TPie.cxx:1334
Float_t fAngle3D
The angle of the pseudo-3d view.
Definition: TPie.h:46
Double_t fHeight
Height of the slice in pixel.
Definition: TPie.h:45
const char * GetPercentFormat()
Definition: TPie.h:77
~TPie()
Destructor.
Definition: TPie.cxx:154
void SetPercentFormat(const char *)
Set the numeric format for the percent value of a slice, default: %3.1f.
Definition: TPie.cxx:1359
void SetLabels(const char *[])
Set the labels for all the slices.
Definition: TPie.cxx:1342
void SetY(Double_t)
Set Y value.
Definition: TPie.cxx:1397
TLegend * GetLegend()
If created before by Paint option or by MakeLegend method return the pointer to the legend,...
Definition: TPie.cxx:703
Bool_t fIs3D
! true if the pseudo-3d is enabled
Definition: TPie.h:44
void Init(Int_t np, Double_t ao, Double_t x, Double_t y, Double_t r)
Common initialization for all constructors.
Definition: TPie.cxx:721
Int_t GetEntryFillColor(Int_t)
Return the color of the slice number "i".
Definition: TPie.cxx:646
void SetLabelsOffset(Float_t)
Set the distance between the label end the external line of the TPie.
Definition: TPie.cxx:1351
virtual void Draw(Option_t *option="l")
Draw the pie chart.
Definition: TPie.cxx:274
void SetEntryRadiusOffset(Int_t, Double_t)
Set the distance, in the direction of the radius of the slice.
Definition: TPie.cxx:1268
Double_t GetAngularOffset()
Definition: TPie.h:62
void SavePrimitive(std::ostream &out, Option_t *opts="")
Save primitive as a C++ statement(s) on output stream out.
Definition: TPie.cxx:1123
void DrawGhost()
This method is for internal use.
Definition: TPie.cxx:297
Float_t GetLabelsOffset()
Definition: TPie.h:74
TLegend * fLegend
!Legend for this piechart
Definition: TPie.h:30
void SetLabelFormat(const char *)
This method is used to customize the label format.
Definition: TPie.cxx:1314
Double_t fRadius
Radius Pie radius.
Definition: TPie.h:35
TPieSlice ** fPieSlices
[fNvals] Slice array of this pie-chart
Definition: TPie.h:43
void SetEntryLineStyle(Int_t, Int_t)
Set the style for the slice's outline. "i" is the slice number.
Definition: TPie.cxx:1236
void MakeSlices(Bool_t force=kFALSE)
Make the slices.
Definition: TPie.cxx:1406
void SetEntryFillColor(Int_t, Int_t)
Set the color for the slice "i".
Definition: TPie.cxx:1252
void SetEntryLineWidth(Int_t, Int_t)
Set the width of the slice's outline. "i" is the slice number.
Definition: TPie.cxx:1244
const char * GetLabelFormat()
Definition: TPie.h:73
void SetCircle(Double_t x=.5, Double_t y=.5, Double_t rad=.4)
Set the coordinates of the circle that describe the pie:
Definition: TPie.cxx:1208
const char * GetEntryLabel(Int_t)
Returns the label of the entry number "i".
Definition: TPie.cxx:638
TString fValueFormat
Vform numeric format for the value.
Definition: TPie.h:39
TString fFractionFormat
Rform numeric format for the fraction of a slice.
Definition: TPie.h:40
virtual void Paint(Option_t *)
Paint a Pie chart in a canvas.
Definition: TPie.cxx:802
void SetAngle3D(Float_t val=30.)
Set the value of for the pseudo 3D view angle, in degree.
Definition: TPie.cxx:1178
Double_t GetEntryRadiusOffset(Int_t)
Return the radial offset's value for the slice number "i".
Definition: TPie.cxx:686
Int_t GetEntryFillStyle(Int_t)
Return the style use to fill the slice number "i".
Definition: TPie.cxx:654
const char * GetValueFormat()
Definition: TPie.h:80
void SetRadius(Double_t)
Set the pie chart's radius' value.
Definition: TPie.cxx:1367
TString fPercentFormat
Pfrom numeric format for the percent of a slice.
Definition: TPie.h:41
TLegend * MakeLegend(Double_t x1=.65, Double_t y1=.65, Double_t x2=.95, Double_t y2=.95, const char *leg_header="")
This method create a legend that explains the contents of the slice for this pie-chart.
Definition: TPie.cxx:768
virtual void ExecuteEvent(Int_t, Int_t, Int_t)
Execute the mouse events.
Definition: TPie.cxx:393
Int_t fNvals
Number of elements.
Definition: TPie.h:42
Int_t DistancetoSlice(Int_t, Int_t)
Returns the slice number at the pixel position (px,py).
Definition: TPie.cxx:188
Int_t GetEntryLineStyle(Int_t)
Return the style used to outline thi "i" slice.
Definition: TPie.cxx:670
Float_t fSum
!Sum for the slice values
Definition: TPie.h:28
Double_t fY
Y coordinate of the pie centre.
Definition: TPie.h:34
void SetEntryLabel(Int_t, const char *text="Slice")
Set slice number "i" label.
Definition: TPie.cxx:1219
Double_t fX
X coordinate of the pie centre.
Definition: TPie.h:33
void SetHeight(Double_t val=.08)
Set the height, in pixel, for the piechart if is drawn using the pseudo-3d mode.
Definition: TPie.cxx:1298
Int_t GetEntryLineWidth(Int_t)
Return the line width used to outline thi "i" slice.
Definition: TPie.cxx:678
Int_t GetEntryLineColor(Int_t)
Return the line color used to outline thi "i" slice.
Definition: TPie.cxx:662
void SetValueFormat(const char *)
Set the numeric format the slices' values.
Definition: TPie.cxx:1381
void SetEntryLineColor(Int_t, Int_t)
Set the color for the slice's outline. "i" is the slice number.
Definition: TPie.cxx:1228
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Evaluate the distance to the chart in gPad.
Definition: TPie.cxx:168
void SetFillColors(Int_t *)
Set the fill colors for all the TPie's slices.
Definition: TPie.cxx:1286
Basic string class.
Definition: TString.h:131
Ssiz_t Length() const
Definition: TString.h:405
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1125
const char * Data() const
Definition: TString.h:364
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:687
TString & Remove(Ssiz_t pos)
Definition: TString.h:668
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:619
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition: TString.h:634
Float_t GetTitleX() const
Definition: TStyle.h:267
Int_t GetOptTitle() const
Definition: TStyle.h:233
Float_t GetTitleY() const
Definition: TStyle.h:268
Style_t GetTitleFont(Option_t *axis="X") const
Return title font.
Definition: TStyle.cxx:1021
Color_t GetTitleFillColor() const
Definition: TStyle.h:258
Style_t GetTitleStyle() const
Definition: TStyle.h:260
Float_t GetLabelOffset(Option_t *axis="X") const
Return label offset.
Definition: TStyle.cxx:945
Width_t GetTitleBorderSize() const
Definition: TStyle.h:262
Int_t GetColorPalette(Int_t i) const
Return color number i in current palette.
Definition: TStyle.cxx:913
Int_t GetTitleAlign()
Definition: TStyle.h:257
Color_t GetTitleTextColor() const
Definition: TStyle.h:259
Float_t GetTitleH() const
Definition: TStyle.h:270
Float_t GetTitleFontSize() const
Definition: TStyle.h:261
Float_t GetTitleW() const
Definition: TStyle.h:269
Base class for several text objects.
Definition: TText.h:23
TText * text
TLine * line
h1 SetFillColor(kGreen)
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
std::string GetName(const std::string &scope_name)
Definition: Cppyy.cxx:146
double dist(Rotation3D const &r1, Rotation3D const &r2)
Definition: 3DDistances.cxx:48
static constexpr double rad
static constexpr double pi
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:212
constexpr Double_t PiOver2()
Definition: TMath.h:52
Double_t ATan2(Double_t, Double_t)
Definition: TMath.h:667
Double_t Sqrt(Double_t x)
Definition: TMath.h:679
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:180
constexpr Double_t PiOver4()
Definition: TMath.h:59
Double_t Cos(Double_t)
Definition: TMath.h:629
constexpr Double_t Pi()
Definition: TMath.h:38
Double_t Sin(Double_t)
Definition: TMath.h:625
Short_t Abs(Short_t d)
Definition: TMathBase.h:120
constexpr Double_t TwoPi()
Definition: TMath.h:45
Definition: first.py:1
TCanvas * style()
Definition: style.C:1
auto * l
Definition: textangle.C:4