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 fN = 0;
78 fP = 0;
79 fLastPoint = -1;
80 fName = "TPolyMarker3D";
81}
82
83////////////////////////////////////////////////////////////////////////////////
84/// 3-D polymarker normal constructor with initialization to 0.
85
87{
88 fName = "TPolyMarker3D";
89 fOption = option;
90 SetMarkerStyle(marker);
92 fLastPoint = -1;
93 if (n <= 0) {
94 fN = 0;
95 fP = 0;
96 return;
97 }
98
99 fN = n;
100 fP = new Float_t [kDimension*fN];
101 for (Int_t i = 0; i < kDimension*fN; i++) fP[i] = 0;
102}
103
104////////////////////////////////////////////////////////////////////////////////
105/// 3-D polymarker constructor. Polymarker is initialized with p.
106
108 Option_t *option)
109{
110 fName = "TPolyMarker3D";
111 SetMarkerStyle(marker);
113 fOption = option;
114 fLastPoint = -1;
115 if (n <= 0) {
116 fN = 0;
117 fP = 0;
118 return;
119 }
120
121 fN = n;
122 fP = new Float_t [kDimension*fN];
123 if (p) {
124 for (Int_t i = 0; i < kDimension*fN; i++)
125 fP[i] = p[i];
126 fLastPoint = fN-1;
127 } else
128 memset(fP,0,kDimension*fN*sizeof(Float_t));
129}
130
131////////////////////////////////////////////////////////////////////////////////
132/// 3-D polymarker constructor. Polymarker is initialized with p
133/// (cast to float).
134
136 Option_t *option)
137{
138 fName = "TPolyMarker3D";
139 SetMarkerStyle(marker);
141 fOption = option;
142 fLastPoint = -1;
143 if (n <= 0) {
144 fN = 0;
145 fP = 0;
146 return;
147 }
148
149 fN = n;
150 fP = new Float_t [kDimension*fN];
151 if (p) {
152 for (Int_t i = 0; i < kDimension*fN; i++)
153 fP[i] = (Float_t) p[i];
154 fLastPoint = fN-1;
155 } else
156 memset(fP,0,kDimension*fN*sizeof(Float_t));
157}
158
159////////////////////////////////////////////////////////////////////////////////
160/// assignment operator
161
163{
164 if(this!=&tp3) {
166 TAttMarker::operator=(tp3);
167 TAtt3D::operator=(tp3);
168 fN=tp3.fN;
169 fP=tp3.fP;
170 fOption=tp3.fOption;
172 fName=tp3.fName;
173 }
174 return *this;
175}
176
177////////////////////////////////////////////////////////////////////////////////
178/// 3-D polymarker destructor.
179
181{
182 fN = 0;
183 if (fP) delete [] fP;
184 fLastPoint = -1;
185}
186
187////////////////////////////////////////////////////////////////////////////////
188/// 3-D polymarker copy ctor.
189
191 TObject(p), TAttMarker(p), TAtt3D(p)
192{
193 fN = 0;
194 fP = 0;
195 fLastPoint = -1;
196 p.Copy(*this);
197}
198
199////////////////////////////////////////////////////////////////////////////////
200/// Copy polymarker to polymarker obj.
201
203{
204 TObject::Copy(obj);
205 ((TPolyMarker3D&)obj).fN = fN;
206 if (fN > 0) {
207 ((TPolyMarker3D&)obj).fP = new Float_t [kDimension*fN];
208 for (Int_t i = 0; i < kDimension*fN; i++) ((TPolyMarker3D&)obj).fP[i] = fP[i];
209 } else {
210 ((TPolyMarker3D&)obj).fP = 0;
211 }
212 ((TPolyMarker3D&)obj).SetMarkerStyle(GetMarkerStyle());
213 ((TPolyMarker3D&)obj).fOption = fOption;
214 ((TPolyMarker3D&)obj).fLastPoint = fLastPoint;
215 ((TPolyMarker3D&)obj).fName = fName;
216}
217
218////////////////////////////////////////////////////////////////////////////////
219/// Compute distance from point px,py to a 3-D polymarker.
220/// Compute the closest distance of approach from point px,py to each segment
221/// of the polymarker.
222/// Returns when the distance found is below DistanceMaximum.
223/// The distance is computed in pixels units.
224
226{
227 const Int_t inaxis = 7;
228 Int_t dist = 9999;
229
230 Int_t puxmin = gPad->XtoAbsPixel(gPad->GetUxmin());
231 Int_t puymin = gPad->YtoAbsPixel(gPad->GetUymin());
232 Int_t puxmax = gPad->XtoAbsPixel(gPad->GetUxmax());
233 Int_t puymax = gPad->YtoAbsPixel(gPad->GetUymax());
234
235 // return if point is not in the user area
236 if (px < puxmin - inaxis) return dist;
237 if (py > puymin + inaxis) return dist;
238 if (px > puxmax + inaxis) return dist;
239 if (py < puymax - inaxis) return dist;
240
241 TView *view = gPad->GetView();
242 if (!view) return dist;
243 Int_t i, dpoint;
244 Float_t xndc[3];
245 Int_t x1,y1;
246 Double_t u,v;
247 for (i=0;i<Size();i++) {
248 view->WCtoNDC(&fP[3*i], xndc);
249 u = (Double_t)xndc[0];
250 v = (Double_t)xndc[1];
251 if (u < gPad->GetUxmin() || u > gPad->GetUxmax()) continue;
252 if (v < gPad->GetUymin() || v > gPad->GetUymax()) continue;
253 x1 = gPad->XtoAbsPixel(u);
254 y1 = gPad->YtoAbsPixel(v);
255 dpoint = Int_t(TMath::Sqrt((((Double_t)px-x1)*((Double_t)px-x1)
256 + ((Double_t)py-y1)*((Double_t)py-y1))));
257 if (dpoint < dist) dist = dpoint;
258 }
259 return dist;
260}
261
262////////////////////////////////////////////////////////////////////////////////
263/// Draws 3-D polymarker with its current attributes.
264
266{
267 AppendPad(option);
268}
269
270////////////////////////////////////////////////////////////////////////////////
271/// Draw this 3-D polymarker with new coordinates. Creates a new
272/// polymarker which will be adopted by the pad in which it is drawn.
273/// Does not change the original polymarker (should be static method).
274
276{
277 TPolyMarker3D *newpolymarker = new TPolyMarker3D();
278 newpolymarker->fN = n;
279 newpolymarker->fP = new Float_t [kDimension*fN];
280 for (Int_t i = 0; i < kDimension*fN; i++) newpolymarker->fP[i] = p[i];
281 newpolymarker->SetMarkerStyle(GetMarkerStyle());
282 newpolymarker->fOption = fOption;
283 newpolymarker->fLastPoint = fLastPoint;
284 newpolymarker->SetBit(kCanDelete);
285 newpolymarker->AppendPad(option);
286}
287
288////////////////////////////////////////////////////////////////////////////////
289/// Execute action corresponding to one event.
290
292{
293 if (!gPad) return;
294 if (gPad->GetView()) gPad->GetView()->ExecuteRotateView(event, px, py);
295}
296
297////////////////////////////////////////////////////////////////////////////////
298/// List this 3-D polymarker.
299
300void TPolyMarker3D::ls(Option_t *option) const
301{
303 std::cout << " TPolyMarker3D N=" << Size() <<" Option="<<option<<std::endl;
304}
305
306////////////////////////////////////////////////////////////////////////////////
307/// Merge polymarkers in the collection in this polymarker
308
310{
311 if (!li) return 0;
312 TIter next(li);
313
314 //first loop to count the number of entries
315 TPolyMarker3D *pm;
316 Int_t npoints = Size();
317 while ((pm = (TPolyMarker3D*)next())) {
318 if (!pm->InheritsFrom(TPolyMarker3D::Class())) {
319 Error("Add","Attempt to add object of class: %s to a %s",pm->ClassName(),this->ClassName());
320 return -1;
321 }
322 npoints += pm->Size();
323 }
324 Int_t currPoint = Size();
325
326 //extend this polymarker to hold npoints
327 SetPoint(npoints-1,0,0,0);
328
329 //merge all polymarkers
330 next.Reset();
331 while ((pm = (TPolyMarker3D*)next())) {
332 Int_t np = pm->Size();
333 Float_t *p = pm->GetP();
334 for (Int_t i = 0; i < np; i++) {
335 SetPoint(currPoint++, p[3*i], p[3*i+1], p[3*i+2]);
336 }
337 }
338 return npoints;
339}
340
341////////////////////////////////////////////////////////////////////////////////
342/// Paint a TPolyMarker3D.
343
345{
346 // No need to continue if there is nothing to paint
347 if (Size() <= 0) return;
348
349 static TBuffer3D buffer(TBuffer3DTypes::kMarker);
350
351 buffer.ClearSectionsValid();
352
353 // Section kCore
354 buffer.fID = this;
355 buffer.fColor = GetMarkerColor();
356 buffer.fTransparency = 0;
357 buffer.fLocalFrame = kFALSE;
359
360 // We fill kCore and kRawSizes on first pass and try with viewer
361 TVirtualViewer3D * viewer3D = gPad->GetViewer3D();
362 if (!viewer3D) return;
363 Int_t reqSections = viewer3D->AddObject(buffer);
364 if (reqSections == TBuffer3D::kNone) {
365 return;
366 }
367
368 if (reqSections & TBuffer3D::kRawSizes) {
369 if (!buffer.SetRawSizes(Size(), 3*Size(), 1, 1, 0, 0)) {
370 return;
371 }
373 }
374
375 if ((reqSections & TBuffer3D::kRaw) && buffer.SectionsValid(TBuffer3D::kRawSizes)) {
376 // Points
377 for (UInt_t i=0; i<3*buffer.NbPnts(); i++) {
378 buffer.fPnts[i] = (Double_t)fP[i];
379 }
380
381 // Transform points - we don't support local->global matrix
382 // so always work in global reference frame
383 if (gGeometry) {
384 Double_t dlocal[3];
385 Double_t dmaster[3];
386 for (UInt_t j=0; j<buffer.NbPnts(); j++) {
387 dlocal[0] = buffer.fPnts[3*j];
388 dlocal[1] = buffer.fPnts[3*j+1];
389 dlocal[2] = buffer.fPnts[3*j+2];
390 gGeometry->Local2Master(&dlocal[0],&dmaster[0]);
391 buffer.fPnts[3*j] = dmaster[0];
392 buffer.fPnts[3*j+1] = dmaster[1];
393 buffer.fPnts[3*j+2] = dmaster[2];
394 }
395 }
396
397 // Basic colors: 0, 1, ... 7
398 Int_t c = (((GetMarkerColor()) %8) -1) * 4;
399 if (c < 0) c = 0;
400
401 // Segments
402 buffer.fSegs[0] = c;
403
405
407 }
408
409 viewer3D->AddObject(buffer);
410}
411
412////////////////////////////////////////////////////////////////////////////////
413/// Paint 3-d histogram h with 3-d polymarkers.
414
416{
417 const Int_t kMaxEntry = 100000;
418 Int_t in, bin, binx, biny, binz;
419
420 TAxis *xaxis = h->GetXaxis();
421 TAxis *yaxis = h->GetYaxis();
422 TAxis *zaxis = h->GetZaxis();
423 Double_t entry = 0;
424 for (binz=zaxis->GetFirst();binz<=zaxis->GetLast();binz++) {
425 for (biny=yaxis->GetFirst();biny<=yaxis->GetLast();biny++) {
426 for (binx=xaxis->GetFirst();binx<=xaxis->GetLast();binx++) {
427 bin = h->GetBin(binx,biny,binz);
428 entry += h->GetBinContent(bin);
429 }
430 }
431 }
432
433 // if histogram has too many entries, rescale it
434 // never draw more than kMaxEntry markers, otherwise this kills
435 // the X server
436 Double_t scale = 1.;
437 if (entry > kMaxEntry) scale = kMaxEntry/Double_t(entry);
438
439 //Create or modify 3-d view object
440 TView *view = gPad->GetView();
441 if (!view) {
442 gPad->Range(-1,-1,1,1);
443 view = TView::CreateView(1,0,0);
444 if (!view) return;
445 }
446 view->SetRange(xaxis->GetBinLowEdge(xaxis->GetFirst()),
447 yaxis->GetBinLowEdge(yaxis->GetFirst()),
448 zaxis->GetBinLowEdge(zaxis->GetFirst()),
449 xaxis->GetBinUpEdge(xaxis->GetLast()),
450 yaxis->GetBinUpEdge(yaxis->GetLast()),
451 zaxis->GetBinUpEdge(zaxis->GetLast()));
452
453 view->PadRange(gPad->GetFrameFillColor());
454
455 if (entry == 0) return;
456 Int_t nmk = Int_t(TMath::Min(Double_t(kMaxEntry),entry));
457 TPolyMarker3D *pm3d = new TPolyMarker3D(nmk);
458 pm3d->SetMarkerStyle(h->GetMarkerStyle());
459 pm3d->SetMarkerColor(h->GetMarkerColor());
460 pm3d->SetMarkerSize(h->GetMarkerSize());
461 gPad->Modified(kTRUE);
462
463 entry = 0;
464 Double_t x,y,z,xw,yw,zw,xp,yp,zp;
465 Int_t ncounts;
466 for (binz=zaxis->GetFirst();binz<=zaxis->GetLast();binz++) {
467 z = zaxis->GetBinLowEdge(binz);
468 zw = zaxis->GetBinWidth(binz);
469 for (biny=yaxis->GetFirst();biny<=yaxis->GetLast();biny++) {
470 y = yaxis->GetBinLowEdge(biny);
471 yw = yaxis->GetBinWidth(biny);
472 for (binx=xaxis->GetFirst();binx<=xaxis->GetLast();binx++) {
473 x = xaxis->GetBinLowEdge(binx);
474 xw = xaxis->GetBinWidth(binx);
475 bin = h->GetBin(binx,biny,binz);
476 ncounts = Int_t(h->GetBinContent(bin)*scale+0.5);
477 for (in=0;in<ncounts;in++) {
478 xp = x + xw*gRandom->Rndm();
479 yp = y + yw*gRandom->Rndm();
480 zp = z + zw*gRandom->Rndm();
481 pm3d->SetPoint(Int_t(entry),xp,yp,zp);
482 entry++;
483 }
484 }
485 }
486 }
487 pm3d->Paint(option);
488 delete pm3d;
489}
490
491////////////////////////////////////////////////////////////////////////////////
492/// Print 3-D polymarker with its attributes on stdout.
493
495{
496 printf("TPolyMarker3D N=%d, Option=%s\n",fN,option);
497 TString opt = option;
498 opt.ToLower();
499 if (opt.Contains("all")) {
500 for (Int_t i=0;i<Size();i++) {
502 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]);
503 }
504 }
505}
506
507////////////////////////////////////////////////////////////////////////////////
508/// Save primitive as a C++ statement(s) on output stream.
509
510void TPolyMarker3D::SavePrimitive(std::ostream &out, Option_t * /*= ""*/)
511{
512 char quote = '"';
513 out<<" "<<std::endl;
514 if (gROOT->ClassSaved(TPolyMarker3D::Class())) {
515 out<<" ";
516 } else {
517 out<<" TPolyMarker3D *";
518 }
519 out<<"pmarker3D = new TPolyMarker3D("<<fN<<","<<GetMarkerStyle()<<","<<quote<<fOption<<quote<<");"<<std::endl;
520 out<<" pmarker3D->SetName("<<quote<<GetName()<<quote<<");"<<std::endl;
521
522 SaveMarkerAttributes(out,"pmarker3D",1,1,1);
523
524 for (Int_t i=0;i<Size();i++) {
525 out<<" pmarker3D->SetPoint("<<i<<","<<fP[3*i]<<","<<fP[3*i+1]<<","<<fP[3*i+2]<<");"<<std::endl;
526 }
527 out<<" pmarker3D->Draw();"<<std::endl;
528}
529
530////////////////////////////////////////////////////////////////////////////////
531/// Change (i.e. set) the name of the TNamed.
532/// WARNING: if the object is a member of a THashTable or THashList container
533/// the container must be Rehash()'ed after SetName(). For example the list
534/// of objects in the current directory is a THashList.
535
537{
538 fName = name;
539 if (gPad && TestBit(kMustCleanup)) gPad->Modified();
540}
541
542////////////////////////////////////////////////////////////////////////////////
543/// Set point following LastPoint to x, y, z.
544/// Returns index of the point (new last point).
545
547{
548 fLastPoint++;
549 SetPoint(fLastPoint, x, y, z);
550 return fLastPoint;
551}
552
553////////////////////////////////////////////////////////////////////////////////
554/// Set point n to x, y, z.
555/// If n is more then the current TPolyMarker3D size (n > fN) then
556/// the polymarker will be resized to contain at least n points.
557
559{
560 if (n < 0) return;
561 if (!fP || n >= fN) {
562 // re-allocate the object
563 Int_t newN = TMath::Max(2*fN,n+1);
564 Float_t *savepoint = new Float_t [kDimension*newN];
565 if (fP && fN){
566 memcpy(savepoint,fP,kDimension*fN*sizeof(Float_t));
567 memset(&savepoint[kDimension*fN],0,(newN-fN)*sizeof(Float_t));
568 delete [] fP;
569 }
570 fP = savepoint;
571 fN = newN;
572 }
573 fP[kDimension*n ] = x;
574 fP[kDimension*n+1] = y;
575 fP[kDimension*n+2] = z;
577}
578
579////////////////////////////////////////////////////////////////////////////////
580/// Re-initialize polymarker with n points from p. If p=0 initialize with 0.
581/// if n <= 0 the current array of points is deleted.
582
584{
585 SetMarkerStyle(marker);
586 fOption = option;
587 if (n <= 0) {
588 fN = 0;
589 fLastPoint = -1;
590 delete [] fP;
591 fP = 0;
592 return;
593 }
594 fN = n;
595 if (fP) delete [] fP;
596 fP = new Float_t [3*fN];
597 if (p) {
598 for (Int_t i = 0; i < fN; i++) {
599 fP[3*i] = p[3*i];
600 fP[3*i+1] = p[3*i+1];
601 fP[3*i+2] = p[3*i+2];
602 }
603 } else {
604 memset(fP,0,kDimension*fN*sizeof(Float_t));
605 }
606 fLastPoint = fN-1;
607}
608
609////////////////////////////////////////////////////////////////////////////////
610/// Re-initialize polymarker with n points from p. If p=0 initialize with 0.
611/// if n <= 0 the current array of points is deleted.
612
614{
615 SetMarkerStyle(marker);
616 fOption = option;
617 if (n <= 0) {
618 fN = 0;
619 fLastPoint = -1;
620 delete [] fP;
621 fP = 0;
622 return;
623 }
624 fN = n;
625 if (fP) delete [] fP;
626 fP = new Float_t [3*fN];
627 if (p) {
628 for (Int_t i = 0; i < fN; i++) {
629 fP[3*i] = (Float_t) p[3*i];
630 fP[3*i+1] = (Float_t) p[3*i+1];
631 fP[3*i+2] = (Float_t) p[3*i+2];
632 }
633 } else {
634 memset(fP,0,kDimension*fN*sizeof(Float_t));
635 }
636 fLastPoint = fN-1;
637}
638
639////////////////////////////////////////////////////////////////////////////////
640/// Stream a 3-D polymarker object.
641
642void TPolyMarker3D::Streamer(TBuffer &b)
643{
644 UInt_t R__s, R__c;
645 if (b.IsReading()) {
646 Version_t R__v = b.ReadVersion(&R__s, &R__c);
647 if (R__v > 2) b.ClassBegin(TPolyMarker3D::IsA());
648 if (R__v > 2) b.ClassMember("TObject");
649 TObject::Streamer(b);
650 if (R__v > 2) b.ClassMember("TAttMarker");
651 TAttMarker::Streamer(b);
652 if (R__v > 2) b.ClassMember("fN","Int_t");
653 b >> fN;
654 if (fN) {
655 if (R__v > 2) b.ClassMember("fP","Float_t", kDimension*fN);
656 fP = new Float_t[kDimension*fN];
657 b.ReadFastArray(fP,kDimension*fN);
658 }
659 fLastPoint = fN-1;
660 if (R__v > 2) b.ClassMember("fOption","TString");
661 fOption.Streamer(b);
662 if (R__v > 2) b.ClassMember("fName","TString");
663 if (R__v > 1) fName.Streamer(b);
664 if (R__v > 2) b.ClassEnd(TPolyMarker3D::IsA());
665 b.CheckByteCount(R__s, R__c, TPolyMarker3D::IsA());
666 } else {
667 R__c = b.WriteVersion(TPolyMarker3D::IsA(), kTRUE);
668 b.ClassBegin(TPolyMarker3D::IsA());
669 b.ClassMember("TObject");
670 TObject::Streamer(b);
671 b.ClassMember("TAttMarker");
672 TAttMarker::Streamer(b);
673 b.ClassMember("fN","Int_t");
674 Int_t size = Size();
675 b << size;
676 if (size) {
677 b.ClassMember("fP","Float_t", kDimension*size);
678 b.WriteFastArray(fP, kDimension*size);
679 }
680 b.ClassMember("fOption","TString");
681 fOption.Streamer(b);
682 b.ClassMember("fName","TString");
683 fName.Streamer(b);
684 b.ClassEnd(TPolyMarker3D::IsA());
685 b.SetByteCount(R__c, kTRUE);
686 }
687}
688
689////////////////////////////////////////////////////////////////////////////////
690/// Fills the parameters x, y, z with the coordinate of the n-th point
691/// n must be between 0 and Size() - 1.
692
694{
695 if (n < 0 || n >= Size()) return;
696 if (!fP) return;
697 x = fP[kDimension*n ];
698 y = fP[kDimension*n+1];
699 z = fP[kDimension*n+2];
700}
701
702////////////////////////////////////////////////////////////////////////////////
703/// Fills the parameters x, y, z with the coordinate of the n-th point
704/// n must be between 0 and Size() - 1.
705
707{
708 if (n < 0 || n >= Size()) return;
709 if (!fP) return;
710 x = (Double_t)fP[kDimension*n ];
711 y = (Double_t)fP[kDimension*n+1];
712 z = (Double_t)fP[kDimension*n+2];
713}
714
715
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define h(i)
Definition RSha256.hxx:106
static const double x1[5]
int Int_t
Definition RtypesCore.h:45
short Version_t
Definition RtypesCore.h:65
short Marker_t
Definition RtypesCore.h:81
unsigned int UInt_t
Definition RtypesCore.h:46
const Bool_t kFALSE
Definition RtypesCore.h:92
double Double_t
Definition RtypesCore.h:59
float Float_t
Definition RtypesCore.h:57
const Bool_t kTRUE
Definition RtypesCore.h:91
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:364
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
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition TAttMarker.h:41
Class to manage histogram axis.
Definition TAxis.h:30
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:113
Bool_t fLocalFrame
Definition TBuffer3D.h:90
Int_t fColor
Definition TBuffer3D.h:88
Short_t fTransparency
Definition TBuffer3D.h:89
TObject * fID
Definition TBuffer3D.h:87
Bool_t SetRawSizes(UInt_t reqPnts, UInt_t reqPntsCapacity, UInt_t reqSegs, UInt_t reqSegsCapacity, UInt_t reqPols, UInt_t reqPolsCapacity)
Set kRaw tessellation section of buffer with supplied sizes.
Double_t * fPnts
Definition TBuffer3D.h:112
Buffer base class used for serializing objects.
Definition TBuffer.h:43
Collection abstract base class.
Definition TCollection.h:63
virtual void Local2Master(Double_t *local, Double_t *master)
Convert one point from local system to master reference system.
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:58
void Reset()
Mother of all ROOT objects.
Definition TObject.h:37
TObject & operator=(const TObject &rhs)
TObject assignment operator.
Definition TObject.h:283
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:187
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:130
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition TObject.cxx:107
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:696
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition TObject.cxx:445
virtual void Copy(TObject &object) const
Copy this to obj.
Definition TObject.cxx:63
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:893
@ kCanDelete
if object in a list can be deleted
Definition TObject.h:58
@ kMustCleanup
if object destructor must call RecursiveRemove()
Definition TObject.h:60
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.
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save primitive as a C++ statement(s) on output stream.
virtual void DrawPolyMarker(Int_t n, Float_t *p, Marker_t marker, Option_t *option="")
Draw this 3-D polymarker with new coordinates.
void SetPoint(Int_t n, Double_t x, Double_t y, Double_t z)
Set point n to x, y, z.
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 ~TPolyMarker3D()
3-D polymarker destructor.
virtual Int_t Size() const
virtual const char * GetName() const
Returns name of object.
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute action corresponding to one event.
virtual void SetName(const char *name)
Change (i.e.
TPolyMarker3D & operator=(const TPolyMarker3D &)
assignment operator
virtual void Paint(Option_t *option="")
Paint a TPolyMarker3D.
virtual Float_t * GetP() const
virtual void Copy(TObject &polymarker) const
Copy polymarker to polymarker obj.
virtual Int_t SetNextPoint(Double_t x, Double_t y, Double_t z)
Set point following LastPoint to x, y, z.
virtual void Print(Option_t *option="") const
Print 3-D polymarker with its attributes on stdout.
virtual void ls(Option_t *option="") const
List this 3-D polymarker.
Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute distance from point px,py to a 3-D polymarker.
virtual void Draw(Option_t *option="")
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:2811
virtual Double_t Rndm()
Machine independent random number generator.
Definition TRandom.cxx:552
Basic string class.
Definition TString.h:136
void ToLower()
Change string to lower-case.
Definition TString.cxx:1145
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:624
See TView3D.
Definition TView.h:25
virtual void WCtoNDC(const Float_t *pw, Float_t *pn)=0
static TView * CreateView(Int_t system=1, const Double_t *rmin=0, const Double_t *rmax=0)
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=0)=0
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
Short_t Max(Short_t a, Short_t b)
Definition TMathBase.h:212
Double_t Sqrt(Double_t x)
Definition TMath.h:691
Short_t Min(Short_t a, Short_t b)
Definition TMathBase.h:180