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