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