Logo ROOT  
Reference Guide
TGLAxis.cxx
Go to the documentation of this file.
1// @(#)root/gl:$Id$
2// Author: Olivier Couet 17/04/2007
3
4/*************************************************************************
5 * Copyright (C) 1995-2006, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11#include "Riostream.h"
12#include "TROOT.h"
13
14#include "TGLIncludes.h"
15#include "TGLUtil.h"
16#include "TGLAxis.h"
17#include "TGLText.h"
18#include "TColor.h"
19#include "TString.h"
20#include "TMath.h"
21#include "THLimitsFinder.h"
22
23/** \class TGLAxis
24\ingroup opengl
25GL Axis.
26
27To draw a 3D axis in a GL window. The labels are drawn using FTGL.
28*/
29
31
32////////////////////////////////////////////////////////////////////////////////
33/// Constructor.
34
35TGLAxis::TGLAxis(): TAttLine(1,1,1), TAttText(20,0.,1,42,0.04)
36{
37 Init();
38}
39
40////////////////////////////////////////////////////////////////////////////////
41/// Default initialization.
42
44{
45 fNDiv = fNDiv1 = fNDiv2 = fNDiv3 = 0;
46 fNTicks1 = fNTicks2 = 0;
47 fTicks1 = 0;
48 fTicks2 = 0;
49 fLabels = 0;
50 fText = 0;
51 fAngle1 = 90.;
52 fAngle2 = 0.;
53 fAngle3 = 0.;
54 fAxisLength = 0.;
55 fWmin = fWmax = 0.;
56 fTickMarksLength = 0.04; // % of fAxisLength
57 fTickMarksOrientation = 2; // can be 0, 1, 2, or 3
58 fLabelsOffset = 0.09; // % of fAxisLength
59 fLabelsSize = 0.06; // % of fAxisLength
60 fGridLength = 0.; // 0 means no grid
61}
62
63////////////////////////////////////////////////////////////////////////////////
64/// Destructor.
65
67{
68 if (fTicks1) delete [] fTicks1;
69 if (fTicks2) delete [] fTicks2;
70 if (fLabels) delete [] fLabels;
71 if (fText) delete fText;
72}
73
74////////////////////////////////////////////////////////////////////////////////
75/// Paint GL Axis.
76///
77/// - p1, p2 : Axis position in the 3D space.
78/// - wmin, wmax : Minimum and maximum values along the axis. wmin < wmax.
79/// - ndiv : Number of axis divisions. It is an integer in the form
80/// "ttsspp" where "tt" is the number of tertiary divisions,
81/// "ss" is the number of secondary divisions and "pp" the
82/// number of primary divisions.
83/// - opt : Options.
84/// "N" - By default the number of divisions is optimized to
85/// get a nice labeling. When option "N" is given, the
86/// number of divisions is not optimized.
87
88void TGLAxis::PaintGLAxis(const Double_t p1[3], const Double_t p2[3],
89 Double_t wmin, Double_t wmax, Int_t ndiv,
90 Option_t *opt)
91{
92 fNDiv = ndiv;
93 if (wmax<=wmin) {
94 fWmax = wmin;
95 fWmin = wmax;
96 } else {
97 fWmax = wmax;
98 fWmin = wmin;
99 }
100
101 // Compute the axis length.
102 Double_t x1 = p1[0], y1 = p1[1], z1 = p1[2];
103 Double_t x2 = p2[0], y2 = p2[1], z2 = p2[2];
105 (y2-y1)*(y2-y1)+
106 (z2-z1)*(z2-z1));
107
108 TicksPositions(opt);
109
110 DoLabels();
111
112 // Paint using GL
113 glPushMatrix();
114
115 // Translation
116 glTranslatef(x1, y1, z1);
117
118 // Rotation in the Z plane
119 Double_t phi=0;
120 Double_t normal[3];
121 normal[0] = 0;
122 normal[1] = 1;
123 normal[2] = 0;
124 if (z1!=z2) {
125 if (y2==y1 && x2==x1) {
126 if (z2<z1) phi = 90;
127 else phi = 270;
128 } else {
129 Double_t p3[3];
130 p3[0] = p2[0]; p3[1] = p2[1]; p3[2] = 0;
131 TMath::Normal2Plane(p1, p2, p3, normal);
132 phi = TMath::ACos(TMath::Abs(z2-z1)/fAxisLength);
133 phi = -(90-(180/TMath::Pi())*phi);
134 }
135 glRotatef(phi, normal[0], normal[1], normal[2]);
136 }
137
138 // Rotation in the XY plane
139 Double_t theta = 0;
140 if (y2!=y1) {
141 if ((x2-x1)>0) {
142 theta = TMath::ATan((y2-y1)/(x2-x1));
143 theta = (180/TMath::Pi())*theta;
144 } else if ((x2-x1)<0) {
145 theta = TMath::ATan((y2-y1)/(x2-x1));
146 theta = 180+(180/TMath::Pi())*theta;
147 } else {
148 if (y2>y1) theta = 90;
149 else theta = 270;
150 }
151 } else {
152 if (x2<x1) theta = 180;
153 }
154 glRotatef(theta, 0., 0., 1.);
155
157
159
161
162 glPopMatrix();
163}
164
165////////////////////////////////////////////////////////////////////////////////
166/// Paint horizontal axis body at position (0,0,0)
167
169{
170 TColor *col;
171 Float_t red = 1.f, green = 1.f, blue = 1.f;
172 if ((col = gROOT->GetColor(GetLineColor())))
173 col->GetRGB(red, green, blue);
174 glColor3d(red, green, blue);
176 glBegin(GL_LINES);
177 glVertex3d(0., 0., 0.);
178 glVertex3d(fAxisLength, 0., 0.);
179 glEnd();
180}
181
182////////////////////////////////////////////////////////////////////////////////
183/// Paint axis tick marks.
184
186{
187 Int_t i;
188
189 // Ticks marks orientation;
190 Double_t yo=0, zo=0;
191 switch (fTickMarksOrientation) {
192 case 0:
193 yo = 0;
194 zo = 1;
195 break;
196 case 1:
197 yo = -1;
198 zo = 0;
199 break;
200 case 2:
201 yo = 0;
202 zo = -1;
203 break;
204 case 3:
205 yo = 1;
206 zo = 0;
207 break;
208 }
209
210 // Paint level 1 tick marks.
211 if (fTicks1) {
212 // Paint the tick marks, if needed.
213 if (fTickMarksLength) {
215 glBegin(GL_LINES);
216 for (i=0; i<fNTicks1 ; i++) {
217 glVertex3f( fTicks1[i], 0, 0);
218 glVertex3f( fTicks1[i], yo*tl, zo*tl);
219 }
220 glEnd();
221 }
222
223 // Paint the grid, if needed, on level 1 tick marks.
224 if (fGridLength) {
225 const UShort_t stipple = 0x8888;
226 glLineStipple(1, stipple);
227 glEnable(GL_LINE_STIPPLE);
228 glBegin(GL_LINES);
229 for (i=0; i<fNTicks1; i++) {
230 glVertex3f( fTicks1[i], 0, 0);
231 glVertex3f( fTicks1[i], -yo*fGridLength, -zo*fGridLength);
232 }
233 glEnd();
234 glDisable(GL_LINE_STIPPLE);
235 }
236 }
237
238 // Paint level 2 tick marks.
239 if (fTicks2) {
240 if (fTickMarksLength) {
242 glBegin(GL_LINES);
243 for (i=0; i<fNTicks2; i++) {
244 glVertex3f( fTicks2[i], 0, 0);
245 glVertex3f( fTicks2[i], yo*tl, zo*tl);
246 }
247 glEnd();
248 }
249 }
250}
251
252////////////////////////////////////////////////////////////////////////////////
253/// Paint axis labels on the main tick marks.
254
256{
257 if (!fLabelsSize) return;
258
259 Double_t x=0,y=0,z=0;
260 Int_t i;
261
262 if (!fText) {
263 fText = new TGLText();
268 }
270
271 switch (fTickMarksOrientation) {
272 case 0:
273 y = 0;
275 break;
276 case 1:
278 z = 0;
279 break;
280 case 2:
281 y = 0;
283 break;
284 case 3:
286 z = 0;
287 break;
288 }
289 for (i=0; i<fNDiv1+1 ; i++) {
290 x = fTicks1[i];
291 fText->PaintGLText(x,y,z,fLabels[i]);
292 }
293}
294
295////////////////////////////////////////////////////////////////////////////////
296/// Compute ticks positions.
297
299{
300 Bool_t optionNoopt /* , optionLog */;
301
302 if (strchr(opt,'N')) optionNoopt = kTRUE; else optionNoopt = kFALSE;
303 // if (strchr(opt,'G')) optionLog = kTRUE; else optionLog = kFALSE;
304
305 // Determine number of tick marks 1, 2 and 3.
306 fNDiv3 = fNDiv/10000;
307 fNDiv2 = (fNDiv-10000*fNDiv3)/100;
308 fNDiv1 = fNDiv%100;
309
310 // Clean the tick marks arrays if they exist.
311 if (fTicks1) {
312 delete [] fTicks1;
313 fTicks1 = 0;
314 }
315 if (fTicks2) {
316 delete [] fTicks2;
317 fTicks2 = 0;
318 }
319
320 // Compute the tick marks positions according to the options.
321 if (optionNoopt) {
323 } else {
325 }
326}
327
328////////////////////////////////////////////////////////////////////////////////
329/// Compute ticks positions. Linear and not optimized.
330
332{
333 Int_t i, j, k;
334 Double_t step1 = fAxisLength/(fNDiv1);
335
336 fNTicks1 = fNDiv1+1;
338
339 // Level 1 tick marks.
340 for (i=0; i<fNTicks1; i++) fTicks1[i] = i*step1;
341
342 // Level 2 tick marks.
343 if (fNDiv2) {
344 Double_t t2;
345 Double_t step2 = step1/fNDiv2;
346 fNTicks2 = fNDiv1*(fNDiv2-1);
348 k = 0;
349 for (i=0; i<fNTicks1-1; i++) {
350 t2 = fTicks1[i]+step2;
351 for (j=0; j<fNDiv2-1; j++) {
352 fTicks2[k] = t2;
353 k++;
354 t2 = t2+step2;
355 }
356 }
357 }
358}
359
360////////////////////////////////////////////////////////////////////////////////
361/// Compute ticks positions. Linear and optimized.
362
364{
365 Int_t i, j, k, nDivOpt;
366 Double_t step1=0, step2=0, wmin2=0, wmax2=0;
367 Double_t wmin = fWmin;
368 Double_t wmax = fWmax;
369
370 // Level 1 tick marks.
372 fWmin, fWmax, nDivOpt,
373 step1, "");
374 fNDiv1 = nDivOpt;
375 fNTicks1 = fNDiv1+1;
377
378 Double_t r = fAxisLength/(wmax-wmin);
379 Double_t w = fWmin;
380 i = 0;
381 while (w<=fWmax) {
382 fTicks1[i] = r*(w-wmin);
383 i++;
384 w = w+step1;
385 }
386
387 // Level 2 tick marks.
388 if (fNDiv2) {
389 Double_t t2;
391 wmin2, wmax2, nDivOpt,
392 step2, "");
393 fNDiv2 = nDivOpt;
394 step2 = TMath::Abs((fTicks1[1]-fTicks1[0])/fNDiv2);
395 Int_t nTickl = (Int_t)(fTicks1[0]/step2);
396 Int_t nTickr = (Int_t)((fAxisLength-fTicks1[fNTicks1-1])/step2);
397 fNTicks2 = fNDiv1*(fNDiv2-1)+nTickl+nTickr;
399 k = 0;
400 for (i=0; i<fNTicks1-1; i++) {
401 t2 = fTicks1[i]+step2;
402 for (j=0; j<fNDiv2-1; j++) {
403 fTicks2[k] = t2;
404 k++;
405 t2 = t2+step2;
406 }
407 }
408 if (nTickl) {
409 t2 = fTicks1[0]-step2;
410 for (i=0; i<nTickl; i++) {
411 fTicks2[k] = t2;
412 k++;
413 t2 = t2-step2;
414 }
415 }
416 if (nTickr) {
417 t2 = fTicks1[fNTicks1-1]+step2;
418 for (i=0; i<nTickr; i++) {
419 fTicks2[k] = t2;
420 k++;
421 t2 = t2+step2;
422 }
423 }
424 }
425}
426
427////////////////////////////////////////////////////////////////////////////////
428/// Do labels.
429
431{
432 if (fLabels) delete [] fLabels;
433 fLabels = new TString[fNTicks1];
434 Int_t i;
435
436 // Make regular labels between fWmin and fWmax.
437 Double_t dw = (fWmax-fWmin)/(fNDiv1);
438 for (i=0; i<fNTicks1; i++) {
439 fLabels[i] = Form("%g",fWmin+i*dw);
440 }
441}
442
443////////////////////////////////////////////////////////////////////////////////
444/// Set labels' angles.
445
447{
448 fAngle1 = a1;
449 fAngle2 = a2;
450 fAngle3 = a3;
451}
#define GL_LINES
Definition: GL_glu.h:284
ROOT::R::TRInterface & r
Definition: Object.C:4
static const double x2[5]
static const double x1[5]
unsigned short UShort_t
Definition: RtypesCore.h:38
int Int_t
Definition: RtypesCore.h:43
const Bool_t kFALSE
Definition: RtypesCore.h:90
double Double_t
Definition: RtypesCore.h:57
float Float_t
Definition: RtypesCore.h:55
const Bool_t kTRUE
Definition: RtypesCore.h:89
const char Option_t
Definition: RtypesCore.h:64
#define ClassImp(name)
Definition: Rtypes.h:361
#define gROOT
Definition: TROOT.h:406
char * Form(const char *fmt,...)
Line Attributes class.
Definition: TAttLine.h:18
virtual Color_t GetLineColor() const
Return the line color.
Definition: TAttLine.h:33
virtual Width_t GetLineWidth() const
Return the line width.
Definition: TAttLine.h:35
Text Attributes class.
Definition: TAttText.h:18
virtual void SetTextAlign(Short_t align=11)
Set the text alignment.
Definition: TAttText.h:41
virtual Short_t GetTextAlign() const
Return the text alignment.
Definition: TAttText.h:32
virtual Font_t GetTextFont() const
Return the text font.
Definition: TAttText.h:35
virtual Color_t GetTextColor() const
Return the text color.
Definition: TAttText.h:34
virtual void SetTextColor(Color_t tcolor=1)
Set the text color.
Definition: TAttText.h:43
virtual void SetTextSize(Float_t tsize=1)
Set the text size.
Definition: TAttText.h:46
The color creation and management class.
Definition: TColor.h:19
virtual void GetRGB(Float_t &r, Float_t &g, Float_t &b) const
Definition: TColor.h:51
GL Axis.
Definition: TGLAxis.h:22
Double_t fLabelsOffset
Definition: TGLAxis.h:41
Double_t * fTicks1
Definition: TGLAxis.h:33
void Init()
Default initialization.
Definition: TGLAxis.cxx:43
void TicksPositionsOpt()
Compute ticks positions. Linear and optimized.
Definition: TGLAxis.cxx:363
void 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.
Definition: TGLAxis.cxx:88
Double_t fTickMarksLength
Definition: TGLAxis.h:39
Double_t fAxisLength
Definition: TGLAxis.h:36
void PaintGLAxisBody()
Paint horizontal axis body at position (0,0,0)
Definition: TGLAxis.cxx:168
void TicksPositions(Option_t *opt="")
Compute ticks positions.
Definition: TGLAxis.cxx:298
Double_t fAngle1
Definition: TGLAxis.h:45
Double_t fGridLength
Definition: TGLAxis.h:43
void SetLabelsAngles(Double_t a1, Double_t a2, Double_t a3)
Set labels' angles.
Definition: TGLAxis.cxx:446
Int_t fNDiv
Definition: TGLAxis.h:27
Int_t fNDiv3
Definition: TGLAxis.h:30
void PaintGLAxisLabels()
Paint axis labels on the main tick marks.
Definition: TGLAxis.cxx:255
Double_t fAngle3
Definition: TGLAxis.h:47
void DoLabels()
Do labels.
Definition: TGLAxis.cxx:430
TGLAxis()
Constructor.
Definition: TGLAxis.cxx:35
Double_t fLabelsSize
Definition: TGLAxis.h:42
void TicksPositionsNoOpt()
Compute ticks positions. Linear and not optimized.
Definition: TGLAxis.cxx:331
Double_t fWmin
Definition: TGLAxis.h:37
Int_t fNDiv2
Definition: TGLAxis.h:29
TGLText * fText
Definition: TGLAxis.h:44
TString * fLabels
Definition: TGLAxis.h:35
Int_t fTickMarksOrientation
Definition: TGLAxis.h:40
Int_t fNDiv1
Definition: TGLAxis.h:28
Double_t fWmax
Definition: TGLAxis.h:38
void PaintGLAxisTickMarks()
Paint axis tick marks.
Definition: TGLAxis.cxx:185
Int_t fNTicks2
Definition: TGLAxis.h:32
Double_t fAngle2
Definition: TGLAxis.h:46
virtual ~TGLAxis()
Destructor.
Definition: TGLAxis.cxx:66
Double_t * fTicks2
Definition: TGLAxis.h:34
Int_t fNTicks1
Definition: TGLAxis.h:31
GL Text.
Definition: TGLText.h:19
void SetGLTextAngles(Double_t a1, Double_t a2, Double_t a3)
Set the text rotation angles.
Definition: TGLText.cxx:167
void PaintGLText(Double_t x, Double_t y, Double_t z, const char *text)
Draw text.
Definition: TGLText.cxx:96
void SetGLTextFont(Font_t fontnumber)
Definition: TGLText.cxx:177
static Float_t LineWidth()
Get the line-width, taking the global scaling into account.
Definition: TGLUtil.cxx:1938
static void Optimize(Double_t A1, Double_t A2, Int_t nold, Double_t &BinLow, Double_t &BinHigh, Int_t &nbins, Double_t &BWID, Option_t *option="")
Static function to compute reasonable axis limits.
Basic string class.
Definition: TString.h:131
Double_t y[n]
Definition: legend1.C:17
Double_t x[n]
Definition: legend1.C:17
T * Normal2Plane(const T v1[3], const T v2[3], const T v3[3], T normal[3])
Calculate a normal vector of a plane.
Definition: TMath.h:1179
Double_t ACos(Double_t)
Definition: TMath.h:658
Double_t ATan(Double_t)
Definition: TMath.h:665
Double_t Sqrt(Double_t x)
Definition: TMath.h:681
constexpr Double_t Pi()
Definition: TMath.h:38
Short_t Abs(Short_t d)
Definition: TMathBase.h:120