Logo ROOT   6.08/07
Reference Guide
TPolyLine3D.cxx
Go to the documentation of this file.
1 // @(#)root/g3d:$Id$
2 // Author: Nenad Buncic 17/08/95
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 "Riostream.h"
13 #include "TROOT.h"
14 #include "TPolyLine3D.h"
15 #include "TVirtualPad.h"
16 #include "TView.h"
17 #include "TVirtualViewer3D.h"
18 #include "TBuffer3D.h"
19 #include "TBuffer3DTypes.h"
20 #include "TGeometry.h"
21 #include "TMath.h"
22 
23 #include <assert.h>
24 
26 
27 /** \class TPolyLine3D
28 \ingroup g3d
29 A 3-dimensional polyline. It has 4 different constructors.
30 
31 First one, without any parameters TPolyLine3D(), we call 'default
32 constructor' and it's used in a case that just an initialisation is
33 needed (i.e. pointer declaration).
34 
35 Example:
36 
37 ~~~ {.cpp}
38  TPolyLine3D *pl1 = new TPolyLine3D;
39 ~~~
40 
41 Second one is 'normal constructor' with, usually, one parameter
42 n (number of points), and it just allocates a space for the points.
43 
44 Example:
45 
46 ~~~ {.cpp}
47  TPolyLine3D pl1(150);
48 ~~~
49 
50 Third one allocates a space for the points, and also makes
51 initialisation from the given array.
52 
53 Example:
54 
55 ~~~ {.cpp}
56  TPolyLine3D pl1(150, pointerToAnArray);
57 ~~~
58 
59 Fourth one is, almost, similar to the constructor above, except
60 initialisation is provided with three independent arrays (array of
61 x coordinates, y coordinates and z coordinates).
62 
63 Example:
64 
65 ~~~ {.cpp}
66  TPolyLine3D pl1(150, xArray, yArray, zArray);
67 ~~~
68 
69 Example:
70 
71 Begin_Macro(source)
72 {
73  TCanvas *c1 = new TCanvas("c1","c1",500,500);
74  TView *view = TView::CreateView(1);
75  view->SetRange(0,0,0,2,2,2);
76  const Int_t n = 500;
77  r = new TRandom();
78  Double_t x, y, z;
79  TPolyLine3D *l = new TPolyLine3D(n);
80  for (Int_t i=0;i<n;i++) {
81  r->Sphere(x, y, z, 1);
82  l->SetPoint(i,x+1,y+1,z+1);
83  }
84  l->Draw();
85 }
86 End_Macro
87 */
88 
89 ////////////////////////////////////////////////////////////////////////////////
90 /// 3-D polyline default constructor.
91 
93 {
94  fN = 0;
95  fP = 0;
96  fLastPoint = -1;
97 }
98 
99 ////////////////////////////////////////////////////////////////////////////////
100 /// 3-D polyline normal constructor with initialization to 0.
101 /// If n < 0 the default size (2 points) is set.
102 
104 {
105  fOption = option;
107  fLastPoint = -1;
108  if (n <= 0) {
109  fN = 0;
110  fP = 0;
111  return;
112  }
113 
114  fN = n;
115  fP = new Float_t[3*fN];
116  for (Int_t i=0; i<3*fN; i++) fP[i] = 0;
117 }
118 
119 ////////////////////////////////////////////////////////////////////////////////
120 /// 3-D polyline normal constructor. Polyline is intialized with p.
121 /// If n < 0 the default size (2 points) is set.
122 
124 {
125  fOption = option;
127  fLastPoint = -1;
128  if (n <= 0) {
129  fN = 0;
130  fP = 0;
131  return;
132  }
133 
134  fN = n;
135  fP = new Float_t[3*fN];
136  for (Int_t i=0; i<3*n; i++) {
137  fP[i] = p[i];
138  }
139  fLastPoint = fN-1;
140 }
141 
142 ////////////////////////////////////////////////////////////////////////////////
143 /// 3-D polyline normal constructor. Polyline is initialized with p
144 /// (cast to float). If n < 0 the default size (2 points) is set.
145 
147 {
148  fOption = option;
150  fLastPoint = -1;
151  if (n <= 0) {
152  fN = 0;
153  fP = 0;
154  return;
155  }
156 
157  fN = n;
158  fP = new Float_t[3*fN];
159  for (Int_t i=0; i<3*n; i++) {
160  fP[i] = (Float_t) p[i];
161  }
162  fLastPoint = fN-1;
163 }
164 
165 ////////////////////////////////////////////////////////////////////////////////
166 /// 3-D polyline normal constructor. Polyline is initialized withe the
167 /// x, y ,z arrays. If n < 0 the default size (2 points) is set.
168 
170 {
171  fOption = option;
173  fLastPoint = -1;
174  if (n <= 0) {
175  fN = 0;
176  fP = 0;
177  return;
178  }
179 
180  fN = n;
181  fP = new Float_t[3*fN];
182  Int_t j = 0;
183  for (Int_t i=0; i<n;i++) {
184  fP[j] = x[i];
185  fP[j+1] = y[i];
186  fP[j+2] = z[i];
187  j += 3;
188  }
189  fLastPoint = fN-1;
190 }
191 
192 ////////////////////////////////////////////////////////////////////////////////
193 /// 3-D polyline normal constructor. Polyline is initialized withe the
194 /// x, y, z arrays (which are cast to float).
195 /// If n < 0 the default size (2 points) is set.
196 
198 {
199  fOption = option;
201  fLastPoint = -1;
202  if (n <= 0) {
203  fN = 0;
204  fP = 0;
205  return;
206  }
207 
208  fN = n;
209  fP = new Float_t[3*fN];
210  Int_t j = 0;
211  for (Int_t i=0; i<n;i++) {
212  fP[j] = (Float_t) x[i];
213  fP[j+1] = (Float_t) y[i];
214  fP[j+2] = (Float_t) z[i];
215  j += 3;
216  }
217  fLastPoint = fN-1;
218 }
219 
220 ////////////////////////////////////////////////////////////////////////////////
221 /// assignment operator
222 
224 {
225  if(this!=&pl) {
226  TObject::operator=(pl);
227  TAttLine::operator=(pl);
228  TAtt3D::operator=(pl);
229  fN=pl.fN;
230  fP=pl.fP;
231  fOption=pl.fOption;
233  }
234  return *this;
235 }
236 
237 ////////////////////////////////////////////////////////////////////////////////
238 /// 3-D polyline destructor.
239 
241 {
242  if (fP) delete [] fP;
243 }
244 
245 ////////////////////////////////////////////////////////////////////////////////
246 /// 3-D polyline copy ctor.
247 
248 TPolyLine3D::TPolyLine3D(const TPolyLine3D &polyline) : TObject(polyline), TAttLine(polyline), TAtt3D(polyline)
249 {
250  fP = 0;
251  fLastPoint = 0;
252  fN = 0;
253  ((TPolyLine3D&)polyline).TPolyLine3D::Copy(*this);
254 }
255 
256 ////////////////////////////////////////////////////////////////////////////////
257 /// Copy polyline to polyline obj.
258 
259 void TPolyLine3D::Copy(TObject &obj) const
260 {
261  TObject::Copy(obj);
262  TAttLine::Copy(((TPolyLine3D&)obj));
263  ((TPolyLine3D&)obj).fN = fN;
264  if (((TPolyLine3D&)obj).fP)
265  delete [] ((TPolyLine3D&)obj).fP;
266  if (fN > 0) {
267  ((TPolyLine3D&)obj).fP = new Float_t[3*fN];
268  for (Int_t i=0; i<3*fN;i++) {((TPolyLine3D&)obj).fP[i] = fP[i];}
269  } else {
270  ((TPolyLine3D&)obj).fP = 0;
271  }
272  ((TPolyLine3D&)obj).fOption = fOption;
273  ((TPolyLine3D&)obj).fLastPoint = fLastPoint;
274 }
275 
276 ////////////////////////////////////////////////////////////////////////////////
277 /// Compute distance from point px,py to a 3-D polyline.
278 /// Compute the closest distance of approach from point px,py to each segment
279 /// of the polyline.
280 /// Returns when the distance found is below DistanceMaximum.
281 /// The distance is computed in pixels units.
282 
284 {
285  const Int_t inaxis = 7;
286  Int_t dist = 9999;
287 
288  Int_t puxmin = gPad->XtoAbsPixel(gPad->GetUxmin());
289  Int_t puymin = gPad->YtoAbsPixel(gPad->GetUymin());
290  Int_t puxmax = gPad->XtoAbsPixel(gPad->GetUxmax());
291  Int_t puymax = gPad->YtoAbsPixel(gPad->GetUymax());
292 
293  // return if point is not in the user area
294  if (px < puxmin - inaxis) return dist;
295  if (py > puymin + inaxis) return dist;
296  if (px > puxmax + inaxis) return dist;
297  if (py < puymax - inaxis) return dist;
298 
299  TView *view = gPad->GetView();
300  if (!view) return dist;
301 
302  Int_t i, dsegment;
303  Double_t x1,y1,x2,y2;
304  Float_t xndc[3];
305  for (i=0;i<Size()-1;i++) {
306  view->WCtoNDC(&fP[3*i], xndc);
307  x1 = xndc[0];
308  y1 = xndc[1];
309  view->WCtoNDC(&fP[3*i+3], xndc);
310  x2 = xndc[0];
311  y2 = xndc[1];
312  dsegment = DistancetoLine(px,py,x1,y1,x2,y2);
313  if (dsegment < dist) dist = dsegment;
314  }
315  return dist;
316 }
317 
318 ////////////////////////////////////////////////////////////////////////////////
319 /// Draw this 3-D polyline with its current attributes.
320 
322 {
323  AppendPad(option);
324 }
325 
326 ////////////////////////////////////////////////////////////////////////////////
327 /// Draw cube outline with 3d polylines.
328 ///
329 /// ~~~ {.cpp}
330 /// xmin = fRmin[0] xmax = fRmax[0]
331 /// ymin = fRmin[1] ymax = fRmax[1]
332 /// zmin = fRmin[2] zmax = fRmax[2]
333 ///
334 /// (xmin,ymax,zmax) +---------+ (xmax,ymax,zmax)
335 /// / /|
336 /// / / |
337 /// / / |
338 ///(xmin,ymin,zmax) +---------+ |
339 /// | | + (xmax,ymax,zmin)
340 /// | | /
341 /// | | /
342 /// | |/
343 /// +---------+
344 /// (xmin,ymin,zmin) (xmax,ymin,zmin)
345 /// ~~~
346 
348 {
349  Double_t xmin = rmin[0]; Double_t xmax = rmax[0];
350  Double_t ymin = rmin[1]; Double_t ymax = rmax[1];
351  Double_t zmin = rmin[2]; Double_t zmax = rmax[2];
352 
353  TPolyLine3D *pl3d = (TPolyLine3D *)outline->First();
354  if (!pl3d) {
355  TView *view = gPad->GetView();
356  if (!view) return;
357  TPolyLine3D *p1 = new TPolyLine3D(4);
358  TPolyLine3D *p2 = new TPolyLine3D(4);
359  TPolyLine3D *p3 = new TPolyLine3D(4);
360  TPolyLine3D *p4 = new TPolyLine3D(4);
361  p1->SetLineColor(view->GetLineColor());
362  p1->SetLineStyle(view->GetLineStyle());
363  p1->SetLineWidth(view->GetLineWidth());
364  p1->Copy(*p2);
365  p1->Copy(*p3);
366  p1->Copy(*p4);
367  outline->Add(p1);
368  outline->Add(p2);
369  outline->Add(p3);
370  outline->Add(p4);
371  }
372 
373  pl3d = (TPolyLine3D *)outline->First();
374 
375  if (pl3d) {
376  pl3d->SetPoint(0, xmin, ymin, zmin);
377  pl3d->SetPoint(1, xmax, ymin, zmin);
378  pl3d->SetPoint(2, xmax, ymax, zmin);
379  pl3d->SetPoint(3, xmin, ymax, zmin);
380  }
381 
382  pl3d = (TPolyLine3D *)outline->After(pl3d);
383 
384  if (pl3d) {
385  pl3d->SetPoint(0, xmax, ymin, zmin);
386  pl3d->SetPoint(1, xmax, ymin, zmax);
387  pl3d->SetPoint(2, xmax, ymax, zmax);
388  pl3d->SetPoint(3, xmax, ymax, zmin);
389  }
390 
391  pl3d = (TPolyLine3D *)outline->After(pl3d);
392 
393  if (pl3d) {
394  pl3d->SetPoint(0, xmax, ymin, zmax);
395  pl3d->SetPoint(1, xmin, ymin, zmax);
396  pl3d->SetPoint(2, xmin, ymax, zmax);
397  pl3d->SetPoint(3, xmax, ymax, zmax);
398  }
399 
400  pl3d = (TPolyLine3D *)outline->After(pl3d);
401 
402  if (pl3d) {
403  pl3d->SetPoint(0, xmin, ymin, zmax);
404  pl3d->SetPoint(1, xmin, ymin, zmin);
405  pl3d->SetPoint(2, xmin, ymax, zmin);
406  pl3d->SetPoint(3, xmin, ymax, zmax);
407  }
408 }
409 
410 ////////////////////////////////////////////////////////////////////////////////
411 /// Draw 3-D polyline with new coordinates. Creates a new polyline which
412 /// will be adopted by the pad in which it is drawn. Does not change the
413 /// original polyline (should be static method).
414 
416 {
417  TPolyLine3D *newpolyline = new TPolyLine3D();
418  Int_t size = 3*Size();
419  newpolyline->fN =n;
420  newpolyline->fP = new Float_t[size];
421  for (Int_t i=0; i<size;i++) { newpolyline->fP[i] = p[i];}
422  TAttLine::Copy(*newpolyline);
423  newpolyline->fOption = fOption;
424  newpolyline->fLastPoint = fLastPoint;
425  newpolyline->SetBit(kCanDelete);
426  newpolyline->AppendPad(option);
427 }
428 
429 ////////////////////////////////////////////////////////////////////////////////
430 /// Execute action corresponding to one event.
431 
433 {
434  if (!gPad) return;
435  if (gPad->GetView()) gPad->GetView()->ExecuteRotateView(event, px, py);
436 }
437 
438 ////////////////////////////////////////////////////////////////////////////////
439 /// List this 3-D polyline.
440 
441 void TPolyLine3D::ls(Option_t *option) const
442 {
444  std::cout <<"PolyLine3D N=" <<fN<<" Option="<<option<<std::endl;
445 }
446 
447 ////////////////////////////////////////////////////////////////////////////////
448 /// Merge polylines in the collection in this polyline
449 
451 {
452  if (!li) return 0;
453  TIter next(li);
454 
455  //first loop to count the number of entries
456  TPolyLine3D *pl;
457  Int_t npoints = 0;
458  while ((pl = (TPolyLine3D*)next())) {
459  if (!pl->InheritsFrom(TPolyLine3D::Class())) {
460  Error("Add","Attempt to add object of class: %s to a %s",pl->ClassName(),this->ClassName());
461  return -1;
462  }
463  npoints += pl->Size();
464  }
465 
466  //extend this polyline to hold npoints
467  SetPoint(npoints-1,0,0,0);
468 
469  //merge all polylines
470  next.Reset();
471  while ((pl = (TPolyLine3D*)next())) {
472  Int_t np = pl->Size();
473  Float_t *p = pl->GetP();
474  for (Int_t i=0;i<np;i++) {
475  SetPoint(i,p[3*i],p[3*i+1],p[3*i+2]);
476  }
477  }
478 
479  return npoints;
480 }
481 
482 ////////////////////////////////////////////////////////////////////////////////
483 /// Paint a TPolyLine3D.
484 
485 void TPolyLine3D::Paint(Option_t * /* option */ )
486 {
487  UInt_t i;
488 
489  // No need to continue if there is nothing to paint
490  if (Size() <= 0) return;
491 
492  static TBuffer3D buffer(TBuffer3DTypes::kLine);
493 
494  // TPolyLine3D can only be described by filling the TBuffer3D 'tesselation'
495  // parts - so there are no 'optional' sections - we just fill everything.
496 
497  buffer.ClearSectionsValid();
498 
499  // Section kCore
500  buffer.fID = this;
501  buffer.fColor = GetLineColor();
502  buffer.fTransparency = 0;
503  buffer.fLocalFrame = kFALSE;
505 
506  // We fill kCore and kRawSizes on first pass and try with viewer
507  TVirtualViewer3D * viewer3D = gPad->GetViewer3D();
508  if (!viewer3D) return;
509  Int_t reqSections = viewer3D->AddObject(buffer);
510  if (reqSections == TBuffer3D::kNone) {
511  return;
512  }
513 
514  if (reqSections & TBuffer3D::kRawSizes) {
515  Int_t nbPnts = Size();
516  Int_t nbSegs = nbPnts-1;
517  if (!buffer.SetRawSizes(nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, 0, 0)) {
518  return;
519  }
520  buffer.SetSectionsValid(TBuffer3D::kRawSizes);
521  }
522 
523  if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
524  // Points
525  for (i=0; i<3*buffer.NbPnts(); i++) {
526  buffer.fPnts[i] = (Double_t)fP[i];
527  }
528 
529  // Transform points
530  if (gGeometry && !buffer.fLocalFrame) {
531  Double_t dlocal[3];
532  Double_t dmaster[3];
533  for (UInt_t j=0; j<buffer.NbPnts(); j++) {
534  dlocal[0] = buffer.fPnts[3*j];
535  dlocal[1] = buffer.fPnts[3*j+1];
536  dlocal[2] = buffer.fPnts[3*j+2];
537  gGeometry->Local2Master(&dlocal[0],&dmaster[0]);
538  buffer.fPnts[3*j] = dmaster[0];
539  buffer.fPnts[3*j+1] = dmaster[1];
540  buffer.fPnts[3*j+2] = dmaster[2];
541  }
542  }
543 
544  // Basic colors: 0, 1, ... 8
545  Int_t c = (((GetLineColor()) %8) -1) * 4;
546  if (c < 0) c = 0;
547 
548  // Segments
549  for (i = 0; i < buffer.NbSegs(); i++) {
550  buffer.fSegs[3*i ] = c;
551  buffer.fSegs[3*i+1] = i;
552  buffer.fSegs[3*i+2] = i+1;
553  }
554 
556 
557  buffer.SetSectionsValid(TBuffer3D::kRaw);
558  }
559 
560  viewer3D->AddObject(buffer);
561 }
562 
563 ////////////////////////////////////////////////////////////////////////////////
564 /// Dump this 3-D polyline with its attributes on stdout.
565 
566 void TPolyLine3D::Print(Option_t *option) const
567 {
568  printf(" TPolyLine3D N=%d, Option=%s\n",fN,option);
569  TString opt = option;
570  opt.ToLower();
571  if (opt.Contains("all")) {
572  for (Int_t i=0;i<Size();i++) {
573  printf(" x[%d]=%g, y[%d]=%g, z[%d]=%g\n",i,fP[3*i],i,fP[3*i+1],i,fP[3*i+2]);
574  }
575  }
576 }
577 
578 ////////////////////////////////////////////////////////////////////////////////
579 /// Save primitive as a C++ statement(s) on output stream.
580 
581 void TPolyLine3D::SavePrimitive(std::ostream &out, Option_t * /*= ""*/)
582 {
583  char quote = '"';
584  out<<" "<<std::endl;
585  if (gROOT->ClassSaved(TPolyLine3D::Class())) {
586  out<<" ";
587  } else {
588  out<<" TPolyLine3D *";
589  }
590  Int_t size=Size();
591  out<<"pline3D = new TPolyLine3D("<<fN<<","<<quote<<fOption<<quote<<");"<<std::endl;
592 
593  SaveLineAttributes(out,"pline3D",1,1,1);
594 
595  if (size > 0) {
596  for (Int_t i=0;i<size;i++)
597  out<<" pline3D->SetPoint("<<i<<","<<fP[3*i]<<","<<fP[3*i+1]<<","<<fP[3*i+2]<<");"<<std::endl;
598  }
599  out<<" pline3D->Draw();"<<std::endl;
600 }
601 
602 ////////////////////////////////////////////////////////////////////////////////
603 /// Set point following LastPoint to x, y, z.
604 /// Returns index of the point (new last point).
605 
607 {
608  fLastPoint++;
609  SetPoint(fLastPoint, x, y, z);
610  return fLastPoint;
611 }
612 
613 ////////////////////////////////////////////////////////////////////////////////
614 /// Set point n to x, y, z.
615 /// If n is more then the current TPolyLine3D size (n > fN) then
616 /// the polyline will be resized to contain at least n points.
617 
619 {
620  if (n < 0) return;
621  if (!fP || n >= fN) {
622  // re-allocate the object
623  Int_t newN = TMath::Max(2*fN,n+1);
624  Float_t *savepoint = new Float_t [3*newN];
625  if (fP && fN){
626  memcpy(savepoint,fP,3*fN*sizeof(Float_t));
627  memset(&savepoint[3*fN],0,(newN-fN)*sizeof(Float_t));
628  delete [] fP;
629  }
630  fP = savepoint;
631  fN = newN;
632  }
633  fP[3*n ] = x;
634  fP[3*n+1] = y;
635  fP[3*n+2] = z;
637 }
638 
639 ////////////////////////////////////////////////////////////////////////////////
640 /// Re-initialize polyline with n points (0,0,0).
641 /// if n <= 0 the current array of points is deleted.
642 
644 {
645  fOption = option;
646  if (n <= 0) {
647  fN = 0;
648  fLastPoint = -1;
649  delete [] fP;
650  fP = 0;
651  return;
652  }
653  fN = n;
654  if (fP) delete [] fP;
655  fP = new Float_t[3*fN];
656  memset(fP,0,3*fN*sizeof(Float_t));
657  fLastPoint = fN-1;
658 }
659 
660 ////////////////////////////////////////////////////////////////////////////////
661 /// Re-initialize polyline with n points from p. If p=0 initialize with 0.
662 /// if n <= 0 the current array of points is deleted.
663 
665 {
666  fOption = option;
667  if (n <= 0) {
668  fN = 0;
669  fLastPoint = -1;
670  delete [] fP;
671  fP = 0;
672  return;
673  }
674  fN = n;
675  if (fP) delete [] fP;
676  fP = new Float_t[3*fN];
677  if (p) {
678  for (Int_t i=0; i<fN;i++) {
679  fP[3*i] = p[3*i];
680  fP[3*i+1] = p[3*i+1];
681  fP[3*i+2] = p[3*i+2];
682  }
683  } else {
684  memset(fP,0,3*fN*sizeof(Float_t));
685  }
686  fLastPoint = fN-1;
687 }
688 
689 ////////////////////////////////////////////////////////////////////////////////
690 /// Re-initialize polyline with n points from p. If p=0 initialize with 0.
691 /// if n <= 0 the current array of points is deleted.
692 
694 {
695  fOption = option;
696  if (n <= 0) {
697  fN = 0;
698  fLastPoint = -1;
699  delete [] fP;
700  fP = 0;
701  return;
702  }
703  fN = n;
704  if (fP) delete [] fP;
705  fP = new Float_t[3*fN];
706  if (p) {
707  for (Int_t i=0; i<fN;i++) {
708  fP[3*i] = (Float_t) p[3*i];
709  fP[3*i+1] = (Float_t) p[3*i+1];
710  fP[3*i+2] = (Float_t) p[3*i+2];
711  }
712  } else {
713  memset(fP,0,3*fN*sizeof(Float_t));
714  }
715  fLastPoint = fN-1;
716 }
717 
718 ////////////////////////////////////////////////////////////////////////////////
719 /// Stream a 3-D polyline object.
720 
721 void TPolyLine3D::Streamer(TBuffer &b)
722 {
723  UInt_t R__s, R__c;
724  if (b.IsReading()) {
725  b.ReadVersion(&R__s, &R__c);
726  TObject::Streamer(b);
727  TAttLine::Streamer(b);
728  b >> fN;
729  if (fN) {
730  fP = new Float_t[3*fN];
731  b.ReadFastArray(fP,3*fN);
732  }
733  fOption.Streamer(b);
734  fLastPoint = fN-1;
735  b.CheckByteCount(R__s, R__c, TPolyLine3D::IsA());
736  } else {
737  R__c = b.WriteVersion(TPolyLine3D::IsA(), kTRUE);
738  TObject::Streamer(b);
739  TAttLine::Streamer(b);
740  Int_t size = Size();
741  b << size;
742  if (size) b.WriteFastArray(fP, 3*size);
743  fOption.Streamer(b);
744  b.SetByteCount(R__c, kTRUE);
745  }
746 }
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition: TAttLine.h:49
Bool_t IsReading() const
Definition: TBuffer.h:83
double dist(Rotation3D const &r1, Rotation3D const &r2)
Definition: 3DDistances.cxx:48
float xmin
Definition: THbookFile.cxx:93
TPolyLine3D()
3-D polyline default constructor.
Definition: TPolyLine3D.cxx:92
static double p3(double t, double a, double b, double c, double d)
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute action corresponding to one event.
Float_t * fP
[3*fN] Array of 3-D coordinates (x,y,z)
Definition: TPolyLine3D.h:45
virtual void WCtoNDC(const Float_t *pw, Float_t *pn)=0
float Float_t
Definition: RtypesCore.h:53
virtual void SetPolyLine(Int_t n, Option_t *option="")
Re-initialize polyline with n points (0,0,0).
return c
const char Option_t
Definition: RtypesCore.h:62
float ymin
Definition: THbookFile.cxx:93
See TView3D.
Definition: TView.h:29
Use this attribute class when an object should have 3D capabilities.
Definition: TAtt3D.h:21
Buffer base class used for serializing objects.
Definition: TBuffer.h:42
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)=0
#define gROOT
Definition: TROOT.h:364
UInt_t NbPnts() const
Definition: TBuffer3D.h:82
void ClearSectionsValid()
Clear any sections marked valid.
Definition: TBuffer3D.cxx:286
Basic string class.
Definition: TString.h:137
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1089
int Int_t
Definition: RtypesCore.h:41
A 3-dimensional polyline.
Definition: TPolyLine3D.h:41
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual Int_t Merge(TCollection *list)
Merge polylines in the collection in this polyline.
virtual UInt_t WriteVersion(const TClass *cl, Bool_t useBcnt=kFALSE)=0
virtual void Modify()
Change current line attributes if necessary.
Definition: TAttLine.cxx:232
const char * Class
Definition: TXMLSetup.cxx:64
virtual Width_t GetLineWidth() const
Return the line width.
Definition: TAttLine.h:41
void Reset()
Definition: TCollection.h:161
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:739
if object in a list can be deleted
Definition: TObject.h:56
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition: TObject.cxx:165
virtual Style_t GetLineStyle() const
Return the line style.
Definition: TAttLine.h:40
static const double x2[5]
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition: TObject.cxx:188
Double_t x[n]
Definition: legend1.C:17
virtual void ls(Option_t *option="") const
List this 3-D polyline.
void Copy(TAttLine &attline) const
Copy this line attributes to a new TAttLine.
Definition: TAttLine.cxx:162
virtual void Copy(TObject &object) const
Copy this to obj.
Definition: TObject.cxx:123
virtual void SaveLineAttributes(std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1, Int_t widdef=1)
Save line attributes as C++ statement(s) on output stream out.
Definition: TAttLine.cxx:260
Float_t * GetP() const
Definition: TPolyLine3D.h:67
virtual ~TPolyLine3D()
3-D polyline destructor.
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute distance from point px,py to a 3-D polyline.
UInt_t NbSegs() const
Definition: TBuffer3D.h:83
static double p2(double t, double a, double b, double c)
virtual void Local2Master(Double_t *local, Double_t *master)
Convert one point from local system to master reference system.
Definition: TGeometry.cxx:407
Double_t * fPnts
Definition: TBuffer3D.h:114
Abstract 3D shapes viewer.
TObject & operator=(const TObject &rhs)
TObject assignment operator.
Definition: TObject.cxx:103
void SetSectionsValid(UInt_t mask)
Definition: TBuffer3D.h:67
A doubly linked list.
Definition: TList.h:47
Bool_t fLocalFrame
Definition: TBuffer3D.h:92
static void DrawOutlineCube(TList *outline, Double_t *rmin, Double_t *rmax)
Draw cube outline with 3d polylines.
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:46
TString fOption
options
Definition: TPolyLine3D.h:46
float ymax
Definition: THbookFile.cxx:93
Bool_t SectionsValid(UInt_t mask) const
Definition: TBuffer3D.h:69
virtual TObject * First() const
Return the first object in the list. Returns 0 when list is empty.
Definition: TList.cxx:557
virtual void Print(Option_t *option="") const
Dump this 3-D polyline with its attributes on stdout.
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:488
Collection abstract base class.
Definition: TCollection.h:48
TPolyLine3D & operator=(const TPolyLine3D &polylin)
assignment operator
virtual Int_t AddObject(const TBuffer3D &buffer, Bool_t *addChildren=0)=0
unsigned int UInt_t
Definition: RtypesCore.h:42
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:925
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:19
static void IndentLevel()
Functions used by ls() to indent an object hierarchy.
Definition: TROOT.cxx:2618
virtual void SetByteCount(UInt_t cntpos, Bool_t packInVersion=kFALSE)=0
virtual TObject * After(const TObject *obj) const
Returns the object after object obj.
Definition: TList.cxx:289
Int_t fN
Number of points.
Definition: TPolyLine3D.h:44
static double p1(double t, double a, double b)
float xmax
Definition: THbookFile.cxx:93
TObject * fID
Definition: TBuffer3D.h:89
virtual void ReadFastArray(Bool_t *b, Int_t n)=0
virtual void WriteFastArray(const Bool_t *b, Int_t n)=0
virtual Color_t GetLineColor() const
Return the line color.
Definition: TAttLine.h:39
virtual void Draw(Option_t *option="")
Draw this 3-D polyline with its current attributes.
virtual void SetPoint(Int_t point, Double_t x, Double_t y, Double_t z)
Set point n to x, y, z.
static const double x1[5]
#define ClassImp(name)
Definition: Rtypes.h:279
double Double_t
Definition: RtypesCore.h:55
Int_t fLastPoint
The index of the last filled point.
Definition: TPolyLine3D.h:47
virtual Int_t SetNextPoint(Double_t x, Double_t y, Double_t z)
Set point following LastPoint to x, y, z.
Int_t fColor
Definition: TBuffer3D.h:90
Double_t y[n]
Definition: legend1.C:17
virtual void Paint(Option_t *option="")
Paint a TPolyLine3D.
Int_t DistancetoLine(Int_t px, Int_t py, Double_t xp1, Double_t yp1, Double_t xp2, Double_t yp2)
Compute distance from point px,py to a line.
Definition: TAttLine.cxx:196
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:567
virtual void Copy(TObject &polyline) const
Copy polyline to polyline obj.
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition: TAttLine.h:48
Mother of all ROOT objects.
Definition: TObject.h:37
you should not use this method at all Int_t Int_t z
Definition: TRolke.cxx:630
Int_t * fSegs
Definition: TBuffer3D.h:115
virtual void Add(TObject *obj)
Definition: TList.h:81
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:202
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
#define gPad
Definition: TVirtualPad.h:289
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save primitive as a C++ statement(s) on output stream.
Short_t fTransparency
Definition: TBuffer3D.h:91
virtual void DrawPolyLine(Int_t n, Float_t *p, Option_t *option="")
Draw 3-D polyline with new coordinates.
const Bool_t kTRUE
Definition: Rtypes.h:91
R__EXTERN TGeometry * gGeometry
Definition: TGeometry.h:162
virtual Int_t Size() const
Definition: TPolyLine3D.h:80
const Int_t n
Definition: legend1.C:16
Line Attributes class.
Definition: TAttLine.h:24
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0