ROOT  6.06/09
Reference Guide
TGraphPolargram.cxx
Go to the documentation of this file.
1 // @(#)root/graf:$Id$
2 // Author: Sebastian Boser, Mathieu Demaret 02/02/06
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 /** \class TGraphPolargram
13 \ingroup BasicGraphics
14 
15 To draw polar axis
16 
17 TGraphPolargram draw the polar axis of the TGraphPolar.
18 
19 Example:
20 
21 Begin_Macro(source)
22 {
23  TCanvas * CPol = new TCanvas("CPol","TGraphPolar Examples",500,500);
24 
25  Double_t rmin=0;
26  Double_t rmax=TMath::Pi()*2;
27  Double_t r[1000];
28  Double_t theta[1000];
29 
30  TF1 * fp1 = new TF1("fplot","cos(x)",rmin,rmax);
31  for (Int_t ipt = 0; ipt < 1000; ipt++) {
32  r[ipt] = ipt*(rmax-rmin)/1000+rmin;
33  theta[ipt] = fp1->Eval(r[ipt]);
34  }
35  TGraphPolar * grP1 = new TGraphPolar(1000,r,theta);
36  grP1->SetLineColor(2);
37  grP1->Draw("AOL");
38 
39  return CPol;
40 }
41 End_Macro
42 */
43 
44 #include "TGraphPolar.h"
45 #include "TGraphPolargram.h"
46 #include "TGaxis.h"
47 #include "THLimitsFinder.h"
48 #include "TVirtualPad.h"
49 #include "TROOT.h"
50 #include "TLatex.h"
51 #include "TEllipse.h"
52 #include "TMath.h"
53 
55 
56 ////////////////////////////////////////////////////////////////////////////////
57 /// TGraphPolargram Constructor.
58 
60  Double_t tmin, Double_t tmax):
61  TNamed(name,"Polargram")
62 {
63  Init();
64  fNdivRad = 508;
65  fNdivPol = 508;
67  fRwrmax = rmax;
68  fRwrmin = rmin;
69  fRwtmin = tmin;
70  fRwtmax = tmax;
71 }
72 
73 ////////////////////////////////////////////////////////////////////////////////
74 /// Short constructor used in the case of a spider plot.
75 
77  TNamed(name,"Polargram")
78 {
79  Init();
80  fNdivRad = 0;
81  fNdivPol = 0;
83  fRwrmax = 1;
84  fRwrmin = 0;
85  fRwtmax = 0;
86  fRwtmin = 0;
87 }
88 
89 ////////////////////////////////////////////////////////////////////////////////
90 /// TGraphPolargram destructor.
91 
93 {
94  if (fPolarLabels != NULL) delete [] fPolarLabels;
95 }
96 
97 ////////////////////////////////////////////////////////////////////////////////
98 /// Set the Polar range.
99 /// \param[in] tmin the start number.
100 /// \param[in] tmax the end number.
101 
103 {
104  if (tmin < tmax) {
105  fRwtmin = tmin;
106  fRwtmax = tmax;
107  }
108  if (gPad) gPad->Modified();
109 }
110 
111 ////////////////////////////////////////////////////////////////////////////////
112 /// Everything within the circle belongs to the TGraphPolargram.
113 
115 {
116  Int_t i;
117  Double_t x = gPad->AbsPixeltoX(px);
118  Double_t y = gPad->AbsPixeltoY(py);
119 
120  // Check if close to a (major) radial line.
121  Double_t rad = TMath::Sqrt(x*x+y*y);
122  Int_t div = (Int_t)rad*(fNdivRad%100);
123  Double_t dr = TMath::Min(TMath::Abs(rad-div*1./(fNdivRad%100)),
124  TMath::Abs(rad-(div+1)*1./(fNdivRad%100)));
125  Int_t drad = gPad->XtoPixel(dr)-gPad->XtoPixel(0);
126 
127  // Check if close to a (major) Polar line.
128  // This is not a proper calculation, but rather fast.
129  Int_t dt = kMaxPixel;
130  for (i=0; i<(fNdivPol%100); i++) {
131  Double_t theta = i*2*TMath::Pi()/(fNdivPol%100);
132 
133  // Attention: px,py in pixel units, line given in user coordinates.
134  Int_t dthis = DistancetoLine(px,py,0.,0.,TMath::Cos(theta),
135  TMath::Sin(theta));
136 
137  // Fails if we are outside box described by the line.
138  // (i.e for all hor/vert lines)
139  if (dthis==9999) {
140 
141  // Outside -> Get distance to endpoint of line.
142  if (rad>1) {
143  dthis = (Int_t)TMath::Sqrt(
144  TMath::Power(px-gPad->XtoPixel(TMath::Cos(theta)),2)+
145  TMath::Power(py-gPad->YtoPixel(TMath::Sin(theta)),2));
146  } else {
147 
148  // Check for horizontal line.
149  if (((TMath::Abs(theta-TMath::Pi())<0.1) &&
150  ((px-gPad->XtoPixel(0))<0)) ||
151  ((TMath::Abs(theta)<0.1) &&
152  ((px-gPad->XtoPixel(0))>0))) {
153  dthis = TMath::Abs(py-gPad->YtoPixel(0.));
154  }
155 
156  //Check for vertical line.
157  if (((TMath::Abs(theta-TMath::PiOver2())<0.1) &&
158  ((py-gPad->YtoPixel(0))>0)) ||
159  ((TMath::Abs(theta-3*TMath::PiOver2())<0.1) &&
160  (py-gPad->YtoPixel(0))<0)) {
161  dthis = TMath::Abs(px-gPad->XtoPixel(0.));
162  }
163  if (dthis==9999) {
164 
165  // Inside, but out of box for nonorthogonal line ->
166  // get distance to start point.
167  dthis = (Int_t)TMath::Sqrt(
168  TMath::Power(px-gPad->XtoPixel(0.),2)+
169  TMath::Power(py-gPad->YtoPixel(0.),2));
170  }
171  }
172  }
173 
174  // Take distance to closes line.
175  dt = TMath::Min(dthis,dt);
176  }
177  return TMath::Min(drad, dt);
178 }
179 
180 ////////////////////////////////////////////////////////////////////////////////
181 /// Draw Polargram.
182 
184 {
185  Paint(options);
186  AppendPad(options);
187 }
188 
189 ////////////////////////////////////////////////////////////////////////////////
190 /// Indicate that there is something to click here.
191 
193 {
194  if (!gPad) return;
195 
196  Int_t kMaxDiff = 20;
197  static Int_t d1, d2, d3, px1, py1, px3, py3;
198  static Bool_t p1, p2, p3, p4, p5, p6, p7, p8;
199  Double_t px2, py2;
200  p2 = p3 = p4 = p5 = p6 = p7 = p8 = kFALSE;
201  if (!gPad->IsEditable()) return;
202  switch (event) {
203  case kMouseMotion:
204  px1 = gPad->XtoAbsPixel(TMath::Cos(GetAngle()));
205  py1 = gPad->YtoAbsPixel(TMath::Sin(GetAngle()));
206  d1 = TMath::Abs(px1 - px) + TMath::Abs(py1-py); //simply take sum of pixels differences
207  p1 = kFALSE;
208  px2 = gPad->XtoAbsPixel(-1);
209  py2 = gPad->YtoAbsPixel(1);
210  d2 = (Int_t)(TMath::Abs(px2 - px) + TMath::Abs(py2 - py)) ;
211  px3 = gPad->XtoAbsPixel(-1);
212  py3 = gPad->YtoAbsPixel(-1);
213  d3 = TMath::Abs(px3 - px) + TMath::Abs(py3 - py) ; //simply take sum of pixels differences
214  // check if point is close to the radial axis
215  if (d1 < kMaxDiff) {
216  gPad->SetCursor(kMove);
217  p1 = kTRUE;
218  }
219  // check if point is close to the left high axis
220  if ( d2 < kMaxDiff) {
221  gPad->SetCursor(kHand);
222  p7 = kTRUE;
223  }
224  // check if point is close to the left down axis
225  if ( d3 < kMaxDiff) {
226  gPad->SetCursor(kHand);
227  p8 = kTRUE;
228  }
229  // check if point is close to a main circle
230  if (!p1 && !p7 ) {
231  p6 = kTRUE;
232  gPad->SetCursor(kHand);
233  }
234  break;
235 
236  case kButton1Down:
237  // Record initial coordinates
238  //px4 = px;
239  //py4 = py;
240 
241  case kButton1Motion:
242  if (p1) {
243  px2 = gPad->AbsPixeltoX(px);
244  py2 = gPad->AbsPixeltoY(py);
245  if ( px2 < 0 && py2 < 0) {p2 = kTRUE;};
246  if ( px2 < 0 && py2 > 0 ) {p3 = kTRUE;};
247  if ( px2 > 0 && py2 > 0 ) {p4 = kTRUE;};
248  if ( px2 > 0 && py2 < 0 ) {p5 = kTRUE;};
249  px2 = TMath::ACos(TMath::Abs(px2));
250  py2 = TMath::ASin(TMath::Abs(py2));
251  if (p2) {
252  fAxisAngle = TMath::Pi()+(px2+py2)/2;
253  p2 = kFALSE;
254  };
255  if (p3) {
256  fAxisAngle = TMath::Pi()-(px2+py2)/2;
257  p3 = kFALSE;
258  };
259  if (p4) {
260  fAxisAngle = (px2+py2)/2;
261  p4 = kFALSE;
262  };
263  if (p5) {
264  fAxisAngle = -(px2+py2)/2;
265  p5 = kFALSE;
266  };
267  }
268  break;
269 
270  case kButton1Up:
271  Paint();
272  }
273 }
274 
275 ////////////////////////////////////////////////////////////////////////////////
276 /// Find the alignement rule to apply for TText::SetTextAlign(Short_t).
277 
279 {
280  Double_t pi = TMath::Pi();
281 
282  while(angle < 0 || angle > 2*pi){
283  if(angle < 0) angle+=2*pi;
284  if(angle > 2*pi) angle-=2*pi;
285  }
287  if(angle > 0 && angle < pi/2) return 11;
288  else if(angle > pi/2 && angle < pi) return 31;
289  else if(angle > pi && angle < 3*pi/2) return 33;
290  else if(angle > 3*pi/2 && angle < 2*pi) return 13;
291  else if(angle == 0 || angle == 2*pi) return 12;
292  else if(angle == pi/2) return 21;
293  else if(angle == pi) return 32;
294  else if(angle == 3*pi/2) return 23;
295  else return 0;
296  }
297  else{
298  if(angle >= 0 && angle <= pi/2) return 12;
299  else if((angle > pi/2 && angle <= pi) || (angle > pi && angle <= 3*pi/2)) return 32;
300  else if(angle > 3*pi/2 && angle <= 2*pi) return 12;
301  else return 0;
302  }
303 }
304 
305 ////////////////////////////////////////////////////////////////////////////////
306 /// Determine the orientation of the polar labels according to their angle.
307 
309 {
310  Double_t pi = TMath::Pi();
311  Double_t convraddeg = 180.0/pi;
312 
313  while(angle < 0 || angle > 2*pi){
314  if(angle < 0) angle+=2*pi;
315  if(angle > 2*pi) angle-=2*pi;
316  }
317 
318  if(angle >= 0 && angle <= pi/2) return angle*convraddeg;
319  else if(angle > pi/2 && angle <= pi) return (angle + pi)*convraddeg;
320  else if(angle > pi && angle <= 3*pi/2) return (angle - pi)*convraddeg;
321  else if(angle > 3*pi/2 && angle <= 2*pi) return angle*convraddeg;
322  else return 0;
323 }
324 
325 ////////////////////////////////////////////////////////////////////////////////
326 /// Initialize some of the fields of TGraphPolargram.
327 
329 {
330  fAxisAngle = 0;
331  fCutRadial = 0;
332  fDegree = kFALSE;
333  fGrad = kFALSE;
334  fLineStyle = 3;
335  fPolarLabelColor = 1;
336  fPolarLabelFont = 62;
337  fPolarOffset = 0.04;
338  fPolarTextSize = 0.04;
339  fRadialOffset = 0.025;
340  fRadian = kTRUE;
341  fRadialLabelColor = 1;
342  fRadialLabelFont = 62;
343  fRadialTextSize = 0.035;
344  fTickpolarSize = 0.02;
345 }
346 
347 ////////////////////////////////////////////////////////////////////////////////
348 /// Paint TGraphPolargram.
349 
351 {
352  Int_t optionpoldiv, optionraddiv;
353  Bool_t optionLabels = kTRUE;
354 
355  TString opt = chopt;
356  opt.ToUpper();
357 
358  if(opt.Contains('P')) optionpoldiv=1; else optionpoldiv=0;
359  if(opt.Contains('R')) optionraddiv=1; else optionraddiv=0;
362  if(!opt.Contains('P') && !opt.Contains('R')) optionpoldiv=optionraddiv=1;
363  if(opt.Contains('N')) optionLabels = kFALSE;
364 
365  if(optionraddiv) PaintRadialDivisions(kTRUE);
367  if(optionpoldiv) PaintPolarDivisions(optionLabels);
368 }
369 
370 ////////////////////////////////////////////////////////////////////////////////
371 /// This is simplified from TEllipse::PaintEllipse.
372 /// Draw this ellipse with new coordinates.
373 
375  Double_t phimin, Double_t phimax, Double_t theta)
376 {
377  Int_t i;
378  const Int_t np = 200; // Number of point to draw circle
379  static Double_t x[np+3], y[np+3];
380 
381  // Set number of points approximatively proportional to the ellipse
382  // circumference.
383 
384  Double_t circ = TMath::Pi()*2*r*(phimax-phimin)/36;
385  Int_t n = (Int_t)(np*circ/((gPad->GetX2()-gPad->GetX1())+
386  (gPad->GetY2()-gPad->GetY1())));
387  if (n < 8) n = 8;
388  if (n > np) n = np;
389  Double_t angle,dx,dy;
390  Double_t dphi = (phimax-phimin)*TMath::Pi()/(180*n);
391  Double_t ct = TMath::Cos(TMath::Pi()*theta/180);
392  Double_t st = TMath::Sin(TMath::Pi()*theta/180);
393  for (i=0; i<=n; i++) {
394  angle = phimin*TMath::Pi()/180 + Double_t(i)*dphi;
395  dx = r*TMath::Cos(angle);
396  dy = r*TMath::Sin(angle);
397  x[i] = x1 + dx*ct - dy*st;
398  y[i] = y1 + dx*st + dy*ct;
399  }
400  gPad->PaintPolyLine(n+1,x,y);
401 }
402 
403 ////////////////////////////////////////////////////////////////////////////////
404 /// Draw Polar divisions.
405 /// Check for editable pad or create default.
406 
408 {
409  Int_t i, j, rnum, rden, first, last;
410  if (!gPad) return ;
411 
412  gPad->RangeAxis(-1,-1,1,1);
413  gPad->Range(-1.25,-1.25,1.25,1.25);
414  Int_t ndivMajor = fNdivPol%100;
415  Int_t ndivMinor = fNdivPol/100;
416 
417  if (!gPad->GetLogy()) {
418  for (i=0; i<ndivMajor; i++) {
419  Double_t txtval = fRwtmin + i*(fRwtmax-fRwtmin)/ndivMajor;
420  Double_t theta = i*2*TMath::Pi()/ndivMajor;
421  Double_t costheta = TMath::Cos(theta);
422  Double_t sintheta = TMath::Sin(theta);
423  Double_t tantheta = TMath::Tan(theta);
424  Double_t costhetas = (1+fPolarOffset)*costheta;
425  Double_t sinthetas = (1+fPolarOffset)*sintheta;
426  Double_t corr = 0.01;
427 
428  TLatex *textangular = new TLatex();
429  textangular->SetTextColor(GetPolarColorLabel());
430  textangular->SetTextFont(GetPolarLabelFont());
431 
432  const char* form = (char *)" ";
433  TGaxis axis;
435  // Polar numbers are aligned with their axis.
436  if(fPolarLabels == NULL && optionLabels){;
437  if (fRadian) {
438  // Radian case.
439  ReduceFraction(2*i, ndivMajor, rnum, rden); // Reduces the fraction.
440  if (rnum == 0) form = Form("%d",rnum);
441  if (rnum == 1 && rden == 1) form = Form("#pi");
442  if (rnum == 1 && rden != 1) form = Form("#frac{#pi}{%d}",rden);
443  if (rnum != 1 && rden == 1 && i !=0) form= Form("%d#pi",rnum);
444  if (rnum != 1 && rden != 1) form = Form("#frac{%d#pi}{%d}",rnum,rden);
445  textangular->SetTextAlign(FindAlign(theta));
446  textangular->PaintLatex(costhetas,
447  sinthetas, FindTextAngle(theta),
448  GetPolarLabelSize(), form);
449  } else {
450  // Any other cases: numbers are aligned with their axis.
451  form = Form("%5.3g",txtval);
452  axis.LabelsLimits(form,first,last);
453  TString s = Form("%s",form);
454  if (first != 0) s.Remove(0, first);
455  textangular->SetTextAlign(FindAlign(theta));
456  textangular->PaintLatex(costhetas,
457  sinthetas, FindTextAngle(theta),
458  GetPolarLabelSize(), s);
459  }
460  } else if (fPolarLabels){
461  // print the specified polar labels
462  textangular->SetTextAlign(FindAlign(theta));
463  textangular->PaintLatex(costhetas,sinthetas,FindTextAngle(theta),
465  }
466  } else {
467  // Polar numbers are shown horizontally.
468  if(fPolarLabels == NULL && optionLabels){
469  if (fRadian) {
470  // Radian case
471  ReduceFraction(2*i, ndivMajor, rnum, rden);
472  if (rnum == 0) form = Form("%d",rnum);
473  if (rnum == 1 && rden == 1) form = Form("#pi");
474  if (rnum == 1 && rden != 1) form = Form("#frac{#pi}{%d}",rden);
475  if (rnum != 1 && rden == 1 && i !=0) form = Form("%d#pi",rnum);
476  if (rnum != 1 && rden != 1) form = Form("#frac{%d#pi}{%d}",rnum,rden);
477  if(theta >= 3*TMath::Pi()/12.0 && theta < 2*TMath::Pi()/3.0) corr=0.04;
478  textangular->SetTextAlign(FindAlign(theta));
479  textangular->PaintLatex(costhetas,corr+sinthetas,0,
480  GetPolarLabelSize(),form);
481  } else {
482  // Any other cases where numbers are shown horizontally.
483  form = Form("%5.3g",txtval);
484  axis.LabelsLimits(form,first,last);
485  TString s = Form("%s",form);
486  if (first != 0) s.Remove(0, first);
487  if(theta >= 3*TMath::Pi()/12.0 && theta < 2*TMath::Pi()/3.0) corr=0.04;
488  textangular->SetTextAlign(FindAlign(theta));
489  textangular->PaintLatex(costhetas, //j'ai efface des offset la
490  corr+sinthetas,0,GetPolarLabelSize(),s);
491  }
492  } else if (fPolarLabels) {
493  // print the specified polar labels
494  textangular->SetTextAlign(FindAlign(theta));
495  textangular->PaintText(costhetas,sinthetas,fPolarLabels[i]);
496  }
497  }
499  //Check if SetTickPolar is activated, and draw tick marks
500  Bool_t issettickpolar = gPad->GetTicky();
501 
502  if (issettickpolar) {
503  if (theta != 0 && theta !=TMath::Pi()) {
504  gPad->PaintLine((sintheta-GetTickpolarSize())/tantheta,sintheta-GetTickpolarSize(),
505  (sintheta+GetTickpolarSize())/tantheta,sintheta+GetTickpolarSize());
506  }
507  if (theta == 0 || theta ==TMath::Pi()) {
508  gPad->PaintLine(1-GetTickpolarSize(),0,1+GetTickpolarSize(),0);
509  gPad->PaintLine(-1+GetTickpolarSize(),0,-1-GetTickpolarSize(),0);
510  }
511  }
514  gPad->PaintLine(0.,0.,costheta,sintheta);
515  delete textangular;
516  // Add minor lines w/o text.
517  Int_t oldLineStyle = GetLineStyle();
518  TAttLine::SetLineStyle(2); //Minor lines always in this style.
519  TAttLine::Modify(); //Changes line attributes apart from style.
520  for (j=1; j<ndivMinor; j++) {
521  Double_t thetamin = theta+j*2*TMath::Pi()/(ndivMajor*ndivMinor);
522  gPad->PaintLine(0.,0.,TMath::Cos(thetamin),TMath::Sin(thetamin));
523  }
524  TAttLine::SetLineStyle(oldLineStyle);
526  }
527  } else {
528  Int_t big = (Int_t)fRwtmax;
529  Int_t test= 1;
530  while (big >= 10) {
531  big = big/10;
532  test++;
533  }
534  for (i=1; i<=test; i++) {
535  Double_t txtval = pow((double)10,(double)(i-1));
536  Double_t theta = (i-1)*2*TMath::Pi()/(double)(test);
537  Double_t costheta = TMath::Cos(theta);
538  Double_t sintheta = TMath::Sin(theta);
539  Double_t tantheta = TMath::Tan(theta);
540  Double_t costhetas = (1+fPolarOffset)*costheta;
541  Double_t sinthetas = (1+fPolarOffset)*sintheta;
542  Double_t corr = 0.01;
543 
544  TLatex *textangular = new TLatex();
545  textangular->SetTextColor(GetPolarColorLabel());
546  textangular->SetTextFont(GetPolarLabelFont());
547 
548  const char* form = (char *)" ";
549  TGaxis axis;
550 
552  if(fPolarLabels==NULL && optionLabels){
553  // Polar numbers are aligned with their axis.
554  form = Form("%5.3g",txtval);
555  axis.LabelsLimits(form,first,last);
556  TString s = Form("%s",form);
557  if (first != 0) s.Remove(0, first);
558  textangular->SetTextAlign(FindAlign(theta));
559  textangular->PaintLatex(costhetas,
560  sinthetas, FindTextAngle(theta), GetPolarLabelSize(), s);
561  }
562  else if (fPolarLabels){
563  // print the specified polar labels
564  textangular->SetTextAlign(FindAlign(theta));
565  textangular->PaintText(costhetas,sinthetas,fPolarLabels[i]);
566  }
567 
568  } else {
569  if(fPolarLabels==NULL && optionLabels){
570  // Polar numbers are shown horizontally.
571  form = Form("%5.3g",txtval);
572  axis.LabelsLimits(form,first,last);
573  TString s = Form("%s",form);
574  if (first != 0) s.Remove(0, first);
575  if(theta >= 3*TMath::Pi()/12.0 && theta < 2*TMath::Pi()/3.0) corr=0.04;
576  textangular->SetTextAlign(FindAlign(theta));
577  textangular->PaintLatex(costhetas,
578  corr+sinthetas,0,GetPolarLabelSize(),s);
579  } else if (fPolarLabels){
580  // print the specified polar labels
581  textangular->SetTextAlign(FindAlign(theta));
582  textangular->PaintText(costhetas,sinthetas,fPolarLabels[i]);
583  }
584  }
585 
587  //Check if SetTickPolar is activated, and draw tick marks
588  Bool_t issettickpolar = gPad->GetTicky();
589  if (issettickpolar) {
590  if (theta != 0 && theta !=TMath::Pi()) {
591  gPad->PaintLine((sintheta-GetTickpolarSize())/tantheta,sintheta-GetTickpolarSize(),
592  (sintheta+GetTickpolarSize())/tantheta,sintheta+GetTickpolarSize());
593  }
594  if (theta == 0 || theta ==TMath::Pi()) {
595  gPad->PaintLine(1-GetTickpolarSize(),0,1+GetTickpolarSize(),0);
596  gPad->PaintLine(-1+GetTickpolarSize(),0,-1-GetTickpolarSize(),0);
597  }
598  }
601  gPad->PaintLine(0.,0.,costheta,sintheta);
602  delete textangular;
603  // Add minor lines w/o text.
604  Int_t oldLineStyle = GetLineStyle();
605  TAttLine::SetLineStyle(2); //Minor lines always in this style.
606  TAttLine::Modify(); //Changes line attributes apart from style.
607  Double_t a=0;
608  Double_t b,c,d;
609  b = TMath::Log(10)*test;
610  d= 2*TMath::Pi()/(double)test;
611  for (j=1; j<9; j++) {
612  a=TMath::Log(j+1)-TMath::Log(j)+a;
613  c=a/b*6.28+d*(i-1);
614  gPad->PaintLine(0.,0.,TMath::Cos(c),TMath::Sin(c));
615  }
616  TAttLine::SetLineStyle(oldLineStyle);
618  }
619  }
620 }
621 
622 ////////////////////////////////////////////////////////////////////////////////
623 /// Paint radial divisions.
624 /// Check for editable pad or create default.
625 
627 {
628  static char chopt[8] = "";
629  Int_t i,j;
630  Int_t ndiv = TMath::Abs(fNdivRad);
631  Int_t ndivMajor = ndiv%100;
632  Int_t ndivMinor = ndiv/100;
633  Int_t ndivmajor = 0;
634  Double_t frwrmin = 0., frwrmax = 0., binWidth = 0;
635 
636  THLimitsFinder::Optimize(fRwrmin,fRwrmax,ndivMajor,frwrmin,
637  frwrmax, ndivmajor,binWidth,"");
638 
639  if (!gPad) return ;
640  if (!gPad->GetLogx()) {
641  gPad->RangeAxis(-1,-1,1,1);
642  gPad->Range(-1.25,-1.25,1.25,1.25);
643  Double_t umin = fRwrmin;
644  Double_t umax = fRwrmax;
645  Double_t rmajmin = (frwrmin-fRwrmin)/(fRwrmax-fRwrmin);
646  Double_t rmajmax = (frwrmax-fRwrmin)/(fRwrmax-fRwrmin);
647  Double_t dist = (rmajmax-rmajmin)/ndivmajor;
648  Int_t ndivminor = 0;
649 
650  chopt[0] = 0;
651  strncat(chopt, "SDH", 3);
652  if (fNdivRad < 0) strncat(chopt, "N",1);
653  if(drawaxis){
654  // Paint axis.
655  TGaxis axis;
660  axis.PaintAxis(0, 0, TMath::Cos(GetAngle()), TMath::Sin(GetAngle()),
661  umin, umax, ndiv, chopt, 0., kFALSE);
662  }
663 
664  // Paint Circles.
665  // First paint main circle.
666  PaintCircle(0.,0.,1,0.,360,0);
667  // Optimised case.
668  if (fNdivRad>0 ) {
669  Double_t frwrmini = 0., frwrmaxi = 0., binWidth2 =0;
670  THLimitsFinder::Optimize(frwrmin,frwrmin+binWidth,ndivMinor,frwrmini,
671  frwrmaxi, ndivminor,binWidth2,"");
672  Double_t dist2 = dist/(ndivminor);
673  // Paint major circles.
674  for (i=1; i<=ndivmajor+2; i++) {
677  PaintCircle(0.,0.,rmajmin,0.,360,0);
678 
679  //Paint minor circles.
682  for (j=1; j<ndivminor+1; j++) {
683  if (rmajmin+j*dist2<=1) PaintCircle(0.,0.,rmajmin+j*dist2,0.,360,0);
684  }
685  rmajmin = (frwrmin-fRwrmin)/(fRwrmax-fRwrmin)+(i-1)*dist;
686  }
687  // Non-optimized case.
688  } else {
689 
690  // Paint major circles.
691  for (i=1; i<=ndivMajor; i++) {
694  Double_t rmaj = i*1./ndivMajor;
695  PaintCircle(0.,0.,rmaj,0.,360,0);
696 
697  // Paint minor circles.
698  for (j=1; j<ndivMinor; j++) {
701  PaintCircle(0.,0.,rmaj- j*1./(ndivMajor*ndivMinor),0.,360,0);
702  }
703  }
704  }
705  } else {
706  // Draw Log scale on radial axis if option activated.
707  Int_t big = (Int_t)fRwrmax;
708  Int_t test= 0;
709  while (big >= 10) {
710  big = big/10;
711  test++;
712  }
713  for (i=1; i<=test; i++) {
716  Double_t ecart;
717  ecart = ((double) i)/ ((double) test);
718  PaintCircle(0.,0.,ecart,0,360,0);
721  Double_t a=0;
722  Double_t b,c,d;
723  b = TMath::Log(10)*test;
724  d = 1/(double)test;
725  for (j=1; j<9; j++) {
726  a = TMath::Log(j+1)-TMath::Log(j)+a;
727  c = a/b+d*(i-1);
728  PaintCircle(0,0.,c,0.,360,0);
729  }
730  }
731  }
734 }
735 
736 ////////////////////////////////////////////////////////////////////////////////
737 /// Reduce fractions.
738 
740 {
741  Int_t a = 0;
742  Int_t b = 0;
743  Int_t i = 0;
744  Int_t j = 0;
745  a = den;
746  b = num;
747  if (b > a) {
748  j = b;
749  } else {
750  j = a;
751  }
752  for (i=j; i > 1; i--) {
753  if ((a % i == 0) && (b % i == 0)) {
754  a = a/i;
755  b = b/i;
756  }
757  }
758  rden = a;
759  rnum = b;
760 }
761 
762 ////////////////////////////////////////////////////////////////////////////////
763 /// Set axis angle.
764 
766 {
767  fAxisAngle = angle/180*TMath::Pi();
768 }
769 
770 ////////////////////////////////////////////////////////////////////////////////
771 /// Set the number of Polar divisions: enter a number ij with 0<i<99 and 0<j<99
772 /// - i sets the major Polar divisions.
773 /// - j sets the minor Polar divisions.
774 
776 {
777  if (ndiv > 0)
778  fNdivPol = ndiv;
779  if (gPad) gPad->Modified();
780 }
781 
782 ////////////////////////////////////////////////////////////////////////////////
783 /// Set the number of radial divisions: enter a number ij with 0<i<99 and 0<j<99
784 /// - i sets the major radial divisions.
785 /// - j sets the minor radial divisions.
786 
788 {
789  fNdivRad = ndiv;
790  if (gPad) gPad->Modified();
791 }
792 
793 ////////////////////////////////////////////////////////////////////////////////
794 /// Set some specified polar labels, used in the case of a spider plot.
795 
797 {
798  if(fPolarLabels == NULL)
800  fPolarLabels[div]=label;
801  if (gPad) gPad->Modified();
802 }
803 
804 ////////////////////////////////////////////////////////////////////////////////
805 /// Set Polar labels color.
806 
808 {
809  fPolarLabelColor = tcolorangular;
810 }
811 
812 ////////////////////////////////////////////////////////////////////////////////
813 
815 {
816  // Set Polar label font.
817 
818  fPolarLabelFont = tfontangular;
819 }
820 
821 ////////////////////////////////////////////////////////////////////////////////
822 /// Set angular labels size.
823 
825 {
826  fPolarTextSize = angularsize;
827 }
828 
829 ////////////////////////////////////////////////////////////////////////////////
830 /// Set the labels offset.
831 
833 {
834  fPolarOffset = angularOffset;
835  if (gPad) gPad->Modified();
836 }
837 
838 ////////////////////////////////////////////////////////////////////////////////
839 /// Set radial labels color.
840 
842 {
843  fRadialLabelColor = tcolorradial;
844 }
845 
846 ////////////////////////////////////////////////////////////////////////////////
847 /// Set radial label font.
848 
850 {
851  fRadialLabelFont = tfontradial;
852 }
853 
854 ////////////////////////////////////////////////////////////////////////////////
855 /// Set radial labels size.
856 
858 {
859  fRadialTextSize = radialsize;
860 }
861 
862 ////////////////////////////////////////////////////////////////////////////////
863 /// Set the labels offset.
864 
866 {
867  fRadialOffset = radialOffset;
868  if (gPad) gPad->Modified();
869 }
870 
871 ////////////////////////////////////////////////////////////////////////////////
872 /// Allows to change range Polar.
873 /// \param[in] tmin the start number.
874 /// \param[in] tmax the end number.
875 
877 {
878  fDegree = kFALSE;
879  fGrad = kFALSE;
880  fRadian = kFALSE;
881 
882  if (tmin < tmax) {
883  fRwtmin = tmin;
884  fRwtmax = tmax;
885  }
886  if (gPad) gPad->Modified();
887 }
888 
889 ////////////////////////////////////////////////////////////////////////////////
890 /// Set the radial range.
891 /// \param[in] rmin radius at center of the circle.
892 /// \param[in] rmax radius at the intersection of the right X axis part and the circle.
893 
895 {
896  if (rmin < rmax) {
897  fRwrmin = rmin;
898  fRwrmax = rmax;
899  }
900  if (gPad) gPad->Modified();
901 }
902 
903 ////////////////////////////////////////////////////////////////////////////////
904 /// Set polar ticks size.
905 
907 {
908  fTickpolarSize = tickpolarsize;
909 }
910 
911 ////////////////////////////////////////////////////////////////////////////////
912 /// The Polar circle is labelled using degree.
913 
915 {
916  fDegree = kTRUE;
917  fGrad = kFALSE;
918  fRadian = kFALSE;
919  ChangeRangePolar(0,360);
920 }
921 
922 ////////////////////////////////////////////////////////////////////////////////
923 /// The Polar circle is labelled using gradian.
924 
926 {
927  fGrad = kTRUE;
928  fRadian = kFALSE;
929  fDegree = kFALSE;
930  ChangeRangePolar(0,200);
931 }
932 
933 ////////////////////////////////////////////////////////////////////////////////
934 /// The Polar circle is labelled using radian.
935 
937 {
938  fRadian = kTRUE;
939  fGrad = kFALSE;
940  fDegree = kFALSE;
942 }
943 
944 ////////////////////////////////////////////////////////////////////////////////
945 /// Set range from 0 to 2*pi
946 
948 {
949  SetRangePolar(0,2*TMath::Pi());
950 }
virtual Style_t GetLineStyle() const
Definition: TAttLine.h:48
void SetToRadian()
The Polar circle is labelled using radian.
double dist(Rotation3D const &r1, Rotation3D const &r2)
Definition: 3DDistances.cxx:48
void SetTickpolarSize(Double_t tickpolarsize=0.02)
Set polar ticks size.
void SetTwoPi()
Set range from 0 to 2*pi.
Color_t fPolarLabelColor
Double_t fPolarOffset
void Init()
Initialize some of the fields of TGraphPolargram.
static double p3(double t, double a, double b, double c, double d)
void SetRadialLabelSize(Double_t radialsize=0.035)
Set radial labels size.
Double_t Log(Double_t x)
Definition: TMath.h:526
const double pi
const char Option_t
Definition: RtypesCore.h:62
void PaintCircle(Double_t x, Double_t y, Double_t r, Double_t phimin, Double_t phimax, Double_t theta)
This is simplified from TEllipse::PaintEllipse.
void SetToDegree()
The Polar circle is labelled using degree.
void Draw(Option_t *options="")
Draw Polargram.
static void Optimize(Double_t A1, Double_t A2, Int_t nold, Double_t &BinLow, Double_t &BinHigh, Int_t &nbins, Double_t &BWID, Option_t *option="")
static function to compute reasonable axis limits
void ToUpper()
Change string to upper case.
Definition: TString.cxx:1101
Double_t fTickpolarSize
Double_t GetPolarLabelSize()
Double_t GetRadialOffset()
Basic string class.
Definition: TString.h:137
ClassImp(TGraphPolargram)
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:170
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
const Bool_t kFALSE
Definition: Rtypes.h:92
Double_t fRadialTextSize
virtual void Modify()
Change current line attributes if necessary.
Definition: TAttLine.cxx:229
short Font_t
Definition: RtypesCore.h:75
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:501
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:732
Int_t FindAlign(Double_t angle)
Find the alignement rule to apply for TText::SetTextAlign(Short_t).
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition: TObject.cxx:164
Double_t fRadialOffset
void ChangeRangePolar(Double_t tmin, Double_t tmax)
Set the Polar range.
virtual void SetTextFont(Font_t tfont=62)
Definition: TAttText.h:59
void SetRadialLabelFont(Font_t tfontradial=62)
Set radial label font.
Double_t x[n]
Definition: legend1.C:17
The TNamed class is the base class for all named ROOT classes.
Definition: TNamed.h:33
void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Indicate that there is something to click here.
To draw Mathematical Formula.
Definition: TLatex.h:33
const Int_t kMaxPixel
Definition: GuiTypes.h:370
static double p2(double t, double a, double b, double c)
double pow(double, double)
void Paint(Option_t *options="")
[fNdivPol] Specified polar labels
void PaintPolarDivisions(Bool_t noLabels)
Draw Polar divisions.
void SetPolarOffset(Double_t PolarOffset=0.04)
Set the labels offset.
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
void SetNdivRadial(Int_t Ndiv=508)
Set the number of radial divisions: enter a number ij with 0
void SetAxisAngle(Double_t angle=0)
Set axis angle.
TGraphPolargram(const char *name, Double_t rmin, Double_t rmax, Double_t tmin, Double_t tmax)
TGraphPolargram Constructor.
void SetLabelSize(Float_t labelsize)
Definition: TGaxis.h:118
short Color_t
Definition: RtypesCore.h:79
virtual void SetTextAlign(Short_t align=11)
Definition: TAttText.h:55
Color_t fRadialLabelColor
Style_t fLineStyle
Definition: TAttLine.h:36
virtual void PaintAxis(Double_t xmin, Double_t ymin, Double_t xmax, Double_t ymax, Double_t &wmin, Double_t &wmax, Int_t &ndiv, Option_t *chopt="", Double_t gridlength=0, Bool_t drawGridOnly=kFALSE)
Control function to draw an axis.
Definition: TGaxis.cxx:705
ROOT::R::TRInterface & r
Definition: Object.C:4
return
Definition: TBase64.cxx:62
void PaintRadialDivisions(Bool_t drawaxis)
Paint radial divisions.
void SetRangeRadial(Double_t rmin, Double_t rmax)
Set the radial range.
Double_t GetTickpolarSize()
void ReduceFraction(Int_t Num, Int_t Denom, Int_t &rnum, Int_t &rden)
Reduce fractions.
Bool_t TestBit(UInt_t f) const
Definition: TObject.h:173
char * Form(const char *fmt,...)
To draw polar axis.
The axis painter class.
Definition: TGaxis.h:39
Font_t GetPolarLabelFont()
Double_t ACos(Double_t)
Definition: TMath.h:445
static double p1(double t, double a, double b)
void SetPolarLabelSize(Double_t angularsize=0.04)
Set angular labels size.
Double_t fPolarTextSize
virtual ~TGraphPolargram()
TGraphPolargram destructor.
Double_t Cos(Double_t)
Definition: TMath.h:424
void SetPolarLabel(Int_t div, const TString &label)
Set some specified polar labels, used in the case of a spider plot.
Double_t Pi()
Definition: TMath.h:44
TString * fPolarLabels
TString & Remove(Ssiz_t pos)
Definition: TString.h:616
void SetLabelOffset(Float_t labeloffset)
Definition: TGaxis.h:117
static const double x1[5]
double Double_t
Definition: RtypesCore.h:55
Double_t y[n]
Definition: legend1.C:17
void SetLabelFont(Int_t labelfont)
Definition: TGaxis.h:116
Int_t DistancetoLine(Int_t px, Int_t py, Double_t xp1, Double_t yp1, Double_t xp2, Double_t yp2)
Compute distance from point px,py to a line.
Definition: TAttLine.cxx:193
virtual void SetLineStyle(Style_t lstyle)
Definition: TAttLine.h:56
#define name(a, b)
Definition: linkTestLib0.cpp:5
Int_t DistancetoPrimitive(Int_t px, Int_t py)
Everything within the circle belongs to the TGraphPolargram.
void SetToGrad()
The Polar circle is labelled using gradian.
Double_t PiOver2()
Definition: TMath.h:46
void SetNdivPolar(Int_t Ndiv=508)
Set the number of Polar divisions: enter a number ij with 0
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:567
Double_t Sin(Double_t)
Definition: TMath.h:421
Double_t ASin(Double_t)
Definition: TMath.h:439
#define NULL
Definition: Rtypes.h:82
void SetLabelColor(Int_t labelcolor)
Definition: TGaxis.h:115
#define gPad
Definition: TVirtualPad.h:288
virtual void SetTextColor(Color_t tcolor=1)
Definition: TAttText.h:57
void ResetBit(UInt_t f)
Definition: TObject.h:172
void SetPolarLabelColor(Color_t tcolorangular=1)
Set Polar labels color.
Color_t GetRadialColorLabel()
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
const Bool_t kTRUE
Definition: Rtypes.h:91
void SetRadialOffset(Double_t RadialOffset=0.025)
Set the labels offset.
Color_t GetPolarColorLabel()
virtual void PaintText(Double_t x, Double_t y, const char *text)
Draw this text with new coordinates.
Definition: TText.cxx:740
const Int_t n
Definition: legend1.C:16
void SetRangePolar(Double_t tmin, Double_t tmax)
Allows to change range Polar.
Double_t Tan(Double_t)
Definition: TMath.h:427
void SetRadialLabelColor(Color_t tcolorradial=1)
Set radial labels color.
Double_t FindTextAngle(Double_t theta)
Determine the orientation of the polar labels according to their angle.
Font_t GetRadialLabelFont()
Double_t GetRadialLabelSize()
void LabelsLimits(const char *label, Int_t &first, Int_t &last)
Internal method to find first and last character of a label.
Definition: TGaxis.cxx:2113
Double_t GetAngle()
void SetPolarLabelFont(Font_t tfontangular=62)