#include "TEveProjectionAxesGL.h"
#include "TEveProjectionAxes.h"
#include "TEveProjectionManager.h"
#include "TGLRnrCtx.h"
#include "TGLIncludes.h"
#include "TMath.h"
ClassImp(TEveProjectionAxesGL);
TEveProjectionAxesGL::TEveProjectionAxesGL() :
TEveTextGL(),
fRange(300),
fLabelSize(0.02),
fLabelOff(0.5),
fTMSize(0.02),
fAxesModel(0),
fProjection(0)
{
fDLCache = kFALSE;
fMultiColor = kTRUE;
}
Bool_t TEveProjectionAxesGL::SetModel(TObject* obj, const Option_t* )
{
if (SetModelCheckClass(obj, TEveProjectionAxes::Class())) {
fAxesModel = dynamic_cast<TEveProjectionAxes*>(obj);
if (fAxesModel->GetManager() == 0)
return kFALSE;
TEveTextGL::fM = fAxesModel;
return kTRUE;
}
return kFALSE;
}
void TEveProjectionAxesGL::SetBBox()
{
SetAxisAlignedBBox(((TEveProjectionAxes*)fExternalObj)->AssertBBox());
}
void TEveProjectionAxesGL::RenderText(const char* txt, Float_t x, Float_t y) const
{
if (fFont.GetMode() < TGLFont::kTexture) {
glRasterPos3f(0, 0, 0);
glBitmap(0, 0, 0, 0, x, y, 0);
fFont.Render(txt);
} else {
glPushMatrix();
glTranslatef(x, y, 0);
fFont.Render(txt);
glPopMatrix();
}
}
void TEveProjectionAxesGL::DrawTickMarks(Float_t y) const
{
if (fFont.GetMode() == TGLFont::kTexture) glDisable(GL_TEXTURE_2D);
glBegin(GL_LINES);
for (std::list<TM_t>::iterator it = fTMList.begin(); it != fTMList.end(); ++it)
{
glVertex2f((*it).first, 0);
glVertex2f((*it).first, y);
}
glEnd();
if (fFont.GetMode() == TGLFont::kTexture) glEnable(GL_TEXTURE_2D);
}
void TEveProjectionAxesGL::DrawHInfo() const
{
Float_t tmH = fTMSize*fRange;
DrawTickMarks(-tmH);
Float_t off = tmH + fLabelOff*tmH;
Float_t llx, lly, llz, urx, ury, urz;
const char* txt;
for (std::list<TM_t>::iterator it = fTMList.begin(); it != fTMList.end(); ++it)
{
glPushMatrix();
glTranslatef((*it).first, -off, 0);
txt = TEveUtil::FormAxisValue((*it).second);
fFont.BBox(txt, llx, lly, llz, urx, ury, urz);
Float_t xd = -0.5f*(urx+llx);
if (txt[0] == '-')
xd -= 0.5f * (urx-llx) / strlen(txt);
RenderText(txt, xd, -ury);
glPopMatrix();
}
fTMList.clear();
}
void TEveProjectionAxesGL::DrawVInfo() const
{
Float_t tmH = fTMSize*fRange;
glPushMatrix();
glRotatef(90, 0, 0, 1);
DrawTickMarks(tmH);
glPopMatrix();
Float_t off = fLabelOff*tmH + tmH;
Float_t llx, lly, llz, urx, ury, urz;
const char* txt;
for (std::list<TM_t>::iterator it = fTMList.begin(); it!= fTMList.end(); ++it)
{
glPushMatrix();
glTranslatef(-off, (*it).first, 0);
txt = TEveUtil::FormAxisValue((*it).second);
fFont.BBox(txt, llx, lly, llz, urx, ury, urz);
RenderText(txt, -urx, -0.38f*(ury+lly));
glPopMatrix();
}
fTMList.clear();
}
void TEveProjectionAxesGL::SplitInterval(Float_t p1, Float_t p2, Int_t ax) const
{
Float_t down = fProjection->GetLimit(ax, kFALSE)*0.95;
p1 = TMath::Max(p1, down);
Float_t up = fProjection->GetLimit(ax, kTRUE)*0.95;
p2 = TMath::Min(p2, up);
if (fAxesModel->GetStepMode() == TEveProjectionAxes::kValue)
{
SplitIntervalByVal(fProjection->GetValForScreenPos(ax, p1), fProjection->GetValForScreenPos(ax, p2), ax);
}
else if (fAxesModel->GetStepMode() == TEveProjectionAxes::kPosition)
{
SplitIntervalByPos(p1, p2, ax);
}
}
void TEveProjectionAxesGL::SplitIntervalByPos(Float_t p1, Float_t p2, Int_t ax) const
{
TEveVector zeroPos;
fProjection->ProjectVector(zeroPos);
Float_t step = (p2-p1)/fAxesModel->GetNumTickMarks();
Float_t p = zeroPos.fX;
while (p > p1) {
fTMList.push_back(TM_t(p , fProjection->GetValForScreenPos(ax, p)));
p -= step;
}
p = zeroPos.fX + step;
while (p < p2) {
fTMList.push_back(TM_t(p , fProjection->GetValForScreenPos(ax, p)));
p += step;
}
}
void TEveProjectionAxesGL::SplitIntervalByVal(Float_t v1, Float_t v2, Int_t ax) const
{
Float_t step = (v2-v1)/fAxesModel->GetNumTickMarks();
Float_t v = 0.f;
while (v > v1) {
fTMList.push_back(TM_t(fProjection->GetScreenVal(ax, v) , v));
v -= step;
}
v = step;
while (v < v2) {
fTMList.push_back(TM_t(fProjection->GetScreenVal(ax, v) , v));
v += step;
}
}
void TEveProjectionAxesGL::DirectDraw(TGLRnrCtx& rnrCtx) const
{
fProjection = fAxesModel->GetManager()->GetProjection();
Float_t* bbox = fAxesModel->GetBBox();
fRange = bbox[1] - bbox[0];
TEveVector zeroPos;
fProjection->ProjectVector(zeroPos);
glBegin(GL_LINES);
glVertex3f(bbox[0], bbox[2], 0.);
glVertex3f(bbox[1], bbox[2], 0.);
glVertex3f(bbox[0], bbox[2], 0.);
glVertex3f(bbox[0], bbox[3], 0.);
glEnd();
Float_t d = 10;
if (fAxesModel->GetDrawCenter()) {
Float_t* c = fProjection->GetProjectedCenter();
TGLUtil::Color3f(1., 0., 0.);
glBegin(GL_LINES);
glVertex3f(c[0] +d, c[1], c[2]); glVertex3f(c[0] - d, c[1] , c[2]);
glVertex3f(c[0] , c[1] +d, c[2]); glVertex3f(c[0] , c[1] -d, c[2]);
glVertex3f(c[0] , c[1], c[2] + d); glVertex3f(c[0] , c[1] , c[2] - d);
glEnd();
}
if (fAxesModel->GetDrawOrigin()) {
TEveVector zero;
fProjection->ProjectVector(zero);
TGLUtil::Color3f(1., 1., 1.);
glBegin(GL_LINES);
glVertex3f(zero[0] + d, zero[1], zero[2]); glVertex3f(zero[0] - d, zero[1] , zero[2]);
glVertex3f(zero[0] , zero[1] + d, zero[2]); glVertex3f(zero[0] , zero[1] - d, zero[2]);
glVertex3f(zero[0] , zero[1], zero[2] + d); glVertex3f(zero[0] , zero[1] , zero[2] - d);
glEnd();
}
SetFont(rnrCtx);
fFont.PreRender(fAxesModel->GetAutoLighting(), fAxesModel->GetLighting());
glPushMatrix();
glTranslatef(zeroPos.fX, bbox[3]*1.1, 0);
Float_t llx, lly, llz, urx, ury, urz;
fFont.BBox(fAxesModel->GetText(), llx, lly, llz, urx, ury, urz);
RenderText(fAxesModel->GetText(), -llx, 0);
glPopMatrix();
glPushMatrix();
glTranslatef(0, bbox[2], 0);
SplitInterval(bbox[0], bbox[1], 0);
DrawHInfo();
glPopMatrix();
glPushMatrix();
glTranslatef(bbox[0], 0, 0);
SplitInterval(bbox[2], bbox[3], 1);
DrawVInfo();
glPopMatrix();
fFont.PostRender();
fProjection = 0;
}
Last change: Wed Jun 25 08:37:48 2008
Last generated: 2008-06-25 08:37
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.