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)
234 pl.TPolyLine3D::Copy(*this);
235 return *this;
236}
237
238////////////////////////////////////////////////////////////////////////////////
239/// 3-D polyline destructor.
240
242{
243 if (fP) delete [] fP;
244}
245
246////////////////////////////////////////////////////////////////////////////////
247/// 3-D polyline copy ctor.
248
249TPolyLine3D::TPolyLine3D(const TPolyLine3D &polyline) : TObject(polyline), TAttLine(polyline), TAtt3D(polyline)
250{
251 fP = 0;
252 fLastPoint = 0;
253 fN = 0;
254 polyline.TPolyLine3D::Copy(*this);
255}
256
257////////////////////////////////////////////////////////////////////////////////
258/// Copy polyline to polyline obj.
259
261{
262 auto &tgt = static_cast<TPolyLine3D &>(obj);
263 TObject::Copy(obj);
264 TAttLine::Copy(tgt);
265 tgt.fN = fN;
266 if (tgt.fP)
267 delete [] tgt.fP;
268 if (fN > 0) {
269 tgt.fP = new Float_t[3*fN];
270 for (Int_t i=0; i<3*fN;i++)
271 tgt.fP[i] = fP[i];
272 } else {
273 tgt.fP = nullptr;
274 }
275 tgt.fOption = fOption;
276 tgt.fLastPoint = fLastPoint;
277}
278
279////////////////////////////////////////////////////////////////////////////////
280/// Compute distance from point px,py to a 3-D polyline.
281/// Compute the closest distance of approach from point px,py to each segment
282/// of the polyline.
283/// Returns when the distance found is below DistanceMaximum.
284/// The distance is computed in pixels units.
285
287{
288 const Int_t inaxis = 7;
289 Int_t dist = 9999;
290
291 Int_t puxmin = gPad->XtoAbsPixel(gPad->GetUxmin());
292 Int_t puymin = gPad->YtoAbsPixel(gPad->GetUymin());
293 Int_t puxmax = gPad->XtoAbsPixel(gPad->GetUxmax());
294 Int_t puymax = gPad->YtoAbsPixel(gPad->GetUymax());
295
296 // return if point is not in the user area
297 if (px < puxmin - inaxis) return dist;
298 if (py > puymin + inaxis) return dist;
299 if (px > puxmax + inaxis) return dist;
300 if (py < puymax - inaxis) return dist;
301
302 TView *view = gPad->GetView();
303 if (!view) return dist;
304
305 Int_t i, dsegment;
306 Double_t x1,y1,x2,y2;
307 Float_t xndc[3];
308 for (i=0;i<Size()-1;i++) {
309 view->WCtoNDC(&fP[3*i], xndc);
310 x1 = xndc[0];
311 y1 = xndc[1];
312 view->WCtoNDC(&fP[3*i+3], xndc);
313 x2 = xndc[0];
314 y2 = xndc[1];
315 dsegment = DistancetoLine(px,py,x1,y1,x2,y2);
316 if (dsegment < dist) dist = dsegment;
317 }
318 return dist;
319}
320
321////////////////////////////////////////////////////////////////////////////////
322/// Draw this 3-D polyline with its current attributes.
323
325{
326 AppendPad(option);
327}
328
329////////////////////////////////////////////////////////////////////////////////
330/// Draw cube outline with 3d polylines.
331///
332/// ~~~ {.cpp}
333/// xmin = fRmin[0] xmax = fRmax[0]
334/// ymin = fRmin[1] ymax = fRmax[1]
335/// zmin = fRmin[2] zmax = fRmax[2]
336///
337/// (xmin,ymax,zmax) +---------+ (xmax,ymax,zmax)
338/// / /|
339/// / / |
340/// / / |
341///(xmin,ymin,zmax) +---------+ |
342/// | | + (xmax,ymax,zmin)
343/// | | /
344/// | | /
345/// | |/
346/// +---------+
347/// (xmin,ymin,zmin) (xmax,ymin,zmin)
348/// ~~~
349
351{
352 Double_t xmin = rmin[0]; Double_t xmax = rmax[0];
353 Double_t ymin = rmin[1]; Double_t ymax = rmax[1];
354 Double_t zmin = rmin[2]; Double_t zmax = rmax[2];
355
356 TPolyLine3D *pl3d = (TPolyLine3D *)outline->First();
357 if (!pl3d) {
358 TView *view = gPad->GetView();
359 if (!view) return;
360 TPolyLine3D *p1 = new TPolyLine3D(4);
361 TPolyLine3D *p2 = new TPolyLine3D(4);
362 TPolyLine3D *p3 = new TPolyLine3D(4);
363 TPolyLine3D *p4 = new TPolyLine3D(4);
364 p1->SetLineColor(view->GetLineColor());
365 p1->SetLineStyle(view->GetLineStyle());
366 p1->SetLineWidth(view->GetLineWidth());
367 p1->Copy(*p2);
368 p1->Copy(*p3);
369 p1->Copy(*p4);
370 outline->Add(p1);
371 outline->Add(p2);
372 outline->Add(p3);
373 outline->Add(p4);
374 }
375
376 pl3d = (TPolyLine3D *)outline->First();
377
378 if (pl3d) {
379 pl3d->SetPoint(0, xmin, ymin, zmin);
380 pl3d->SetPoint(1, xmax, ymin, zmin);
381 pl3d->SetPoint(2, xmax, ymax, zmin);
382 pl3d->SetPoint(3, xmin, ymax, zmin);
383 }
384
385 pl3d = (TPolyLine3D *)outline->After(pl3d);
386
387 if (pl3d) {
388 pl3d->SetPoint(0, xmax, ymin, zmin);
389 pl3d->SetPoint(1, xmax, ymin, zmax);
390 pl3d->SetPoint(2, xmax, ymax, zmax);
391 pl3d->SetPoint(3, xmax, ymax, zmin);
392 }
393
394 pl3d = (TPolyLine3D *)outline->After(pl3d);
395
396 if (pl3d) {
397 pl3d->SetPoint(0, xmax, ymin, zmax);
398 pl3d->SetPoint(1, xmin, ymin, zmax);
399 pl3d->SetPoint(2, xmin, ymax, zmax);
400 pl3d->SetPoint(3, xmax, ymax, zmax);
401 }
402
403 pl3d = (TPolyLine3D *)outline->After(pl3d);
404
405 if (pl3d) {
406 pl3d->SetPoint(0, xmin, ymin, zmax);
407 pl3d->SetPoint(1, xmin, ymin, zmin);
408 pl3d->SetPoint(2, xmin, ymax, zmin);
409 pl3d->SetPoint(3, xmin, ymax, zmax);
410 }
411}
412
413////////////////////////////////////////////////////////////////////////////////
414/// Draw 3-D polyline with new coordinates. Creates a new polyline which
415/// will be adopted by the pad in which it is drawn. Does not change the
416/// original polyline (should be static method).
417
419{
420 TPolyLine3D *newpolyline = new TPolyLine3D();
421 Int_t size = 3*Size();
422 newpolyline->fN =n;
423 newpolyline->fP = new Float_t[size];
424 for (Int_t i=0; i<size;i++) { newpolyline->fP[i] = p[i];}
425 TAttLine::Copy(*newpolyline);
426 newpolyline->fOption = fOption;
427 newpolyline->fLastPoint = fLastPoint;
428 newpolyline->SetBit(kCanDelete);
429 newpolyline->AppendPad(option);
430}
431
432////////////////////////////////////////////////////////////////////////////////
433/// Execute action corresponding to one event.
434
436{
437 if (!gPad) return;
438 if (gPad->GetView()) gPad->GetView()->ExecuteRotateView(event, px, py);
439}
440
441////////////////////////////////////////////////////////////////////////////////
442/// List this 3-D polyline.
443
444void TPolyLine3D::ls(Option_t *option) const
445{
447 std::cout <<"PolyLine3D N=" <<fN<<" Option="<<option<<std::endl;
448}
449
450////////////////////////////////////////////////////////////////////////////////
451/// Merge polylines in the collection in this polyline
452
454{
455 if (!li) return 0;
456 TIter next(li);
457
458 //first loop to count the number of entries
459 TPolyLine3D *pl;
460 Int_t npoints = 0;
461 while ((pl = (TPolyLine3D*)next())) {
462 if (!pl->InheritsFrom(TPolyLine3D::Class())) {
463 Error("Add","Attempt to add object of class: %s to a %s",pl->ClassName(),this->ClassName());
464 return -1;
465 }
466 npoints += pl->Size();
467 }
468
469 //extend this polyline to hold npoints
470 SetPoint(npoints-1,0,0,0);
471
472 //merge all polylines
473 next.Reset();
474 while ((pl = (TPolyLine3D*)next())) {
475 Int_t np = pl->Size();
476 Float_t *p = pl->GetP();
477 for (Int_t i=0;i<np;i++) {
478 SetPoint(i,p[3*i],p[3*i+1],p[3*i+2]);
479 }
480 }
481
482 return npoints;
483}
484
485////////////////////////////////////////////////////////////////////////////////
486/// Paint a TPolyLine3D.
487
488void TPolyLine3D::Paint(Option_t * /* option */ )
489{
490 UInt_t i;
491
492 // No need to continue if there is nothing to paint
493 if (Size() <= 0) return;
494
495 static TBuffer3D buffer(TBuffer3DTypes::kLine);
496
497 // TPolyLine3D can only be described by filling the TBuffer3D 'tesselation'
498 // parts - so there are no 'optional' sections - we just fill everything.
499
500 buffer.ClearSectionsValid();
501
502 // Section kCore
503 buffer.fID = this;
504 buffer.fColor = GetLineColor();
505 buffer.fTransparency = 0;
506 buffer.fLocalFrame = kFALSE;
508
509 // We fill kCore and kRawSizes on first pass and try with viewer
510 TVirtualViewer3D * viewer3D = gPad->GetViewer3D();
511 if (!viewer3D) return;
512 Int_t reqSections = viewer3D->AddObject(buffer);
513 if (reqSections == TBuffer3D::kNone) {
514 return;
515 }
516
517 if (reqSections & TBuffer3D::kRawSizes) {
518 Int_t nbPnts = Size();
519 Int_t nbSegs = nbPnts-1;
520 if (!buffer.SetRawSizes(nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, 0, 0)) {
521 return;
522 }
524 }
525
526 if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
527 // Points
528 for (i=0; i<3*buffer.NbPnts(); i++) {
529 buffer.fPnts[i] = (Double_t)fP[i];
530 }
531
532 // Transform points
533 if (gGeometry && !buffer.fLocalFrame) {
534 Double_t dlocal[3];
535 Double_t dmaster[3];
536 for (UInt_t j=0; j<buffer.NbPnts(); j++) {
537 dlocal[0] = buffer.fPnts[3*j];
538 dlocal[1] = buffer.fPnts[3*j+1];
539 dlocal[2] = buffer.fPnts[3*j+2];
540 gGeometry->Local2Master(&dlocal[0],&dmaster[0]);
541 buffer.fPnts[3*j] = dmaster[0];
542 buffer.fPnts[3*j+1] = dmaster[1];
543 buffer.fPnts[3*j+2] = dmaster[2];
544 }
545 }
546
547 // Basic colors: 0, 1, ... 8
548 Int_t c = (((GetLineColor()) %8) -1) * 4;
549 if (c < 0) c = 0;
550
551 // Segments
552 for (i = 0; i < buffer.NbSegs(); i++) {
553 buffer.fSegs[3*i ] = c;
554 buffer.fSegs[3*i+1] = i;
555 buffer.fSegs[3*i+2] = i+1;
556 }
557
559
561 }
562
563 viewer3D->AddObject(buffer);
564}
565
566////////////////////////////////////////////////////////////////////////////////
567/// Dump this 3-D polyline with its attributes on stdout.
568
569void TPolyLine3D::Print(Option_t *option) const
570{
571 printf(" TPolyLine3D N=%d, Option=%s\n",fN,option);
572 TString opt = option;
573 opt.ToLower();
574 if (opt.Contains("all")) {
575 for (Int_t i=0;i<Size();i++) {
576 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]);
577 }
578 }
579}
580
581////////////////////////////////////////////////////////////////////////////////
582/// Save primitive as a C++ statement(s) on output stream.
583
584void TPolyLine3D::SavePrimitive(std::ostream &out, Option_t * /*= ""*/)
585{
586 char quote = '"';
587 out<<" "<<std::endl;
588 if (gROOT->ClassSaved(TPolyLine3D::Class())) {
589 out<<" ";
590 } else {
591 out<<" TPolyLine3D *";
592 }
593 Int_t size=Size();
594 out<<"pline3D = new TPolyLine3D("<<fN<<","<<quote<<fOption<<quote<<");"<<std::endl;
595
596 SaveLineAttributes(out,"pline3D",1,1,1);
597
598 if (size > 0) {
599 for (Int_t i=0;i<size;i++)
600 out<<" pline3D->SetPoint("<<i<<","<<fP[3*i]<<","<<fP[3*i+1]<<","<<fP[3*i+2]<<");"<<std::endl;
601 }
602 out<<" pline3D->Draw();"<<std::endl;
603}
604
605////////////////////////////////////////////////////////////////////////////////
606/// Set point following LastPoint to x, y, z.
607/// Returns index of the point (new last point).
608
610{
611 fLastPoint++;
612 SetPoint(fLastPoint, x, y, z);
613 return fLastPoint;
614}
615
616////////////////////////////////////////////////////////////////////////////////
617/// Set point n to x, y, z.
618/// If n is more then the current TPolyLine3D size (n > fN) then
619/// the polyline will be resized to contain at least n points.
620
622{
623 if (n < 0) return;
624 if (!fP || n >= fN) {
625 // re-allocate the object
626 Int_t newN = TMath::Max(2*fN,n+1);
627 Float_t *savepoint = new Float_t [3*newN];
628 if (fP && fN){
629 memcpy(savepoint,fP,3*fN*sizeof(Float_t));
630 memset(&savepoint[3*fN],0,(newN-fN)*sizeof(Float_t));
631 delete [] fP;
632 }
633 fP = savepoint;
634 fN = newN;
635 }
636 fP[3*n ] = x;
637 fP[3*n+1] = y;
638 fP[3*n+2] = z;
640}
641
642////////////////////////////////////////////////////////////////////////////////
643/// Re-initialize polyline with n points (0,0,0).
644/// if n <= 0 the current array of points is deleted.
645
647{
648 fOption = option;
649 if (n <= 0) {
650 fN = 0;
651 fLastPoint = -1;
652 delete [] fP;
653 fP = 0;
654 return;
655 }
656 fN = n;
657 if (fP) delete [] fP;
658 fP = new Float_t[3*fN];
659 memset(fP,0,3*fN*sizeof(Float_t));
660 fLastPoint = fN-1;
661}
662
663////////////////////////////////////////////////////////////////////////////////
664/// Re-initialize polyline with n points from p. If p=0 initialize with 0.
665/// if n <= 0 the current array of points is deleted.
666
668{
669 fOption = option;
670 if (n <= 0) {
671 fN = 0;
672 fLastPoint = -1;
673 delete [] fP;
674 fP = 0;
675 return;
676 }
677 fN = n;
678 if (fP) delete [] fP;
679 fP = new Float_t[3*fN];
680 if (p) {
681 for (Int_t i=0; i<fN;i++) {
682 fP[3*i] = p[3*i];
683 fP[3*i+1] = p[3*i+1];
684 fP[3*i+2] = p[3*i+2];
685 }
686 } else {
687 memset(fP,0,3*fN*sizeof(Float_t));
688 }
689 fLastPoint = fN-1;
690}
691
692////////////////////////////////////////////////////////////////////////////////
693/// Re-initialize polyline with n points from p. If p=0 initialize with 0.
694/// if n <= 0 the current array of points is deleted.
695
697{
698 fOption = option;
699 if (n <= 0) {
700 fN = 0;
701 fLastPoint = -1;
702 delete [] fP;
703 fP = 0;
704 return;
705 }
706 fN = n;
707 if (fP) delete [] fP;
708 fP = new Float_t[3*fN];
709 if (p) {
710 for (Int_t i=0; i<fN;i++) {
711 fP[3*i] = (Float_t) p[3*i];
712 fP[3*i+1] = (Float_t) p[3*i+1];
713 fP[3*i+2] = (Float_t) p[3*i+2];
714 }
715 } else {
716 memset(fP,0,3*fN*sizeof(Float_t));
717 }
718 fLastPoint = fN-1;
719}
720
721////////////////////////////////////////////////////////////////////////////////
722/// Stream a 3-D polyline object.
723
724void TPolyLine3D::Streamer(TBuffer &b)
725{
726 UInt_t R__s, R__c;
727 if (b.IsReading()) {
728 b.ReadVersion(&R__s, &R__c);
729 b.ClassBegin(TPolyLine3D::IsA());
730 b.ClassMember("TObject");
731 TObject::Streamer(b);
732 b.ClassMember("TAttLine");
733 TAttLine::Streamer(b);
734 b.ClassMember("fN", "Int_t");
735 b >> fN;
736 if (fN) {
737 fP = new Float_t[3*fN];
738 b.ClassMember("fP", "Float_t", 3 * fN);
739 b.ReadFastArray(fP, 3 * fN);
740 }
741 b.ClassMember("fOption", "TString");
742 fOption.Streamer(b);
743 fLastPoint = fN-1;
744 b.ClassEnd(TPolyLine3D::IsA());
745 b.CheckByteCount(R__s, R__c, TPolyLine3D::IsA());
746 } else {
747 R__c = b.WriteVersion(TPolyLine3D::IsA(), kTRUE);
748 b.ClassBegin(TPolyLine3D::IsA());
749 b.ClassMember("TObject");
750 TObject::Streamer(b);
751 b.ClassMember("TAttLine");
752 TAttLine::Streamer(b);
753 b.ClassMember("fN", "Int_t");
754 Int_t size = Size();
755 b << size;
756 if (size) {
757 b.ClassMember("fP", "Float_t", 3 * size);
758 b.WriteFastArray(fP, 3 * size);
759 }
760 b.ClassMember("fOption", "TString");
761 fOption.Streamer(b);
762 b.ClassEnd(TPolyLine3D::IsA());
763 b.SetByteCount(R__c, kTRUE);
764 }
765}
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
static const double x2[5]
static const double x1[5]
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
int Int_t
Definition RtypesCore.h:45
unsigned int UInt_t
Definition RtypesCore.h:46
const Bool_t kFALSE
Definition RtypesCore.h:101
double Double_t
Definition RtypesCore.h:59
float Float_t
Definition RtypesCore.h:57
const Bool_t kTRUE
Definition RtypesCore.h:100
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:404
#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:245
void Copy(TAttLine &attline) const
Copy this line attributes to a new TAttLine.
Definition TAttLine.cxx:175
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:209
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:273
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:65
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:38
virtual void Add(TObject *obj)
Definition TList.h:81
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:41
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:200
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition TObject.cxx:177
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:766
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:515
virtual void Copy(TObject &object) const
Copy this to obj.
Definition TObject.cxx:133
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:963
@ kCanDelete
if object in a list can be deleted
Definition TObject.h:62
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:2819
Basic string class.
Definition TString.h:136
void ToLower()
Change string to lower-case.
Definition TString.cxx:1150
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:208