TGLAxis.cxx

Go to the documentation of this file.
00001 // @(#)root/gl:$Id: TGLAxis.cxx 35580 2010-09-22 13:27:00Z couet $
00002 // Author:  Olivier Couet  17/04/2007
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2006, 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 #include "Riostream.h"
00012 #include "TROOT.h"
00013 
00014 #include "TGLIncludes.h"
00015 #include "TGLUtil.h"
00016 #include "TGLAxis.h"
00017 #include "TGLText.h"
00018 #include "TColor.h"
00019 #include "TString.h"
00020 #include "TMath.h"
00021 #include "THLimitsFinder.h"
00022 
00023 //______________________________________________________________________________
00024 /* Begin_Html
00025 <center><h2>GL Axis</h2></center>
00026 To draw a 3D axis in a GL window. The labels are drawn using FTGL.
00027 End_Html */
00028 
00029 ClassImp(TGLAxis)
00030 
00031 //______________________________________________________________________________
00032 TGLAxis::TGLAxis(): TAttLine(1,1,1), TAttText(20,0.,1,42,0.04)
00033 {
00034    // Constructor.
00035 
00036    Init();
00037 }
00038 
00039 
00040 //______________________________________________________________________________
00041 void TGLAxis::Init()
00042 {
00043    // Default initialization.
00044 
00045    fNDiv = fNDiv1 = fNDiv2 = fNDiv3 = 0;
00046    fNTicks1 = fNTicks2 = 0;
00047    fTicks1          = 0;
00048    fTicks2          = 0;
00049    fLabels          = 0;
00050    fText            = 0;
00051    fAngle1          = 90.;
00052    fAngle2          = 0.;
00053    fAngle3          = 0.;
00054    fAxisLength      = 0.;
00055    fWmin = fWmax    = 0.;
00056    fTickMarksLength = 0.04; // % of fAxisLength
00057    fTickMarksOrientation = 2; // can be 0, 1, 2, or 3
00058    fLabelsOffset    = 0.09; // % of fAxisLength
00059    fLabelsSize      = 0.06; // % of fAxisLength
00060    fGridLength      = 0.;   // 0 means no grid
00061 }
00062 
00063 
00064 //______________________________________________________________________________
00065 TGLAxis::~TGLAxis()
00066 {
00067    // Destructor.
00068 
00069    if (fTicks1) delete [] fTicks1;
00070    if (fTicks2) delete [] fTicks2;
00071    if (fLabels) delete [] fLabels;
00072    if (fText)   delete fText;
00073 }
00074 
00075 
00076 //______________________________________________________________________________
00077 void TGLAxis::PaintGLAxis(const Double_t p1[3], const Double_t p2[3],
00078                           Double_t wmin,  Double_t wmax, Int_t ndiv,
00079                           Option_t *opt)
00080 {
00081    // Paint GL Axis.
00082    //
00083    // p1, p2     : Axis position in the 3D space.
00084    // wmin, wmax : Minimum and maximum values along the axis. wmin < wmax.
00085    // ndiv       : Number of axis divisions. It is an integer in the form
00086    //              "ttsspp" where "tt" is the number of tertiary divisions,
00087    //              "ss" is the number of secondary divisions and "pp" the
00088    //              number of primary divisions.
00089    // opt        : Options.
00090    //              "N" - By default the number of divisions is optimized to
00091    //                    get a nice labeling. When option "N" is given, the
00092    //                    number of divisions is not optimized.
00093 
00094    fNDiv = ndiv;
00095    if (wmax<=wmin) {
00096       fWmax = wmin;
00097       fWmin = wmax;
00098    } else {
00099       fWmax = wmax;
00100       fWmin = wmin;
00101    }
00102 
00103    // Compute the axis length.
00104    Double_t x1 = p1[0], y1 = p1[1], z1 = p1[2];
00105    Double_t x2 = p2[0], y2 = p2[1], z2 = p2[2];
00106    fAxisLength = TMath::Sqrt((x2-x1)*(x2-x1)+
00107                              (y2-y1)*(y2-y1)+
00108                              (z2-z1)*(z2-z1));
00109 
00110    TicksPositions(opt);
00111 
00112    DoLabels();
00113 
00114    // Paint using GL
00115    glPushMatrix();
00116 
00117    // Translation
00118    glTranslatef(x1, y1, z1);
00119 
00120    // Rotation in the Z plane
00121    Double_t phi=0;
00122    Double_t normal[3];
00123    normal[0] = 0;
00124    normal[1] = 1;
00125    normal[2] = 0;
00126    if (z1!=z2) {
00127       if (y2==y1 && x2==x1) {
00128          if (z2<z1) phi = 90;
00129          else       phi = 270;
00130       } else {
00131          Double_t p3[3];
00132          p3[0] = p2[0]; p3[1] = p2[1]; p3[2] = 0;
00133          TMath::Normal2Plane(p1, p2, p3, normal);
00134          phi = TMath::ACos(TMath::Abs(z2-z1)/fAxisLength);
00135          phi = -(90-(180/TMath::Pi())*phi);
00136       }
00137       glRotatef(phi, normal[0], normal[1], normal[2]);
00138    }
00139 
00140    // Rotation in the XY plane
00141    Double_t theta = 0;
00142    if (y2!=y1) {
00143       if ((x2-x1)>0) {
00144          theta = TMath::ATan((y2-y1)/(x2-x1));
00145          theta = (180/TMath::Pi())*theta;
00146       } else if ((x2-x1)<0) {
00147          theta = TMath::ATan((y2-y1)/(x2-x1));
00148          theta = 180+(180/TMath::Pi())*theta;
00149       } else {
00150          if (y2>y1) theta = 90;
00151          else       theta = 270;
00152       }
00153    } else {
00154       if (x2<x1) theta = 180;
00155    }
00156    glRotatef(theta, 0., 0., 1.);
00157 
00158    PaintGLAxisBody();
00159 
00160    PaintGLAxisTickMarks();
00161 
00162    PaintGLAxisLabels();
00163 
00164    glPopMatrix();
00165 }
00166 
00167 
00168 //______________________________________________________________________________
00169 void TGLAxis::PaintGLAxisBody()
00170 {
00171    // Paint horizontal axis body at position (0,0,0)
00172 
00173    TColor *col;
00174    Float_t red, green, blue;
00175    col = gROOT->GetColor(GetLineColor());
00176    col->GetRGB(red, green, blue);
00177    glColor3d(red, green, blue);
00178    TGLUtil::LineWidth(GetLineWidth());
00179    glBegin(GL_LINES);
00180    glVertex3d(0., 0., 0.);
00181    glVertex3d(fAxisLength, 0., 0.);
00182    glEnd();
00183 }
00184 
00185 
00186 //______________________________________________________________________________
00187 void TGLAxis::PaintGLAxisTickMarks()
00188 {
00189    // Paint axis tick marks.
00190 
00191    Int_t i;
00192 
00193    // Ticks marks orientation;
00194    Double_t yo=0, zo=0;
00195    switch (fTickMarksOrientation) {
00196       case 0:
00197          yo = 0;
00198          zo = 1;
00199       break;
00200       case 1:
00201          yo = -1;
00202          zo = 0;
00203       break;
00204       case 2:
00205          yo = 0;
00206          zo = -1;
00207       break;
00208       case 3:
00209          yo = 1;
00210          zo = 0;
00211       break;
00212    }
00213 
00214    // Paint level 1 tick marks.
00215    if (fTicks1) {
00216       // Paint the tick marks, if needed.
00217       if (fTickMarksLength) {
00218          Double_t tl = fTickMarksLength*fAxisLength;
00219          glBegin(GL_LINES);
00220          for (i=0; i<fNTicks1 ; i++) {
00221             glVertex3f( fTicks1[i], 0, 0);
00222             glVertex3f( fTicks1[i], yo*tl, zo*tl);
00223          }
00224          glEnd();
00225       }
00226 
00227       // Paint the grid, if needed, on level 1 tick marks.
00228       if (fGridLength) {
00229          const UShort_t stipple = 0x8888;
00230          glLineStipple(1, stipple);
00231          glEnable(GL_LINE_STIPPLE);
00232          glBegin(GL_LINES);
00233          for (i=0; i<fNTicks1; i++) {
00234             glVertex3f( fTicks1[i], 0, 0);
00235             glVertex3f( fTicks1[i], -yo*fGridLength, -zo*fGridLength);
00236          }
00237          glEnd();
00238          glDisable(GL_LINE_STIPPLE);
00239       }
00240    }
00241 
00242    // Paint level 2 tick marks.
00243    if (fTicks2) {
00244       if (fTickMarksLength) {
00245          Double_t tl = 0.5*fTickMarksLength*fAxisLength;
00246          glBegin(GL_LINES);
00247          for (i=0; i<fNTicks2; i++) {
00248             glVertex3f( fTicks2[i], 0, 0);
00249             glVertex3f( fTicks2[i], yo*tl, zo*tl);
00250          }
00251          glEnd();
00252       }
00253    }
00254 }
00255 
00256 
00257 //______________________________________________________________________________
00258 void TGLAxis::PaintGLAxisLabels()
00259 {
00260    // Paint axis labels on the main tick marks.
00261 
00262    if (!fLabelsSize) return;
00263 
00264    Double_t x=0,y=0,z=0;
00265    Int_t i;
00266 
00267    if (!fText) {
00268       fText = new TGLText();
00269       fText->SetTextColor(GetTextColor());
00270       fText->SetGLTextFont(GetTextFont());
00271       fText->SetTextSize(fLabelsSize*fAxisLength);
00272       fText->SetTextAlign(GetTextAlign());
00273    }
00274    fText->SetGLTextAngles(fAngle1, fAngle2, fAngle3);
00275 
00276    switch (fTickMarksOrientation) {
00277       case 0:
00278          y = 0;
00279          z = fLabelsOffset*fAxisLength;
00280       break;
00281       case 1:
00282          y = -fLabelsOffset*fAxisLength;
00283          z = 0;
00284       break;
00285       case 2:
00286          y = 0;
00287          z = -fLabelsOffset*fAxisLength;
00288       break;
00289       case 3:
00290          y = fLabelsOffset*fAxisLength;
00291          z = 0;
00292       break;
00293    }
00294    for (i=0; i<fNDiv1+1 ; i++) {
00295       x = fTicks1[i];
00296       fText->PaintGLText(x,y,z,fLabels[i]);
00297    }
00298 }
00299 
00300 
00301 //______________________________________________________________________________
00302 void TGLAxis::TicksPositions(Option_t *opt)
00303 {
00304    // Compute ticks positions.
00305 
00306    Bool_t optionNoopt, optionLog;
00307 
00308    if (strchr(opt,'N')) optionNoopt = kTRUE;  else optionNoopt = kFALSE;
00309    if (strchr(opt,'G')) optionLog   = kTRUE;  else optionLog   = kFALSE;
00310 
00311    // Determine number of tick marks 1, 2 and 3.
00312    fNDiv3 = fNDiv/10000;
00313    fNDiv2 = (fNDiv-10000*fNDiv3)/100;
00314    fNDiv1 = fNDiv%100;
00315 
00316    // Clean the tick marks arrays if they exist.
00317    if (fTicks1) {
00318       delete [] fTicks1;
00319       fTicks1 = 0;
00320    }
00321    if (fTicks2) {
00322       delete [] fTicks2;
00323       fTicks2 = 0;
00324    }
00325 
00326    // Compute the tick marks positions according to the options.
00327    if (optionNoopt) {
00328       TicksPositionsNoOpt();
00329    } else {
00330       TicksPositionsOpt();
00331    }
00332 }
00333 
00334 
00335 //______________________________________________________________________________
00336 void TGLAxis::TicksPositionsNoOpt()
00337 {
00338    // Compute ticks positions. Linear and not optimized.
00339 
00340    Int_t i, j, k;
00341    Double_t step1 = fAxisLength/(fNDiv1);
00342 
00343    fNTicks1 = fNDiv1+1;
00344    fTicks1  = new Double_t[fNTicks1];
00345 
00346    // Level 1 tick marks.
00347    for (i=0; i<fNTicks1; i++) fTicks1[i] = i*step1;
00348 
00349    // Level 2 tick marks.
00350    if (fNDiv2) {
00351       Double_t t2;
00352       Double_t step2 = step1/fNDiv2;
00353       fNTicks2       = fNDiv1*(fNDiv2-1);
00354       fTicks2        = new Double_t[fNTicks2];
00355       k = 0;
00356       for (i=0; i<fNTicks1-1; i++) {
00357          t2 = fTicks1[i]+step2;
00358          for (j=0; j<fNDiv2-1; j++) {
00359             fTicks2[k] = t2;
00360             k++;
00361             t2 = t2+step2;
00362          }
00363       }
00364    }
00365 }
00366 
00367 
00368 //______________________________________________________________________________
00369 void TGLAxis::TicksPositionsOpt()
00370 {
00371    // Compute ticks positions. Linear and optimized.
00372 
00373    Int_t i, j, k, nDivOpt;
00374    Double_t step1=0, step2=0, wmin2, wmax2;
00375    Double_t wmin = fWmin;
00376    Double_t wmax = fWmax;
00377 
00378    // Level 1 tick marks.
00379    THLimitsFinder::Optimize(wmin,  wmax, fNDiv1,
00380                             fWmin, fWmax, nDivOpt,
00381                             step1, "");
00382    fNDiv1   = nDivOpt;
00383    fNTicks1 = fNDiv1+1;
00384    fTicks1  = new Double_t[fNTicks1];
00385 
00386    Double_t r = fAxisLength/(wmax-wmin);
00387    Double_t w = fWmin;
00388    i = 0;
00389    while (w<=fWmax) {
00390       fTicks1[i] = r*(w-wmin);
00391       i++;
00392       w = w+step1;
00393    }
00394 
00395    // Level 2 tick marks.
00396    if (fNDiv2) {
00397       Double_t t2;
00398       THLimitsFinder::Optimize(fWmin, fWmin+step1, fNDiv2,
00399                                wmin2, wmax2, nDivOpt,
00400                                step2, "");
00401       fNDiv2       = nDivOpt;
00402       step2        = TMath::Abs((fTicks1[1]-fTicks1[0])/fNDiv2);
00403       Int_t nTickl = (Int_t)(fTicks1[0]/step2);
00404       Int_t nTickr = (Int_t)((fAxisLength-fTicks1[fNTicks1-1])/step2);
00405       fNTicks2     = fNDiv1*(fNDiv2-1)+nTickl+nTickr;
00406       fTicks2      = new Double_t[fNTicks2];
00407       k = 0;
00408       for (i=0; i<fNTicks1-1; i++) {
00409          t2 = fTicks1[i]+step2;
00410          for (j=0; j<fNDiv2-1; j++) {
00411             fTicks2[k] = t2;
00412             k++;
00413             t2 = t2+step2;
00414          }
00415       }
00416       if (nTickl) {
00417          t2 = fTicks1[0]-step2;
00418          for (i=0; i<nTickl; i++) {
00419             fTicks2[k] = t2;
00420             k++;
00421             t2 = t2-step2;
00422          }
00423       }
00424       if (nTickr) {
00425          t2 = fTicks1[fNTicks1-1]+step2;
00426          for (i=0; i<nTickr; i++) {
00427             fTicks2[k] = t2;
00428             k++;
00429             t2 = t2+step2;
00430          }
00431       }
00432    }
00433 }
00434 
00435 
00436 //______________________________________________________________________________
00437 void TGLAxis::DoLabels()
00438 {
00439    // Do labels.
00440 
00441    if (fLabels) delete [] fLabels;
00442    fLabels = new TString[fNTicks1];
00443    Int_t i;
00444 
00445    // Make regular labels between fWmin and fWmax.
00446    Double_t dw = (fWmax-fWmin)/(fNDiv1);
00447    for (i=0; i<fNTicks1; i++) {
00448       fLabels[i] = Form("%g",fWmin+i*dw);
00449    }
00450 }
00451 
00452 
00453 //______________________________________________________________________________
00454 void TGLAxis::SetLabelsAngles(Double_t a1, Double_t a2, Double_t a3)
00455 {
00456    // Set labels' angles.
00457 
00458    fAngle1 = a1;
00459    fAngle2 = a2;
00460    fAngle3 = a3;
00461 }

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