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