Logo ROOT   6.12/07
Reference Guide
TPolyMarker.cxx
Go to the documentation of this file.
1 // @(#)root/hist:$Id$
2 // Author: Rene Brun 12/12/94
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 "TROOT.h"
14 #include "TVirtualPad.h"
15 #include "TPolyMarker.h"
16 #include "TClass.h"
17 #include "TMath.h"
18 
20 
21 
22 /** \class TPolyMarker
23  \ingroup Hist
24 A PolyMarker is defined by an array on N points in a 2-D space.
25 At each point x[i], y[i] a marker is drawn.
26 Marker attributes are managed by TAttMarker.
27 See TMarker for the list of possible marker types.
28 */
29 
30 ////////////////////////////////////////////////////////////////////////////////
31 /// Default constructor.
32 
34 {
35  fN = 0;
36  fX = fY = 0;
37  fLastPoint = -1;
38 }
39 
40 ////////////////////////////////////////////////////////////////////////////////
41 /// Constructor.
42 
44  :TObject(), TAttMarker()
45 {
46  fOption = option;
48  fLastPoint = -1;
49  if (n <= 0) {
50  fN = 0;
51  fLastPoint = -1;
52  fX = fY = 0;
53  return;
54  }
55  fN = n;
56  fX = new Double_t [fN];
57  fY = new Double_t [fN];
58 }
59 
60 ////////////////////////////////////////////////////////////////////////////////
61 /// Constructor.
62 
64  :TObject(), TAttMarker()
65 {
66  fOption = option;
68  fLastPoint = -1;
69  if (n <= 0) {
70  fN = 0;
71  fLastPoint = -1;
72  fX = fY = 0;
73  return;
74  }
75  fN = n;
76  fX = new Double_t [fN];
77  fY = new Double_t [fN];
78  if (!x || !y) return;
79  for (Int_t i=0; i<fN;i++) { fX[i] = x[i]; fY[i] = y[i]; }
80  fLastPoint = fN-1;
81 }
82 
83 ////////////////////////////////////////////////////////////////////////////////
84 /// Constructor.
85 
87  :TObject(), TAttMarker()
88 {
89  fOption = option;
91  fLastPoint = -1;
92  if (n <= 0) {
93  fN = 0;
94  fLastPoint = -1;
95  fX = fY = 0;
96  return;
97  }
98  fN = n;
99  fX = new Double_t [fN];
100  fY = new Double_t [fN];
101  if (!x || !y) return;
102  for (Int_t i=0; i<fN;i++) { fX[i] = x[i]; fY[i] = y[i]; }
103  fLastPoint = fN-1;
104 }
105 
106 ////////////////////////////////////////////////////////////////////////////////
107 ///assignment operator
108 
110 {
111  if(this!=&pm) {
112  TObject::operator=(pm);
114  fN=pm.fN;
116  // delete first previous existing fX and fY
117  if (fX) delete [] fX;
118  if (fY) delete [] fY;
119  fX=pm.fX;
120  fY=pm.fY;
121  fOption=pm.fOption;
122  }
123  return *this;
124 }
125 
126 ////////////////////////////////////////////////////////////////////////////////
127 /// Destructor.
128 
130 {
131  if (fX) delete [] fX;
132  if (fY) delete [] fY;
133  fLastPoint = -1;
134 }
135 
136 ////////////////////////////////////////////////////////////////////////////////
137 /// Copy constructor.
138 
139 TPolyMarker::TPolyMarker(const TPolyMarker &polymarker) : TObject(polymarker), TAttMarker(polymarker)
140 {
141  fN = 0;
142  fX = fY = 0;
143  fLastPoint = -1;
144  ((TPolyMarker&)polymarker).Copy(*this);
145 }
146 
147 ////////////////////////////////////////////////////////////////////////////////
148 
149 void TPolyMarker::Copy(TObject &obj) const
150 {
151  // Copy.
152 
153  TObject::Copy(obj);
154  TAttMarker::Copy(((TPolyMarker&)obj));
155  ((TPolyMarker&)obj).fN = fN;
156  // delete first previous existing fX and fY
157  if (((TPolyMarker&)obj).fX) delete [] (((TPolyMarker&)obj).fX);
158  if (((TPolyMarker&)obj).fY) delete [] (((TPolyMarker&)obj).fY);
159  if (fN > 0) {
160  ((TPolyMarker&)obj).fX = new Double_t [fN];
161  ((TPolyMarker&)obj).fY = new Double_t [fN];
162  for (Int_t i=0; i<fN;i++) { ((TPolyMarker&)obj).fX[i] = fX[i], ((TPolyMarker&)obj).fY[i] = fY[i]; }
163  } else {
164  ((TPolyMarker&)obj).fX = 0;
165  ((TPolyMarker&)obj).fY = 0;
166  }
167  ((TPolyMarker&)obj).fOption = fOption;
168  ((TPolyMarker&)obj).fLastPoint = fLastPoint;
169 }
170 
171 ////////////////////////////////////////////////////////////////////////////////
172 /// Compute distance from point px,py to a polymarker.
173 ///
174 /// Compute the closest distance of approach from point px,py to each point
175 /// of the polymarker.
176 /// Returns when the distance found is below DistanceMaximum.
177 /// The distance is computed in pixels units.
178 
180 {
181  const Int_t big = 9999;
182 
183  // check if point is near one of the points
184  Int_t i, pxp, pyp, d;
185  Int_t distance = big;
186 
187  for (i=0;i<Size();i++) {
188  pxp = gPad->XtoAbsPixel(gPad->XtoPad(fX[i]));
189  pyp = gPad->YtoAbsPixel(gPad->YtoPad(fY[i]));
190  d = TMath::Abs(pxp-px) + TMath::Abs(pyp-py);
191  if (d < distance) distance = d;
192  }
193  return distance;
194 }
195 
196 ////////////////////////////////////////////////////////////////////////////////
197 /// Draw.
198 
200 {
201  AppendPad(option);
202 }
203 
204 ////////////////////////////////////////////////////////////////////////////////
205 /// Draw polymarker.
206 
208 {
209  TPolyMarker *newpolymarker = new TPolyMarker(n,x,y);
210  TAttMarker::Copy(*newpolymarker);
211  newpolymarker->fOption = fOption;
212  newpolymarker->SetBit(kCanDelete);
213  newpolymarker->AppendPad();
214 }
215 
216 ////////////////////////////////////////////////////////////////////////////////
217 /// Execute action corresponding to one event.
218 ///
219 /// This member function must be implemented to realize the action
220 /// corresponding to the mouse click on the object in the window
221 
223 {
224 }
225 
226 ////////////////////////////////////////////////////////////////////////////////
227 /// ls.
228 
230 {
232  printf("TPolyMarker N=%d\n",fN);
233 }
234 
235 ////////////////////////////////////////////////////////////////////////////////
236 /// Merge polymarkers in the collection in this polymarker.
237 
239 {
240  if (!li) return 0;
241  TIter next(li);
242 
243  //first loop to count the number of entries
244  TPolyMarker *pm;
245  Int_t npoints = 0;
246  while ((pm = (TPolyMarker*)next())) {
247  if (!pm->InheritsFrom(TPolyMarker::Class())) {
248  Error("Add","Attempt to add object of class: %s to a %s",pm->ClassName(),this->ClassName());
249  return -1;
250  }
251  npoints += pm->Size();
252  }
253 
254  //extend this polymarker to hold npoints
255  SetPoint(npoints-1,0,0);
256 
257  //merge all polymarkers
258  next.Reset();
259  while ((pm = (TPolyMarker*)next())) {
260  Int_t np = pm->Size();
261  Double_t *x = pm->GetX();
262  Double_t *y = pm->GetY();
263  for (Int_t i=0;i<np;i++) {
264  SetPoint(i,x[i],y[i]);
265  }
266  }
267 
268  return npoints;
269 }
270 
271 ////////////////////////////////////////////////////////////////////////////////
272 /// Paint.
273 
275 {
276  PaintPolyMarker(fLastPoint+1, fX, fY, option);
277 }
278 
279 ////////////////////////////////////////////////////////////////////////////////
280 /// Paint polymarker.
281 
283 {
284  if (n <= 0) return;
285  TAttMarker::Modify(); //Change marker attributes only if necessary
286  Double_t *xx = x;
287  Double_t *yy = y;
288  if (gPad->GetLogx()) {
289  xx = new Double_t[n];
290  for (Int_t ix=0;ix<n;ix++) xx[ix] = gPad->XtoPad(x[ix]);
291  }
292  if (gPad->GetLogy()) {
293  yy = new Double_t[n];
294  for (Int_t iy=0;iy<n;iy++) yy[iy] = gPad->YtoPad(y[iy]);
295  }
296  gPad->PaintPolyMarker(n,xx,yy,option);
297  if (x != xx) delete [] xx;
298  if (y != yy) delete [] yy;
299 }
300 
301 ////////////////////////////////////////////////////////////////////////////////
302 /// Print polymarker.
303 
305 {
306  printf("TPolyMarker N=%d\n",fN);
307 }
308 
309 ////////////////////////////////////////////////////////////////////////////////
310 /// Save primitive as a C++ statement(s) on output stream out.
311 
312 void TPolyMarker::SavePrimitive(std::ostream &out, Option_t *option /*= ""*/)
313 {
314  char quote = '"';
315  out<<" "<<std::endl;
316  out<<" Double_t *dum = 0;"<<std::endl;
317  if (gROOT->ClassSaved(TPolyMarker::Class())) {
318  out<<" ";
319  } else {
320  out<<" TPolyMarker *";
321  }
322  out<<"pmarker = new TPolyMarker("<<fN<<",dum,dum,"<<quote<<fOption<<quote<<");"<<std::endl;
323 
324  SaveMarkerAttributes(out,"pmarker",1,1,1);
325 
326  for (Int_t i=0;i<Size();i++) {
327  out<<" pmarker->SetPoint("<<i<<","<<fX[i]<<","<<fY[i]<<");"<<std::endl;
328  }
329  out<<" pmarker->Draw("
330  <<quote<<option<<quote<<");"<<std::endl;
331 }
332 
333 ////////////////////////////////////////////////////////////////////////////////
334 /// Set point following LastPoint to x, y.
335 /// Returns index of the point (new last point).
336 
338 {
339  fLastPoint++;
340  SetPoint(fLastPoint, x, y);
341  return fLastPoint;
342 }
343 
344 ////////////////////////////////////////////////////////////////////////////////
345 /// Set point number n.
346 /// if n is greater than the current size, the arrays are automatically
347 /// extended
348 
350 {
351  if (n < 0) return;
352  if (!fX || !fY || n >= fN) {
353  // re-allocate the object
354  Int_t newN = TMath::Max(2*fN,n+1);
355  Double_t *savex = new Double_t [newN];
356  Double_t *savey = new Double_t [newN];
357  if (fX && fN){
358  memcpy(savex,fX,fN*sizeof(Double_t));
359  memset(&savex[fN],0,(newN-fN)*sizeof(Double_t));
360  delete [] fX;
361  }
362  if (fY && fN){
363  memcpy(savey,fY,fN*sizeof(Double_t));
364  memset(&savey[fN],0,(newN-fN)*sizeof(Double_t));
365  delete [] fY;
366  }
367  fX = savex;
368  fY = savey;
369  fN = newN;
370  }
371  fX[n] = x;
372  fY[n] = y;
374 }
375 
376 ////////////////////////////////////////////////////////////////////////////////
377 /// If n <= 0 the current arrays of points are deleted.
378 
380 {
381  if (n <= 0) {
382  fN = 0;
383  fLastPoint = -1;
384  delete [] fX;
385  delete [] fY;
386  fX = fY = 0;
387  return;
388  }
389  SetPoint(n-1,0,0);
390 }
391 
392 ////////////////////////////////////////////////////////////////////////////////
393 /// If n <= 0 the current arrays of points are deleted.
394 
396 {
397  if (n <= 0) {
398  fN = 0;
399  fLastPoint = -1;
400  delete [] fX;
401  delete [] fY;
402  fX = fY = 0;
403  return;
404  }
405  fN =n;
406  if (fX) delete [] fX;
407  if (fY) delete [] fY;
408  fX = new Double_t[fN];
409  fY = new Double_t[fN];
410  for (Int_t i=0; i<fN;i++) {
411  if (x) fX[i] = (Double_t)x[i];
412  if (y) fY[i] = (Double_t)y[i];
413  }
414  fOption = option;
415  fLastPoint = fN-1;
416 }
417 
418 ////////////////////////////////////////////////////////////////////////////////
419 /// If n <= 0 the current arrays of points are deleted.
420 
422 {
423  if (n <= 0) {
424  fN = 0;
425  fLastPoint = -1;
426  delete [] fX;
427  delete [] fY;
428  fX = fY = 0;
429  return;
430  }
431  fN =n;
432  if (fX) delete [] fX;
433  if (fY) delete [] fY;
434  fX = new Double_t[fN];
435  fY = new Double_t[fN];
436  for (Int_t i=0; i<fN;i++) {
437  if (x) fX[i] = x[i];
438  if (y) fY[i] = y[i];
439  }
440  fOption = option;
441  fLastPoint = fN-1;
442 }
443 
444 ////////////////////////////////////////////////////////////////////////////////
445 /// Stream a class object.
446 
447 void TPolyMarker::Streamer(TBuffer &R__b)
448 {
449  if (R__b.IsReading()) {
450  UInt_t R__s, R__c;
451  Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
452  if (R__v > 1) {
453  R__b.ReadClassBuffer(TPolyMarker::Class(), this, R__v, R__s, R__c);
454  return;
455  }
456  //====process old versions before automatic schema evolution
457  TObject::Streamer(R__b);
458  TAttMarker::Streamer(R__b);
459  R__b >> fN;
460  fX = new Double_t[fN];
461  fY = new Double_t[fN];
462  Int_t i;
463  Float_t xold,yold;
464  for (i=0;i<fN;i++) {R__b >> xold; fX[i] = xold;}
465  for (i=0;i<fN;i++) {R__b >> yold; fY[i] = yold;}
466  fOption.Streamer(R__b);
467  R__b.CheckByteCount(R__s, R__c, TPolyMarker::IsA());
468  //====end of old versions
469 
470  } else {
472  }
473 }
virtual void Draw(Option_t *option="")
Draw.
Bool_t IsReading() const
Definition: TBuffer.h:83
Double_t * GetX() const
Definition: TPolyMarker.h:56
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute action corresponding to one event.
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
short Version_t
Definition: RtypesCore.h:61
virtual void DrawPolyMarker(Int_t n, Double_t *x, Double_t *y, Option_t *option="")
Draw polymarker.
float Float_t
Definition: RtypesCore.h:53
const char Option_t
Definition: RtypesCore.h:62
virtual void PaintPolyMarker(Int_t n, Double_t *x, Double_t *y, Option_t *option="")
Paint polymarker.
TString fOption
Definition: TPolyMarker.h:37
virtual ~TPolyMarker()
Destructor.
Buffer base class used for serializing objects.
Definition: TBuffer.h:40
TPolyMarker()
Default constructor.
Definition: TPolyMarker.cxx:33
virtual Int_t Merge(TCollection *list)
Merge polymarkers in the collection in this polymarker.
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)=0
#define gROOT
Definition: TROOT.h:402
int Int_t
Definition: RtypesCore.h:41
virtual Int_t Size() const
Definition: TPolyMarker.h:69
void Copy(TAttMarker &attmarker) const
Copy this marker attributes to a new TAttMarker.
Definition: TAttMarker.cxx:209
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save primitive as a C++ statement(s) on output stream out.
Short_t Abs(Short_t d)
Definition: TMathBase.h:108
TPolyMarker & operator=(const TPolyMarker &)
assignment operator
void Reset()
Definition: TCollection.h:250
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:694
if object in a list can be deleted
Definition: TObject.h:58
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition: TObject.cxx:105
Marker Attributes class.
Definition: TAttMarker.h:19
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition: TObject.cxx:128
Double_t x[n]
Definition: legend1.C:17
void Class()
Definition: Class.C:29
virtual void Copy(TObject &object) const
Copy this to obj.
Definition: TObject.cxx:61
TObject & operator=(const TObject &rhs)
TObject assignment operator.
Definition: TObject.h:271
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:244
virtual Int_t SetNextPoint(Double_t x, Double_t y)
Set point following LastPoint to x, y.
virtual void Paint(Option_t *option="")
Paint.
virtual void SetPolyMarker(Int_t n)
If n <= 0 the current arrays of points are deleted.
virtual void ls(Option_t *option="") const
ls.
virtual void Modify()
Change current marker attributes if necessary.
Definition: TAttMarker.cxx:219
Double_t * fX
Definition: TPolyMarker.h:35
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:443
Collection abstract base class.
Definition: TCollection.h:63
unsigned int UInt_t
Definition: RtypesCore.h:42
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:880
static void IndentLevel()
Functions used by ls() to indent an object hierarchy.
Definition: TROOT.cxx:2746
virtual void SetPoint(Int_t point, Double_t x, Double_t y)
Set point number n.
virtual void Print(Option_t *option="") const
Print polymarker.
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
A PolyMarker is defined by an array on N points in a 2-D space.
Definition: TPolyMarker.h:31
#define ClassImp(name)
Definition: Rtypes.h:359
double Double_t
Definition: RtypesCore.h:55
Double_t y[n]
Definition: legend1.C:17
Binding & operator=(OUT(*fun)(void))
Mother of all ROOT objects.
Definition: TObject.h:37
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute distance from point px,py to a polymarker.
Double_t * fY
Definition: TPolyMarker.h:36
virtual void Copy(TObject &polymarker) const
Copy this to obj.
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:200
#define gPad
Definition: TVirtualPad.h:285
Double_t * GetY() const
Definition: TPolyMarker.h:57
Int_t fLastPoint
Definition: TPolyMarker.h:34
const Int_t n
Definition: legend1.C:16
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0