// @(#)root/gl:$Id$
// Author:  Olivier Couet  17/04/2007

/*************************************************************************
 * Copyright (C) 1995-2006, Rene Brun and Fons Rademakers.               *
 * All rights reserved.                                                  *
 *                                                                       *
 * For the licensing terms see $ROOTSYS/LICENSE.                         *
 * For the list of contributors see $ROOTSYS/README/CREDITS.             *
 *************************************************************************/
#include "Riostream.h"
#include "TROOT.h"

#include "TGLIncludes.h"
#include "TGLUtil.h"
#include "TGLAxis.h"
#include "TGLText.h"
#include "TColor.h"
#include "TString.h"
#include "TMath.h"
#include "THLimitsFinder.h"

//______________________________________________________________________________
/* Begin_Html
<center><h2>GL Axis</h2></center>
To draw a 3D axis in a GL window. The labels are drawn using FTGL.
End_Html */

ClassImp(TGLAxis)

//______________________________________________________________________________
TGLAxis::TGLAxis(): TAttLine(1,1,1), TAttText(20,0.,1,42,0.04)
{
   // Constructor.

   Init();
}


//______________________________________________________________________________
void TGLAxis::Init()
{
   // Default initialization.

   fNDiv = fNDiv1 = fNDiv2 = fNDiv3 = 0;
   fNTicks1 = fNTicks2 = 0;
   fTicks1          = 0;
   fTicks2          = 0;
   fLabels          = 0;
   fText            = 0;
   fAngle1          = 90.;
   fAngle2          = 0.;
   fAngle3          = 0.;
   fAxisLength      = 0.;
   fWmin = fWmax    = 0.;
   fTickMarksLength = 0.04; // % of fAxisLength
   fTickMarksOrientation = 2; // can be 0, 1, 2, or 3
   fLabelsOffset    = 0.09; // % of fAxisLength
   fLabelsSize      = 0.06; // % of fAxisLength
   fGridLength      = 0.;   // 0 means no grid
}


//______________________________________________________________________________
TGLAxis::~TGLAxis()
{
   // Destructor.

   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)
{
   // Paint GL Axis.
   //
   // p1, p2     : Axis position in the 3D space.
   // wmin, wmax : Minimum and maximum values along the axis. wmin < wmax.
   // ndiv       : Number of axis divisions. It is an integer in the form
   //              "ttsspp" where "tt" is the number of tertiary divisions,
   //              "ss" is the number of secondary divisions and "pp" the
   //              number of primary divisions.
   // opt        : Options.
   //              "N" - By default the number of divisions is optimized to
   //                    get a nice labeling. When option "N" is given, the
   //                    number of divisions is not optimized.

   fNDiv = ndiv;
   if (wmax<=wmin) {
      fWmax = wmin;
      fWmin = wmax;
   } else {
      fWmax = wmax;
      fWmin = wmin;
   }

   // Compute the axis length.
   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();

   // Paint using GL
   glPushMatrix();

   // Translation
   glTranslatef(x1, y1, z1);

   // Rotation in the Z plane
   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]);
   }

   // Rotation in the XY plane
   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()
{
   // Paint horizontal axis body at position (0,0,0)

   TColor *col;
   Float_t red = 1.f, green = 1.f, blue = 1.f;
   if ((col = gROOT->GetColor(GetLineColor())))
      col->GetRGB(red, green, blue);
   glColor3d(red, green, blue);
   TGLUtil::LineWidth(GetLineWidth());
   glBegin(GL_LINES);
   glVertex3d(0., 0., 0.);
   glVertex3d(fAxisLength, 0., 0.);
   glEnd();
}


//______________________________________________________________________________
void TGLAxis::PaintGLAxisTickMarks()
{
   // Paint axis tick marks.

   Int_t i;

   // Ticks marks orientation;
   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;
   }

   // Paint level 1 tick marks.
   if (fTicks1) {
      // Paint the tick marks, if needed.
      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();
      }

      // Paint the grid, if needed, on level 1 tick marks.
      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);
      }
   }

   // Paint level 2 tick marks.
   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()
{
   // Paint axis labels on the main tick marks.

   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)
{
   // Compute ticks positions.

   Bool_t optionNoopt /* , optionLog */;

   if (strchr(opt,'N')) optionNoopt = kTRUE;  else optionNoopt = kFALSE;
   // if (strchr(opt,'G')) optionLog   = kTRUE;  else optionLog   = kFALSE;

   // Determine number of tick marks 1, 2 and 3.
   fNDiv3 = fNDiv/10000;
   fNDiv2 = (fNDiv-10000*fNDiv3)/100;
   fNDiv1 = fNDiv%100;

   // Clean the tick marks arrays if they exist.
   if (fTicks1) {
      delete [] fTicks1;
      fTicks1 = 0;
   }
   if (fTicks2) {
      delete [] fTicks2;
      fTicks2 = 0;
   }

   // Compute the tick marks positions according to the options.
   if (optionNoopt) {
      TicksPositionsNoOpt();
   } else {
      TicksPositionsOpt();
   }
}


