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