//Begin_Html
/*
<img src="gif/graphpolar.gif">
*/
//End_Html
#include "TGraphPolar.h"
ClassImp(TGraphPolar);
TGraphPolar::TGraphPolar(Int_t n,
const Double_t* x, const Double_t* y, const Double_t *ex, const Double_t* ey)
: TGraphErrors(n,x,y,ex,ey)
{
fXpol = 0;
fYpol = 0;
fPolargram = 0;
fOptionAxis = kFALSE;
SetEditable(kFALSE);
};
TGraphPolar::~TGraphPolar()
{
if (fXpol) delete fXpol;
if (fYpol) delete fYpol;
}
Int_t TGraphPolar::DistancetoPrimitive(Int_t px, Int_t py)
{
Double_t* xold = fX; fX = fXpol;
Double_t* yold = fY; fY = fYpol;
Int_t dist = TGraphErrors::DistancetoPrimitive(px,py);
fX = xold;
fY = yold;
return dist;
}
void TGraphPolar::ExecuteEvent(Int_t , Int_t , Int_t )
{
gPad->SetCursor(kHand);
}
void TGraphPolar::SetMaxRadial(Double_t maximum)
{
if (fPolargram) fPolargram->SetRangeRadial(fPolargram->GetRMin(),maximum);
}
void TGraphPolar::SetMinRadial(Double_t minimum)
{
if (fPolargram) fPolargram->SetRangeRadial(minimum, fPolargram->GetRMax());
}
void TGraphPolar::SetMaxPolar(Double_t maximum)
{
if (fPolargram) fPolargram->SetRangePolar(fPolargram->GetTMin(),maximum);
}
void TGraphPolar::SetMinPolar(Double_t minimum)
{
if (fPolargram) fPolargram->SetRangePolar(minimum, fPolargram->GetTMax());
}
void TGraphPolar::Draw(Option_t* options)
{
TString opt = options;
opt.ToUpper();
opt.ReplaceAll("SAME","");
if (opt.Contains("A")) fOptionAxis = kTRUE;
opt.ReplaceAll("A","");
AppendPad(opt);
}
void TGraphPolar::Paint(Option_t* options)
{
Int_t ipt, i;
if (fNpoints<1) return;
TString opt = options;
opt.ToUpper();
Double_t rwrmin, rwrmax, rwtmin, rwtmax;
if (gPad) {
if (fPolargram)
if (!gPad->FindObject(fPolargram->GetName())) fPolargram=0;
if (!fPolargram) {
TListIter padObjIter(gPad->GetListOfPrimitives());
while (TObject* AnyObj = padObjIter.Next()) {
if (TString(AnyObj->ClassName()).CompareTo("TGraphPolargram",
TString::kExact)==0)
fPolargram = (TGraphPolargram*)AnyObj;
}
}
}
if (!fPolargram) {
rwrmin = fY[0]; rwrmax = fY[fNpoints-1];
rwtmin = fX[0]; rwtmax = fX[fNpoints-1];
for (ipt = 0; ipt < fNpoints; ipt++) {
if (fEX) {
if (fX[ipt] -fEX[ipt] < rwtmin) rwtmin = fX[ipt]-fEX[ipt];
if (fX[ipt] +fEX[ipt] > rwtmax) rwtmax = fX[ipt]+fEX[ipt];
} else {
if (fX[ipt] < rwtmin) rwtmin=fX[ipt];
if (fX[ipt] > rwtmax) rwtmax=fX[ipt];
}
if (fEY) {
if (fY[ipt] -fEY[ipt] < rwrmin) rwrmin = fY[ipt]-fEY[ipt];
if (fY[ipt] +fEY[ipt] > rwrmax) rwrmax = fY[ipt]+fEY[ipt];
} else {
if (fY[ipt] < rwrmin) rwrmin=fY[ipt];
if (fY[ipt] > rwrmax) rwrmax=fY[ipt];
}
}
if (rwrmin == rwrmax) rwrmax += 1.;
if (rwtmin == rwtmax) rwtmax += 1.;
Double_t dr = (rwrmax-rwrmin);
Double_t dt = (rwtmax-rwtmin);
rwrmax += 0.1*dr;
rwrmin -= 0.1*dr;
rwtmax += dt/fNpoints;
} else {
rwrmin = fPolargram->GetRMin();
rwrmax = fPolargram->GetRMax();
rwtmin = fPolargram->GetTMin();
rwtmax = fPolargram->GetTMax();
}
if ((!fPolargram)||fOptionAxis) {
fPolargram = new TGraphPolargram("Polargram",rwrmin,rwrmax,rwtmin,rwtmax);
fPolargram->Draw(opt);
fOptionAxis = kFALSE;
}
if (fXpol) delete fXpol;
if (fYpol) delete fYpol;
fXpol = new Double_t[fNpoints];
fYpol = new Double_t[fNpoints];
Double_t radiusNDC = rwrmax-rwrmin;
Double_t thetaNDC = (rwtmax-rwtmin)/(2*kPi);
if (opt.Contains("E")) {
if (fEY) {
for (i=0; i<fNpoints; i++) {
Double_t eymin, eymax, exmin,exmax;
exmin = (fY[i]-fEY[i]-rwrmin)/radiusNDC*
TMath::Cos((fX[i]-rwtmin)/thetaNDC);
eymin = (fY[i]-fEY[i]-rwrmin)/radiusNDC*
TMath::Sin((fX[i]-rwtmin)/thetaNDC);
exmax = (fY[i]+fEY[i]-rwrmin)/radiusNDC*
TMath::Cos((fX[i]-rwtmin)/thetaNDC);
eymax = (fY[i]+fEY[i]-rwrmin)/radiusNDC*
TMath::Sin((fX[i]-rwtmin)/thetaNDC);
TAttLine::Modify();
gPad->PaintLine(exmin,eymin,exmax,eymax);
}
}
if (fEX) {
for (i=0; i<fNpoints; i++) {
Double_t rad = (fY[i]-rwrmin)/radiusNDC;
Double_t phimin = (fX[i]-fEX[i]-rwtmin)/thetaNDC*180/kPi;
Double_t phimax = (fX[i]+fEX[i]-rwtmin)/thetaNDC*180/kPi;
TAttLine::Modify();
fPolargram->PaintCircle(0,0,rad,phimin,phimax,0);
}
}
}
for (i=0; i<fNpoints; i++) {
fXpol[i] = (fY[i]-rwrmin)/radiusNDC*TMath::Cos((fX[i]-rwtmin)/thetaNDC);
fYpol[i] = (fY[i]-rwrmin)/radiusNDC*TMath::Sin((fX[i]-rwtmin)/thetaNDC);
}
TGraph::PaintGraph(fNpoints, fXpol, fYpol,opt);
}
ClassImp(TGraphPolargram);
TGraphPolargram::TGraphPolargram(const char* name, Double_t rmin, Double_t rmax,
Double_t tmin, Double_t tmax):
TNamed(name,"Polargram")
{
fRwrmin = rmin;
fRwrmax = rmax;
fRwtmin = tmin;
fRwtmax = tmax;
fNdivRad = 502;
fNdivPol = 304;
fLabelOffset = 0.03;
}
TGraphPolargram::~TGraphPolargram()
{
}
Int_t TGraphPolargram::DistancetoPrimitive(Int_t px, Int_t py)
{
Int_t i;
Double_t x = gPad->AbsPixeltoX(px);
Double_t y = gPad->AbsPixeltoY(py);
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);
Int_t dt = kMaxPixel;
for (i=0; i<(fNdivPol%100); i++) {
Double_t theta = i*2*kPi/(fNdivPol%100);
Int_t dthis = DistancetoLine(px,py,0.,0.,TMath::Cos(theta),
TMath::Sin(theta));
if (dthis==9999) {
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 {
if (((TMath::Abs(theta-kPi)<0.1)&&((px-gPad->XtoPixel(0))<0))
||((TMath::Abs(theta)<0.1)&&((px-gPad->XtoPixel(0))>0))) {
dthis = abs(py-gPad->YtoPixel(0.));
}
if (((TMath::Abs(theta-kPi/2)<0.1)&&((py-gPad->YtoPixel(0))>0))
||((TMath::Abs(theta-3*kPi/2)<0.1)&&(py-gPad->YtoPixel(0))<0)) {
dthis = abs(px-gPad->XtoPixel(0.));
}
if (dthis==9999) {
dthis = (Int_t)TMath::Sqrt(
TMath::Power(px-gPad->XtoPixel(0.),2)+
TMath::Power(py-gPad->YtoPixel(0.),2));
}
}
}
dt = TMath::Min(dthis,dt);
}
return TMath::Min(drad, dt);
}
void TGraphPolargram::ExecuteEvent(Int_t , Int_t , Int_t )
{
gPad->SetCursor(kHand);
}
void TGraphPolargram::SetNdivRadial(Int_t ndiv)
{
if (ndiv > 0) fNdivRad = ndiv;
if (gPad) gPad->Modified();
}
void TGraphPolargram::SetNdivPolar(Int_t ndiv)
{
if (ndiv > 0) fNdivPol = ndiv;
if (gPad) gPad->Modified();
}
void TGraphPolargram::SetLabelOffset(Double_t labelOffset)
{
fLabelOffset = labelOffset;
if (gPad) gPad->Modified();
}
void TGraphPolargram::SetRangeRadial(Double_t rmin, Double_t rmax)
{
if (rmin < rmax) {
fRwrmin = rmin;
fRwrmax = rmax;
}
if (gPad) gPad->Modified();
}
void TGraphPolargram::SetRangePolar(Double_t tmin, Double_t tmax)
{
if (tmin < tmax) {
fRwtmin = tmin;
fRwtmax = tmax;
}
if (gPad) gPad->Modified();
}
void TGraphPolargram::Draw(Option_t* options)
{
TGraphPolargram::Paint(options);
AppendPad(options);
}
void TGraphPolargram::PaintCircle(Double_t x1, Double_t y1, Double_t r,
Double_t phimin, Double_t phimax, Double_t theta)
{
Int_t i;
const Int_t np = 200;
static Double_t x[np+3], y[np+3];
Double_t circ = kPi*2*r*(phimax-phimin)/360;
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)*kPi/(180*n);
Double_t ct = TMath::Cos(kPi*theta/180);
Double_t st = TMath::Sin(kPi*theta/180);
for (i=0;i<=n;i++) {
angle = phimin*kPi/180 + Double_t(i)*dphi;
dx = r*TMath::Cos(angle);
dy = r*TMath::Sin(angle);
x[i] = gPad->XtoPad(x1 + dx*ct - dy*st);
y[i] = gPad->YtoPad(y1 + dx*st + dy*ct);
}
gPad->PaintPolyLine(n+1,x,y);
}
void TGraphPolargram::Paint(Option_t* )
{
Int_t i,j;
if (!gPad) return ;
if (!gPad->IsEditable()) gROOT->GetMakeDefCanvas();
gPad->RangeAxis(-1,-1,1,1);
gPad->Range(-1.25,-1.25,1.25,1.25);
Int_t ndivMajor = fNdivRad%100;
Int_t ndivMinor = fNdivRad/100;
for (i=1; i<=ndivMajor; i++) {
TAttLine::Modify();
Double_t rmaj = i*1./ndivMajor;
PaintCircle(0.,0.,rmaj,0.,360,0);
Double_t txtval = fRwrmin+i*(fRwrmax-fRwrmin)/ndivMajor;
SetTextAlign(23);
TAttText::Modify();
gPad->PaintText(rmaj,-fLabelOffset,Form("%5.3g",txtval));
Int_t oldLineStyle = GetLineStyle();
TAttLine::SetLineStyle(2);
TAttLine::Modify();
for (j=1; j<ndivMinor; j++) {
PaintCircle(0.,0.,rmaj- j*1./(ndivMajor*ndivMinor),0.,360,0);
}
SetLineStyle(oldLineStyle);
}
ndivMajor = fNdivPol%100;
ndivMinor = fNdivPol/100;
for (i=0; i<ndivMajor; i++) {
Double_t txtval = fRwtmin + i*(fRwtmax-fRwtmin)/ndivMajor;
Double_t theta = i*2*kPi/ndivMajor;
TAttLine::Modify();
gPad->PaintLine(0.,0.,TMath::Cos(theta),TMath::Sin(theta));
if ((theta < kPi/2)||(theta > 3*kPi/2)) SetTextAlign(12);
else SetTextAlign(32);
if (TMath::Abs(theta - kPi/2)<0.2) SetTextAlign(21);
if (TMath::Abs(theta - 3*kPi/2)<0.2) SetTextAlign(23);
TAttText::Modify();
gPad->PaintText((1+fLabelOffset)*TMath::Cos(theta),
(1+fLabelOffset)*TMath::Sin(theta),
Form("%5.3g",txtval));
Int_t oldLineStyle = GetLineStyle();
TAttLine::SetLineStyle(2);
TAttLine::Modify();
for (j=1; j<ndivMinor; j++) {
Double_t thetamin = theta+j*2*kPi/(ndivMajor*ndivMinor);
gPad->PaintLine(0.,0.,TMath::Cos(thetamin),TMath::Sin(thetamin));
}
SetLineStyle(oldLineStyle);
}
}
ROOT page - Class index - Class Hierarchy - Top of the page
This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.