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 auto 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}
103
104////////////////////////////////////////////////////////////////////////////////
105/// 3-D polyline normal constructor with initialization to 0.
106/// If n < 0 the default size (2 points) is set.
107
109{
110 fOption = option;
112 if (n <= 0)
113 return;
114
115 fN = n;
116 fP = new Float_t[3*fN];
117 for (Int_t i=0; i<3*fN; i++) fP[i] = 0;
118}
119
120////////////////////////////////////////////////////////////////////////////////
121/// 3-D polyline normal constructor. Polyline is intialized with p.
122/// If n < 0 the default size (2 points) is set.
123
125{
126 fOption = option;
128 if (n <= 0)
129 return;
130
131 fN = n;
132 fP = new Float_t[3*fN];
133 for (Int_t i=0; i<3*n; i++) {
134 fP[i] = p[i];
135 }
136 fLastPoint = fN-1;
137}
138
139////////////////////////////////////////////////////////////////////////////////
140/// 3-D polyline normal constructor. Polyline is initialized with p
141/// (cast to float). If n < 0 the default size (2 points) is set.
142
144{
145 fOption = option;
147 if (n <= 0)
148 return;
149
150 fN = n;
151 fP = new Float_t[3*fN];
152 for (Int_t i=0; i<3*n; i++) {
153 fP[i] = (Float_t) p[i];
154 }
155 fLastPoint = fN-1;
156}
157
158////////////////////////////////////////////////////////////////////////////////
159/// 3-D polyline normal constructor. Polyline is initialized withe the
160/// x, y ,z arrays. If n < 0 the default size (2 points) is set.
161
163{
164 fOption = option;
166 if (n <= 0)
167 return;
168
169 fN = n;
170 fP = new Float_t[3*fN];
171 Int_t j = 0;
172 for (Int_t i=0; i<n;i++) {
173 fP[j] = x[i];
174 fP[j+1] = y[i];
175 fP[j+2] = z[i];
176 j += 3;
177 }
178 fLastPoint = fN-1;
179}
180
181////////////////////////////////////////////////////////////////////////////////
182/// 3-D polyline normal constructor. Polyline is initialized withe the
183/// x, y, z arrays (which are cast to float).
184/// If n < 0 the default size (2 points) is set.
185
187{
188 fOption = option;
190 if (n <= 0)
191 return;
192
193 fN = n;
194 fP = new Float_t[3*fN];
195 Int_t j = 0;
196 for (Int_t i=0; i<n;i++) {
197 fP[j] = (Float_t) x[i];
198 fP[j+1] = (Float_t) y[i];
199 fP[j+2] = (Float_t) z[i];
200 j += 3;
201 }
202 fLastPoint = fN-1;
203}
204
205////////////////////////////////////////////////////////////////////////////////
206/// assignment operator
207
209{
210 if(this != &pl)
211 pl.TPolyLine3D::Copy(*this);
212 return *this;
213}
214
215////////////////////////////////////////////////////////////////////////////////
216/// 3-D polyline destructor.
217
219{
220 if (fP) delete [] fP;
221}
222
223////////////////////////////////////////////////////////////////////////////////
224/// 3-D polyline copy ctor.
225
226TPolyLine3D::TPolyLine3D(const TPolyLine3D &polyline) : TObject(polyline), TAttLine(polyline), TAtt3D(polyline)
227{
228 polyline.TPolyLine3D::Copy(*this);
229}
230
231////////////////////////////////////////////////////////////////////////////////
232/// Copy polyline to polyline obj.
233
235{
236 auto &tgt = static_cast<TPolyLine3D &>(obj);
237 TObject::Copy(obj);
238 TAttLine::Copy(tgt);
239 tgt.fN = fN;
240 if (tgt.fP)
241 delete [] tgt.fP;
242 if (fN > 0) {
243 tgt.fP = new Float_t[3*fN];
244 for (Int_t i=0; i<3*fN;i++)
245 tgt.fP[i] = fP[i];
246 } else {
247 tgt.fP = nullptr;
248 }
249 tgt.fOption = fOption;
250 tgt.fLastPoint = fLastPoint;
251}
252
253////////////////////////////////////////////////////////////////////////////////
254/// Compute distance from point px,py to a 3-D polyline.
255/// Compute the closest distance of approach from point px,py to each segment
256/// of the polyline.
257/// Returns when the distance found is below DistanceMaximum.
258/// The distance is computed in pixels units.
259
261{
262 const Int_t inaxis = 7;
263 Int_t dist = 9999;
264
265 Int_t puxmin = gPad->XtoAbsPixel(gPad->GetUxmin());
266 Int_t puymin = gPad->YtoAbsPixel(gPad->GetUymin());
267 Int_t puxmax = gPad->XtoAbsPixel(gPad->GetUxmax());
268 Int_t puymax = gPad->YtoAbsPixel(gPad->GetUymax());
269
270 // return if point is not in the user area
271 if (px < puxmin - inaxis) return dist;
272 if (py > puymin + inaxis) return dist;
273 if (px > puxmax + inaxis) return dist;
274 if (py < puymax - inaxis) return dist;
275
276 TView *view = gPad->GetView();
277 if (!view) return dist;
278
279 Int_t i, dsegment;
281 Float_t xndc[3];
282 for (i=0;i<Size()-1;i++) {
283 view->WCtoNDC(&fP[3*i], xndc);
284 x1 = xndc[0];
285 y1 = xndc[1];
286 view->WCtoNDC(&fP[3*i+3], xndc);
287 x2 = xndc[0];
288 y2 = xndc[1];
289 dsegment = DistancetoLine(px,py,x1,y1,x2,y2);
290 if (dsegment < dist) dist = dsegment;
291 }
292 return dist;
293}
294
295////////////////////////////////////////////////////////////////////////////////
296/// Draw this 3-D polyline with its current attributes.
297
299{
301}
302
303////////////////////////////////////////////////////////////////////////////////
304/// Draw cube outline with 3d polylines.
305///
306/// ~~~ {.cpp}
307/// xmin = fRmin[0] xmax = fRmax[0]
308/// ymin = fRmin[1] ymax = fRmax[1]
309/// zmin = fRmin[2] zmax = fRmax[2]
310///
311/// (xmin,ymax,zmax) +---------+ (xmax,ymax,zmax)
312/// / /|
313/// / / |
314/// / / |
315///(xmin,ymin,zmax) +---------+ |
316/// | | + (xmax,ymax,zmin)
317/// | | /
318/// | | /
319/// | |/
320/// +---------+
321/// (xmin,ymin,zmin) (xmax,ymin,zmin)
322/// ~~~
323
325{
326 Double_t xmin = rmin[0]; Double_t xmax = rmax[0];
327 Double_t ymin = rmin[1]; Double_t ymax = rmax[1];
328 Double_t zmin = rmin[2]; Double_t zmax = rmax[2];
329
330 TPolyLine3D *pl3d = (TPolyLine3D *)outline->First();
331 if (!pl3d) {
332 TView *view = gPad->GetView();
333 if (!view) return;
334 TPolyLine3D *p1 = new TPolyLine3D(4);
335 TPolyLine3D *p2 = new TPolyLine3D(4);
336 TPolyLine3D *p3 = new TPolyLine3D(4);
337 TPolyLine3D *p4 = new TPolyLine3D(4);
338 p1->SetLineColor(view->GetLineColor());
339 p1->SetLineStyle(view->GetLineStyle());
340 p1->SetLineWidth(view->GetLineWidth());
341 p1->Copy(*p2);
342 p1->Copy(*p3);
343 p1->Copy(*p4);
344 outline->Add(p1);
345 outline->Add(p2);
346 outline->Add(p3);
347 outline->Add(p4);
348 }
349
350 pl3d = (TPolyLine3D *)outline->First();
351
352 if (pl3d) {
353 pl3d->SetPoint(0, xmin, ymin, zmin);
354 pl3d->SetPoint(1, xmax, ymin, zmin);
355 pl3d->SetPoint(2, xmax, ymax, zmin);
356 pl3d->SetPoint(3, xmin, ymax, zmin);
357 }
358
359 pl3d = (TPolyLine3D *)outline->After(pl3d);
360
361 if (pl3d) {
362 pl3d->SetPoint(0, xmax, ymin, zmin);
363 pl3d->SetPoint(1, xmax, ymin, zmax);
364 pl3d->SetPoint(2, xmax, ymax, zmax);
365 pl3d->SetPoint(3, xmax, ymax, zmin);
366 }
367
368 pl3d = (TPolyLine3D *)outline->After(pl3d);
369
370 if (pl3d) {
371 pl3d->SetPoint(0, xmax, ymin, zmax);
372 pl3d->SetPoint(1, xmin, ymin, zmax);
373 pl3d->SetPoint(2, xmin, ymax, zmax);
374 pl3d->SetPoint(3, xmax, ymax, zmax);
375 }
376
377 pl3d = (TPolyLine3D *)outline->After(pl3d);
378
379 if (pl3d) {
380 pl3d->SetPoint(0, xmin, ymin, zmax);
381 pl3d->SetPoint(1, xmin, ymin, zmin);
382 pl3d->SetPoint(2, xmin, ymax, zmin);
383 pl3d->SetPoint(3, xmin, ymax, zmax);
384 }
385}
386
387////////////////////////////////////////////////////////////////////////////////
388/// Draw 3-D polyline with new coordinates. Creates a new polyline which
389/// will be adopted by the pad in which it is drawn. Does not change the
390/// original polyline (should be static method).
391
393{
394 TPolyLine3D *newpolyline = new TPolyLine3D();
395 Int_t size = 3*Size();
396 newpolyline->fN =n;
397 newpolyline->fP = new Float_t[size];
398 for (Int_t i=0; i<size;i++) { newpolyline->fP[i] = p[i];}
399 TAttLine::Copy(*newpolyline);
400 newpolyline->fOption = fOption;
401 newpolyline->fLastPoint = fLastPoint;
402 newpolyline->SetBit(kCanDelete);
403 newpolyline->AppendPad(option);
404}
405
406////////////////////////////////////////////////////////////////////////////////
407/// Execute action corresponding to one event.
408
410{
411 if (!gPad) return;
412 if (gPad->GetView()) gPad->GetView()->ExecuteRotateView(event, px, py);
413}
414
415////////////////////////////////////////////////////////////////////////////////
416/// List this 3-D polyline.
417
419{
421 std::cout <<"PolyLine3D N=" <<fN<<" Option="<<option<<std::endl;
422}
423
424////////////////////////////////////////////////////////////////////////////////
425/// Merge polylines in the collection in this polyline
426
428{
429 if (!li) return 0;
430 TIter next(li);
431
432 //first loop to count the number of entries
433 TPolyLine3D *pl;
434 Int_t npoints = 0;
435 while ((pl = (TPolyLine3D*)next())) {
436 if (!pl->InheritsFrom(TPolyLine3D::Class())) {
437 Error("Add","Attempt to add object of class: %s to a %s",pl->ClassName(),this->ClassName());
438 return -1;
439 }
440 npoints += pl->Size();
441 }
442
443 //extend this polyline to hold npoints
444 SetPoint(npoints-1,0,0,0);
445
446 //merge all polylines
447 next.Reset();
448 while ((pl = (TPolyLine3D*)next())) {
449 Int_t np = pl->Size();
450 Float_t *p = pl->GetP();
451 for (Int_t i=0;i<np;i++) {
452 SetPoint(i,p[3*i],p[3*i+1],p[3*i+2]);
453 }
454 }
455
456 return npoints;
457}
458
459////////////////////////////////////////////////////////////////////////////////
460/// Paint a TPolyLine3D.
461
462void TPolyLine3D::Paint(Option_t * /* option */ )
463{
464 UInt_t i;
465
466 // No need to continue if there is nothing to paint
467 if (Size() <= 0) return;
468
469 static TBuffer3D buffer(TBuffer3DTypes::kLine);
470
471 // TPolyLine3D can only be described by filling the TBuffer3D 'tesselation'
472 // parts - so there are no 'optional' sections - we just fill everything.
473
474 buffer.ClearSectionsValid();
475
476 // Section kCore
477 buffer.fID = this;
478 buffer.fColor = GetLineColor();
479 buffer.fTransparency = 0;
480 buffer.fLocalFrame = kFALSE;
482
483 // We fill kCore and kRawSizes on first pass and try with viewer
484 TVirtualViewer3D * viewer3D = gPad->GetViewer3D();
485 if (!viewer3D) return;
486 Int_t reqSections = viewer3D->AddObject(buffer);
487 if (reqSections == TBuffer3D::kNone) {
488 return;
489 }
490
491 if (reqSections & TBuffer3D::kRawSizes) {
492 Int_t nbPnts = Size();
493 Int_t nbSegs = nbPnts-1;
494 if (!buffer.SetRawSizes(nbPnts, 3*nbPnts, nbSegs, 3*nbSegs, 0, 0)) {
495 return;
496 }
498 }
499
500 if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
501 // Points
502 for (i=0; i<3*buffer.NbPnts(); i++) {
503 buffer.fPnts[i] = (Double_t)fP[i];
504 }
505
506 // Transform points
507 if (gGeometry && !buffer.fLocalFrame) {
508 Double_t dlocal[3];
509 Double_t dmaster[3];
510 for (UInt_t j=0; j<buffer.NbPnts(); j++) {
511 dlocal[0] = buffer.fPnts[3*j];
512 dlocal[1] = buffer.fPnts[3*j+1];
513 dlocal[2] = buffer.fPnts[3*j+2];
514 gGeometry->Local2Master(&dlocal[0],&dmaster[0]);
515 buffer.fPnts[3*j] = dmaster[0];
516 buffer.fPnts[3*j+1] = dmaster[1];
517 buffer.fPnts[3*j+2] = dmaster[2];
518 }
519 }
520
521 // Basic colors: 0, 1, ... 8
522 Int_t c = (((GetLineColor()) %8) -1) * 4;
523 if (c < 0) c = 0;
524
525 // Segments
526 for (i = 0; i < buffer.NbSegs(); i++) {
527 buffer.fSegs[3*i ] = c;
528 buffer.fSegs[3*i+1] = i;
529 buffer.fSegs[3*i+2] = i+1;
530 }
531
533
535 }
536
537 viewer3D->AddObject(buffer);
538}
539
540////////////////////////////////////////////////////////////////////////////////
541/// Dump this 3-D polyline with its attributes on stdout.
542
544{
545 printf(" TPolyLine3D N=%d, Option=%s\n",fN,option);
546 TString opt = option;
547 opt.ToLower();
548 if (opt.Contains("all")) {
549 for (Int_t i=0;i<Size();i++) {
550 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]);
551 }
552 }
553}
554
555////////////////////////////////////////////////////////////////////////////////
556/// Save primitive as a C++ statement(s) on output stream.
557
558void TPolyLine3D::SavePrimitive(std::ostream &out, Option_t * /*= ""*/)
559{
560 char quote = '"';
561 out<<" "<<std::endl;
562 if (gROOT->ClassSaved(TPolyLine3D::Class())) {
563 out<<" ";
564 } else {
565 out<<" TPolyLine3D *";
566 }
567 Int_t size=Size();
568 out<<"pline3D = new TPolyLine3D("<<fN<<","<<quote<<fOption<<quote<<");"<<std::endl;
569
570 SaveLineAttributes(out,"pline3D",1,1,1);
571
572 if (size > 0) {
573 for (Int_t i=0;i<size;i++)
574 out<<" pline3D->SetPoint("<<i<<","<<fP[3*i]<<","<<fP[3*i+1]<<","<<fP[3*i+2]<<");"<<std::endl;
575 }
576 out<<" pline3D->Draw();"<<std::endl;
577}
578
579////////////////////////////////////////////////////////////////////////////////
580/// Set point following LastPoint to x, y, z.
581/// Returns index of the point (new last point).
582
584{
585 fLastPoint++;
586 SetPoint(fLastPoint, x, y, z);
587 return fLastPoint;
588}
589
590////////////////////////////////////////////////////////////////////////////////
591/// Set point n to x, y, z.
592/// If n is more then the current TPolyLine3D size (n > fN) then
593/// the polyline will be resized to contain at least n points.
594
596{
597 if (n < 0) return;
598 if (!fP || n >= fN) {
599 // re-allocate the object
600 Int_t newN = TMath::Max(2*fN,n+1);
601 Float_t *savepoint = new Float_t [3*newN];
602 if (fP && fN){
603 memcpy(savepoint,fP,3*fN*sizeof(Float_t));
604 memset(&savepoint[3*fN],0,(newN-fN)*sizeof(Float_t));
605 delete [] fP;
606 }
607 fP = savepoint;
608 fN = newN;
609 }
610 fP[3*n ] = x;
611 fP[3*n+1] = y;
612 fP[3*n+2] = z;
614}
615
616////////////////////////////////////////////////////////////////////////////////
617/// Re-initialize polyline with n points (0,0,0).
618/// if n <= 0 the current array of points is deleted.
619
621{
622 fOption = option;
623 if (n <= 0) {
624 fN = 0;
625 fLastPoint = -1;
626 delete [] fP;
627 fP = nullptr;
628 return;
629 }
630 fN = n;
631 if (fP) delete [] fP;
632 fP = new Float_t[3*fN];
633 memset(fP,0,3*fN*sizeof(Float_t));
634 fLastPoint = fN-1;
635}
636
637////////////////////////////////////////////////////////////////////////////////
638/// Re-initialize polyline with n points from p. If p=0 initialize with 0.
639/// if n <= 0 the current array of points is deleted.
640
642{
643 fOption = option;
644 if (n <= 0) {
645 fN = 0;
646 fLastPoint = -1;
647 delete [] fP;
648 fP = nullptr;
649 return;
650 }
651 fN = n;
652 if (fP) delete [] fP;
653 fP = new Float_t[3*fN];
654 if (p) {
655 for (Int_t i=0; i<fN;i++) {
656 fP[3*i] = p[3*i];
657 fP[3*i+1] = p[3*i+1];
658 fP[3*i+2] = p[3*i+2];
659 }
660 } else {
661 memset(fP,0,3*fN*sizeof(Float_t));
662 }
663 fLastPoint = fN-1;
664}
665
666////////////////////////////////////////////////////////////////////////////////
667/// Re-initialize polyline with n points from p. If p=0 initialize with 0.
668/// if n <= 0 the current array of points is deleted.
669
671{
672 fOption = option;
673 if (n <= 0) {
674 fN = 0;
675 fLastPoint = -1;
676 delete [] fP;
677 fP = nullptr;
678 return;
679 }
680 fN = n;
681 if (fP) delete [] fP;
682 fP = new Float_t[3*fN];
683 if (p) {
684 for (Int_t i=0; i<fN;i++) {
685 fP[3*i] = (Float_t) p[3*i];
686 fP[3*i+1] = (Float_t) p[3*i+1];
687 fP[3*i+2] = (Float_t) p[3*i+2];
688 }
689 } else {
690 memset(fP,0,3*fN*sizeof(Float_t));
691 }
692 fLastPoint = fN-1;
693}
694
695////////////////////////////////////////////////////////////////////////////////
696/// Stream a 3-D polyline object.
697
699{
700 UInt_t R__s, R__c;
701 if (b.IsReading()) {
702 b.ReadVersion(&R__s, &R__c);
703 b.ClassBegin(TPolyLine3D::IsA());
704 b.ClassMember("TObject");
706 b.ClassMember("TAttLine");
708 b.ClassMember("fN", "Int_t");
709 b >> fN;
710 if (fN) {
711 fP = new Float_t[3*fN];
712 b.ClassMember("fP", "Float_t", 3 * fN);
713 b.ReadFastArray(fP, 3 * fN);
714 }
715 b.ClassMember("fOption", "TString");
717 fLastPoint = fN-1;
718 b.ClassEnd(TPolyLine3D::IsA());
719 b.CheckByteCount(R__s, R__c, TPolyLine3D::IsA());
720 } else {
721 R__c = b.WriteVersion(TPolyLine3D::IsA(), kTRUE);
722 b.ClassBegin(TPolyLine3D::IsA());
723 b.ClassMember("TObject");
725 b.ClassMember("TAttLine");
727 b.ClassMember("fN", "Int_t");
728 Int_t size = Size();
729 b << size;
730 if (size) {
731 b.ClassMember("fP", "Float_t", 3 * size);
732 b.WriteFastArray(fP, 3 * size);
733 }
734 b.ClassMember("fOption", "TString");
736 b.ClassEnd(TPolyLine3D::IsA());
737 b.SetByteCount(R__c, kTRUE);
738 }
739}
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
float Float_t
Definition RtypesCore.h:57
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
double Double_t
Definition RtypesCore.h:59
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:382
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t np
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t TPoint TPoint const char y1
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 void Streamer(TBuffer &)
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:247
void Copy(TAttLine &attline) const
Copy this line attributes to a new TAttLine.
Definition TAttLine.cxx:177
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:211
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:275
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:114
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:113
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
TObject * After(const TObject *obj) const override
Returns the object after object obj.
Definition TList.cxx:328
void Add(TObject *obj) override
Definition TList.h:83
TObject * First() const override
Return the first object in the list. Returns 0 when list is empty.
Definition TList.cxx:657
Mother of all ROOT objects.
Definition TObject.h:41
virtual void Streamer(TBuffer &)
Stream an object of class TObject.
Definition TObject.cxx:894
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:213
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition TObject.cxx:190
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:786
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:530
virtual void Copy(TObject &object) const
Copy this to obj.
Definition TObject.cxx:146
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:993
@ 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.
void Paint(Option_t *option="") override
Paint a TPolyLine3D.
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 Int_t Merge(TCollection *list)
Merge polylines in the collection in this 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
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save primitive as a C++ statement(s) on output stream.
~TPolyLine3D() override
3-D polyline destructor.
Float_t * GetP() const
Definition TPolyLine3D.h:58
void ExecuteEvent(Int_t event, Int_t px, Int_t py) override
Execute action corresponding to one event.
virtual Int_t Size() const
Definition TPolyLine3D.h:71
TString fOption
options
Definition TPolyLine3D.h:37
void ls(Option_t *option="") const override
List this 3-D polyline.
void Draw(Option_t *option="") override
Draw this 3-D polyline with its current attributes.
virtual void SetPolyLine(Int_t n, Option_t *option="")
Re-initialize polyline with n points (0,0,0).
static TClass * Class()
void Copy(TObject &polyline) const override
Copy polyline to polyline obj.
void Streamer(TBuffer &) override
Stream a 3-D polyline object.
virtual void DrawPolyLine(Int_t n, Float_t *p, Option_t *option="")
Draw 3-D polyline with new coordinates.
void Print(Option_t *option="") const override
Dump this 3-D polyline with its attributes on stdout.
TClass * IsA() const override
Definition TPolyLine3D.h:75
Int_t DistancetoPrimitive(Int_t px, Int_t py) override
Compute distance from point px,py to a 3-D polyline.
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:2895
Basic string class.
Definition TString.h:139
void ToLower()
Change string to lower-case.
Definition TString.cxx:1182
virtual void Streamer(TBuffer &)
Stream a string object.
Definition TString.cxx:1412
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:632
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=nullptr)=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)
Returns the largest of a and b.
Definition TMathBase.h:250