TEveProjectionAxesGL.cxx

Go to the documentation of this file.
00001 // @(#)root/eve:$Id: TEveProjectionAxesGL.cxx 36384 2010-10-20 14:26:41Z matevz $
00002 // Author: Alja Mrak-Tadel 2007
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2007, Rene Brun and Fons Rademakers.               *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 #include "TEveProjectionAxesGL.h"
00013 #include "TEveProjectionAxes.h"
00014 #include "TEveProjectionManager.h"
00015 #include "THLimitsFinder.h"
00016 
00017 #include "TGLIncludes.h"
00018 #include "TGLRnrCtx.h"
00019 #include "TGLFontManager.h"
00020 #include "TGLCamera.h"
00021 
00022 #include "TMath.h"
00023 
00024 //______________________________________________________________________________
00025 //
00026 // OpenGL renderer class for TEveProjectionAxes.
00027 //
00028 
00029 ClassImp(TEveProjectionAxesGL);
00030 
00031 //______________________________________________________________________________
00032 TEveProjectionAxesGL::TEveProjectionAxesGL() :
00033    TGLObject(),
00034    fM(0),
00035    fProjection(0)
00036 {
00037    // Constructor.
00038 
00039    fDLCache    = kFALSE; // Disable display list.
00040 }
00041 
00042 //______________________________________________________________________________
00043 Bool_t TEveProjectionAxesGL::SetModel(TObject* obj, const Option_t* /*opt*/)
00044 {
00045    // Set model object.
00046    // Virtual from TGLObject.
00047 
00048    fM = SetModelDynCast<TEveProjectionAxes>(obj);
00049    fAxisPainter.SetAttAxis(fM);
00050    return fM->GetManager() ? kTRUE : kFALSE;
00051 }
00052 
00053 //______________________________________________________________________________
00054 void TEveProjectionAxesGL::SetBBox()
00055 {
00056    // Fill the bounding-box data of the logical-shape.
00057    // Virtual from TGLObject.
00058 
00059    SetAxisAlignedBBox(((TEveProjectionAxes*)fExternalObj)->AssertBBox());
00060 }
00061 
00062 //______________________________________________________________________________
00063 void TEveProjectionAxesGL::FilterOverlappingLabels(Int_t idx, Float_t ref) const
00064 {
00065    TGLAxisPainter::LabVec_t &orig = fAxisPainter.RefLabVec();
00066    if (orig.size() == 0) return;
00067 
00068    Float_t center = fM->GetManager()->GetCenter()[idx];
00069 
00070    // Get index of label closest to the distortion center.
00071    // Needed to keep simetry around center.
00072    Int_t minIdx = 0;
00073    Int_t cnt = 0;
00074    Float_t currD = 0;
00075    Float_t minD = TMath::Abs(orig[0].first -center);
00076    for (TGLAxisPainter::LabVec_t::iterator it = orig.begin(); it != orig.end(); ++it)
00077    {
00078       currD = TMath::Abs((*it).first - center);
00079       if (minD > currD)
00080       {
00081          minD = currD;
00082          minIdx = cnt;
00083       }
00084       cnt++;
00085    }
00086 
00087    // Minimum allowed distance 4* font size.
00088    TGLAxisPainter::LabVec_t  filtered;
00089    filtered.push_back(orig[minIdx]);
00090    Int_t size = orig.size();
00091    Float_t minDist = 4*fM->GetLabelSize()*ref;
00092    Float_t pos = 0;
00093 
00094    // Go from center to minimum.
00095    if (minIdx > 0)
00096    {
00097       pos =  orig[minIdx].first;
00098       for (Int_t i=minIdx-1; i>=0; --i)
00099       {
00100          if (TMath::Abs(pos - orig[i].first) > minDist)
00101          {
00102             filtered.push_back(orig[i]);
00103             pos = orig[i].first;
00104          }
00105       }
00106    }
00107 
00108    // Go from center to maximum.
00109    if (minIdx < (size -1))
00110    {
00111       pos =  orig[minIdx].first;
00112       for (Int_t i=minIdx+1; i<size; ++i)
00113       {
00114          if (TMath::Abs(orig[i].first - pos) > minDist)
00115          {
00116             filtered.push_back(orig[i]);
00117             pos = orig[i].first;
00118          }
00119       }
00120    }
00121 
00122    // Set labels list and text format.
00123    if (filtered.size() >= 2)
00124    {
00125       if ( minIdx > 0 )
00126          fAxisPainter.SetTextFormat(orig.front().second, orig.back().second,  orig[minIdx].second - orig[minIdx-1].second);
00127       else
00128          fAxisPainter.SetTextFormat(orig.front().second, orig.back().second,  orig[minIdx+1].second - orig[minIdx].second);
00129 
00130       fAxisPainter.RefLabVec().swap(filtered);
00131    }
00132    else
00133    {
00134       fAxisPainter.SetTextFormat(orig.front().second, orig.back().second,  orig.back().second - orig.front().second);
00135    }
00136 }
00137 
00138 //______________________________________________________________________________
00139 void TEveProjectionAxesGL::SplitInterval(Float_t p1, Float_t p2, Int_t ax) const
00140 {
00141    // Build an array of tick-mark position-value pairs.
00142 
00143    fAxisPainter.RefLabVec().clear();
00144    fAxisPainter.RefTMVec().clear();
00145 
00146    // Get list of label position-value pairs.
00147 
00148 
00149    // Minimum/maximum are defined at the front/back element of list.
00150    fAxisPainter.RefTMVec().push_back(TGLAxisPainter::TM_t(p1, -1));
00151 
00152    if (fM->GetLabMode() == TEveProjectionAxes::kValue)
00153    {
00154       SplitIntervalByVal(p1, p2, ax);
00155    }
00156    else if (fM->GetLabMode() == TEveProjectionAxes::kPosition)
00157    {
00158       SplitIntervalByPos(p1, p2, ax);
00159    }
00160 
00161 
00162    FilterOverlappingLabels(0, p2 -p1);
00163 
00164    // Minimum/maximum are defined at the front/back element of list.
00165    fAxisPainter.RefTMVec().push_back(TGLAxisPainter::TM_t(p2, -1));
00166 }
00167 
00168 //______________________________________________________________________________
00169 void TEveProjectionAxesGL::SplitIntervalByPos(Float_t p1, Float_t p2, Int_t ax) const
00170 {
00171    // Add tick-marks at equidistant position.
00172 
00173    // Limits.
00174    Int_t n1a = TMath::FloorNint(fM->GetNdivisions() / 100);
00175    Int_t n2a = fM->GetNdivisions() - n1a * 100;
00176    Int_t bn1, bn2;
00177    Double_t bw1, bw2; // bin with first second order
00178    Double_t bl1, bh1, bl2, bh2; // bin low, high first second order
00179    THLimitsFinder::Optimize(p1, p2, n1a, bl1, bh1, bn1, bw1);
00180    THLimitsFinder::Optimize(bl1, bl1+bw1, n2a, bl2, bh2, bn2, bw2);
00181 
00182    Int_t n1=TMath::CeilNint(p1/bw1);
00183    Int_t n2=TMath::FloorNint(p2/bw1);
00184 
00185    TGLAxisPainter::LabVec_t &labVec = fAxisPainter.RefLabVec();
00186    TGLAxisPainter::TMVec_t  &tmVec =  fAxisPainter.RefTMVec();
00187 
00188    Float_t p = n1*bw1;
00189    Float_t pMinor = p;
00190    for (Int_t l=n1; l<=n2; l++)
00191    {
00192       // Labels.
00193       labVec.push_back( TGLAxisPainter::Lab_t(p , fProjection->GetValForScreenPos(ax, p)));
00194 
00195       // Tick-marks.
00196       tmVec.push_back(TGLAxisPainter::TM_t(p, 0));
00197       pMinor = p+bw2;
00198       for (Int_t i=1; i<bn2; i++)
00199       {
00200          if (pMinor > p2)  break;
00201          tmVec.push_back( TGLAxisPainter::TM_t(pMinor, 1));
00202          pMinor += bw2;
00203       }
00204       p += bw1;
00205    }
00206 
00207    // Complete second order tick-marks.
00208    pMinor = n1*bw1 -bw2;
00209    while ( pMinor > p1)
00210    {
00211       tmVec.push_back(TGLAxisPainter::TM_t(pMinor, 1));
00212       pMinor -=bw2;
00213    }
00214 }
00215 
00216 //______________________________________________________________________________
00217 void TEveProjectionAxesGL::SplitIntervalByVal(Float_t p1, Float_t p2, Int_t ax) const
00218 {
00219    // Add tick-marks on fixed value step.
00220 
00221    Float_t v1 = fProjection->GetValForScreenPos(ax, p1);
00222    Float_t v2 = fProjection->GetValForScreenPos(ax, p2);
00223 
00224    TGLAxisPainter::LabVec_t &labVec =  fAxisPainter.RefLabVec();
00225    TGLAxisPainter::TMVec_t  &tmVec  =  fAxisPainter.RefTMVec();
00226 
00227    // Limits
00228    Int_t n1a = TMath::FloorNint(fM->GetNdivisions() / 100);
00229    Int_t n2a = fM->GetNdivisions() - n1a * 100;
00230    Int_t bn1, bn2;
00231    Double_t bw1, bw2;           // bin width first / second order
00232    Double_t bl1, bh1, bl2, bh2; // bin low, high first / second order
00233    THLimitsFinder::Optimize(v1, v2, n1a, bl1, bh1, bn1, bw1);
00234    THLimitsFinder::Optimize(bl1, bl1+bw1, n2a, bl2, bh2, bn2, bw2);
00235 
00236    Float_t pFirst, pSecond; // position of first, second order of tickmarks
00237    Float_t v = bl1;
00238    // step
00239    for (Int_t l=0; l<=bn1; l++)
00240    {
00241       // Labels.
00242       pFirst = fProjection->GetScreenVal(ax, v);
00243       labVec.push_back(TGLAxisPainter::Lab_t(pFirst , v));
00244       tmVec.push_back(TGLAxisPainter::TM_t(pFirst, 0));
00245 
00246       // Tickmarks.
00247       for (Int_t k=1; k<bn2; k++)
00248       {
00249          pSecond = fProjection->GetScreenVal(ax, v+k*bw2);
00250          if (pSecond > p2)  break;
00251          tmVec.push_back(TGLAxisPainter::TM_t(pSecond, 1));
00252       }
00253       v += bw1;
00254    }
00255 
00256    // Complete second order tick-marks.
00257    v = bl1 -bw2;
00258    while ( v > v1)
00259    {
00260       pSecond = fProjection->GetScreenVal(ax, v);
00261       if (pSecond < p1)  break;
00262       tmVec.push_back(TGLAxisPainter::TM_t(pSecond, 1));
00263       v -= bw2;
00264    }
00265 }
00266 
00267 //______________________________________________________________________________
00268 void TEveProjectionAxesGL::GetRange(Int_t ax, Float_t frustMin, Float_t frustMax, Float_t& min, Float_t& max) const
00269 {
00270    // Get range from bounding box of projection manager
00271 
00272 
00273    // Compare frustum range with bbox and take larger.
00274 
00275    Float_t frng = (frustMax -frustMin)*0.4;
00276    Float_t c = 0.5*(frustMax +frustMin);
00277    min = c - frng;
00278    max = c + frng;
00279 
00280    // Check projection  limits.
00281    // Set limit factor in case of divergence.
00282    Float_t dLim = fProjection->GetLimit(ax, 0);
00283    Float_t uLim = fProjection->GetLimit(ax, 1);
00284    if (min < dLim) min = dLim*0.98;
00285    if (max > uLim) max   = uLim*0.98;
00286 }
00287 
00288 //______________________________________________________________________________
00289 void TEveProjectionAxesGL::Draw(TGLRnrCtx& rnrCtx) const
00290 {
00291    // Draw function for TEveProjectionAxesGL. Skips line-pass of outline mode.
00292 
00293    if (rnrCtx.IsDrawPassOutlineLine())
00294       return;
00295 
00296    TGLObject::Draw(rnrCtx);
00297 }
00298 
00299 //______________________________________________________________________________
00300 void TEveProjectionAxesGL::DirectDraw(TGLRnrCtx& rnrCtx) const
00301 {
00302    // Actual rendering code.
00303    // Virtual from TGLLogicalShape.
00304 
00305    if (rnrCtx.Selection() || rnrCtx.Highlight()) return;
00306 
00307    glPushAttrib(GL_ENABLE_BIT | GL_POLYGON_BIT);
00308 
00309    glDisable(GL_LIGHTING);
00310    glColorMaterial(GL_FRONT_AND_BACK, GL_DIFFUSE);
00311    glEnable(GL_COLOR_MATERIAL);
00312    glPolygonMode(GL_FRONT_AND_BACK, GL_FILL);
00313    glDisable(GL_CULL_FACE);
00314 
00315    // Draw on front-clipping plane.
00316    Float_t old_depth_range[2];
00317    glGetFloatv(GL_DEPTH_RANGE, old_depth_range);
00318    glDepthRange(0, 0.001);
00319 
00320    // Frustum size.
00321    TGLCamera &camera = rnrCtx.RefCamera();
00322    Float_t l = -camera.FrustumPlane(TGLCamera::kLeft).D();
00323    Float_t r =  camera.FrustumPlane(TGLCamera::kRight).D();
00324    Float_t t =  camera.FrustumPlane(TGLCamera::kTop).D();
00325    Float_t b = -camera.FrustumPlane(TGLCamera::kBottom).D();
00326 
00327    if (fM->fUseColorSet)
00328    {
00329        TGLUtil::Color(rnrCtx.ColorSet().Markup());
00330        fAxisPainter.SetUseAxisColors(kFALSE);
00331    }
00332 
00333    fProjection = fM->GetManager()->GetProjection();
00334    glDisable(GL_LIGHTING);
00335    // Projection center and origin marker.
00336    {
00337       Float_t d = ((r-l) > (b-t)) ? (b-t) : (r-l);
00338       d *= 0.02f;
00339       if (fM->GetDrawCenter())
00340       {
00341          Float_t* c = fProjection->GetProjectedCenter();
00342          TGLUtil::LineWidth(1);
00343          glBegin(GL_LINES);
00344          glVertex3f(c[0] + d, c[1], c[2]); glVertex3f(c[0] - d, c[1], c[2]);
00345          glVertex3f(c[0], c[1] + d, c[2]); glVertex3f(c[0], c[1] - d, c[2]);
00346          glVertex3f(c[0], c[1], c[2] + d); glVertex3f(c[0], c[1], c[2] - d);
00347          glEnd();
00348       }
00349       if (fM->GetDrawOrigin())
00350       {
00351          TEveVector zero;
00352          fProjection->ProjectVector(zero, 0);
00353          TGLUtil::LineWidth(1);
00354          glBegin(GL_LINES);
00355          glVertex3f(zero[0] + d, zero[1], zero[2]); glVertex3f(zero[0] - d, zero[1], zero[2]);
00356          glVertex3f(zero[0], zero[1] + d, zero[2]); glVertex3f(zero[0], zero[1] - d, zero[2]);
00357          glVertex3f(zero[0], zero[1], zero[2] + d); glVertex3f(zero[0], zero[1], zero[2] - d);
00358          glEnd();
00359       }
00360    }
00361 
00362    //
00363    // Axes.
00364    {
00365       using namespace TMath;
00366       GLint   vp[4];
00367       glGetIntegerv(GL_VIEWPORT, vp);
00368       Float_t refLength =  TMath::Sqrt((TMath::Power(vp[2]-vp[0], 2) + TMath::Power(vp[3]-vp[1], 2)));
00369       Float_t tickLength = TMath::Sqrt((TMath::Power(r-l, 2) + TMath::Power(t-b, 2)));
00370       fAxisPainter.SetFontMode(TGLFont::kPixmap);
00371       fAxisPainter.SetLabelFont(rnrCtx, TGLFontManager::GetFontNameFromId(fM->GetLabelFont()),  TMath::CeilNint(refLength*0.02), tickLength*fM->GetLabelSize());
00372 
00373       Float_t min, max;
00374       // X-axis.
00375       if (fM->fAxesMode == TEveProjectionAxes::kAll ||
00376           fM->fAxesMode == TEveProjectionAxes::kHorizontal)
00377       {
00378          GetRange(0, l, r, min, max);
00379          SplitInterval(min, max, 0);
00380 
00381          fAxisPainter.RefDir().Set(1, 0, 0);
00382          fAxisPainter.RefTMOff(0).Set(0, tickLength, 0);
00383 
00384          // Bottom.
00385          glPushMatrix();
00386          glTranslatef( 0, b, 0);
00387          fAxisPainter.SetLabelAlign(TGLFont::kCenterH, TGLFont::kTop);
00388          fAxisPainter.RnrLabels();
00389          fAxisPainter.RnrLines();
00390          glPopMatrix();
00391 
00392          // Top.
00393          glPushMatrix();
00394          glTranslatef( 0, t, 0);
00395          fAxisPainter.SetLabelAlign(TGLFont::kCenterH, TGLFont::kBottom);
00396          fAxisPainter.RefTMOff(0).Negate();
00397          fAxisPainter.RnrLabels();
00398          fAxisPainter.RnrLines();
00399          glPopMatrix();
00400       }
00401 
00402       // Y-axis.
00403       if (fM->fAxesMode == TEveProjectionAxes::kAll ||
00404           fM->fAxesMode == TEveProjectionAxes::kVertical)
00405       {
00406          GetRange(1, b, t, min, max);
00407          SplitInterval(min, max, 1);
00408 
00409          fAxisPainter.RefDir().Set(0, 1, 0);
00410          fAxisPainter.RefTMOff(0).Set(tickLength, 0 , 0);
00411 
00412          // Left.
00413          glPushMatrix();
00414          glTranslatef(l, 0, 0);
00415          fAxisPainter.SetLabelAlign(TGLFont::kLeft, TGLFont::kCenterV);
00416          fAxisPainter.RnrLabels();
00417          fAxisPainter.RnrLines();
00418          glPopMatrix();
00419 
00420          // Right.
00421          glPushMatrix();
00422          glTranslatef(r, 0, 0);
00423          fAxisPainter.SetLabelAlign(TGLFont::kRight, TGLFont::kCenterV);
00424          fAxisPainter.RefTMOff(0).Negate();
00425          fAxisPainter.RnrLabels();
00426          fAxisPainter.RnrLines();
00427          glPopMatrix();
00428       }
00429    }
00430    glDepthRange(old_depth_range[0], old_depth_range[1]);
00431 
00432    glPopAttrib();
00433 }

Generated on Tue Jul 5 14:16:13 2011 for ROOT_528-00b_version by  doxygen 1.5.1