//______________________________________________________________________________
void TGLAxis::TicksPositionsNoOpt()
{
   // Compute ticks positions. Linear and not optimized.

   Int_t i, j, k;
   Double_t step1 = fAxisLength/(fNDiv1);

   fNTicks1 = fNDiv1+1;
   fTicks1  = new Double_t[fNTicks1];

   // Level 1 tick marks.
   for (i=0; i<fNTicks1; i++) fTicks1[i] = i*step1;

   // Level 2 tick marks.
   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()
{
   // Compute ticks positions. Linear and optimized.

   Int_t i, j, k, nDivOpt;
   Double_t step1=0, step2=0, wmin2=0, wmax2=0;
   Double_t wmin = fWmin;
   Double_t wmax = fWmax;

   // Level 1 tick marks.
   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;
   }

   // Level 2 tick marks.
   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()
{
   // Do labels.

   if (fLabels) delete [] fLabels;
   fLabels = new TString[fNTicks1];
   Int_t i;

   // Make regular labels between fWmin and fWmax.
   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)
{
   // Set labels' angles.

   fAngle1 = a1;
   fAngle2 = a2;
   fAngle3 = a3;
}
 TGLAxis.cxx:1
 TGLAxis.cxx:2
 TGLAxis.cxx:3
 TGLAxis.cxx:4
 TGLAxis.cxx:5
 TGLAxis.cxx:6
 TGLAxis.cxx:7
 TGLAxis.cxx:8
 TGLAxis.cxx:9
 TGLAxis.cxx:10
 TGLAxis.cxx:11
 TGLAxis.cxx:12
 TGLAxis.cxx:13
 TGLAxis.cxx:14
 TGLAxis.cxx:15
 TGLAxis.cxx:16
 TGLAxis.cxx:17
 TGLAxis.cxx:18
 TGLAxis.cxx:19
 TGLAxis.cxx:20
 TGLAxis.cxx:21
 TGLAxis.cxx:22
 TGLAxis.cxx:23
 TGLAxis.cxx:24
 TGLAxis.cxx:25
 TGLAxis.cxx:26
 TGLAxis.cxx:27
 TGLAxis.cxx:28
 TGLAxis.cxx:29
 TGLAxis.cxx:30
 TGLAxis.cxx:31
 TGLAxis.cxx:32
 TGLAxis.cxx:33
 TGLAxis.cxx:34
 TGLAxis.cxx:35
 TGLAxis.cxx:36
 TGLAxis.cxx:37
 TGLAxis.cxx:38
 TGLAxis.cxx:39
 TGLAxis.cxx:40
 TGLAxis.cxx:41
 TGLAxis.cxx:42
 TGLAxis.cxx:43
 TGLAxis.cxx:44
 TGLAxis.cxx:45
 TGLAxis.cxx:46
 TGLAxis.cxx:47
 TGLAxis.cxx:48
 TGLAxis.cxx:49
 TGLAxis.cxx:50
 TGLAxis.cxx:51
 TGLAxis.cxx:52
 TGLAxis.cxx:53
 TGLAxis.cxx:54
 TGLAxis.cxx:55
 TGLAxis.cxx:56
 TGLAxis.cxx:57
 TGLAxis.cxx:58
 TGLAxis.cxx:59
 TGLAxis.cxx:60
 TGLAxis.cxx:61
 TGLAxis.cxx:62
 TGLAxis.cxx:63
 TGLAxis.cxx:64
 TGLAxis.cxx:65
 TGLAxis.cxx:66
 TGLAxis.cxx:67
 TGLAxis.cxx:68
 TGLAxis.cxx:69
 TGLAxis.cxx:70
 TGLAxis.cxx:71
 TGLAxis.cxx:72
 TGLAxis.cxx:73
 TGLAxis.cxx:74
 TGLAxis.cxx:75
 TGLAxis.cxx:76
 TGLAxis.cxx:77
 TGLAxis.cxx:78
 TGLAxis.cxx:79
 TGLAxis.cxx:80
 TGLAxis.cxx:81
 TGLAxis.cxx:82
 TGLAxis.cxx:83
 TGLAxis.cxx:84
 TGLAxis.cxx:85
 TGLAxis.cxx:86
 TGLAxis.cxx:87
 TGLAxis.cxx:88
 TGLAxis.cxx:89
 TGLAxis.cxx:90
 TGLAxis.cxx:91
 TGLAxis.cxx:92
 TGLAxis.cxx:93
 TGLAxis.cxx:94
 TGLAxis.cxx:95
 TGLAxis.cxx:96
 TGLAxis.cxx:97
 TGLAxis.cxx:98
 TGLAxis.cxx:99
 TGLAxis.cxx:100
 TGLAxis.cxx:101
 TGLAxis.cxx:102
 TGLAxis.cxx:103
 TGLAxis.cxx:104
 TGLAxis.cxx:105
 TGLAxis.cxx:106
 TGLAxis.cxx:107
 TGLAxis.cxx:108
 TGLAxis.cxx:109
 TGLAxis.cxx:110
 TGLAxis.cxx:111
 TGLAxis.cxx:112
 TGLAxis.cxx:113
 TGLAxis.cxx:114
 TGLAxis.cxx:115
 TGLAxis.cxx:116
 TGLAxis.cxx:117
 TGLAxis.cxx:118
 TGLAxis.cxx:119
 TGLAxis.cxx:120
 TGLAxis.cxx:121
 TGLAxis.cxx:122
 TGLAxis.cxx:123
 TGLAxis.cxx:124
 TGLAxis.cxx:125
 TGLAxis.cxx:126
 TGLAxis.cxx:127
 TGLAxis.cxx:128
 TGLAxis.cxx:129
 TGLAxis.cxx:130
 TGLAxis.cxx:131
 TGLAxis.cxx:132
 TGLAxis.cxx:133
 TGLAxis.cxx:134
 TGLAxis.cxx:135
 TGLAxis.cxx:136
 TGLAxis.cxx:137
 TGLAxis.cxx:138
 TGLAxis.cxx:139
 TGLAxis.cxx:140
 TGLAxis.cxx:141
 TGLAxis.cxx:142
 TGLAxis.cxx:143
 TGLAxis.cxx:144
 TGLAxis.cxx:145
 TGLAxis.cxx:146
 TGLAxis.cxx:147
 TGLAxis.cxx:148
 TGLAxis.cxx:149
 TGLAxis.cxx:150
 TGLAxis.cxx:151
 TGLAxis.cxx:152
 TGLAxis.cxx:153
 TGLAxis.cxx:154
 TGLAxis.cxx:155
 TGLAxis.cxx:156
 TGLAxis.cxx:157
 TGLAxis.cxx:158
 TGLAxis.cxx:159
 TGLAxis.cxx:160
 TGLAxis.cxx:161
 TGLAxis.cxx:162
 TGLAxis.cxx:163
 TGLAxis.cxx:164
 TGLAxis.cxx:165
 TGLAxis.cxx:166
 TGLAxis.cxx:167
 TGLAxis.cxx:168
 TGLAxis.cxx:169
 TGLAxis.cxx:170
 TGLAxis.cxx:171
 TGLAxis.cxx:172
 TGLAxis.cxx:173
 TGLAxis.cxx:174
 TGLAxis.cxx:175
 TGLAxis.cxx:176
 TGLAxis.cxx:177
 TGLAxis.cxx:178
 TGLAxis.cxx:179
 TGLAxis.cxx:180
 TGLAxis.cxx:181
 TGLAxis.cxx:182
 TGLAxis.cxx:183
 TGLAxis.cxx:184
 TGLAxis.cxx:185
 TGLAxis.cxx:186
 TGLAxis.cxx:187
 TGLAxis.cxx:188
 TGLAxis.cxx:189
 TGLAxis.cxx:190
 TGLAxis.cxx:191
 TGLAxis.cxx:192
 TGLAxis.cxx:193
 TGLAxis.cxx:194
 TGLAxis.cxx:195
 TGLAxis.cxx:196
 TGLAxis.cxx:197
 TGLAxis.cxx:198
 TGLAxis.cxx:199
 TGLAxis.cxx:200
 TGLAxis.cxx:201
 TGLAxis.cxx:202
 TGLAxis.cxx:203
 TGLAxis.cxx:204
 TGLAxis.cxx:205
 TGLAxis.cxx:206
 TGLAxis.cxx:207
 TGLAxis.cxx:208
 TGLAxis.cxx:209
 TGLAxis.cxx:210
 TGLAxis.cxx:211
 TGLAxis.cxx:212
 TGLAxis.cxx:213
 TGLAxis.cxx:214
 TGLAxis.cxx:215
 TGLAxis.cxx:216
 TGLAxis.cxx:217
 TGLAxis.cxx:218
 TGLAxis.cxx:219
 TGLAxis.cxx:220
 TGLAxis.cxx:221
 TGLAxis.cxx:222
 TGLAxis.cxx:223
 TGLAxis.cxx:224
 TGLAxis.cxx:225
 TGLAxis.cxx:226
 TGLAxis.cxx:227
 TGLAxis.cxx:228
 TGLAxis.cxx:229
 TGLAxis.cxx:230
 TGLAxis.cxx:231
 TGLAxis.cxx:232
 TGLAxis.cxx:233
 TGLAxis.cxx:234
 TGLAxis.cxx:235
 TGLAxis.cxx:236
 TGLAxis.cxx:237
 TGLAxis.cxx:238
 TGLAxis.cxx:239
 TGLAxis.cxx:240
 TGLAxis.cxx:241
 TGLAxis.cxx:242
 TGLAxis.cxx:243
 TGLAxis.cxx:244
 TGLAxis.cxx:245
 TGLAxis.cxx:246
 TGLAxis.cxx:247
 TGLAxis.cxx:248
 TGLAxis.cxx:249
 TGLAxis.cxx:250
 TGLAxis.cxx:251
 TGLAxis.cxx:252
 TGLAxis.cxx:253
 TGLAxis.cxx:254
 TGLAxis.cxx:255
 TGLAxis.cxx:256
 TGLAxis.cxx:257
 TGLAxis.cxx:258
 TGLAxis.cxx:259
 TGLAxis.cxx:260
 TGLAxis.cxx:261
 TGLAxis.cxx:262
 TGLAxis.cxx:263
 TGLAxis.cxx:264
 TGLAxis.cxx:265
 TGLAxis.cxx:266
 TGLAxis.cxx:267
 TGLAxis.cxx:268
 TGLAxis.cxx:269
 TGLAxis.cxx:270
 TGLAxis.cxx:271
 TGLAxis.cxx:272
 TGLAxis.cxx:273
 TGLAxis.cxx:274
 TGLAxis.cxx:275
 TGLAxis.cxx:276
 TGLAxis.cxx:277
 TGLAxis.cxx:278
 TGLAxis.cxx:279
 TGLAxis.cxx:280
 TGLAxis.cxx:281
 TGLAxis.cxx:282
 TGLAxis.cxx:283
 TGLAxis.cxx:284
 TGLAxis.cxx:285
 TGLAxis.cxx:286
 TGLAxis.cxx:287
 TGLAxis.cxx:288
 TGLAxis.cxx:289
 TGLAxis.cxx:290
 TGLAxis.cxx:291
 TGLAxis.cxx:292
 TGLAxis.cxx:293
 TGLAxis.cxx:294
 TGLAxis.cxx:295
 TGLAxis.cxx:296
 TGLAxis.cxx:297
 TGLAxis.cxx:298
 TGLAxis.cxx:299
 TGLAxis.cxx:300
 TGLAxis.cxx:301
 TGLAxis.cxx:302
 TGLAxis.cxx:303
 TGLAxis.cxx:304
 TGLAxis.cxx:305
 TGLAxis.cxx:306
 TGLAxis.cxx:307
 TGLAxis.cxx:308
 TGLAxis.cxx:309
 TGLAxis.cxx:310
 TGLAxis.cxx:311
 TGLAxis.cxx:312
 TGLAxis.cxx:313
 TGLAxis.cxx:314
 TGLAxis.cxx:315
 TGLAxis.cxx:316
 TGLAxis.cxx:317
 TGLAxis.cxx:318
 TGLAxis.cxx:319
 TGLAxis.cxx:320
 TGLAxis.cxx:321
 TGLAxis.cxx:322
 TGLAxis.cxx:323
 TGLAxis.cxx:324
 TGLAxis.cxx:325
 TGLAxis.cxx:326
 TGLAxis.cxx:327
 TGLAxis.cxx:328
 TGLAxis.cxx:329
 TGLAxis.cxx:330
 TGLAxis.cxx:331
 TGLAxis.cxx:332
 TGLAxis.cxx:333
 TGLAxis.cxx:334
 TGLAxis.cxx:335
 TGLAxis.cxx:336
 TGLAxis.cxx:337
 TGLAxis.cxx:338
 TGLAxis.cxx:339
 TGLAxis.cxx:340
 TGLAxis.cxx:341
 TGLAxis.cxx:342
 TGLAxis.cxx:343
 TGLAxis.cxx:344
 TGLAxis.cxx:345
 TGLAxis.cxx:346
 TGLAxis.cxx:347
 TGLAxis.cxx:348
 TGLAxis.cxx:349
 TGLAxis.cxx:350
 TGLAxis.cxx:351
 TGLAxis.cxx:352
 TGLAxis.cxx:353
 TGLAxis.cxx:354
 TGLAxis.cxx:355
 TGLAxis.cxx:356
 TGLAxis.cxx:357
 TGLAxis.cxx:358
 TGLAxis.cxx:359
 TGLAxis.cxx:360
 TGLAxis.cxx:361
 TGLAxis.cxx:362
 TGLAxis.cxx:363
 TGLAxis.cxx:364
 TGLAxis.cxx:365
 TGLAxis.cxx:366
 TGLAxis.cxx:367
 TGLAxis.cxx:368
 TGLAxis.cxx:369
 TGLAxis.cxx:370
 TGLAxis.cxx:371
 TGLAxis.cxx:372
 TGLAxis.cxx:373
 TGLAxis.cxx:374
 TGLAxis.cxx:375
 TGLAxis.cxx:376
 TGLAxis.cxx:377
 TGLAxis.cxx:378
 TGLAxis.cxx:379
 TGLAxis.cxx:380
 TGLAxis.cxx:381
 TGLAxis.cxx:382
 TGLAxis.cxx:383
 TGLAxis.cxx:384
 TGLAxis.cxx:385
 TGLAxis.cxx:386
 TGLAxis.cxx:387
 TGLAxis.cxx:388
 TGLAxis.cxx:389
 TGLAxis.cxx:390
 TGLAxis.cxx:391
 TGLAxis.cxx:392
 TGLAxis.cxx:393
 TGLAxis.cxx:394
 TGLAxis.cxx:395
 TGLAxis.cxx:396
 TGLAxis.cxx:397
 TGLAxis.cxx:398
 TGLAxis.cxx:399
 TGLAxis.cxx:400
 TGLAxis.cxx:401
 TGLAxis.cxx:402
 TGLAxis.cxx:403
 TGLAxis.cxx:404
 TGLAxis.cxx:405
 TGLAxis.cxx:406
 TGLAxis.cxx:407
 TGLAxis.cxx:408
 TGLAxis.cxx:409
 TGLAxis.cxx:410
 TGLAxis.cxx:411
 TGLAxis.cxx:412
 TGLAxis.cxx:413
 TGLAxis.cxx:414
 TGLAxis.cxx:415
 TGLAxis.cxx:416
 TGLAxis.cxx:417
 TGLAxis.cxx:418
 TGLAxis.cxx:419
 TGLAxis.cxx:420
 TGLAxis.cxx:421
 TGLAxis.cxx:422
 TGLAxis.cxx:423
 TGLAxis.cxx:424
 TGLAxis.cxx:425
 TGLAxis.cxx:426
 TGLAxis.cxx:427
 TGLAxis.cxx:428
 TGLAxis.cxx:429
 TGLAxis.cxx:430
 TGLAxis.cxx:431
 TGLAxis.cxx:432
 TGLAxis.cxx:433
 TGLAxis.cxx:434
 TGLAxis.cxx:435
 TGLAxis.cxx:436
 TGLAxis.cxx:437
 TGLAxis.cxx:438
 TGLAxis.cxx:439
 TGLAxis.cxx:440
 TGLAxis.cxx:441
 TGLAxis.cxx:442
 TGLAxis.cxx:443
 TGLAxis.cxx:444
 TGLAxis.cxx:445
 TGLAxis.cxx:446
 TGLAxis.cxx:447
 TGLAxis.cxx:448
 TGLAxis.cxx:449
 TGLAxis.cxx:450
 TGLAxis.cxx:451
 TGLAxis.cxx:452
 TGLAxis.cxx:453
 TGLAxis.cxx:454
 TGLAxis.cxx:455
 TGLAxis.cxx:456
 TGLAxis.cxx:457
 TGLAxis.cxx:458
 TGLAxis.cxx:459
 TGLAxis.cxx:460
 TGLAxis.cxx:461