// @(#)root/graf:$Id$ // Author: Sebastian Boser, Mathieu Demaret 02/02/06 /************************************************************************* * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. * * All rights reserved. * * * * For the licensing terms see $ROOTSYS/LICENSE. * * For the list of contributors see $ROOTSYS/README/CREDITS. * *************************************************************************/ //______________________________________________________________________________ /* Begin_Html
Example:
End_Html
Begin_Macro(source)
{
TCanvas * CPol = new TCanvas("CPol","TGraphPolar Examples",500,500);
Double_t rmin=0;
Double_t rmax=TMath::Pi()*2;
Double_t r[1000];
Double_t theta[1000];
TF1 * fp1 = new TF1("fplot","cos(x)",rmin,rmax);
for (Int_t ipt = 0; ipt < 1000; ipt++) {
r[ipt] = ipt*(rmax-rmin)/1000+rmin;
theta[ipt] = fp1->Eval(r[ipt]);
}
TGraphPolar * grP1 = new TGraphPolar(1000,r,theta);
grP1->SetLineColor(2);
grP1->Draw("AOL");
return CPol;
}
End_Macro */
#include "TGraphPolar.h"
#include "TGraphPolargram.h"
#include "TGaxis.h"
#include "THLimitsFinder.h"
#include "TVirtualPad.h"
#include "TROOT.h"
#include "TLatex.h"
#include "TEllipse.h"
#include "TMath.h"
ClassImp(TGraphPolargram);
//______________________________________________________________________________
TGraphPolargram::TGraphPolargram(const char* name, Double_t rmin, Double_t rmax,
Double_t tmin, Double_t tmax):
TNamed(name,"Polargram")
{
// TGraphPolargram Constructor.
Init();
fNdivRad = 508;
fNdivPol = 508;
fPolarLabels = NULL;
fRwrmax = rmax;
fRwrmin = rmin;
fRwtmin = tmin;
fRwtmax = tmax;
}
//______________________________________________________________________________
TGraphPolargram::TGraphPolargram(const char* name):
TNamed(name,"Polargram")
{
// Short constructor used in the case of a spider plot.
Init();
fNdivRad = 0;
fNdivPol = 0;
fPolarLabels = NULL;
fRwrmax = 1;
fRwrmin = 0;
fRwtmax = 0;
fRwtmin = 0;
}
//______________________________________________________________________________
TGraphPolargram::~TGraphPolargram()
{
// TGraphPolargram destructor.
if (fPolarLabels != NULL) delete [] fPolarLabels;
}
//______________________________________________________________________________
void TGraphPolargram::ChangeRangePolar(Double_t tmin, Double_t tmax)
{
// Set the Polar range.
// tmin is the start number.
// tmax is the end number.
if (tmin < tmax) {
fRwtmin = tmin;
fRwtmax = tmax;
}
if (gPad) gPad->Modified();
}
//______________________________________________________________________________
Int_t TGraphPolargram::DistancetoPrimitive(Int_t px, Int_t py)
{
// Everything within the circle belongs to the TGraphPolargram.
Int_t i;
Double_t x = gPad->AbsPixeltoX(px);
Double_t y = gPad->AbsPixeltoY(py);
// Check if close to a (major) radial line.
Double_t rad = TMath::Sqrt(x*x+y*y);
Int_t div = (Int_t)rad*(fNdivRad%100);
Double_t dr = TMath::Min(TMath::Abs(rad-div*1./(fNdivRad%100)),
TMath::Abs(rad-(div+1)*1./(fNdivRad%100)));
Int_t drad = gPad->XtoPixel(dr)-gPad->XtoPixel(0);
// Check if close to a (major) Polar line.
// This is not a proper calculation, but rather fast.
Int_t dt = kMaxPixel;
for (i=0; i<(fNdivPol%100); i++) {
Double_t theta = i*2*TMath::Pi()/(fNdivPol%100);
// Attention: px,py in pixel units, line given in user coordinates.
Int_t dthis = DistancetoLine(px,py,0.,0.,TMath::Cos(theta),
TMath::Sin(theta));
// Fails if we are outside box discribed by the line.
// (i.e for all hor/vert lines)
if (dthis==9999) {
// Outside -> Get distance to endpoint of line.
if (rad>1) {
dthis = (Int_t)TMath::Sqrt(
TMath::Power(px-gPad->XtoPixel(TMath::Cos(theta)),2)+
TMath::Power(py-gPad->YtoPixel(TMath::Sin(theta)),2));
} else {
// Check for horizontal line.
if (((TMath::Abs(theta-TMath::Pi())<0.1) &&
((px-gPad->XtoPixel(0))<0)) ||
((TMath::Abs(theta)<0.1) &&
((px-gPad->XtoPixel(0))>0))) {
dthis = TMath::Abs(py-gPad->YtoPixel(0.));
}
//Check for vertical line.
if (((TMath::Abs(theta-TMath::PiOver2())<0.1) &&
((py-gPad->YtoPixel(0))>0)) ||
((TMath::Abs(theta-3*TMath::PiOver2())<0.1) &&
(py-gPad->YtoPixel(0))<0)) {
dthis = TMath::Abs(px-gPad->XtoPixel(0.));
}
if (dthis==9999) {
// Inside, but out of box for nonorthogonal line ->
// get distance to start point.
dthis = (Int_t)TMath::Sqrt(
TMath::Power(px-gPad->XtoPixel(0.),2)+
TMath::Power(py-gPad->YtoPixel(0.),2));
}
}
}
// Take distance to closes line.
dt = TMath::Min(dthis,dt);
}
return TMath::Min(drad, dt);
}
//______________________________________________________________________________
void TGraphPolargram::Draw(Option_t* options)
{
// Draw Polargram.
Paint(options);
AppendPad(options);
}
//______________________________________________________________________________
void TGraphPolargram::ExecuteEvent(Int_t event, Int_t px, Int_t py)
{
// Indicate that there is something to click here.
Int_t kMaxDiff = 20;
static Int_t d1, d2, d3, px1, py1, px3, py3;
static Bool_t p1, p2, p3, p4, p5, p6, p7, p8;
Double_t px2, py2;
p2 = p3 = p4 = p5 = p6 = p7 = p8 = kFALSE;
if (!gPad->IsEditable()) return;
switch (event) {
case kMouseMotion:
px1 = gPad->XtoAbsPixel(TMath::Cos(GetAngle()));
py1 = gPad->YtoAbsPixel(TMath::Sin(GetAngle()));
d1 = TMath::Abs(px1 - px) + TMath::Abs(py1-py); //simply take sum of pixels differences
p1 = kFALSE;
px2 = gPad->XtoAbsPixel(-1);
py2 = gPad->YtoAbsPixel(1);
d2 = (Int_t)(TMath::Abs(px2 - px) + TMath::Abs(py2 - py)) ;
px3 = gPad->XtoAbsPixel(-1);
py3 = gPad->YtoAbsPixel(-1);
d3 = TMath::Abs(px3 - px) + TMath::Abs(py3 - py) ; //simply take sum of pixels differences
// check if point is close to the radial axis
if (d1 < kMaxDiff) {
gPad->SetCursor(kMove);
p1 = kTRUE;
}
// check if point is close to the left high axis
if ( d2 < kMaxDiff) {
gPad->SetCursor(kHand);
p7 = kTRUE;
}
// check if point is close to the left down axis
if ( d3 < kMaxDiff) {
gPad->SetCursor(kHand);
p8 = kTRUE;
}
// check if point is close to a main circle
if (!p1 && !p7 ) {
p6 = kTRUE;
gPad->SetCursor(kHand);
}
break;
case kButton1Down:
// Record initial coordinates
//px4 = px;
//py4 = py;
case kButton1Motion:
if (p1) {
px2 = gPad->AbsPixeltoX(px);
py2 = gPad->AbsPixeltoY(py);
if ( px2 < 0 && py2 < 0) {p2 = kTRUE;};
if ( px2 < 0 && py2 > 0 ) {p3 = kTRUE;};
if ( px2 > 0 && py2 > 0 ) {p4 = kTRUE;};
if ( px2 > 0 && py2 < 0 ) {p5 = kTRUE;};
px2 = TMath::ACos(TMath::Abs(px2));
py2 = TMath::ASin(TMath::Abs(py2));
if (p2) {
fAxisAngle = TMath::Pi()+(px2+py2)/2;
p2 = kFALSE;
};
if (p3) {
fAxisAngle = TMath::Pi()-(px2+py2)/2;
p3 = kFALSE;
};
if (p4) {
fAxisAngle = (px2+py2)/2;
p4 = kFALSE;
};
if (p5) {
fAxisAngle = -(px2+py2)/2;
p5 = kFALSE;
};
}
break;
case kButton1Up:
Paint();
}
}
//______________________________________________________________________________
Int_t TGraphPolargram::FindAlign(Double_t angle)
{
// Find the alignement rule to apply for TText::SetTextAlign(Short_t).
Double_t pi = TMath::Pi();
while(angle < 0 || angle > 2*pi){
if(angle < 0) angle+=2*pi;
if(angle > 2*pi) angle-=2*pi;
}
if(!TestBit(TGraphPolargram::kLabelOrtho)){
if(angle > 0 && angle < pi/2) return 11;
else if(angle > pi/2 && angle < pi) return 31;
else if(angle > pi && angle < 3*pi/2) return 33;
else if(angle > 3*pi/2 && angle < 2*pi) return 13;
else if(angle == 0 || angle == 2*pi) return 12;
else if(angle == pi/2) return 21;
else if(angle == pi) return 32;
else if(angle == 3*pi/2) return 23;
else return 0;
}
else{
if(angle >= 0 && angle <= pi/2) return 12;
else if((angle > pi/2 && angle <= pi) || (angle > pi && angle <= 3*pi/2)) return 32;
else if(angle > 3*pi/2 && angle <= 2*pi) return 12;
else return 0;
}
}
//______________________________________________________________________________
Double_t TGraphPolargram::FindTextAngle(Double_t angle)
{
// Determine the orientation of the polar labels according to their angle.
Double_t pi = TMath::Pi();
Double_t convraddeg = 180.0/pi;
while(angle < 0 || angle > 2*pi){
if(angle < 0) angle+=2*pi;
if(angle > 2*pi) angle-=2*pi;
}
if(angle >= 0 && angle <= pi/2) return angle*convraddeg;
else if(angle > pi/2 && angle <= pi) return (angle + pi)*convraddeg;
else if(angle > pi && angle <= 3*pi/2) return (angle - pi)*convraddeg;
else if(angle > 3*pi/2 && angle <= 2*pi) return angle*convraddeg;
else return 0;
}
//______________________________________________________________________________
void TGraphPolargram::Init()
{
// Initiallize some of the fields of TGraphPolargram.
fAxisAngle = 0;
fCutRadial = 0;
fDegree = kFALSE;
fGrad = kFALSE;
fLineStyle = 3;
fPolarLabelColor = 1;
fPolarLabelFont = 62;
fPolarOffset = 0.04;
fPolarTextSize = 0.04;
fRadialOffset = 0.025;
fRadian = kTRUE;
fRadialLabelColor = 1;
fRadialLabelFont = 62;
fRadialTextSize = 0.035;
fTickpolarSize = 0.02;
}
//______________________________________________________________________________
void TGraphPolargram::Paint(Option_t * chopt)
{
// Paint TGraphPolargram.
Int_t optionpoldiv, optionraddiv;
Bool_t optionLabels = kTRUE;
TString opt = chopt;
opt.ToUpper();
if(opt.Contains('P')) optionpoldiv=1; else optionpoldiv=0;
if(opt.Contains('R')) optionraddiv=1; else optionraddiv=0;
if(opt.Contains('O')) SetBit(TGraphPolargram::kLabelOrtho);
else ResetBit(TGraphPolargram::kLabelOrtho);
if(!opt.Contains('P') && !opt.Contains('R')) optionpoldiv=optionraddiv=1;
if(opt.Contains('N')) optionLabels = kFALSE;
if(optionraddiv) PaintRadialDivisions(kTRUE);
else PaintRadialDivisions(kFALSE);
if(optionpoldiv) PaintPolarDivisions(optionLabels);
}
//______________________________________________________________________________
void TGraphPolargram::PaintCircle(Double_t x1, Double_t y1, Double_t r,
Double_t phimin, Double_t phimax, Double_t theta)
{
// This is simplifed from TEllipse::PaintEllipse.
// Draw this ellipse with new coordinates.
Int_t i;
const Int_t np = 200; // Number of point to draw circle
static Double_t x[np+3], y[np+3];
// Set number of points approximatively proportional to the ellipse
// circumference.
Double_t circ = TMath::Pi()*2*r*(phimax-phimin)/36;
Int_t n = (Int_t)(np*circ/((gPad->GetX2()-gPad->GetX1())+
(gPad->GetY2()-gPad->GetY1())));
if (n < 8) n = 8;
if (n > np) n = np;
Double_t angle,dx,dy;
Double_t dphi = (phimax-phimin)*TMath::Pi()/(180*n);
Double_t ct = TMath::Cos(TMath::Pi()*theta/180);
Double_t st = TMath::Sin(TMath::Pi()*theta/180);
for (i=0; i<=n; i++) {
angle = phimin*TMath::Pi()/180 + Double_t(i)*dphi;
dx = r*TMath::Cos(angle);
dy = r*TMath::Sin(angle);
x[i] = x1 + dx*ct - dy*st;
y[i] = y1 + dx*st + dy*ct;
}
gPad->PaintPolyLine(n+1,x,y);
}
//______________________________________________________________________________
void TGraphPolargram::PaintPolarDivisions(Bool_t optionLabels)
{
// Draw Polar divisions.
// Check for editable pad or create default.
Int_t i, j, rnum, rden, first, last;
if (!gPad) return ;
gPad->RangeAxis(-1,-1,1,1);
gPad->Range(-1.25,-1.25,1.25,1.25);
Int_t ndivMajor = fNdivPol%100;
Int_t ndivMinor = fNdivPol/100;
if (!gPad->GetLogy()) {
for (i=0; i