Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGLCylinder.cxx
Go to the documentation of this file.
1// @(#)root/gl:$Id$
2// Author: Timur Pocheptsov 03/08/2004
3// NOTE: This code moved from obsoleted TGLSceneObject.h / .cxx - see these
4// attic files for previous CVS history
5
6/*************************************************************************
7 * Copyright (C) 1995-2006, Rene Brun and Fons Rademakers. *
8 * All rights reserved. *
9 * *
10 * For the licensing terms see $ROOTSYS/LICENSE. *
11 * For the list of contributors see $ROOTSYS/README/CREDITS. *
12 *************************************************************************/
13
14#include "TGLCylinder.h"
15#include "TGLRnrCtx.h"
16#include "TGLIncludes.h"
17
18#include "TBuffer3D.h"
19#include "TBuffer3DTypes.h"
20#include "TMath.h"
21
22// For debug tracing
23#include "TClass.h"
24#include "TError.h"
25
28
30{
31protected:
32 // active LOD (level of detail) - quality
34
37
38 //normals for top and bottom (for cuts)
41
42 void GetNormal(const TGLVertex3 &vertex, TGLVector3 &normal)const;
45
46public:
49 virtual ~TGLMesh() { }
50 virtual void Draw() const = 0;
51};
52
53//segment contains 3 quad strips:
54//one for inner and outer sides, two for top and bottom
55class TubeSegMesh : public TGLMesh {
56private:
57 // Allocate space for highest quality (LOD) meshes
60
61public:
63 Double_t phi1, Double_t phi2, const TGLVector3 &l = gLowNormalDefault,
65
66 void Draw() const override;
67};
68
69//four quad strips:
70//outer, inner, top, bottom
71class TubeMesh : public TGLMesh
72{
73private:
74 // Allocate space for highest quality (LOD) meshes
77
78public:
81
82 void Draw() const override;
83};
84
85//One quad mesh and 2 triangle funs
86class TCylinderMesh : public TGLMesh {
87private:
88 // Allocate space for highest quality (LOD) meshes
91
92public:
95
96 void Draw() const override;
97};
98
99//One quad mesh and 2 triangle fans
101{
102private:
103 // Allocate space for highest quality (LOD) meshes
106
107public:
110 void Draw() const override;
111};
112
114 const TGLVector3 &l, const TGLVector3 &h) :
115 fLOD(LOD),
116 fRmin1(r1), fRmax1(r2), fRmin2(r3), fRmax2(r4),
117 fDz(dz), fNlow(l), fNhigh(h)
118{
119 // constructor
120}
121
122////////////////////////////////////////////////////////////////////////////////
123/// get normal
124
126{
127 if( fDz < 1.e-10 ) {
128 n[0] = 0.;
129 n[1] = 0.;
130 n[2] = 1.;
131 }
132 Double_t z = (fRmax1 - fRmax2) / (2 * fDz);
133 Double_t mag = TMath::Sqrt(v[0] * v[0] + v[1] * v[1] + z * z);
134 if( mag > 1.e-10 ) {
135 n[0] = v[0] / mag;
136 n[1] = v[1] / mag;
137 n[2] = z / mag;
138 } else {
139 n[0] = v[0];
140 n[1] = v[1];
141 n[2] = z;
142 }
143}
144
145////////////////////////////////////////////////////////////////////////////////
146/// get Z coordinate
147
149{
150 Double_t newz = 0;
151 if (z < 0) newz = -fDz - (x * fNlow[0] + y * fNlow[1]) / fNlow[2];
152 else newz = fDz - (x * fNhigh[0] + y * fNhigh[1]) / fNhigh[2];
153
154 return newz;
155}
156
157////////////////////////////////////////////////////////////////////////////////
158/// make vertex
159
161{
162 static TGLVertex3 vert(0., 0., 0.);
163 vert[0] = x;
164 vert[1] = y;
165 vert[2] = GetZcoord(x, y, z);
166
167 return vert;
168}
169
170////////////////////////////////////////////////////////////////////////////////
171
173 Double_t phi1, Double_t phi2,
174 const TGLVector3 &l, const TGLVector3 &h)
175 :TGLMesh(LOD, r1, r2, r3, r4, dz, l, h), fMesh(), fNorm()
176
177{
178 // constructor
179 const Double_t delta = (phi2 - phi1) / LOD;
180 Double_t currAngle = phi1;
181
182 Bool_t even = kTRUE;
183 Double_t c = TMath::Cos(currAngle);
184 Double_t s = TMath::Sin(currAngle);
185 const Int_t topShift = (fLOD + 1) * 4 + 8;
186 const Int_t botShift = (fLOD + 1) * 6 + 8;
187 Int_t j = 4 * (fLOD + 1) + 2;
188
189 //defining all three strips here, first strip is non-closed here
190 for (Int_t i = 0, e = (fLOD + 1) * 2; i < e; ++i) {
191 if (even) {
192 fMesh[i] = MakeVertex(fRmax2 * c, fRmax2 * s, fDz);
193 fMesh[j] = MakeVertex(fRmin2 * c, fRmin2 * s, fDz);
194 fMesh[i + topShift] = MakeVertex(fRmin2 * c, fRmin2 * s, fDz);
195 fMesh[i + botShift] = MakeVertex(fRmax1 * c, fRmax1 * s, - fDz);
196 GetNormal(fMesh[j], fNorm[j]);
197 fNorm[j].Negate();
198 even = kFALSE;
199 } else {
200 fMesh[i] = MakeVertex(fRmax1 * c, fRmax1 * s, - fDz);
201 fMesh[j + 1] = MakeVertex(fRmin1 * c, fRmin1 * s, -fDz);
202 fMesh[i + topShift] = MakeVertex(fRmax2 * c, fRmax2 * s, fDz);
203 fMesh[i + botShift] = MakeVertex(fRmin1 * c, fRmin1 * s, - fDz);
204 GetNormal(fMesh[j + 1], fNorm[j + 1]);
205 fNorm[j + 1].Negate();
206 even = kTRUE;
207 currAngle += delta;
208 c = TMath::Cos(currAngle);
209 s = TMath::Sin(currAngle);
210 j -= 2;
211 }
212
213 GetNormal(fMesh[i], fNorm[i]);
214 fNorm[i + topShift] = fNhigh;
215 fNorm[i + botShift] = fNlow;
216 }
217
218 //closing first strip
219 Int_t ind = 2 * (fLOD + 1);
220 TGLVector3 norm(0., 0., 0.);
221
222 fMesh[ind] = fMesh[ind - 2];
223 fMesh[ind + 1] = fMesh[ind - 1];
224 fMesh[ind + 2] = fMesh[ind + 4];
225 fMesh[ind + 3] = fMesh[ind + 5];
226 TMath::Normal2Plane(fMesh[ind].CArr(), fMesh[ind + 1].CArr(), fMesh[ind + 2].CArr(),
227 norm.Arr());
228 fNorm[ind] = norm;
229 fNorm[ind + 1] = norm;
230 fNorm[ind + 2] = norm;
231 fNorm[ind + 3] = norm;
232
233 ind = topShift - 4;
234 fMesh[ind] = fMesh[ind - 2];
235 fMesh[ind + 1] = fMesh[ind - 1];
236 fMesh[ind + 2] = fMesh[0];
237 fMesh[ind + 3] = fMesh[1];
238 TMath::Normal2Plane(fMesh[ind].CArr(), fMesh[ind + 1].CArr(), fMesh[ind + 2].CArr(),
239 norm.Arr());
240 fNorm[ind] = norm;
241 fNorm[ind + 1] = norm;
242 fNorm[ind + 2] = norm;
243 fNorm[ind + 3] = norm;
244}
245
246////////////////////////////////////////////////////////////////////////////////
247///Tube segment is drawn as three quad strips
248///1. enabling vertex arrays
249
251{
252 glEnableClientState(GL_VERTEX_ARRAY);
253 glEnableClientState(GL_NORMAL_ARRAY);
254 //2. setting arrays
255 glVertexPointer(3, GL_DOUBLE, sizeof(TGLVertex3), fMesh[0].CArr());
256 glNormalPointer(GL_DOUBLE, sizeof(TGLVector3), fNorm[0].CArr());
257 //3. draw first strip
258 glDrawArrays(GL_QUAD_STRIP, 0, 4 * (fLOD + 1) + 8);
259 //4. draw top and bottom strips
260 glDrawArrays(GL_QUAD_STRIP, 4 * (fLOD + 1) + 8, 2 * (fLOD + 1));
261 glDrawArrays(GL_QUAD_STRIP, 6 * (fLOD + 1) + 8, 2 * (fLOD + 1));
262
263 glDisableClientState(GL_VERTEX_ARRAY);
264 glDisableClientState(GL_NORMAL_ARRAY);
265}
266
267////////////////////////////////////////////////////////////////////////////////
268/// constructor
269
271 const TGLVector3 &l, const TGLVector3 &h)
272 :TGLMesh(LOD, r1, r2, r3, r4, z, l, h), fMesh(), fNorm()
273{
274 const Double_t delta = TMath::TwoPi() / fLOD;
275 Double_t currAngle = 0.;
276
277 Bool_t even = kTRUE;
278 Double_t c = TMath::Cos(currAngle);
279 Double_t s = TMath::Sin(currAngle);
280
281 const Int_t topShift = (fLOD + 1) * 4;
282 const Int_t botShift = (fLOD + 1) * 6;
283 Int_t j = 4 * (fLOD + 1) - 2;
284
285 //defining all four strips here
286 for (Int_t i = 0, e = (fLOD + 1) * 2; i < e; ++i) {
287 if (even) {
288 fMesh[i] = MakeVertex(fRmax2 * c, fRmax2 * s, fDz);
289 fMesh[j] = MakeVertex(fRmin2 * c, fRmin2 * s, fDz);
290 fMesh[i + topShift] = MakeVertex(fRmin2 * c, fRmin2 * s, fDz);
291 fMesh[i + botShift] = MakeVertex(fRmax1 * c, fRmax1 * s, - fDz);
292 GetNormal(fMesh[j], fNorm[j]);
293 fNorm[j].Negate();
294 even = kFALSE;
295 } else {
296 fMesh[i] = MakeVertex(fRmax1 * c, fRmax1 * s, - fDz);
297 fMesh[j + 1] = MakeVertex(fRmin1 * c, fRmin1 * s, -fDz);
298 fMesh[i + topShift] = MakeVertex(fRmax2 * c, fRmax2 * s, fDz);
299 fMesh[i + botShift] = MakeVertex(fRmin1 * c, fRmin1 * s, - fDz);
300 GetNormal(fMesh[j + 1], fNorm[j + 1]);
301 fNorm[j + 1].Negate();
302 even = kTRUE;
303 currAngle += delta;
304 c = TMath::Cos(currAngle);
305 s = TMath::Sin(currAngle);
306 j -= 2;
307 }
308
309 GetNormal(fMesh[i], fNorm[i]);
310 fNorm[i + topShift] = fNhigh;
311 fNorm[i + botShift] = fNlow;
312 }
313}
314
315////////////////////////////////////////////////////////////////////////////////
316///Tube is drawn as four quad strips
317
318void TubeMesh::Draw() const
319{
320 glEnableClientState(GL_VERTEX_ARRAY);
321 glEnableClientState(GL_NORMAL_ARRAY);
322
323 glVertexPointer(3, GL_DOUBLE, sizeof(TGLVertex3), fMesh[0].CArr());
324 glNormalPointer(GL_DOUBLE, sizeof(TGLVector3), fNorm[0].CArr());
325 //draw outer and inner strips
326 glDrawArrays(GL_QUAD_STRIP, 0, 2 * (fLOD + 1));
327 glDrawArrays(GL_QUAD_STRIP, 2 * (fLOD + 1), 2 * (fLOD + 1));
328 //draw top and bottom strips
329 glDrawArrays(GL_QUAD_STRIP, 4 * (fLOD + 1), 2 * (fLOD + 1));
330 glDrawArrays(GL_QUAD_STRIP, 6 * (fLOD + 1), 2 * (fLOD + 1));
331 //5. disabling vertex arrays
332 glDisableClientState(GL_VERTEX_ARRAY);
333 glDisableClientState(GL_NORMAL_ARRAY);
334}
335
336////////////////////////////////////////////////////////////////////////////////
337/// constructor
338
340 const TGLVector3 &l, const TGLVector3 &h)
341 :TGLMesh(LOD, 0., r1, 0., r2, dz, l, h), fMesh(), fNorm()
342{
343 const Double_t delta = TMath::TwoPi() / fLOD;
344 Double_t currAngle = 0.;
345
346 Bool_t even = kTRUE;
347 Double_t c = TMath::Cos(currAngle);
348 Double_t s = TMath::Sin(currAngle);
349
350 //central point of top fan
351 Int_t topShift = (fLOD + 1) * 2;
352 fMesh[topShift][0] = fMesh[topShift][1] = 0., fMesh[topShift][2] = fDz;
353 fNorm[topShift] = fNhigh;
354 ++topShift;
355
356 //central point of bottom fun
357 Int_t botShift = topShift + 2 * (fLOD + 1);
358 fMesh[botShift][0] = fMesh[botShift][1] = 0., fMesh[botShift][2] = -fDz;
359 fNorm[botShift] = fNlow;
360 ++botShift;
361
362 //defining 1 strip and 2 fans
363 for (Int_t i = 0, e = (fLOD + 1) * 2, j = 0; i < e; ++i) {
364 if (even) {
365 fMesh[i] = MakeVertex(fRmax2 * c, fRmax2 * s, fDz);
366 fMesh[j + topShift] = MakeVertex(fRmin2 * c, fRmin2 * s, fDz);
367 fMesh[j + botShift] = MakeVertex(fRmax1 * c, fRmax1 * s, - fDz);
368 even = kFALSE;
369 } else {
370 fMesh[i] = MakeVertex(fRmax1 * c, fRmax1 * s, - fDz);
371 even = kTRUE;
372 currAngle += delta;
373 c = TMath::Cos(currAngle);
374 s = TMath::Sin(currAngle);
375 ++j;
376 }
377
378 GetNormal(fMesh[i], fNorm[i]);
379 fNorm[i + topShift] = fNhigh;
380 fNorm[i + botShift] = fNlow;
381 }
382}
383
384////////////////////////////////////////////////////////////////////////////////
385/// draw cylinder mesh
386
388{
389 glEnableClientState(GL_VERTEX_ARRAY);
390 glEnableClientState(GL_NORMAL_ARRAY);
391
392 glVertexPointer(3, GL_DOUBLE, sizeof(TGLVertex3), fMesh[0].CArr());
393 glNormalPointer(GL_DOUBLE, sizeof(TGLVector3), fNorm[0].CArr());
394
395 //draw quad strip
396 glDrawArrays(GL_QUAD_STRIP, 0, 2 * (fLOD + 1));
397 //draw top and bottom funs
398 glDrawArrays(GL_TRIANGLE_FAN, 2 * (fLOD + 1), fLOD + 2);
399 glDrawArrays(GL_TRIANGLE_FAN, 3 * (fLOD + 1) + 1, fLOD + 2);
400
401 glDisableClientState(GL_VERTEX_ARRAY);
402 glDisableClientState(GL_NORMAL_ARRAY);
403}
404
405////////////////////////////////////////////////////////////////////////////////
406///One quad mesh and two fans
407
409 Double_t phi2, const TGLVector3 &l,
410 const TGLVector3 &h)
411 :TGLMesh(LOD, 0., r1, 0., r2, dz, l, h), fMesh(), fNorm()
412{
413 Double_t delta = (phi2 - phi1) / fLOD;
414 Double_t currAngle = phi1;
415
416 Bool_t even = kTRUE;
417 Double_t c = TMath::Cos(currAngle);
418 Double_t s = TMath::Sin(currAngle);
419
420 const TGLVertex3 vTop(0., 0., fDz);
421 const TGLVertex3 vBot(0., 0., - fDz);
422
423 //center of top fan
424 Int_t topShift = (fLOD + 1) * 2 + 8;
425 fMesh[topShift] = vTop;
426 fNorm[topShift] = fNhigh;
427 ++topShift;
428
429 //center of bottom fan
430 Int_t botShift = topShift + fLOD + 1;
431 fMesh[botShift] = vBot;
432 fNorm[botShift] = fNlow;
433 ++botShift;
434
435 //defining strip and two fans
436 //strip is not closed here
437 Int_t i = 0;
438 for (Int_t e = (fLOD + 1) * 2, j = 0; i < e; ++i) {
439 if (even) {
440 fMesh[i] = MakeVertex(fRmax2 * c, fRmax2 * s, fDz);
441 fMesh[j + topShift] = MakeVertex(fRmax2 * c, fRmax2 * s, fDz);
442 fMesh[j + botShift] = MakeVertex(fRmax1 * c, fRmax1 * s, - fDz);
443 even = kFALSE;
444 fNorm[j + topShift] = fNhigh;
445 fNorm[j + botShift] = fNlow;
446 } else {
447 fMesh[i] = MakeVertex(fRmax1 * c, fRmax1 * s, - fDz);
448 even = kTRUE;
449 currAngle += delta;
450 c = TMath::Cos(currAngle);
451 s = TMath::Sin(currAngle);
452 ++j;
453 }
454
455 GetNormal(fMesh[i], fNorm[i]);
456 }
457
458 //closing first strip
459 Int_t ind = 2 * (fLOD + 1);
460 TGLVector3 norm(0., 0., 0.);
461
462 fMesh[ind] = fMesh[ind - 2];
463 fMesh[ind + 1] = fMesh[ind - 1];
464 fMesh[ind + 2] = vTop;
465 fMesh[ind + 3] = vBot;
466 TMath::Normal2Plane(fMesh[ind].CArr(), fMesh[ind + 1].CArr(), fMesh[ind + 2].CArr(),
467 norm.Arr());
468 fNorm[ind] = norm;
469 fNorm[ind + 1] = norm;
470 fNorm[ind + 2] = norm;
471 fNorm[ind + 3] = norm;
472
473 ind += 4;
474 fMesh[ind] = vTop;
475 fMesh[ind + 1] = vBot;
476 fMesh[ind + 2] = fMesh[0];
477 fMesh[ind + 3] = fMesh[1];
478 TMath::Normal2Plane(fMesh[ind].CArr(), fMesh[ind + 1].CArr(), fMesh[ind + 2].CArr(),
479 norm.Arr());
480 fNorm[ind] = norm;
481 fNorm[ind + 1] = norm;
482 fNorm[ind + 2] = norm;
483 fNorm[ind + 3] = norm;
484}
485
486////////////////////////////////////////////////////////////////////////////////
487///Cylinder segment is drawn as one quad strip and
488///two triangle fans
489///1. enabling vertex arrays
490
492{
493 glEnableClientState(GL_VERTEX_ARRAY);
494 glEnableClientState(GL_NORMAL_ARRAY);
495 //2. setting arrays
496 glVertexPointer(3, GL_DOUBLE, sizeof(TGLVertex3), fMesh[0].CArr());
497 glNormalPointer(GL_DOUBLE, sizeof(TGLVector3), fNorm[0].CArr());
498 //3. draw quad strip
499 glDrawArrays(GL_QUAD_STRIP, 0, 2 * (fLOD + 1) + 8);
500 //4. draw top and bottom funs
501 glDrawArrays(GL_TRIANGLE_FAN, 2 * (fLOD + 1) + 8, fLOD + 2);
502 // glDrawArrays(GL_TRIANGLE_FAN, 3 * (fLOD + 1) + 9, fLOD + 2);
503 //5. disabling vertex arrays
504 glDisableClientState(GL_VERTEX_ARRAY);
505 glDisableClientState(GL_NORMAL_ARRAY);
506}
507
508
509/** \class TGLCylinder
510\ingroup opengl
511Implements a native ROOT-GL cylinder that can be rendered at
512different levels of detail.
513*/
514
516
517////////////////////////////////////////////////////////////////////////////////
518/// Copy out relevant parts of buffer - we create and delete mesh
519/// parts on demand in DirectDraw() and they are DL cached
520
522 TGLLogicalShape(buffer)
523{
524 fDLSize = 14;
525
526 fR1 = buffer.fRadiusInner;
527 fR2 = buffer.fRadiusOuter;
528 fR3 = buffer.fRadiusInner;
529 fR4 = buffer.fRadiusOuter;
530 fDz = buffer.fHalfLength;
531
534
535 switch (buffer.Type())
536 {
537 default:
539 {
541 fPhi1 = 0;
542 fPhi2 = 360;
543 break;
544 }
545
548 {
549 fSegMesh = kTRUE;
550
551 const TBuffer3DTubeSeg * segBuffer = dynamic_cast<const TBuffer3DTubeSeg *>(&buffer);
552 if (!segBuffer) {
553 Error("TGLCylinder::TGLCylinder", "cannot cast TBuffer3D");
554 return;
555 }
556
557 fPhi1 = segBuffer->fPhiMin;
558 fPhi2 = segBuffer->fPhiMax;
559 if (fPhi2 < fPhi1) fPhi2 += 360.;
562
563 if (buffer.Type() == TBuffer3DTypes::kCutTube) {
564 const TBuffer3DCutTube * cutBuffer = dynamic_cast<const TBuffer3DCutTube *>(&buffer);
565 if (!cutBuffer) {
566 Error("TGLCylinder::TGLCylinder", "cannot cast TBuffer3D");
567 return;
568 }
569
570 for (UInt_t i =0; i < 3; i++) {
571 fLowPlaneNorm[i] = cutBuffer->fLowPlaneNorm[i];
572 fHighPlaneNorm[i] = cutBuffer->fHighPlaneNorm[i];
573 }
574 }
575 break;
576 }
577 }
578}
579
580////////////////////////////////////////////////////////////////////////////////
581///destructor
582
584{
585}
586
587////////////////////////////////////////////////////////////////////////////////
588/// Return display-list offset for given LOD.
589/// Calculation based on what is done in virtual QuantizeShapeLOD below.
590
592{
593 UInt_t off = 0;
594 if (lod >= 100) off = 0;
595 else if (lod < 10) off = lod / 2;
596 else off = lod / 10 + 4;
597 return off;
598}
599
600////////////////////////////////////////////////////////////////////////////////
601/// Factor in scene/viewer LOD and quantize.
602
604{
605 Int_t lod = ((Int_t)shapeLOD * (Int_t)combiLOD) / 100;
606
607 if (lod >= 100)
608 {
609 lod = 100;
610 }
611 else if (lod > 10)
612 { // Round LOD above 10 to nearest 10
613 Double_t quant = 0.1 * ((static_cast<Double_t>(lod)) + 0.5);
614 lod = 10 * static_cast<Int_t>(quant);
615 }
616 else
617 { // Round LOD below 10 to nearest 2
618 Double_t quant = 0.5 * ((static_cast<Double_t>(lod)) + 0.5);
619 lod = 2 * static_cast<Int_t>(quant);
620 }
621 return static_cast<Short_t>(lod);
622}
623
624////////////////////////////////////////////////////////////////////////////////
625/// Debug tracing
626
628{
629 if (gDebug > 4) {
630 Info("TGLCylinder::DirectDraw", "this %zd (class %s) LOD %d",
631 (size_t)this, IsA()->GetName(), rnrCtx.ShapeLOD());
632 }
633
634 // As we are now support display list caching we can create, draw and
635 // delete mesh parts of suitable LOD (quality) here - it will be cached
636 // into a display list by base-class TGLLogicalShape::Draw(),
637 // against our id and the LOD value. So this will only occur once
638 // for a certain cylinder/LOD combination
639 std::vector<TGLMesh *> meshParts;
640
641 // Create mesh parts
642 if (!fSegMesh) {
643 meshParts.push_back(new TubeMesh (rnrCtx.ShapeLOD(), fR1, fR2, fR3, fR4,
645 } else {
646 meshParts.push_back(new TubeSegMesh(rnrCtx.ShapeLOD(), fR1, fR2, fR3, fR4,
647 fDz, fPhi1, fPhi2,
649 }
650
651 // Draw mesh parts
652 for (UInt_t i = 0; i < meshParts.size(); ++i) meshParts[i]->Draw();
653
654 // Delete mesh parts
655 for (UInt_t i = 0; i < meshParts.size(); ++i) {
656 delete meshParts[i];
657 meshParts[i] = nullptr;//not to have invalid pointer for pseudo-destructor call :)
658 }
659}
#define GL_TRIANGLE_FAN
Definition GL_glu.h:289
#define GL_QUAD_STRIP
Definition GL_glu.h:291
#define c(i)
Definition RSha256.hxx:101
#define h(i)
Definition RSha256.hxx:106
#define e(i)
Definition RSha256.hxx:103
int Int_t
Definition RtypesCore.h:45
short Short_t
Definition RtypesCore.h:39
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
#define ClassImp(name)
Definition Rtypes.h:382
void Info(const char *location, const char *msgfmt,...)
Use this function for informational messages.
Definition TError.cxx:218
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:185
TGLVector3 gHighNormalDefault(0., 0., 1.)
TGLVector3 gLowNormalDefault(0., 0., -1.)
Int_t gDebug
Definition TROOT.cxx:597
Cut tube segment description class - see TBuffer3DTypes for producer classes.
Definition TBuffer3D.h:212
Double_t fLowPlaneNorm[3]
Definition TBuffer3D.h:224
Double_t fHighPlaneNorm[3]
Definition TBuffer3D.h:225
Tube segment description class - see TBuffer3DTypes for producer classes.
Definition TBuffer3D.h:185
Double_t fPhiMax
Definition TBuffer3D.h:204
Double_t fPhiMin
Definition TBuffer3D.h:203
Complete tube description class - see TBuffer3DTypes for producer classes.
Definition TBuffer3D.h:157
Double_t fRadiusInner
Definition TBuffer3D.h:175
Double_t fRadiusOuter
Definition TBuffer3D.h:176
Double_t fHalfLength
Definition TBuffer3D.h:177
Int_t Type() const
Definition TBuffer3D.h:85
TGLVertex3 fMesh[(TGLRnrCtx::kLODHigh+1) *4+2]
TGLVector3 fNorm[(TGLRnrCtx::kLODHigh+1) *4+2]
TCylinderMesh(UInt_t LOD, Double_t r1, Double_t r2, Double_t dz, const TGLVector3 &l=gLowNormalDefault, const TGLVector3 &h=gHighNormalDefault)
constructor
void Draw() const override
draw cylinder mesh
TGLVector3 fNorm[(TGLRnrCtx::kLODHigh+1) *4+10]
TGLVertex3 fMesh[(TGLRnrCtx::kLODHigh+1) *4+10]
TCylinderSegMesh(UInt_t LOD, Double_t r1, Double_t r2, Double_t dz, Double_t phi1, Double_t phi2, const TGLVector3 &l=gLowNormalDefault, const TGLVector3 &h=gHighNormalDefault)
One quad mesh and two fans.
void Draw() const override
Cylinder segment is drawn as one quad strip and two triangle fans.
Implements a native ROOT-GL cylinder that can be rendered at different levels of detail.
Definition TGLCylinder.h:22
TGLVector3 fHighPlaneNorm
Definition TGLCylinder.h:28
void DirectDraw(TGLRnrCtx &rnrCtx) const override
Debug tracing.
Double_t fR2
Definition TGLCylinder.h:24
Double_t fR1
Definition TGLCylinder.h:24
TGLCylinder(const TBuffer3DTube &buffer)
Copy out relevant parts of buffer - we create and delete mesh parts on demand in DirectDraw() and the...
Double_t fDz
Definition TGLCylinder.h:25
TClass * IsA() const override
Definition TGLCylinder.h:50
~TGLCylinder() override
destructor
Short_t QuantizeShapeLOD(Short_t shapeLOD, Short_t combiLOD) const override
Factor in scene/viewer LOD and quantize.
Double_t fR4
Definition TGLCylinder.h:24
Double_t fR3
Definition TGLCylinder.h:24
TGLVector3 fLowPlaneNorm
Definition TGLCylinder.h:28
Double_t fPhi1
Definition TGLCylinder.h:26
Bool_t fSegMesh
Definition TGLCylinder.h:29
UInt_t DLOffset(Short_t lod) const override
Return display-list offset for given LOD.
Double_t fPhi2
Definition TGLCylinder.h:26
Abstract logical shape - a GL 'drawable' - base for all shapes - faceset sphere etc.
Int_t fDLSize
display-list id base
UInt_t fLOD
virtual ~TGLMesh()
const TGLVertex3 & MakeVertex(Double_t x, Double_t y, Double_t z) const
make vertex
Double_t fDz
Double_t fRmin2
Double_t fRmin1
virtual void Draw() const =0
TGLMesh(UInt_t LOD, Double_t r1, Double_t r2, Double_t r3, Double_t r4, Double_t dz, const TGLVector3 &l=gLowNormalDefault, const TGLVector3 &h=gHighNormalDefault)
Double_t GetZcoord(Double_t x, Double_t y, Double_t z) const
get Z coordinate
void GetNormal(const TGLVertex3 &vertex, TGLVector3 &normal) const
get normal
Double_t fRmax2
TGLVector3 fNlow
TGLVector3 fNhigh
Double_t fRmax1
The TGLRnrCtx class aggregates data for a given redering context as needed by various parts of the RO...
Definition TGLRnrCtx.h:41
Short_t ShapeLOD() const
Definition TGLRnrCtx.h:177
3 component (x/y/z) vector class.
Definition TGLUtil.h:248
3 component (x/y/z) vertex class.
Definition TGLUtil.h:84
Double_t * Arr()
Definition TGLUtil.h:127
void Negate()
Definition TGLUtil.h:141
TGLVertex3 fMesh[(TGLRnrCtx::kLODHigh+1) *8]
TubeMesh(UInt_t LOD, Double_t r1, Double_t r2, Double_t r3, Double_t r4, Double_t dz, const TGLVector3 &l=gLowNormalDefault, const TGLVector3 &h=gHighNormalDefault)
constructor
TGLVector3 fNorm[(TGLRnrCtx::kLODHigh+1) *8]
void Draw() const override
Tube is drawn as four quad strips.
TGLVertex3 fMesh[(TGLRnrCtx::kLODHigh+1) *8+8]
TGLVector3 fNorm[(TGLRnrCtx::kLODHigh+1) *8+8]
TubeSegMesh(UInt_t LOD, Double_t r1, Double_t r2, Double_t r3, Double_t r4, Double_t dz, Double_t phi1, Double_t phi2, const TGLVector3 &l=gLowNormalDefault, const TGLVector3 &h=gHighNormalDefault)
void Draw() const override
Tube segment is drawn as three quad strips.
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
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:1212
constexpr Double_t DegToRad()
Conversion from degree to radian: .
Definition TMath.h:79
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:662
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Definition TMath.h:594
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:588
constexpr Double_t TwoPi()
Definition TMath.h:44
th1 Draw()
TLine l
Definition textangle.C:4