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