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