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