Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TXTRU.cxx
Go to the documentation of this file.
1// @@(#)root/g3d:$Id$
2// Author: Robert Hatcher (rhatcher@fnal.gov) 2000.09.06
3
4/*************************************************************************
5 * Copyright (C) 1995-2000, 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 "TXTRU.h"
13
14#include "TBuffer3D.h"
15#include "TBuffer3DTypes.h"
16#include "TGeometry.h"
17#include "TMath.h"
18
19#include <iostream>
20#include <iomanip>
21
22
23/** \class TXTRU
24\ingroup g3d
25A poly-extrusion.
26
27\image html g3d_xtru.png
28
29XTRU is a poly-extrusion with fixed outline shape in x-y,
30a sequence of z extents (segments) and two end faces perpendicular
31to the z axis. The x-y outline is defined by an ordered list of
32points; the overall scale of the outline scales linearly between
33z points and the center can have an x-y offset specified
34at each segment end.
35
36A TXTRU has the following parameters:
37
38 - name: name of the shape
39 - title: shape's title
40 - material: (see TMaterial)
41 - nxy: number of x-y vertex points constituting the outline --
42 this number should be at least 3
43 - nz: number of planes perpendicular to the z axis where
44 the scaling dimension of the section is given --
45 this number should be at least 2
46 - Xvtx: array [nxy] of X coordinates of vertices
47 - Yvtx: array [nxy] of Y coordinates of vertices
48 - z: array [nz] of z plane positions
49 - scale: array [nz] of scale factors
50 - x0: array [nz] of x offsets
51 - y0: array [nz] of y offsets
52
53All XTRU shapes are correctly rendered in wire mode but can encounter
54difficulty when rendered as a solid with hidden surfaces. These
55exceptions occur if the outline shape is not a convex polygon.
56Both the X3D and OpenGL renderers expect polygons to be convex.
57The OpenGL spec specifies that points defining a polygon using the
58GL_POLYGON primitive may be rendered as the convex hull of that set.
59
60Solid rendering under X3D can also give unexpected artifacts if
61the combination of x-y-z offsets and scales for the segments are
62chosen in such a manner that they represent a concave shape when
63sliced along a plane parallel to the z axis.
64
65Choosing sets of point that represent a malformed polygon is
66not supported, but testing for such a condition is not implemented
67and thus it is left to the user to avoid this mistake.
68
69\image html g3d_polytype.png
70*/
71
72////////////////////////////////////////////////////////////////////////////////
73/// TXTRU shape - default constructor
74
76{
77}
78
79////////////////////////////////////////////////////////////////////////////////
80/// TXTRU shape - normal constructor
81///
82/// Parameters of Nxy positions must be entered via TXTRU::DefineVertex
83/// Parameters of Nz positions must be entered via TXTRU::DefineSection
84
85TXTRU::TXTRU(const char *name, const char *title, const char *material,
87 : TShape (name,title,material)
88{
89 // start in a known state even if "Error" is encountered
90 if ( nxy < 3 ) {
91 Error(name,"number of x-y points for %s must be at least three!",name);
92 return;
93 }
94 if ( nz < 2 ) {
95 Error(name,"number of z points for %s must be at least two!",name);
96 return;
97 }
98
99 // allocate space for Nxy vertex points
100 fNxy = nxy;
101 fNxyAlloc = nxy;
102 fXvtx = new Float_t [fNxyAlloc];
103 fYvtx = new Float_t [fNxyAlloc];
104 // zero out the vertex points
105 for (Int_t i = 0; i < fNxyAlloc; i++) {
106 fXvtx[i] = 0.;
107 fYvtx[i] = 0.;
108 }
109
110 // allocate space for Nz sections
111 fNz = nz;
112 fNzAlloc = nz;
113 fZ = new Float_t [fNzAlloc];
114 fScale = new Float_t [fNzAlloc];
115 fX0 = new Float_t [fNzAlloc];
116 fY0 = new Float_t [fNzAlloc];
117 // zero out the z points
118 for (Int_t j = 0; j < fNzAlloc; j++) {
119 fZ[j] = 0.;
120 fScale[j] = 0.;
121 fX0[j] = 0.;
122 fY0[j] = 0.;
123 }
124
125}
126
127////////////////////////////////////////////////////////////////////////////////
128/// TXTRU copy constructor
129
131{
132 xtru.TXTRU::Copy(*this);
133}
134
135////////////////////////////////////////////////////////////////////////////////
136/// TXTRU destructor deallocates arrays
137
139{
140 if (fXvtx) delete [] fXvtx;
141 if (fYvtx) delete [] fYvtx;
142 fXvtx = nullptr;
143 fYvtx = nullptr;
144 fNxy = 0;
145 fNxyAlloc = 0;
146
147 if (fZ) delete [] fZ;
148 if (fScale) delete [] fScale;
149 if (fX0) delete [] fX0;
150 if (fY0) delete [] fY0;
151 fZ = nullptr;
152 fScale = nullptr;
153 fX0 = nullptr;
154 fY0 = nullptr;
155 fNz = 0;
156 fNzAlloc = 0;
157
160}
161
162////////////////////////////////////////////////////////////////////////////////
163/// Deep assignment operator
164
166{
167 // protect against self-assignment
168 if (this != &rhs)
169 rhs.TXTRU::Copy(*this);
170
171 return *this;
172}
173
174////////////////////////////////////////////////////////////////////////////////
175/// TXTRU Copy method
176
177void TXTRU::Copy(TObject &obj) const
178{
179 auto &tgt = static_cast<TXTRU &>(obj);
180
181 delete [] tgt.fXvtx; tgt.fXvtx = nullptr;
182 delete [] tgt.fYvtx; tgt.fYvtx = nullptr;
183 delete [] tgt.fZ; tgt.fZ = nullptr;
184 delete [] tgt.fScale; tgt.fScale = nullptr;
185 delete [] tgt.fX0; tgt.fX0 = nullptr;
186 delete [] tgt.fY0; tgt.fY0 = nullptr;
187
188 TObject::Copy(obj);
189 tgt.fNxy = fNxy;
190 tgt.fNxyAlloc = fNxyAlloc;
191 tgt.fXvtx = new Float_t [fNxyAlloc];
192 tgt.fYvtx = new Float_t [fNxyAlloc];
193 for (Int_t i = 0; i < fNxyAlloc; i++) {
194 tgt.fXvtx[i] = fXvtx[i];
195 tgt.fYvtx[i] = fYvtx[i];
196 }
197
198 tgt.fNz = fNz;
199 tgt.fNzAlloc = fNzAlloc;
200 tgt.fZ = new Float_t [fNzAlloc];
201 tgt.fScale = new Float_t [fNzAlloc];
202 tgt.fX0 = new Float_t [fNzAlloc];
203 tgt.fY0 = new Float_t [fNzAlloc];
204 for (Int_t j = 0; j < fNzAlloc; j++) {
205 tgt.fZ[j] = fZ[j];
206 tgt.fScale[j] = fScale[j];
207 tgt.fX0[j] = fX0[j];
208 tgt.fY0[j] = fY0[j];
209 }
210
211 tgt.fPolygonShape = fPolygonShape;
212 tgt.fZOrdering = fZOrdering;
213}
214
215////////////////////////////////////////////////////////////////////////////////
216/// Set z section iz information
217/// expand size of array if necessary
218
220{
221 if (iz < 0) return;
222
223 // setting a new section makes things unverified
225
226 if (iz >= fNzAlloc) {
227 // re-allocate the z positions/scales
228 Int_t newNalloc = iz + 1;
229 Float_t *newZ = new Float_t [newNalloc];
230 Float_t *newS = new Float_t [newNalloc];
231 Float_t *newX = new Float_t [newNalloc];
232 Float_t *newY = new Float_t [newNalloc];
233 Int_t i = 0;
234 for (i = 0; i < newNalloc; i++) {
235 if (i<fNz) {
236 // copy the old points
237 newZ[i] = fZ[i];
238 newS[i] = fScale[i];
239 newX[i] = fX0[i];
240 newY[i] = fY0[i];
241 } else {
242 // zero out the new points
243 newZ[i] = 0;
244 newS[i] = 0;
245 newX[i] = 0;
246 newY[i] = 0;
247 }
248 }
249 delete [] fZ;
250 delete [] fScale;
251 delete [] fX0;
252 delete [] fY0;
253 fZ = newZ;
254 fScale = newS;
255 fX0 = newX;
256 fY0 = newY;
258 }
259
260 // filled z "iz" means indices 0...iz have values -> iz+1 entries
261 fNz = TMath::Max(iz+1,fNz);
262
263 fZ[iz] = z;
264 fScale[iz] = scale;
265 fX0[iz] = x0;
266 fY0[iz] = y0;
267}
268
269////////////////////////////////////////////////////////////////////////////////
270/// Set vertex point ipt to (x,y)
271/// expand size of array if necessary
272
274 if (ipt < 0) return;
275
276 // setting a new vertex makes things unverified
278
279 if (ipt >= fNxyAlloc) {
280 // re-allocate the outline points
281 Int_t newNalloc = ipt + 1;
282 Float_t *newX = new Float_t [newNalloc];
283 Float_t *newY = new Float_t [newNalloc];
284 Int_t i = 0;
285 for (i = 0; i < newNalloc; i++) {
286 if (i<fNxy) {
287 // copy the old points
288 newX[i] = fXvtx[i];
289 newY[i] = fYvtx[i];
290 } else {
291 // zero out the new points
292 newX[i] = 0;
293 newY[i] = 0;
294 }
295 }
296 delete [] fXvtx;
297 delete [] fYvtx;
298 fXvtx = newX;
299 fYvtx = newY;
301 }
302
303 // filled point "ipt" means indices 0...ipt have values -> ipt+1 entries
304 fNxy = TMath::Max(ipt+1,fNxy);
305
306 fXvtx[ipt] = x;
307 fYvtx[ipt] = y;
308}
309
310////////////////////////////////////////////////////////////////////////////////
311/// Compute the distance from point px,py to a TXTRU
312/// by calculating the closest approach to each corner
313
315{
316 Int_t numPoints = fNz*fNxy;
317 return ShapeDistancetoPrimitive(numPoints,px,py);
318}
319
320////////////////////////////////////////////////////////////////////////////////
321/// Return x coordinate of a vertex point
322
324 if ((n < 0) || (n >= fNxy)) {
325 Error(fName,"no such point %d [of %d]",n,fNxy);
326 return 0.0;
327 }
328 return fXvtx[n];
329}
330
331////////////////////////////////////////////////////////////////////////////////
332/// Return y coordinate of a vertex point
333
335 if ((n < 0) || (n >= fNxy)) {
336 Error(fName,"no such point %d [of %d]",n,fNxy);
337 return 0.0;
338 }
339 return fYvtx[n];
340}
341
342////////////////////////////////////////////////////////////////////////////////
343/// Return x0 shift of a z section
344
346 if ((n < 0) || (n >= fNz)) {
347 Error(fName,"no such section %d [of %d]",n,fNz);
348 return 0.0;
349 }
350 return fX0[n];
351}
352
353////////////////////////////////////////////////////////////////////////////////
354/// Return y0 shift of a z section
355
357 if ((n < 0) || (n >= fNz)) {
358 Error(fName,"no such section %d [of %d]",n,fNz);
359 return 0.0;
360 }
361 return fY0[n];
362}
363
364////////////////////////////////////////////////////////////////////////////////
365/// Return scale factor for a z section
366
368 if ((n < 0) || (n >= fNz)) {
369 Error(fName,"no such section %d [of %d]",n,fNz);
370 return 0.0;
371 }
372 return fScale[n];
373}
374
375////////////////////////////////////////////////////////////////////////////////
376/// Return z of a z section
377
379 if ((n < 0) || (n >= fNz)) {
380 Error(fName,"no such section %d [of %d]",n,fNz);
381 return 0.0;
382 }
383 return fZ[n];
384}
385
386////////////////////////////////////////////////////////////////////////////////
387/// Dump the info of this TXTRU shape
388/// Option:
389/// - "xy" to get x-y information
390/// - "z" to get z information
391/// - "alloc" to show full allocated arrays (not just used values)
392
394{
395 TString opt = option;
396 opt.ToLower();
397
398 printf("TXTRU %s Nxy=%d [of %d] Nz=%d [of %d] Option=%s\n",
400
401 const char *shape = nullptr;
402 const char *zorder = nullptr;
403
404 switch (fPolygonShape) {
405 case kUncheckedXY: shape = "Unchecked "; break;
406 case kMalformedXY: shape = "Malformed "; break;
407 case kConvexCCW: shape = "Convex CCW "; break;
408 case kConvexCW: shape = "Convex CW "; break;
409 case kConcaveCCW: shape = "Concave CCW"; break;
410 case kConcaveCW: shape = "Concave CW "; break;
411 }
412
413 switch (fZOrdering) {
414 case kUncheckedZ: zorder = "Unchecked Z"; break;
415 case kMalformedZ: zorder = "Malformed Z"; break;
416 case kConvexIncZ: zorder = "Convex Increasing Z"; break;
417 case kConvexDecZ: zorder = "Convex Decreasing Z"; break;
418 case kConcaveIncZ: zorder = "Concave Increasing Z"; break;
419 case kConcaveDecZ: zorder = "Concave Decreasing Z"; break;
420 }
421
422 printf(" XY shape '%s', '%s'\n",shape,zorder);
423
424 Int_t nxy, nz;
425
426 if (opt.Contains("alloc")) {
427 nxy = fNxy;
428 nz = fNz;
429 } else {
430 nxy = fNxyAlloc;
431 nz = fNzAlloc;
432 }
433
434 const char *name = nullptr;
435 Float_t *p=nullptr;
436 Int_t nlimit = 0;
437 Bool_t print_vtx = opt.Contains("xy");
438 Bool_t print_z = opt.Contains("z");
439
440 Int_t ixyz=0;
441 for (ixyz=0; ixyz<6; ixyz++) {
442 switch (ixyz) {
443 case 0: p = fXvtx; name = "x"; nlimit = nxy; break;
444 case 1: p = fYvtx; name = "y"; nlimit = nxy; break;
445 case 2: p = fZ; name = "z"; nlimit = nz; break;
446 case 3: p = fScale; name = "scale"; nlimit = nz; break;
447 case 4: p = fX0; name = "x0"; nlimit = nz; break;
448 case 5: p = fY0; name = "y0"; nlimit = nz; break;
449 }
450 if (ixyz<=1 && !print_vtx) continue;
451 if (ixyz>=2 && !print_z) continue;
452
453 printf(" Float_t %s[] = \n { %10g",name,*p++);
454 Int_t i=1;
455 for (i=1;i<nlimit;i++) {
456 printf(", %10g",*p++);
457 if (i%6==5) printf("\n ");
458 }
459 printf(" };\n");
460 }
461
462}
463
464////////////////////////////////////////////////////////////////////////////////
465/// Create TXTRU points in buffer
466/// order as expected by other methods (counterclockwise xy, increasing z)
467
469{
470 if (points) {
471 Int_t ipt, ixy, iz, ioff;
472 Float_t x, y;
473
474 // put xy in counterclockwise order
477
478 // put z
481
482 ipt = 0; // point number
483 Int_t i=0;
484 for (i=0; i<fNz; i++) { // loop over sections
485 iz = (reversez) ? fNz-1 - i : i;
486 Int_t j=0;
487 for (j=0; j<fNxy; j++) { // loop over points in section
488 ixy = (iscw) ? fNxy-1 - j : j;
489 ioff = ipt*3; // 3 words per point (x,y,z)
490 x = fXvtx[ixy];
491 y = fYvtx[ixy];
492 points[ioff ] = x*fScale[iz] + fX0[iz];
493 points[ioff+1] = y*fScale[iz] + fY0[iz];
494 points[ioff+2] = fZ[iz];
495 ipt++;
496 }
497 }
498 }
499}
500
501////////////////////////////////////////////////////////////////////////////////
502/// Return total X3D needed by TNode::ls (when called with option "x")
503
504void TXTRU::Sizeof3D() const
505{
506 gSize3D.numPoints += fNz*fNxy;
507 gSize3D.numSegs += (2*fNz-1)*fNxy;
508 gSize3D.numPolys += (fNz-1)*fNxy+2;
509}
510
511////////////////////////////////////////////////////////////////////////////////
512/// (Dis)Enable the splitting of concave polygon outlines into
513/// multiple convex polygons. This would make for better rendering
514/// in solid mode, but introduces extra, potentially confusing, lines
515/// in wireframe mode.
516///
517/// *** Not yet implemented ***
518
520{
521 fSplitConcave = split;
522
523 // Not implemented yet
524 if (split) {
526 std::cout << TNamed::GetName()
527 << " TXTRU::SplitConcavePolygon is not yet implemented" << std::endl;
528 }
529
530}
531
532////////////////////////////////////////////////////////////////////////////////
533/// Truncate the vertex list
534
536 if ((npts < 0) || (npts > fNxy)) {
537 Error(fName,"truncate to %d impossible on %d points",npts,fNxy);
538 return;
539 }
540 fNxy = npts;
541 return;
542}
543
544////////////////////////////////////////////////////////////////////////////////
545/// Truncate the z section list
546
548 if ((nz < 0) || (nz > fNz)) {
549 Error(fName,"truncate to %d impossible on %d points",nz,fNz);
550 return;
551 }
552 fNz = nz;
553 return;
554}
555
556////////////////////////////////////////////////////////////////////////////////
557/// Determine ordering over which to process points, segments, surfaces
558/// so that they render correctly. Generally this has to do
559/// with getting outward normals in the hidden/solid surface case.
560
562{
564
565 // Check on polygon's shape
566 // Convex vs. Concave and ClockWise vs. Counter-ClockWise
567 plus = minus = zero = 0;
568 Int_t ixy=0;
569 for (ixy=0; ixy<fNxy; ixy++) {
570 // calculate the cross product of the two segments that
571 // meet at vertex "ixy"
572 // concave polygons have a mixture of + and - values
573 Int_t ixyprev = (ixy + fNxy - 1)%fNxy;
574 Int_t ixynext = (ixy + fNxy + 1)%fNxy;
575
580
582
583 if (xprod > 0) {
584 plus += xprod;
585 } else if (xprod < 0) {
586 minus -= xprod;
587 } else {
588 zero++;
589 }
590 }
591
592 if (fNxy<3) {
593 // no check yet written for checking that the segments don't cross
595 } else {
596 if (plus==0 || minus==0) {
597 // convex polygons have all of one sign
598 if (plus>minus) {
600 } else {
602 }
603 } else {
604 // concave
605 if (plus>minus) {
607 } else {
609 }
610 }
611 }
612
613 // Check on z ordering
614 // Convex vs. Concave and increasing or decreasing in z
615 plus = minus = zero = 0;
617 Int_t iz=0;
618 for (iz=0; iz<fNz; iz++) {
619 // calculate the cross product of the two segments that
620 // meet at vertex "iz"
621 // concave polygons have a mixture of + and - values
622 Int_t izprev = (iz + fNz - 1)%fNz;
623 Int_t iznext = (iz + fNz + 1)%fNz;
624
625 Float_t dzprev = fZ[iz] - fZ[izprev];
627 Float_t dznext = fZ[iznext] - fZ[iz];
629
630 // special cases for end faces
631 if (iz==0) {
632 dzprev = 0;
633 dsprev = fScale[0];
634 } else if (iz==fNz-1) {
635 dznext = 0;
636 dsnext = -fScale[iz];
637 }
638
640
641 if (xprod > 0) {
642 plus += xprod;
643 } else if (xprod < 0) {
644 minus -= xprod;
645 } else {
646 zero++;
647 }
648 // also check for scale factors that change sign...
649 if (fScale[iz]*fScale[iznext] < 0) scaleSignChange = kTRUE;
650 }
651
652 if (fNz<1 || scaleSignChange) {
653 // no check yet written for checking that the segments don't cross
655 } else {
656 if (plus==0 || minus==0) {
657 // convex polygons have all of one sign
658 if (plus>minus) {
660 } else {
662 }
663 } else {
664 // concave
665 if (plus>minus) {
667 } else {
669 }
670 }
671 }
672}
673
674////////////////////////////////////////////////////////////////////////////////
675/// Dump the vertex points for visual inspection
676
677void TXTRU::DumpPoints(int npoints, float *pointbuff) const
678{
679 std::cout << "TXTRU::DumpPoints - " << npoints << " points" << std::endl;
680 int ioff = 0;
681 float x,y,z;
682 int ipt=0;
683 for (ipt=0; ipt<npoints; ipt++) {
684 x = pointbuff[ioff++];
685 y = pointbuff[ioff++];
686 z = pointbuff[ioff++];
687 printf(" [%4d] %6.1f %6.1f %6.1f \n",ipt,x,y,z);
688 }
689}
690
691////////////////////////////////////////////////////////////////////////////////
692/// Dump the segment info for visual inspection
693
695{
696 std::cout << "TXTRU::DumpSegments - " << nsegments << " segments" << std::endl;
697 int ioff = 0;
698 int icol, p1, p2;
699 int iseg=0;
700 for (iseg=0; iseg<nsegments; iseg++) {
701 icol = segbuff[ioff++];
702 p1 = segbuff[ioff++];
703 p2 = segbuff[ioff++];
704 printf(" [%4d] %3d (%4d,%4d)\n",iseg,icol,p1,p2);
705 }
706}
707
708////////////////////////////////////////////////////////////////////////////////
709/// Dump the derived polygon info for visual inspection
710
712{
713 std::cout << "TXTRU::DumpPolygons - " << npolygons << " polygons" << std::endl;
714 int ioff = 0;
715 int icol, nseg, iseg;
716 int ipoly=0;
717 for (ipoly=0; ipoly<npolygons; ipoly++) {
718 icol = polybuff[ioff++];
719 nseg = polybuff[ioff++];
720#ifndef R__MACOSX
721 std::cout << " [" << std::setw(4) << ipoly << "] icol " << std::setw(3) << icol
722 << " nseg " << std::setw(3) << nseg << " (";
723#else
724 printf(" [%d4] icol %d3 nseg %d3 (", ipoly, icol, nseg);
725#endif
726 for (iseg=0; iseg<nseg-1; iseg++) {
727 std::cout << polybuff[ioff++] << ",";
728 }
729 std::cout << polybuff[ioff++] << ")" << std::endl;
730 }
731 std::cout << " buffer size " << bufsize << " last used " << --ioff << std::endl;
732}
733
734////////////////////////////////////////////////////////////////////////////////
735/// Get buffer 3d.
736
738{
739 static TBuffer3D buffer(TBuffer3DTypes::kGeneric);
740
742
744 // Check that the polygon is well formed
745 // convex vs. concave, z ordered monotonically
746
749 const_cast<TXTRU *>(this)->CheckOrdering();
750 }
752 Int_t nbSegs = fNxy*(2*fNz-1);
753 Int_t nbPols = fNxy*(fNz-1)+2;
754 if (buffer.SetRawSizes(nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, nbPols, 6*(nbPols-2)+2*(2+fNxy))) {
756 }
757 }
759 // Points
760 SetPoints(buffer.fPnts);
761 if (!buffer.fLocalFrame) {
762 TransformPoints(buffer.fPnts, buffer.NbPnts());
763 }
764
766
767 Int_t i,j, k;
769 indx = indx2 = 0;
770
771 // Segments
772 for (i=0; i<fNz; i++) {
773 // loop Z planes
774 indx2 = i*fNxy;
775 // loop polygon segments
776 for (j=0; j<fNxy; j++) {
777 k = (j+1)%fNxy;
778 buffer.fSegs[indx++] = c;
779 buffer.fSegs[indx++] = indx2+j;
780 buffer.fSegs[indx++] = indx2+k;
781 }
782 } // total: fNz*fNxy polygon segments
783 for (i=0; i<fNz-1; i++) {
784 // loop Z planes
785 indx2 = i*fNxy;
786 // loop polygon segments
787 for (j=0; j<fNxy; j++) {
788 k = j + fNxy;
789 buffer.fSegs[indx++] = c;
790 buffer.fSegs[indx++] = indx2+j;
791 buffer.fSegs[indx++] = indx2+k;
792 }
793 } // total (fNz-1)*fNxy lateral segments
794
795 // Polygons
796 indx = 0;
797
798 // fill lateral polygons
799 for (i=0; i<fNz-1; i++) {
800 indx2 = i*fNxy;
801 for (j=0; j<fNxy; j++) {
802 k = (j+1)%fNxy;
803 buffer.fPols[indx++] = c+j%3;
804 buffer.fPols[indx++] = 4;
805 buffer.fPols[indx++] = indx2+j;
806 buffer.fPols[indx++] = fNz*fNxy+indx2+k;
807 buffer.fPols[indx++] = indx2+fNxy+j;
808 buffer.fPols[indx++] = fNz*fNxy+indx2+j;
809 }
810 } // total (fNz-1)*fNxy polys
811 buffer.fPols[indx++] = c+2;
812 buffer.fPols[indx++] = fNxy;
813 indx2 = 0;
814 for (j = fNxy - 1; j >= 0; --j) {
815 buffer.fPols[indx++] = indx2+j;
816 }
817
818 buffer.fPols[indx++] = c;
819 buffer.fPols[indx++] = fNxy;
820 indx2 = (fNz-1)*fNxy;
821
822 for (j=0; j<fNxy; j++) {
823 buffer.fPols[indx++] = indx2+j;
824 }
825
827 }
828
829 return buffer;
830}
#define c(i)
Definition RSha256.hxx:101
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.
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t points
char name[80]
Definition TGX11.cxx:110
#define gSize3D
Definition X3DBuffer.h:40
Generic 3D primitive description class.
Definition TBuffer3D.h:18
Int_t * fPols
Definition TBuffer3D.h:115
UInt_t NbPnts() const
Definition TBuffer3D.h:80
void SetSectionsValid(UInt_t mask)
Definition TBuffer3D.h:65
Int_t * fSegs
Definition TBuffer3D.h:114
Bool_t fLocalFrame
Definition TBuffer3D.h:90
Bool_t SetRawSizes(UInt_t reqPnts, UInt_t reqPntsCapacity, UInt_t reqSegs, UInt_t reqSegsCapacity, UInt_t reqPols, UInt_t reqPolsCapacity)
Set kRaw tessellation section of buffer with supplied sizes.
Double_t * fPnts
Definition TBuffer3D.h:113
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
TString fName
Definition TNamed.h:32
Mother of all ROOT objects.
Definition TObject.h:41
virtual void Copy(TObject &object) const
Copy this to obj.
Definition TObject.cxx:159
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1071
This is the base class for all geometry shapes.
Definition TShape.h:35
Int_t GetBasicColor() const
Get basic color.
Definition TShape.cxx:240
Int_t ShapeDistancetoPrimitive(Int_t numPoints, Int_t px, Int_t py)
Distance to primitive.
Definition TShape.cxx:116
virtual void FillBuffer3D(TBuffer3D &buffer, Int_t reqSections) const
We have to set kRawSize (unless already done) to allocate buffer space before kRaw can be filled.
Definition TShape.cxx:210
void TransformPoints(Double_t *points, UInt_t NbPnts) const
Transform points (LocalToMaster)
Definition TShape.cxx:189
Basic string class.
Definition TString.h:138
void ToLower()
Change string to lower-case.
Definition TString.cxx:1189
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:640
A poly-extrusion.
Definition TXTRU.h:22
Float_t * fScale
Definition TXTRU.h:68
Float_t * fXvtx
Definition TXTRU.h:65
void CheckOrdering()
Determine ordering over which to process points, segments, surfaces so that they render correctly.
Definition TXTRU.cxx:561
virtual void DefineVertex(Int_t pointNum, Float_t x, Float_t y)
Set vertex point ipt to (x,y) expand size of array if necessary.
Definition TXTRU.cxx:273
void SplitConcavePolygon(Bool_t split=kTRUE)
(Dis)Enable the splitting of concave polygon outlines into multiple convex polygons.
Definition TXTRU.cxx:519
virtual Float_t GetSectionX0(Int_t secNum) const
Return x0 shift of a z section.
Definition TXTRU.cxx:345
EZChecked fZOrdering
Definition TXTRU.h:80
void DumpSegments(int nsegments, int *segbuff) const
Dump the segment info for visual inspection.
Definition TXTRU.cxx:694
TXTRU & operator=(const TXTRU &rhs)
Deep assignment operator.
Definition TXTRU.cxx:165
Float_t * fY0
Definition TXTRU.h:70
@ kConvexDecZ
Definition TXTRU.h:76
@ kConcaveDecZ
Definition TXTRU.h:77
@ kMalformedZ
Definition TXTRU.h:75
@ kConvexIncZ
Definition TXTRU.h:76
@ kUncheckedZ
Definition TXTRU.h:75
@ kConcaveIncZ
Definition TXTRU.h:77
Int_t DistancetoPrimitive(Int_t px, Int_t py) override
Compute the distance from point px,py to a TXTRU by calculating the closest approach to each corner.
Definition TXTRU.cxx:314
Float_t * fYvtx
Definition TXTRU.h:66
virtual Float_t GetOutlinePointX(Int_t pointNum) const
Return x coordinate of a vertex point.
Definition TXTRU.cxx:323
Bool_t fSplitConcave
Definition TXTRU.h:85
virtual Float_t GetSectionScale(Int_t secNum) const
Return scale factor for a z section.
Definition TXTRU.cxx:367
Int_t fNxyAlloc
Definition TXTRU.h:62
virtual Float_t GetSectionY0(Int_t secNum) const
Return y0 shift of a z section.
Definition TXTRU.cxx:356
EXYChecked fPolygonShape
Definition TXTRU.h:79
void Copy(TObject &xtru) const override
TXTRU Copy method.
Definition TXTRU.cxx:177
void DumpPoints(int npoints, float *pointbuff) const
Dump the vertex points for visual inspection.
Definition TXTRU.cxx:677
void DumpPolygons(int npolygons, int *polybuff, int bufsize) const
Dump the derived polygon info for visual inspection.
Definition TXTRU.cxx:711
virtual void TruncateNxy(Int_t npts)
Truncate the vertex list.
Definition TXTRU.cxx:535
void Print(Option_t *option="") const override
Dump the info of this TXTRU shape Option:
Definition TXTRU.cxx:393
Int_t fNzAlloc
Definition TXTRU.h:64
virtual Float_t GetSectionZ(Int_t secNum) const
Return z of a z section.
Definition TXTRU.cxx:378
const TBuffer3D & GetBuffer3D(Int_t) const override
Get buffer 3d.
Definition TXTRU.cxx:737
virtual void DefineSection(Int_t secNum, Float_t z, Float_t scale=1., Float_t x0=0., Float_t y0=0.)
Set z section iz information expand size of array if necessary.
Definition TXTRU.cxx:219
Int_t fNxy
Definition TXTRU.h:61
@ kConcaveCW
Definition TXTRU.h:74
@ kMalformedXY
Definition TXTRU.h:72
@ kConvexCCW
Definition TXTRU.h:73
@ kConcaveCCW
Definition TXTRU.h:74
@ kConvexCW
Definition TXTRU.h:73
@ kUncheckedXY
Definition TXTRU.h:72
Float_t * fX0
Definition TXTRU.h:69
void SetPoints(Double_t *points) const override
Create TXTRU points in buffer order as expected by other methods (counterclockwise xy,...
Definition TXTRU.cxx:468
Int_t fNz
Definition TXTRU.h:63
Float_t * fZ
Definition TXTRU.h:67
void Sizeof3D() const override
Return total X3D needed by TNode::ls (when called with option "x")
Definition TXTRU.cxx:504
TXTRU()
TXTRU shape - default constructor.
Definition TXTRU.cxx:75
~TXTRU() override
TXTRU destructor deallocates arrays.
Definition TXTRU.cxx:138
virtual void TruncateNz(Int_t npts)
Truncate the z section list.
Definition TXTRU.cxx:547
virtual Float_t GetOutlinePointY(Int_t pointNum) const
Return y coordinate of a vertex point.
Definition TXTRU.cxx:334
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:251