#include "Riostream.h"
#include "TROOT.h"
#include "TGLIncludes.h"
#include "TGLAxis.h"
#include "TGLText.h"
#include "TColor.h"
#include "TString.h"
#include "TMath.h"
#include "THLimitsFinder.h"
ClassImp(TGLAxis)
TGLAxis::TGLAxis(): TAttLine(1,1,1), TAttText(20,0.,1,42,0.04)
{
Init();
}
void TGLAxis::Init()
{
fTicks1 = 0;
fTicks2 = 0;
fLabels = 0;
fText = 0;
fAngle1 = 90.;
fAngle2 = 0.;
fAngle3 = 0.;
fAxisLength = 0.;
fTickMarksLength = 0.04;
fTickMarksOrientation = 2;
fLabelsOffset = 0.09;
fLabelsSize = 0.06;
fGridLength = 0.;
}
TGLAxis::~TGLAxis()
{
if (fTicks1) delete [] fTicks1;
if (fTicks2) delete [] fTicks2;
if (fLabels) delete [] fLabels;
if (fText) delete fText;
}
void TGLAxis::PaintGLAxis(const Double_t p1[3], const Double_t p2[3],
Double_t wmin, Double_t wmax, Int_t ndiv,
Option_t *opt)
{
fNDiv = ndiv;
if (wmax<=wmin) {
fWmax = wmin;
fWmin = wmax;
} else {
fWmax = wmax;
fWmin = wmin;
}
Double_t x1 = p1[0], y1 = p1[1], z1 = p1[2];
Double_t x2 = p2[0], y2 = p2[1], z2 = p2[2];
fAxisLength = TMath::Sqrt((x2-x1)*(x2-x1)+
(y2-y1)*(y2-y1)+
(z2-z1)*(z2-z1));
TicksPositions(opt);
DoLabels();
glPushMatrix();
glTranslatef(x1, y1, z1);
Double_t phi=0;
Double_t normal[3];
normal[0] = 0;
normal[1] = 1;
normal[2] = 0;
if (z1!=z2) {
if (y2==y1 && x2==x1) {
if (z2<z1) phi = 90;
else phi = 270;
} else {
Double_t p3[3];
p3[0] = p2[0]; p3[1] = p2[1]; p3[2] = 0;
TMath::Normal2Plane(p1, p2, p3, normal);
phi = TMath::ACos(TMath::Abs(z2-z1)/fAxisLength);
phi = -(90-(180/TMath::Pi())*phi);
}
glRotatef(phi, normal[0], normal[1], normal[2]);
}
Double_t theta = 0;
if (y2!=y1) {
if ((x2-x1)>0) {
theta = TMath::ATan((y2-y1)/(x2-x1));
theta = (180/TMath::Pi())*theta;
} else if ((x2-x1)<0) {
theta = TMath::ATan((y2-y1)/(x2-x1));
theta = 180+(180/TMath::Pi())*theta;
} else {
if (y2>y1) theta = 90;
else theta = 270;
}
} else {
if (x2<x1) theta = 180;
}
glRotatef(theta, 0., 0., 1.);
PaintGLAxisBody();
PaintGLAxisTickMarks();
PaintGLAxisLabels();
glPopMatrix();
}
void TGLAxis::PaintGLAxisBody()
{
TColor *col;
Float_t red, green, blue;
col = gROOT->GetColor(GetLineColor());
col->GetRGB(red, green, blue);
glColor3d(red, green, blue);
glLineWidth(GetLineWidth());
glBegin(GL_LINES);
glVertex3d(0., 0., 0.);
glVertex3d(fAxisLength, 0., 0.);
glEnd();
}
void TGLAxis::PaintGLAxisTickMarks()
{
Int_t i;
Double_t yo=0, zo=0;
switch (fTickMarksOrientation) {
case 0:
yo = 0;
zo = 1;
break;
case 1:
yo = -1;
zo = 0;
break;
case 2:
yo = 0;
zo = -1;
break;
case 3:
yo = 1;
zo = 0;
break;
}
if (fTicks1) {
if (fTickMarksLength) {
Double_t tl = fTickMarksLength*fAxisLength;
glBegin(GL_LINES);
for (i=0; i<fNTicks1 ; i++) {
glVertex3f( fTicks1[i], 0, 0);
glVertex3f( fTicks1[i], yo*tl, zo*tl);
}
glEnd();
}
if (fGridLength) {
const UShort_t stipple = 0x8888;
glLineStipple(1, stipple);
glEnable(GL_LINE_STIPPLE);
glBegin(GL_LINES);
for (i=0; i<fNTicks1; i++) {
glVertex3f( fTicks1[i], 0, 0);
glVertex3f( fTicks1[i], -yo*fGridLength, -zo*fGridLength);
}
glEnd();
glDisable(GL_LINE_STIPPLE);
}
}
if (fTicks2) {
if (fTickMarksLength) {
Double_t tl = 0.5*fTickMarksLength*fAxisLength;
glBegin(GL_LINES);
for (i=0; i<fNTicks2; i++) {
glVertex3f( fTicks2[i], 0, 0);
glVertex3f( fTicks2[i], yo*tl, zo*tl);
}
glEnd();
}
}
}
void TGLAxis::PaintGLAxisLabels()
{
if (!fLabelsSize) return;
Double_t x=0,y=0,z=0;
Int_t i;
if (!fText) {
fText = new TGLText();
fText->SetTextColor(GetTextColor());
fText->SetGLTextFont(GetTextFont());
fText->SetTextSize(fLabelsSize*fAxisLength);
fText->SetTextAlign(GetTextAlign());
}
fText->SetGLTextAngles(fAngle1, fAngle2, fAngle3);
switch (fTickMarksOrientation) {
case 0:
y = 0;
z = fLabelsOffset*fAxisLength;
break;
case 1:
y = -fLabelsOffset*fAxisLength;
z = 0;
break;
case 2:
y = 0;
z = -fLabelsOffset*fAxisLength;
break;
case 3:
y = fLabelsOffset*fAxisLength;
z = 0;
break;
}
for (i=0; i<fNDiv1+1 ; i++) {
x = fTicks1[i];
fText->PaintGLText(x,y,z,fLabels[i]);
}
}
void TGLAxis::TicksPositions(Option_t *opt)
{
Bool_t optionNoopt, optionLog;
if (strchr(opt,'N')) optionNoopt = kTRUE; else optionNoopt = kFALSE;
if (strchr(opt,'G')) optionLog = kTRUE; else optionLog = kFALSE;
fNDiv3 = fNDiv/10000;
fNDiv2 = (fNDiv-10000*fNDiv3)/100;
fNDiv1 = fNDiv%100;
if (fTicks1) {
delete [] fTicks1;
fTicks1 = 0;
}
if (fTicks2) {
delete [] fTicks2;
fTicks2 = 0;
}
if (optionNoopt) {
TicksPositionsNoOpt();
} else {
TicksPositionsOpt();
}
}
void TGLAxis::TicksPositionsNoOpt()
{
Int_t i, j, k;
Double_t step1 = fAxisLength/(fNDiv1);
fNTicks1 = fNDiv1+1;
fTicks1 = new Double_t[fNTicks1];
for (i=0; i<fNTicks1; i++) fTicks1[i] = i*step1;
if (fNDiv2) {
Double_t t2;
Double_t step2 = step1/fNDiv2;
fNTicks2 = fNDiv1*(fNDiv2-1);
fTicks2 = new Double_t[fNTicks2];
k = 0;
for (i=0; i<fNTicks1-1; i++) {
t2 = fTicks1[i]+step2;
for (j=0; j<fNDiv2-1; j++) {
fTicks2[k] = t2;
k++;
t2 = t2+step2;
}
}
}
}
void TGLAxis::TicksPositionsOpt()
{
Int_t i, j, k, nDivOpt;
Double_t step1, step2, wmin2, wmax2;
Double_t wmin = fWmin;
Double_t wmax = fWmax;
THLimitsFinder::Optimize(wmin, wmax, fNDiv1,
fWmin, fWmax, nDivOpt,
step1, "");
fNDiv1 = nDivOpt;
fNTicks1 = fNDiv1+1;
fTicks1 = new Double_t[fNTicks1];
Double_t r = fAxisLength/(wmax-wmin);
Double_t w = fWmin;
i = 0;
while (w<=fWmax) {
fTicks1[i] = r*(w-wmin);
i++;
w = w+step1;
}
if (fNDiv2) {
Double_t t2;
THLimitsFinder::Optimize(fWmin, fWmin+step1, fNDiv2,
wmin2, wmax2, nDivOpt,
step2, "");
fNDiv2 = nDivOpt;
step2 = TMath::Abs((fTicks1[1]-fTicks1[0])/fNDiv2);
Int_t nTickl = (Int_t)(fTicks1[0]/step2);
Int_t nTickr = (Int_t)((fAxisLength-fTicks1[fNTicks1-1])/step2);
fNTicks2 = fNDiv1*(fNDiv2-1)+nTickl+nTickr;
fTicks2 = new Double_t[fNTicks2];
k = 0;
for (i=0; i<fNTicks1-1; i++) {
t2 = fTicks1[i]+step2;
for (j=0; j<fNDiv2-1; j++) {
fTicks2[k] = t2;
k++;
t2 = t2+step2;
}
}
if (nTickl) {
t2 = fTicks1[0]-step2;
for (i=0; i<nTickl; i++) {
fTicks2[k] = t2;
k++;
t2 = t2-step2;
}
}
if (nTickr) {
t2 = fTicks1[fNTicks1-1]+step2;
for (i=0; i<nTickr; i++) {
fTicks2[k] = t2;
k++;
t2 = t2+step2;
}
}
}
}
void TGLAxis::DoLabels()
{
if (fLabels) delete [] fLabels;
fLabels = new TString[fNTicks1];
Int_t i;
Double_t dw = (fWmax-fWmin)/(fNDiv1);
for (i=0; i<fNTicks1; i++) {
fLabels[i] = Form("%g",fWmin+i*dw);
}
}
void TGLAxis::SetLabelsAngles(Double_t a1, Double_t a2, Double_t a3)
{
fAngle1 = a1;
fAngle2 = a2;
fAngle3 = a3;
}
Last change: Wed Jun 25 08:40:43 2008
Last generated: 2008-06-25 08:40
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.