Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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
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
30
31////////////////////////////////////////////////////////////////////////////////
32/// Constructor.
33
34TGLAxis::TGLAxis(): TAttLine(1,1,1), TAttText(20,0.,1,42,0.04)
35{
36 Init();
37}
38
39////////////////////////////////////////////////////////////////////////////////
40/// Default initialization.
41
43{
44 fNDiv = fNDiv1 = fNDiv2 = fNDiv3 = 0;
45 fNTicks1 = fNTicks2 = 0;
46 fTicks1 = nullptr;
47 fTicks2 = nullptr;
48 fLabels = nullptr;
49 fText = nullptr;
50 fAngle1 = 90.;
51 fAngle2 = 0.;
52 fAngle3 = 0.;
53 fAxisLength = 0.;
54 fWmin = fWmax = 0.;
55 fTickMarksLength = 0.04; // % of fAxisLength
56 fTickMarksOrientation = 2; // can be 0, 1, 2, or 3
57 fLabelsOffset = 0.09; // % of fAxisLength
58 fLabelsSize = 0.06; // % of fAxisLength
59 fGridLength = 0.; // 0 means no grid
60}
61
62////////////////////////////////////////////////////////////////////////////////
63/// Destructor.
64
66{
67 if (fTicks1) delete [] fTicks1;
68 if (fTicks2) delete [] fTicks2;
69 if (fLabels) delete [] fLabels;
70 if (fText) delete fText;
71}
72
73////////////////////////////////////////////////////////////////////////////////
74/// Paint GL Axis.
75///
76/// - p1, p2 : Axis position in the 3D space.
77/// - wmin, wmax : Minimum and maximum values along the axis. wmin < wmax.
78/// - ndiv : Number of axis divisions. It is an integer in the form
79/// "ttsspp" where "tt" is the number of tertiary divisions,
80/// "ss" is the number of secondary divisions and "pp" the
81/// number of primary divisions.
82/// - opt : Options.
83/// "N" - By default the number of divisions is optimized to
84/// get a nice labeling. When option "N" is given, the
85/// number of divisions is not optimized.
86
87void TGLAxis::PaintGLAxis(const Double_t p1[3], const Double_t p2[3],
89 Option_t *opt)
90{
91 fNDiv = ndiv;
92 if (wmax<=wmin) {
93 fWmax = wmin;
94 fWmin = wmax;
95 } else {
96 fWmax = wmax;
97 fWmin = wmin;
98 }
99
100 // Compute the axis length.
101 Double_t x1 = p1[0], y1 = p1[1], z1 = p1[2];
102 Double_t x2 = p2[0], y2 = p2[1], z2 = p2[2];
104 (y2-y1)*(y2-y1)+
105 (z2-z1)*(z2-z1));
106
107 TicksPositions(opt);
108
109 DoLabels();
110
111 // Paint using GL
112 glPushMatrix();
113
114 // Translation
115 glTranslatef(x1, y1, z1);
116
117 // Rotation in the Z plane
118 Double_t phi=0;
119 Double_t normal[3];
120 normal[0] = 0;
121 normal[1] = 1;
122 normal[2] = 0;
123 if (z1!=z2) {
124 if (y2==y1 && x2==x1) {
125 if (z2<z1) phi = 90;
126 else phi = 270;
127 } else {
128 Double_t p3[3];
129 p3[0] = p2[0]; p3[1] = p2[1]; p3[2] = 0;
132 phi = -(90-(180/TMath::Pi())*phi);
133 }
134 glRotatef(phi, normal[0], normal[1], normal[2]);
135 }
136
137 // Rotation in the XY plane
138 Double_t theta = 0;
139 if (y2!=y1) {
140 if ((x2-x1)>0) {
141 theta = TMath::ATan((y2-y1)/(x2-x1));
142 theta = (180/TMath::Pi())*theta;
143 } else if ((x2-x1)<0) {
144 theta = TMath::ATan((y2-y1)/(x2-x1));
145 theta = 180+(180/TMath::Pi())*theta;
146 } else {
147 if (y2>y1) theta = 90;
148 else theta = 270;
149 }
150 } else {
151 if (x2<x1) theta = 180;
152 }
153 glRotatef(theta, 0., 0., 1.);
154
156
158
160
161 glPopMatrix();
162}
163
164////////////////////////////////////////////////////////////////////////////////
165/// Paint horizontal axis body at position (0,0,0)
166
168{
169 TColor *col;
170 Float_t red = 1.f, green = 1.f, blue = 1.f;
171 if ((col = gROOT->GetColor(GetLineColor())))
172 col->GetRGB(red, green, blue);
175 glBegin(GL_LINES);
176 glVertex3d(0., 0., 0.);
177 glVertex3d(fAxisLength, 0., 0.);
178 glEnd();
179}
180
181////////////////////////////////////////////////////////////////////////////////
182/// Paint axis tick marks.
183
185{
186 Int_t i;
187
188 // Ticks marks orientation;
189 Double_t yo=0, zo=0;
190 switch (fTickMarksOrientation) {
191 case 0:
192 yo = 0;
193 zo = 1;
194 break;
195 case 1:
196 yo = -1;
197 zo = 0;
198 break;
199 case 2:
200 yo = 0;
201 zo = -1;
202 break;
203 case 3:
204 yo = 1;
205 zo = 0;
206 break;
207 }
208
209 // Paint level 1 tick marks.
210 if (fTicks1) {
211 // Paint the tick marks, if needed.
212 if (fTickMarksLength) {
214 glBegin(GL_LINES);
215 for (i=0; i<fNTicks1 ; i++) {
216 glVertex3f( fTicks1[i], 0, 0);
217 glVertex3f( fTicks1[i], yo*tl, zo*tl);
218 }
219 glEnd();
220 }
221
222 // Paint the grid, if needed, on level 1 tick marks.
223 if (fGridLength) {
224 const UShort_t stipple = 0x8888;
227 glBegin(GL_LINES);
228 for (i=0; i<fNTicks1; i++) {
229 glVertex3f( fTicks1[i], 0, 0);
231 }
232 glEnd();
234 }
235 }
236
237 // Paint level 2 tick marks.
238 if (fTicks2) {
239 if (fTickMarksLength) {
241 glBegin(GL_LINES);
242 for (i=0; i<fNTicks2; i++) {
243 glVertex3f( fTicks2[i], 0, 0);
244 glVertex3f( fTicks2[i], yo*tl, zo*tl);
245 }
246 glEnd();
247 }
248 }
249}
250
251////////////////////////////////////////////////////////////////////////////////
252/// Paint axis labels on the main tick marks.
253
255{
256 if (!fLabelsSize) return;
257
258 Double_t x=0,y=0,z=0;
259 Int_t i;
260
261 if (!fText) {
262 fText = new TGLText();
267 }
269
270 switch (fTickMarksOrientation) {
271 case 0:
272 y = 0;
274 break;
275 case 1:
277 z = 0;
278 break;
279 case 2:
280 y = 0;
282 break;
283 case 3:
285 z = 0;
286 break;
287 }
288 for (i=0; i<fNDiv1+1 ; i++) {
289 x = fTicks1[i];
290 fText->PaintGLText(x,y,z,fLabels[i]);
291 }
292}
293
294////////////////////////////////////////////////////////////////////////////////
295/// Compute ticks positions.
296
298{
299 Bool_t optionNoopt /* , optionLog */;
300
301 if (strchr(opt,'N')) optionNoopt = kTRUE; else optionNoopt = kFALSE;
302 // if (strchr(opt,'G')) optionLog = kTRUE; else optionLog = kFALSE;
303
304 // Determine number of tick marks 1, 2 and 3.
305 fNDiv3 = fNDiv/10000;
306 fNDiv2 = (fNDiv-10000*fNDiv3)/100;
307 fNDiv1 = fNDiv%100;
308
309 // Clean the tick marks arrays if they exist.
310 if (fTicks1) {
311 delete [] fTicks1;
312 fTicks1 = nullptr;
313 }
314 if (fTicks2) {
315 delete [] fTicks2;
316 fTicks2 = nullptr;
317 }
318
319 // Compute the tick marks positions according to the options.
320 if (optionNoopt) {
322 } else {
324 }
325}
326
327////////////////////////////////////////////////////////////////////////////////
328/// Compute ticks positions. Linear and not optimized.
329
331{
332 Int_t i, j, k;
334
335 fNTicks1 = fNDiv1+1;
337
338 // Level 1 tick marks.
339 for (i=0; i<fNTicks1; i++) fTicks1[i] = i*step1;
340
341 // Level 2 tick marks.
342 if (fNDiv2) {
343 Double_t t2;
345 fNTicks2 = fNDiv1*(fNDiv2-1);
347 k = 0;
348 for (i=0; i<fNTicks1-1; i++) {
349 t2 = fTicks1[i]+step2;
350 for (j=0; j<fNDiv2-1; j++) {
351 fTicks2[k] = t2;
352 k++;
353 t2 = t2+step2;
354 }
355 }
356 }
357}
358
359////////////////////////////////////////////////////////////////////////////////
360/// Compute ticks positions. Linear and optimized.
361
363{
364 Int_t i, j, k, nDivOpt;
365 Double_t step1=0, step2=0, wmin2=0, wmax2=0;
368
369 // Level 1 tick marks.
372 step1, "");
373 fNDiv1 = nDivOpt;
374 fNTicks1 = fNDiv1+1;
376
378 Double_t w = fWmin;
379 i = 0;
380 while (w<=fWmax) {
381 fTicks1[i] = r*(w-wmin);
382 i++;
383 w = w+step1;
384 }
385
386 // Level 2 tick marks.
387 if (fNDiv2) {
388 Double_t t2;
391 step2, "");
392 fNDiv2 = nDivOpt;
398 k = 0;
399 for (i=0; i<fNTicks1-1; i++) {
400 t2 = fTicks1[i]+step2;
401 for (j=0; j<fNDiv2-1; j++) {
402 fTicks2[k] = t2;
403 k++;
404 t2 = t2+step2;
405 }
406 }
407 if (nTickl) {
408 t2 = fTicks1[0]-step2;
409 for (i=0; i<nTickl; i++) {
410 fTicks2[k] = t2;
411 k++;
412 t2 = t2-step2;
413 }
414 }
415 if (nTickr) {
417 for (i=0; i<nTickr; i++) {
418 fTicks2[k] = t2;
419 k++;
420 t2 = t2+step2;
421 }
422 }
423 }
424}
425
426////////////////////////////////////////////////////////////////////////////////
427/// Do labels.
428
430{
431 if (fLabels) delete [] fLabels;
432 fLabels = new TString[fNTicks1];
433 Int_t i;
434
435 // Make regular labels between fWmin and fWmax.
437 for (i=0; i<fNTicks1; i++) {
438 fLabels[i] = Form("%g",fWmin+i*dw);
439 }
440}
441
442////////////////////////////////////////////////////////////////////////////////
443/// Set labels' angles.
444
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
float Float_t
Float 4 bytes (float)
Definition RtypesCore.h:71
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t wmin
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t wmax
Option_t Option_t TPoint TPoint const char y1
#define gROOT
Definition TROOT.h:411
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2495
Line Attributes class.
Definition TAttLine.h:20
virtual Color_t GetLineColor() const
Return the line color.
Definition TAttLine.h:35
virtual Width_t GetLineWidth() const
Return the line width.
Definition TAttLine.h:37
Text Attributes class.
Definition TAttText.h:20
virtual void SetTextAlign(Short_t align=11)
Set the text alignment.
Definition TAttText.h:44
virtual Short_t GetTextAlign() const
Return the text alignment.
Definition TAttText.h:34
virtual Font_t GetTextFont() const
Return the text font.
Definition TAttText.h:37
virtual Color_t GetTextColor() const
Return the text color.
Definition TAttText.h:36
virtual void SetTextColor(Color_t tcolor=1)
Set the text color.
Definition TAttText.h:46
virtual void SetTextSize(Float_t tsize=1)
Set the text size.
Definition TAttText.h:49
The color creation and management class.
Definition TColor.h:22
virtual void GetRGB(Float_t &r, Float_t &g, Float_t &b) const
Definition TColor.h:55
static Int_t GetColor(const char *hexcolor)
Static method returning color number for color specified by hex color string of form: "#rrggbb",...
Definition TColor.cxx:1926
Double_t fLabelsOffset
Definition TGLAxis.h:41
Double_t * fTicks1
Definition TGLAxis.h:33
void Init()
Default initialization.
Definition TGLAxis.cxx:42
~TGLAxis() override
Destructor.
Definition TGLAxis.cxx:65
void TicksPositionsOpt()
Compute ticks positions. Linear and optimized.
Definition TGLAxis.cxx:362
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:87
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:167
void TicksPositions(Option_t *opt="")
Compute ticks positions.
Definition TGLAxis.cxx:297
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:445
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:254
Double_t fAngle3
Definition TGLAxis.h:47
void DoLabels()
Do labels.
Definition TGLAxis.cxx:429
TGLAxis()
Constructor.
Definition TGLAxis.cxx:34
Double_t fLabelsSize
Definition TGLAxis.h:42
void TicksPositionsNoOpt()
Compute ticks positions. Linear and not optimized.
Definition TGLAxis.cxx:330
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:184
Int_t fNTicks2
Definition TGLAxis.h:32
Double_t fAngle2
Definition TGLAxis.h:46
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:164
void PaintGLText(Double_t x, Double_t y, Double_t z, const char *text)
Draw text.
Definition TGLText.cxx:93
void SetGLTextFont(Font_t fontnumber)
Definition TGLText.cxx:174
static Float_t LineWidth()
Get the line-width, taking the global scaling into account.
Definition TGLUtil.cxx:1934
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:138
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])
Calculates a normal vector of a plane.
Definition TMath.h:1299
Double_t ACos(Double_t)
Returns the principal value of the arc cosine of x, expressed in radians.
Definition TMath.h:643
Double_t ATan(Double_t)
Returns the principal value of the arc tangent of x, expressed in radians.
Definition TMath.h:651
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:673
constexpr Double_t Pi()
Definition TMath.h:40
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:124