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