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