Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TPolyMarker3D.cxx
Go to the documentation of this file.
1// @(#)root/g3d:$Id$
2// Author: Nenad Buncic 21/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 "TView.h"
13#include "TPolyMarker3D.h"
14#include "TVirtualPad.h"
15#include "TRandom.h"
16#include "TBuffer.h"
17#include "TBuffer3D.h"
18#include "TBuffer3DTypes.h"
19#include "TVirtualViewer3D.h"
20#include "TGeometry.h"
21#include "TH1.h"
22#include "TROOT.h"
23#include "TMath.h"
24
25#include <cassert>
26#include <iostream>
27
29
30constexpr Int_t kDimension = 3;
31
32/** \class TPolyMarker3D
33\ingroup g3d
34A 3D polymarker.
35
36It has three constructors.
37
38First one, without any parameters TPolyMarker3D(), we call 'default
39constructor' and it's used in a case that just an initialisation is
40needed (i.e. pointer declaration).
41
42Example:
43
44~~~ {.cpp}
45 TPolyMarker3D *pm = new TPolyMarker3D;
46~~~
47
48Second one, takes, usually, two parameters, n (number of points) and
49marker (marker style). Third parameter is optional.
50
51Example:
52
53~~~ {.cpp}
54 TPolyMarker3D (150, 1);
55~~~
56
57Third one takes, usually, three parameters, n (number of points), *p
58(pointer to an array of 3D points), and marker (marker style). Fourth
59parameter is optional.
60
61Example:
62
63~~~ {.cpp}
64 Float_t *ptr = new Float_t [150*3];
65 ... ... ...
66 ... ... ...
67 ... ... ...
68 TPolyMarker3D (150, ptr, 1);
69~~~
70*/
71
72////////////////////////////////////////////////////////////////////////////////
73/// 3-D polymarker default constructor.
74
76{
77 fName = "TPolyMarker3D";
78}
79
80////////////////////////////////////////////////////////////////////////////////
81/// 3-D polymarker normal constructor with initialization to 0.
82
84{
85 fName = "TPolyMarker3D";
87 SetMarkerStyle(marker);
89 if (n <= 0)
90 return;
91
92 fN = n;
93 fP = new Float_t [kDimension*fN];
94 for (Int_t i = 0; i < kDimension*fN; i++) fP[i] = 0;
95}
96
97////////////////////////////////////////////////////////////////////////////////
98/// 3-D polymarker constructor. Polymarker is initialized with p.
99
102{
103 fName = "TPolyMarker3D";
104 SetMarkerStyle(marker);
106 fOption = option;
107 if (n <= 0)
108 return;
109
110 fN = n;
111 fP = new Float_t [kDimension*fN];
112 if (p) {
113 for (Int_t i = 0; i < kDimension*fN; i++)
114 fP[i] = p[i];
115 fLastPoint = fN-1;
116 } else
117 memset(fP,0,kDimension*fN*sizeof(Float_t));
118}
119
120////////////////////////////////////////////////////////////////////////////////
121/// 3-D polymarker constructor. Polymarker is initialized with p
122/// (cast to float).
123
126{
127 fName = "TPolyMarker3D";
128 SetMarkerStyle(marker);
130 fOption = option;
131 if (n <= 0)
132 return;
133
134 fN = n;
135 fP = new Float_t [kDimension*fN];
136 if (p) {
137 for (Int_t i = 0; i < kDimension*fN; i++)
138 fP[i] = (Float_t) p[i];
139 fLastPoint = fN-1;
140 } else
141 memset(fP,0,kDimension*fN*sizeof(Float_t));
142}
143
144////////////////////////////////////////////////////////////////////////////////
145/// assignment operator
146
148{
149 if(this != &tp3)
150 tp3.TPolyMarker3D::Copy(*this);
151 return *this;
152}
153
154////////////////////////////////////////////////////////////////////////////////
155/// 3-D polymarker destructor.
156
158{
159 fN = 0;
160 if (fP) delete [] fP;
161 fLastPoint = -1;
162}
163
164////////////////////////////////////////////////////////////////////////////////
165/// 3-D polymarker copy ctor.
166
168{
169 p.TPolyMarker3D::Copy(*this);
170}
171
172////////////////////////////////////////////////////////////////////////////////
173/// Copy polymarker to polymarker obj.
174
176{
177 auto &tgt = static_cast<TPolyMarker3D &>(obj);
178 TObject::Copy(obj);
179 TAttMarker::Copy(tgt);
180 tgt.fN = fN;
181 if (tgt.fP)
182 delete [] tgt.fP;
183 if (fN > 0) {
184 tgt.fP = new Float_t [kDimension*fN];
185 for (Int_t i = 0; i < kDimension*fN; i++)
186 tgt.fP[i] = fP[i];
187 } else {
188 tgt.fP = nullptr;
189 }
190 tgt.fOption = fOption;
191 tgt.fLastPoint = fLastPoint;
192 tgt.fName = fName;
193}
194
195////////////////////////////////////////////////////////////////////////////////
196/// Compute distance from point px,py to a 3-D polymarker.
197/// Compute the closest distance of approach from point px,py to each segment
198/// of the polymarker.
199/// Returns when the distance found is below DistanceMaximum.
200/// The distance is computed in pixels units.
201
203{
204 const Int_t inaxis = 7;
205 Int_t dist = 9999;
206
207 Int_t puxmin = gPad->XtoAbsPixel(gPad->GetUxmin());
208 Int_t puymin = gPad->YtoAbsPixel(gPad->GetUymin());
209 Int_t puxmax = gPad->XtoAbsPixel(gPad->GetUxmax());
210 Int_t puymax = gPad->YtoAbsPixel(gPad->GetUymax());
211
212 // return if point is not in the user area
213 if (px < puxmin - inaxis) return dist;
214 if (py > puymin + inaxis) return dist;
215 if (px > puxmax + inaxis) return dist;
216 if (py < puymax - inaxis) return dist;
217
218 TView *view = gPad->GetView();
219 if (!view) return dist;
220 Int_t i, dpoint;
221 Float_t xndc[3];
222 Int_t x1,y1;
223 Double_t u,v;
224 for (i=0;i<Size();i++) {
225 view->WCtoNDC(&fP[3*i], xndc);
226 u = (Double_t)xndc[0];
227 v = (Double_t)xndc[1];
228 if (u < gPad->GetUxmin() || u > gPad->GetUxmax()) continue;
229 if (v < gPad->GetUymin() || v > gPad->GetUymax()) continue;
230 x1 = gPad->XtoAbsPixel(u);
231 y1 = gPad->YtoAbsPixel(v);
232 dpoint = Int_t(TMath::Sqrt((((Double_t)px-x1)*((Double_t)px-x1)
233 + ((Double_t)py-y1)*((Double_t)py-y1))));
234 if (dpoint < dist) dist = dpoint;
235 }
236 return dist;
237}
238
239////////////////////////////////////////////////////////////////////////////////
240/// Draws 3-D polymarker with its current attributes.
241
243{
245}
246
247////////////////////////////////////////////////////////////////////////////////
248/// Draw this 3-D polymarker with new coordinates. Creates a new
249/// polymarker which will be adopted by the pad in which it is drawn.
250/// Does not change the original polymarker (should be static method).
251
253{
254 TPolyMarker3D *newpolymarker = new TPolyMarker3D();
255 newpolymarker->fN = n;
256 newpolymarker->fP = new Float_t [kDimension*fN];
257 for (Int_t i = 0; i < kDimension*fN; i++) newpolymarker->fP[i] = p[i];
258 newpolymarker->SetMarkerStyle(GetMarkerStyle());
259 newpolymarker->fOption = fOption;
260 newpolymarker->fLastPoint = fLastPoint;
261 newpolymarker->SetBit(kCanDelete);
262 newpolymarker->AppendPad(option);
263}
264
265////////////////////////////////////////////////////////////////////////////////
266/// Execute action corresponding to one event.
267
269{
270 if (!gPad) return;
271 if (gPad->GetView()) gPad->GetView()->ExecuteRotateView(event, px, py);
272}
273
274////////////////////////////////////////////////////////////////////////////////
275/// List this 3-D polymarker.
276
278{
280 std::cout << " TPolyMarker3D N=" << Size() <<" Option="<<option<<std::endl;
281}
282
283////////////////////////////////////////////////////////////////////////////////
284/// Merge polymarkers in the collection in this polymarker
285
287{
288 if (!li) return 0;
289 TIter next(li);
290
291 //first loop to count the number of entries
292 TPolyMarker3D *pm;
293 Int_t npoints = Size();
294 while ((pm = (TPolyMarker3D*)next())) {
296 Error("Add","Attempt to add object of class: %s to a %s",pm->ClassName(),this->ClassName());
297 return -1;
298 }
299 npoints += pm->Size();
300 }
301 Int_t currPoint = Size();
302
303 //extend this polymarker to hold npoints
304 SetPoint(npoints-1,0,0,0);
305
306 //merge all polymarkers
307 next.Reset();
308 while ((pm = (TPolyMarker3D*)next())) {
309 Int_t np = pm->Size();
310 Float_t *p = pm->GetP();
311 for (Int_t i = 0; i < np; i++) {
312 SetPoint(currPoint++, p[3*i], p[3*i+1], p[3*i+2]);
313 }
314 }
315 return npoints;
316}
317
318////////////////////////////////////////////////////////////////////////////////
319/// Paint a TPolyMarker3D.
320
322{
323 // No need to continue if there is nothing to paint
324 if (Size() <= 0) return;
325
326 static TBuffer3D buffer(TBuffer3DTypes::kMarker);
327
328 buffer.ClearSectionsValid();
329
330 // Section kCore
331 buffer.fID = this;
332 buffer.fColor = GetMarkerColor();
333 buffer.fTransparency = 0;
334 buffer.fLocalFrame = kFALSE;
336
337 // We fill kCore and kRawSizes on first pass and try with viewer
338 TVirtualViewer3D * viewer3D = gPad->GetViewer3D();
339 if (!viewer3D) return;
340 Int_t reqSections = viewer3D->AddObject(buffer);
341 if (reqSections == TBuffer3D::kNone) {
342 return;
343 }
344
345 if (reqSections & TBuffer3D::kRawSizes) {
346 if (!buffer.SetRawSizes(Size(), 3*Size(), 1, 1, 0, 0)) {
347 return;
348 }
350 }
351
352 if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
353 // Points
354 for (UInt_t i=0; i<3*buffer.NbPnts(); i++) {
355 buffer.fPnts[i] = (Double_t)fP[i];
356 }
357
358 // Transform points - we don't support local->global matrix
359 // so always work in global reference frame
360 if (gGeometry) {
361 Double_t dlocal[3];
362 Double_t dmaster[3];
363 for (UInt_t j=0; j<buffer.NbPnts(); j++) {
364 dlocal[0] = buffer.fPnts[3*j];
365 dlocal[1] = buffer.fPnts[3*j+1];
366 dlocal[2] = buffer.fPnts[3*j+2];
367 gGeometry->Local2Master(&dlocal[0],&dmaster[0]);
368 buffer.fPnts[3*j] = dmaster[0];
369 buffer.fPnts[3*j+1] = dmaster[1];
370 buffer.fPnts[3*j+2] = dmaster[2];
371 }
372 }
373
374 // Basic colors: 0, 1, ... 7
375 Int_t c = (((GetMarkerColor()) %8) -1) * 4;
376 if (c < 0) c = 0;
377
378 // Segments
379 buffer.fSegs[0] = c;
380
382
384 }
385
386 viewer3D->AddObject(buffer);
387}
388
389////////////////////////////////////////////////////////////////////////////////
390/// Paint 3-d histogram h with 3-d polymarkers.
391
393{
394 const Int_t kMaxEntry = 100000;
395 Int_t in, bin, binx, biny, binz;
396
397 TAxis *xaxis = h->GetXaxis();
398 TAxis *yaxis = h->GetYaxis();
399 TAxis *zaxis = h->GetZaxis();
400 Double_t entry = 0;
401 for (binz=zaxis->GetFirst();binz<=zaxis->GetLast();binz++) {
402 for (biny=yaxis->GetFirst();biny<=yaxis->GetLast();biny++) {
403 for (binx=xaxis->GetFirst();binx<=xaxis->GetLast();binx++) {
404 bin = h->GetBin(binx,biny,binz);
405 entry += h->GetBinContent(bin);
406 }
407 }
408 }
409
410 // if histogram has too many entries, rescale it
411 // never draw more than kMaxEntry markers, otherwise this kills
412 // the X server
413 Double_t scale = 1.;
414 if (entry > kMaxEntry) scale = kMaxEntry/Double_t(entry);
415
416 //Create or modify 3-d view object
417 TView *view = gPad->GetView();
418 if (!view) {
419 gPad->Range(-1,-1,1,1);
420 view = TView::CreateView(1,nullptr,nullptr);
421 if (!view) return;
422 }
423 view->SetRange(xaxis->GetBinLowEdge(xaxis->GetFirst()),
424 yaxis->GetBinLowEdge(yaxis->GetFirst()),
425 zaxis->GetBinLowEdge(zaxis->GetFirst()),
426 xaxis->GetBinUpEdge(xaxis->GetLast()),
427 yaxis->GetBinUpEdge(yaxis->GetLast()),
428 zaxis->GetBinUpEdge(zaxis->GetLast()));
429
430 view->PadRange(gPad->GetFrameFillColor());
431
432 if (entry == 0) return;
433 Int_t nmk = Int_t(TMath::Min(Double_t(kMaxEntry),entry));
434 TPolyMarker3D *pm3d = new TPolyMarker3D(nmk);
435 pm3d->SetMarkerStyle(h->GetMarkerStyle());
436 pm3d->SetMarkerColor(h->GetMarkerColor());
437 pm3d->SetMarkerSize(h->GetMarkerSize());
438 gPad->Modified(kTRUE);
439
440 entry = 0;
441 Double_t x,y,z,xw,yw,zw,xp,yp,zp;
442 Int_t ncounts;
443 for (binz=zaxis->GetFirst();binz<=zaxis->GetLast();binz++) {
444 z = zaxis->GetBinLowEdge(binz);
445 zw = zaxis->GetBinWidth(binz);
446 for (biny=yaxis->GetFirst();biny<=yaxis->GetLast();biny++) {
447 y = yaxis->GetBinLowEdge(biny);
448 yw = yaxis->GetBinWidth(biny);
449 for (binx=xaxis->GetFirst();binx<=xaxis->GetLast();binx++) {
450 x = xaxis->GetBinLowEdge(binx);
451 xw = xaxis->GetBinWidth(binx);
452 bin = h->GetBin(binx,biny,binz);
453 ncounts = Int_t(h->GetBinContent(bin)*scale+0.5);
454 for (in=0;in<ncounts;in++) {
455 xp = x + xw*gRandom->Rndm();
456 yp = y + yw*gRandom->Rndm();
457 zp = z + zw*gRandom->Rndm();
458 pm3d->SetPoint(Int_t(entry),xp,yp,zp);
459 entry++;
460 }
461 }
462 }
463 }
464 pm3d->Paint(option);
465 delete pm3d;
466}
467
468////////////////////////////////////////////////////////////////////////////////
469/// Print 3-D polymarker with its attributes on stdout.
470
472{
473 printf("TPolyMarker3D N=%d, Option=%s\n",fN,option);
474 TString opt = option;
475 opt.ToLower();
476 if (opt.Contains("all")) {
477 for (Int_t i=0;i<Size();i++) {
479 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]);
480 }
481 }
482}
483
484////////////////////////////////////////////////////////////////////////////////
485/// Save primitive as a C++ statement(s) on output stream.
486
487void TPolyMarker3D::SavePrimitive(std::ostream &out, Option_t * /*= ""*/)
488{
489 char quote = '"';
490 out<<" "<<std::endl;
491 if (gROOT->ClassSaved(TPolyMarker3D::Class())) {
492 out<<" ";
493 } else {
494 out<<" TPolyMarker3D *";
495 }
496 out<<"pmarker3D = new TPolyMarker3D("<<fN<<","<<GetMarkerStyle()<<","<<quote<<fOption<<quote<<");"<<std::endl;
497 out<<" pmarker3D->SetName("<<quote<<GetName()<<quote<<");"<<std::endl;
498
499 SaveMarkerAttributes(out,"pmarker3D",1,1,1);
500
501 for (Int_t i=0;i<Size();i++) {
502 out<<" pmarker3D->SetPoint("<<i<<","<<fP[3*i]<<","<<fP[3*i+1]<<","<<fP[3*i+2]<<");"<<std::endl;
503 }
504 out<<" pmarker3D->Draw();"<<std::endl;
505}
506
507////////////////////////////////////////////////////////////////////////////////
508/// Change (i.e. set) the name of the TNamed.
509/// WARNING: if the object is a member of a THashTable or THashList container
510/// the container must be Rehash()'ed after SetName(). For example the list
511/// of objects in the current directory is a THashList.
512
514{
515 fName = name;
516 if (gPad && TestBit(kMustCleanup)) gPad->Modified();
517}
518
519////////////////////////////////////////////////////////////////////////////////
520/// Set point following LastPoint to x, y, z.
521/// Returns index of the point (new last point).
522
524{
525 fLastPoint++;
526 SetPoint(fLastPoint, x, y, z);
527 return fLastPoint;
528}
529
530////////////////////////////////////////////////////////////////////////////////
531/// Set point n to x, y, z.
532/// If n is more then the current TPolyMarker3D size (n > fN) then
533/// the polymarker will be resized to contain at least n points.
534
536{
537 if (n < 0) return;
538 if (!fP || n >= fN) {
539 // re-allocate the object
540 Int_t newN = TMath::Max(2*fN,n+1);
541 Float_t *savepoint = new Float_t [kDimension*newN];
542 if (fP && fN){
543 memcpy(savepoint,fP,kDimension*fN*sizeof(Float_t));
544 memset(&savepoint[kDimension*fN],0,(newN-fN)*sizeof(Float_t));
545 delete [] fP;
546 }
547 fP = savepoint;
548 fN = newN;
549 }
550 fP[kDimension*n ] = x;
551 fP[kDimension*n+1] = y;
552 fP[kDimension*n+2] = z;
554}
555
556////////////////////////////////////////////////////////////////////////////////
557/// Re-initialize polymarker with n points from p. If p=0 initialize with 0.
558/// if n <= 0 the current array of points is deleted.
559
561{
562 SetMarkerStyle(marker);
563 fOption = option;
564 if (n <= 0) {
565 fN = 0;
566 fLastPoint = -1;
567 delete [] fP;
568 fP = nullptr;
569 return;
570 }
571 fN = n;
572 if (fP) delete [] fP;
573 fP = new Float_t [3*fN];
574 if (p) {
575 for (Int_t i = 0; i < fN; i++) {
576 fP[3*i] = p[3*i];
577 fP[3*i+1] = p[3*i+1];
578 fP[3*i+2] = p[3*i+2];
579 }
580 } else {
581 memset(fP,0,kDimension*fN*sizeof(Float_t));
582 }
583 fLastPoint = fN-1;
584}
585
586////////////////////////////////////////////////////////////////////////////////
587/// Re-initialize polymarker with n points from p. If p=0 initialize with 0.
588/// if n <= 0 the current array of points is deleted.
589
591{
592 SetMarkerStyle(marker);
593 fOption = option;
594 if (n <= 0) {
595 fN = 0;
596 fLastPoint = -1;
597 delete [] fP;
598 fP = nullptr;
599 return;
600 }
601 fN = n;
602 if (fP) delete [] fP;
603 fP = new Float_t [3*fN];
604 if (p) {
605 for (Int_t i = 0; i < fN; i++) {
606 fP[3*i] = (Float_t) p[3*i];
607 fP[3*i+1] = (Float_t) p[3*i+1];
608 fP[3*i+2] = (Float_t) p[3*i+2];
609 }
610 } else {
611 memset(fP,0,kDimension*fN*sizeof(Float_t));
612 }
613 fLastPoint = fN-1;
614}
615
616////////////////////////////////////////////////////////////////////////////////
617/// Stream a 3-D polymarker object.
618
620{
621 UInt_t R__s, R__c;
622 if (b.IsReading()) {
623 Version_t R__v = b.ReadVersion(&R__s, &R__c);
624 if (R__v > 2) b.ClassBegin(TPolyMarker3D::IsA());
625 if (R__v > 2) b.ClassMember("TObject");
627 if (R__v > 2) b.ClassMember("TAttMarker");
629 if (R__v > 2) b.ClassMember("fN","Int_t");
630 b >> fN;
631 if (fN) {
632 if (R__v > 2) b.ClassMember("fP","Float_t", kDimension*fN);
633 fP = new Float_t[kDimension*fN];
634 b.ReadFastArray(fP,kDimension*fN);
635 }
636 fLastPoint = fN-1;
637 if (R__v > 2) b.ClassMember("fOption","TString");
639 if (R__v > 2) b.ClassMember("fName","TString");
640 if (R__v > 1) fName.Streamer(b);
641 if (R__v > 2) b.ClassEnd(TPolyMarker3D::IsA());
642 b.CheckByteCount(R__s, R__c, TPolyMarker3D::IsA());
643 } else {
644 R__c = b.WriteVersion(TPolyMarker3D::IsA(), kTRUE);
645 b.ClassBegin(TPolyMarker3D::IsA());
646 b.ClassMember("TObject");
648 b.ClassMember("TAttMarker");
650 b.ClassMember("fN","Int_t");
651 Int_t size = Size();
652 b << size;
653 if (size) {
654 b.ClassMember("fP","Float_t", kDimension*size);
655 b.WriteFastArray(fP, kDimension*size);
656 }
657 b.ClassMember("fOption","TString");
659 b.ClassMember("fName","TString");
661 b.ClassEnd(TPolyMarker3D::IsA());
662 b.SetByteCount(R__c, kTRUE);
663 }
664}
665
666////////////////////////////////////////////////////////////////////////////////
667/// Fills the parameters x, y, z with the coordinate of the n-th point
668/// n must be between 0 and Size() - 1.
669
671{
672 if (n < 0 || n >= Size()) return;
673 if (!fP) return;
674 x = fP[kDimension*n ];
675 y = fP[kDimension*n+1];
676 z = fP[kDimension*n+2];
677}
678
679////////////////////////////////////////////////////////////////////////////////
680/// Fills the parameters x, y, z with the coordinate of the n-th point
681/// n must be between 0 and Size() - 1.
682
684{
685 if (n < 0 || n >= Size()) return;
686 if (!fP) return;
687 x = (Double_t)fP[kDimension*n ];
688 y = (Double_t)fP[kDimension*n+1];
689 z = (Double_t)fP[kDimension*n+2];
690}
691
692
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define h(i)
Definition RSha256.hxx:106
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
int Int_t
Definition RtypesCore.h:45
short Version_t
Definition RtypesCore.h:65
short Marker_t
Definition RtypesCore.h:90
float Float_t
Definition RtypesCore.h:57
constexpr Bool_t kFALSE
Definition RtypesCore.h:101
double Double_t
Definition RtypesCore.h:59
constexpr Bool_t kTRUE
Definition RtypesCore.h:100
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:377
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 x1
Option_t Option_t SetMarkerStyle
Option_t Option_t TPoint TPoint const char y1
char name[80]
Definition TGX11.cxx:110
R__EXTERN TGeometry * gGeometry
Definition TGeometry.h:158
constexpr Int_t kDimension
#define gROOT
Definition TROOT.h:406
R__EXTERN TRandom * gRandom
Definition TRandom.h:62
#define gPad
int currPoint
Definition X3DBuffer.c:17
Use this attribute class when an object should have 3D capabilities.
Definition TAtt3D.h:19
Marker Attributes class.
Definition TAttMarker.h:19
virtual void SaveMarkerAttributes(std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1, Int_t sizdef=1)
Save line attributes as C++ statement(s) on output stream out.
virtual void Modify()
Change current marker attributes if necessary.
virtual Style_t GetMarkerStyle() const
Return the marker style.
Definition TAttMarker.h:32
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition TAttMarker.h:38
virtual Color_t GetMarkerColor() const
Return the marker color.
Definition TAttMarker.h:31
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition TAttMarker.h:40
void Copy(TAttMarker &attmarker) const
Copy this marker attributes to a new TAttMarker.
virtual void Streamer(TBuffer &)
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition TAttMarker.h:45
Class to manage histogram axis.
Definition TAxis.h:31
virtual Double_t GetBinLowEdge(Int_t bin) const
Return low edge of bin.
Definition TAxis.cxx:518
Int_t GetLast() const
Return last bin on the axis i.e.
Definition TAxis.cxx:469
virtual Double_t GetBinWidth(Int_t bin) const
Return bin width.
Definition TAxis.cxx:540
virtual Double_t GetBinUpEdge(Int_t bin) const
Return up edge of bin.
Definition TAxis.cxx:528
Int_t GetFirst() const
Return first bin on the axis i.e.
Definition TAxis.cxx:458
Generic 3D primitive description class.
Definition TBuffer3D.h:18
UInt_t NbPnts() const
Definition TBuffer3D.h:80
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.
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:59
void Reset()
Mother of all ROOT objects.
Definition TObject.h:41
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:199
virtual void Streamer(TBuffer &)
Stream an object of class TObject.
Definition TObject.cxx:888
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:207
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition TObject.cxx:184
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:780
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:525
virtual void Copy(TObject &object) const
Copy this to obj.
Definition TObject.cxx:140
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:987
@ kCanDelete
if object in a list can be deleted
Definition TObject.h:62
@ kMustCleanup
if object destructor must call RecursiveRemove()
Definition TObject.h:64
A 3D polymarker.
static void PaintH3(TH1 *h, Option_t *option)
Paint 3-d histogram h with 3-d polymarkers.
virtual Int_t Merge(TCollection *list)
Merge polymarkers in the collection in this polymarker.
void ExecuteEvent(Int_t event, Int_t px, Int_t py) override
Execute action corresponding to one event.
static TClass * Class()
virtual void DrawPolyMarker(Int_t n, Float_t *p, Marker_t marker, Option_t *option="")
Draw this 3-D polymarker with new coordinates.
~TPolyMarker3D() override
3-D polymarker destructor.
void SetPoint(Int_t n, Double_t x, Double_t y, Double_t z)
Set point n to x, y, z.
const char * GetName() const override
Returns name of object.
virtual void GetPoint(Int_t n, Float_t &x, Float_t &y, Float_t &z) const
Fills the parameters x, y, z with the coordinate of the n-th point n must be between 0 and Size() - 1...
virtual void SetPolyMarker(Int_t n, Float_t *p, Marker_t marker, Option_t *option="")
Re-initialize polymarker with n points from p.
virtual Int_t Size() const
virtual void SetName(const char *name)
Change (i.e.
TPolyMarker3D & operator=(const TPolyMarker3D &)
assignment operator
void Copy(TObject &polymarker) const override
Copy polymarker to polymarker obj.
void Streamer(TBuffer &) override
Stream a 3-D polymarker object.
void ls(Option_t *option="") const override
List this 3-D polymarker.
virtual Float_t * GetP() const
virtual Int_t SetNextPoint(Double_t x, Double_t y, Double_t z)
Set point following LastPoint to x, y, z.
TClass * IsA() const override
void Paint(Option_t *option="") override
Paint a TPolyMarker3D.
void Print(Option_t *option="") const override
Print 3-D polymarker with its attributes on stdout.
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save primitive as a C++ statement(s) on output stream.
Int_t DistancetoPrimitive(Int_t px, Int_t py) override
Compute distance from point px,py to a 3-D polymarker.
void Draw(Option_t *option="") override
Draws 3-D polymarker with its current attributes.
TPolyMarker3D()
3-D polymarker default constructor.
static void IndentLevel()
Functions used by ls() to indent an object hierarchy.
Definition TROOT.cxx:2891
Double_t Rndm() override
Machine independent random number generator.
Definition TRandom.cxx:559
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
static TView * CreateView(Int_t system=1, const Double_t *rmin=nullptr, const Double_t *rmax=nullptr)
Create a concrete default 3-d view via the plug-in manager.
Definition TView.cxx:27
virtual void PadRange(Int_t rback)=0
virtual void SetRange(const Double_t *min, const Double_t *max)=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
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:662
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:198