Logo ROOT   6.10/09
Reference Guide
TPolyLine.cxx
Go to the documentation of this file.
1 // @(#)root/graf:$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 "TMath.h"
15 #include "TVirtualPad.h"
16 #include "TPolyLine.h"
17 #include "TClass.h"
18 
20 
21 /** \class TPolyLine
22 \ingroup BasicGraphics
23 
24 Defined by an array on N points in a 2-D space.
25 
26 One can draw the contour of the polyline or/and its fill area.
27 Example:
28 Begin_Macro(source)
29 {
30  Double_t x[5] = {.2,.7,.6,.25,.2};
31  Double_t y[5] = {.5,.1,.9,.7,.5};
32  TPolyLine *pline = new TPolyLine(5,x,y);
33  pline->SetFillColor(38);
34  pline->SetLineColor(2);
35  pline->SetLineWidth(4);
36  pline->Draw("f");
37  pline->Draw();
38 }
39 End_Macro
40 */
41 
42 ////////////////////////////////////////////////////////////////////////////////
43 /// PolyLine default constructor.
44 
46 {
47  fN = 0;
48  fX = 0;
49  fY = 0;
50  fLastPoint = -1;
51 }
52 
53 ////////////////////////////////////////////////////////////////////////////////
54 /// PolyLine normal constructor without initialisation.
55 /// Allocates n points. The option string is ignored.
56 
58  :TObject(), TAttLine(), TAttFill()
59 {
60  fOption = option;
61  fLastPoint = -1;
62  if (n <= 0) {
63  fN = 0;
64  fLastPoint = -1;
65  fX = fY = 0;
66  return;
67  }
68  fN = n;
69  fX = new Double_t[fN];
70  fY = new Double_t[fN];
71 }
72 
73 ////////////////////////////////////////////////////////////////////////////////
74 /// PolyLine normal constructor (single precision).
75 /// Makes n points with (x, y) coordinates from x and y.
76 /// The option string is ignored.
77 
79  :TObject(), TAttLine(), TAttFill()
80 {
81  fOption = option;
82  fLastPoint = -1;
83  if (n <= 0) {
84  fN = 0;
85  fLastPoint = -1;
86  fX = fY = 0;
87  return;
88  }
89  fN = n;
90  fX = new Double_t[fN];
91  fY = new Double_t[fN];
92  if (!x || !y) return;
93  for (Int_t i=0; i<fN;i++) { fX[i] = x[i]; fY[i] = y[i];}
94  fLastPoint = fN-1;
95 }
96 
97 ////////////////////////////////////////////////////////////////////////////////
98 /// PolyLine normal constructor (double precision).
99 /// Makes n points with (x, y) coordinates from x and y.
100 /// The option string is ignored.
101 
103  :TObject(), TAttLine(), TAttFill()
104 {
105  fOption = option;
106  fLastPoint = -1;
107  if (n <= 0) {
108  fN = 0;
109  fLastPoint = -1;
110  fX = fY = 0;
111  return;
112  }
113  fN = n;
114  fX = new Double_t[fN];
115  fY = new Double_t[fN];
116  if (!x || !y) return;
117  for (Int_t i=0; i<fN;i++) { fX[i] = x[i]; fY[i] = y[i];}
118  fLastPoint = fN-1;
119 }
120 
121 ////////////////////////////////////////////////////////////////////////////////
122 ///assignment operator
123 
125 {
126  if(this!=&pl) {
127  TObject::operator=(pl);
130  fN=pl.fN;
132  fX=pl.fX;
133  fY=pl.fY;
134  fOption=pl.fOption;
135  }
136  return *this;
137 }
138 
139 ////////////////////////////////////////////////////////////////////////////////
140 /// PolyLine default destructor.
141 
143 {
144  if (fX) delete [] fX;
145  if (fY) delete [] fY;
146 }
147 
148 ////////////////////////////////////////////////////////////////////////////////
149 /// PolyLine copy constructor.
150 
151 TPolyLine::TPolyLine(const TPolyLine &polyline) : TObject(polyline), TAttLine(polyline), TAttFill(polyline)
152 {
153  fN = 0;
154  fX = 0;
155  fY = 0;
156  fLastPoint = -1;
157  ((TPolyLine&)polyline).Copy(*this);
158 }
159 
160 ////////////////////////////////////////////////////////////////////////////////
161 /// Copy this polyline to polyline.
162 
163 void TPolyLine::Copy(TObject &obj) const
164 {
165  TObject::Copy(obj);
166  TAttLine::Copy(((TPolyLine&)obj));
167  TAttFill::Copy(((TPolyLine&)obj));
168  ((TPolyLine&)obj).fN = fN;
169  if (fN > 0) {
170  ((TPolyLine&)obj).fX = new Double_t[fN];
171  ((TPolyLine&)obj).fY = new Double_t[fN];
172  for (Int_t i=0; i<fN;i++) {((TPolyLine&)obj).fX[i] = fX[i]; ((TPolyLine&)obj).fY[i] = fY[i];}
173  } else {
174  ((TPolyLine&)obj).fX = 0;
175  ((TPolyLine&)obj).fY = 0;
176  }
177  ((TPolyLine&)obj).fOption = fOption;
178  ((TPolyLine&)obj).fLastPoint = fLastPoint;
179 }
180 
181 ////////////////////////////////////////////////////////////////////////////////
182 /// Returns closest distance in pixels from point (px, py) to a polyline.
183 ///
184 /// First looks for distances to the points of the polyline. Stops search
185 /// and returns if a vertex of the polyline is found to be closer than 10
186 /// pixels. Thus the return value may depend on the ordering of points
187 /// in the polyline.
188 ///
189 /// Then looks for distances to the lines of the polyline. There is no
190 /// arbitrary cutoff; any distance may be found.
191 ///
192 /// Finally checks whether (px, py) is inside a closed and filled polyline.
193 /// (Must be EXACTLY closed. "Filled" means fill color and fill style are
194 /// both non-zero.) If so, returns zero.
195 ///
196 /// Returns 9999 if the polyline has no points.
197 
199 {
200  const Int_t big = 9999;
201  const Int_t kMaxDiff = 10;
202 
203  // check if point is near one of the points
204  Int_t i, pxp, pyp, d;
205  Int_t distance = big;
206  if (Size() <= 0) return distance;
207 
208  for (i=0;i<Size();i++) {
209  pxp = gPad->XtoAbsPixel(gPad->XtoPad(fX[i]));
210  pyp = gPad->YtoAbsPixel(gPad->YtoPad(fY[i]));
211  d = TMath::Abs(pxp-px) + TMath::Abs(pyp-py);
212  if (d < distance) distance = d;
213  }
214  if (distance < kMaxDiff) return distance;
215 
216  // check if point is near one of the connecting lines
217  for (i=0;i<Size()-1;i++) {
218  d = DistancetoLine(px, py, gPad->XtoPad(fX[i]), gPad->YtoPad(fY[i]), gPad->XtoPad(fX[i+1]), gPad->YtoPad(fY[i+1]));
219  if (d < distance) distance = d;
220  }
221 
222  // in case of a closed and filled polyline, check if we are inside
223  if (fFillColor && fFillStyle && fX[0] == fX[fLastPoint] && fY[0] == fY[fLastPoint]) {
224  if (TMath::IsInside(gPad->AbsPixeltoX(px),gPad->AbsPixeltoY(py),fLastPoint+1,fX,fY)) distance = 0;
225  }
226  return distance;
227 }
228 
229 ////////////////////////////////////////////////////////////////////////////////
230 /// Draw this polyline with its current attributes.
231 
233 {
234  AppendPad(option);
235 }
236 
237 ////////////////////////////////////////////////////////////////////////////////
238 /// Draw this polyline with new coordinates.
239 
241 {
242  TPolyLine *newpolyline = new TPolyLine(n,x,y);
243  TAttLine::Copy(*newpolyline);
244  TAttFill::Copy(*newpolyline);
245  newpolyline->fOption = fOption;
246  newpolyline->SetBit(kCanDelete);
247  newpolyline->AppendPad(option);
248 }
249 
250 ////////////////////////////////////////////////////////////////////////////////
251 /// Execute action corresponding to one event.
252 ///
253 /// This member function is called when a polyline is clicked with the locator
254 ///
255 /// If Left button clicked on one of the line end points, this point
256 /// follows the cursor until button is released.
257 ///
258 /// if Middle button clicked, the line is moved parallel to itself
259 /// until the button is released.
260 
262 {
263  if (!gPad) return;
264 
265  Int_t i, d;
266  Double_t xmin, xmax, ymin, ymax, dx, dy, dxr, dyr;
267  const Int_t kMaxDiff = 10;
268  static Bool_t middle;
269  static Int_t ipoint, pxp, pyp;
270  static Int_t px1,px2,py1,py2;
271  static Int_t pxold, pyold, px1old, py1old, px2old, py2old;
272  static Int_t dpx, dpy;
273  static Int_t *x=0, *y=0;
274  Bool_t opaque = gPad->OpaqueMoving();
275 
276  if (!gPad->IsEditable()) return;
277 
278  Int_t np = Size();
279 
280  switch (event) {
281 
282  case kButton1Down:
283  gVirtualX->SetLineColor(-1);
284  TAttLine::Modify(); //Change line attributes only if necessary
285  px1 = gPad->XtoAbsPixel(gPad->GetX1());
286  py1 = gPad->YtoAbsPixel(gPad->GetY1());
287  px2 = gPad->XtoAbsPixel(gPad->GetX2());
288  py2 = gPad->YtoAbsPixel(gPad->GetY2());
289  ipoint = -1;
290 
291 
292  if (x || y) break;
293  x = new Int_t[np+1];
294  y = new Int_t[np+1];
295  for (i=0;i<np;i++) {
296  pxp = gPad->XtoAbsPixel(gPad->XtoPad(fX[i]));
297  pyp = gPad->YtoAbsPixel(gPad->YtoPad(fY[i]));
298  if (!opaque) {
299  gVirtualX->DrawLine(pxp-4, pyp-4, pxp+4, pyp-4);
300  gVirtualX->DrawLine(pxp+4, pyp-4, pxp+4, pyp+4);
301  gVirtualX->DrawLine(pxp+4, pyp+4, pxp-4, pyp+4);
302  gVirtualX->DrawLine(pxp-4, pyp+4, pxp-4, pyp-4);
303  }
304  x[i] = pxp;
305  y[i] = pyp;
306  d = TMath::Abs(pxp-px) + TMath::Abs(pyp-py);
307  if (d < kMaxDiff) ipoint =i;
308  }
309  dpx = 0;
310  dpy = 0;
311  pxold = px;
312  pyold = py;
313  if (ipoint < 0) return;
314  if (ipoint == 0) {
315  px1old = 0;
316  py1old = 0;
317  px2old = gPad->XtoAbsPixel(fX[1]);
318  py2old = gPad->YtoAbsPixel(fY[1]);
319  } else if (ipoint == fN-1) {
320  px1old = gPad->XtoAbsPixel(gPad->XtoPad(fX[fN-2]));
321  py1old = gPad->YtoAbsPixel(gPad->YtoPad(fY[fN-2]));
322  px2old = 0;
323  py2old = 0;
324  } else {
325  px1old = gPad->XtoAbsPixel(gPad->XtoPad(fX[ipoint-1]));
326  py1old = gPad->YtoAbsPixel(gPad->YtoPad(fY[ipoint-1]));
327  px2old = gPad->XtoAbsPixel(gPad->XtoPad(fX[ipoint+1]));
328  py2old = gPad->YtoAbsPixel(gPad->YtoPad(fY[ipoint+1]));
329  }
330  pxold = gPad->XtoAbsPixel(gPad->XtoPad(fX[ipoint]));
331  pyold = gPad->YtoAbsPixel(gPad->YtoPad(fY[ipoint]));
332 
333  break;
334 
335 
336  case kMouseMotion:
337 
338  middle = kTRUE;
339  for (i=0;i<np;i++) {
340  pxp = gPad->XtoAbsPixel(gPad->XtoPad(fX[i]));
341  pyp = gPad->YtoAbsPixel(gPad->YtoPad(fY[i]));
342  d = TMath::Abs(pxp-px) + TMath::Abs(pyp-py);
343  if (d < kMaxDiff) middle = kFALSE;
344  }
345 
346 
347  // check if point is close to an axis
348  if (middle) gPad->SetCursor(kMove);
349  else gPad->SetCursor(kHand);
350  break;
351 
352  case kButton1Motion:
353  if (!opaque) {
354  if (middle) {
355  for(i=0;i<np-1;i++) {
356  gVirtualX->DrawLine(x[i]+dpx, y[i]+dpy, x[i+1]+dpx, y[i+1]+dpy);
357  pxp = x[i]+dpx;
358  pyp = y[i]+dpy;
359  gVirtualX->DrawLine(pxp-4, pyp-4, pxp+4, pyp-4);
360  gVirtualX->DrawLine(pxp+4, pyp-4, pxp+4, pyp+4);
361  gVirtualX->DrawLine(pxp+4, pyp+4, pxp-4, pyp+4);
362  gVirtualX->DrawLine(pxp-4, pyp+4, pxp-4, pyp-4);
363  }
364  pxp = x[np-1]+dpx;
365  pyp = y[np-1]+dpy;
366  gVirtualX->DrawLine(pxp-4, pyp-4, pxp+4, pyp-4);
367  gVirtualX->DrawLine(pxp+4, pyp-4, pxp+4, pyp+4);
368  gVirtualX->DrawLine(pxp+4, pyp+4, pxp-4, pyp+4);
369  gVirtualX->DrawLine(pxp-4, pyp+4, pxp-4, pyp-4);
370  dpx += px - pxold;
371  dpy += py - pyold;
372  pxold = px;
373  pyold = py;
374  for(i=0;i<np-1;i++) {
375  gVirtualX->DrawLine(x[i]+dpx, y[i]+dpy, x[i+1]+dpx, y[i+1]+dpy);
376  pxp = x[i]+dpx;
377  pyp = y[i]+dpy;
378  gVirtualX->DrawLine(pxp-4, pyp-4, pxp+4, pyp-4);
379  gVirtualX->DrawLine(pxp+4, pyp-4, pxp+4, pyp+4);
380  gVirtualX->DrawLine(pxp+4, pyp+4, pxp-4, pyp+4);
381  gVirtualX->DrawLine(pxp-4, pyp+4, pxp-4, pyp-4);
382  }
383  pxp = x[np-1]+dpx;
384  pyp = y[np-1]+dpy;
385  gVirtualX->DrawLine(pxp-4, pyp-4, pxp+4, pyp-4);
386  gVirtualX->DrawLine(pxp+4, pyp-4, pxp+4, pyp+4);
387  gVirtualX->DrawLine(pxp+4, pyp+4, pxp-4, pyp+4);
388  gVirtualX->DrawLine(pxp-4, pyp+4, pxp-4, pyp-4);
389  } else {
390  if (px1old) gVirtualX->DrawLine(px1old, py1old, pxold, pyold);
391  if (px2old) gVirtualX->DrawLine(pxold, pyold, px2old, py2old);
392  gVirtualX->DrawLine(pxold-4, pyold-4, pxold+4, pyold-4);
393  gVirtualX->DrawLine(pxold+4, pyold-4, pxold+4, pyold+4);
394  gVirtualX->DrawLine(pxold+4, pyold+4, pxold-4, pyold+4);
395  gVirtualX->DrawLine(pxold-4, pyold+4, pxold-4, pyold-4);
396  pxold = px;
397  pxold = TMath::Max(pxold, px1);
398  pxold = TMath::Min(pxold, px2);
399  pyold = py;
400  pyold = TMath::Max(pyold, py2);
401  pyold = TMath::Min(pyold, py1);
402  if (px1old) gVirtualX->DrawLine(px1old, py1old, pxold, pyold);
403  if (px2old) gVirtualX->DrawLine(pxold, pyold, px2old, py2old);
404  gVirtualX->DrawLine(pxold-4, pyold-4, pxold+4, pyold-4);
405  gVirtualX->DrawLine(pxold+4, pyold-4, pxold+4, pyold+4);
406  gVirtualX->DrawLine(pxold+4, pyold+4, pxold-4, pyold+4);
407  gVirtualX->DrawLine(pxold-4, pyold+4, pxold-4, pyold-4);
408  }
409  } else {
410  if (middle) {
411  for(i=0;i<np-1;i++) {
412  pxp = x[i]+dpx;
413  pyp = y[i]+dpy;
414  }
415  pxp = x[np-1]+dpx;
416  pyp = y[np-1]+dpy;
417  dpx += px - pxold;
418  dpy += py - pyold;
419  pxold = px;
420  pyold = py;
421  } else {
422  pxold = px;
423  pxold = TMath::Max(pxold, px1);
424  pxold = TMath::Min(pxold, px2);
425  pyold = py;
426  pyold = TMath::Max(pyold, py2);
427  pyold = TMath::Min(pyold, py1);
428  }
429  if (x && y) {
430  if (middle) {
431  for(i=0;i<np;i++) {
432  fX[i] = gPad->PadtoX(gPad->AbsPixeltoX(x[i]+dpx));
433  fY[i] = gPad->PadtoY(gPad->AbsPixeltoY(y[i]+dpy));
434  }
435  } else {
436  fX[ipoint] = gPad->PadtoX(gPad->AbsPixeltoX(pxold));
437  fY[ipoint] = gPad->PadtoY(gPad->AbsPixeltoY(pyold));
438  }
439  }
440  gPad->Modified(kTRUE);
441  }
442  break;
443 
444  case kButton1Up:
445 
446  // Compute x,y range
447  xmin = gPad->GetUxmin();
448  xmax = gPad->GetUxmax();
449  ymin = gPad->GetUymin();
450  ymax = gPad->GetUymax();
451  dx = xmax-xmin;
452  dy = ymax-ymin;
453  dxr = dx/(1 - gPad->GetLeftMargin() - gPad->GetRightMargin());
454  dyr = dy/(1 - gPad->GetBottomMargin() - gPad->GetTopMargin());
455 
456  // Range() could change the size of the pad pixmap and therefore should
457  // be called before the other paint routines
458  gPad->Range(xmin - dxr*gPad->GetLeftMargin(),
459  ymin - dyr*gPad->GetBottomMargin(),
460  xmax + dxr*gPad->GetRightMargin(),
461  ymax + dyr*gPad->GetTopMargin());
462  gPad->RangeAxis(xmin, ymin, xmax, ymax);
463 
464  if (x && y) {
465  if (middle) {
466  for(i=0;i<np;i++) {
467  fX[i] = gPad->PadtoX(gPad->AbsPixeltoX(x[i]+dpx));
468  fY[i] = gPad->PadtoY(gPad->AbsPixeltoY(y[i]+dpy));
469  }
470  } else {
471  fX[ipoint] = gPad->PadtoX(gPad->AbsPixeltoX(pxold));
472  fY[ipoint] = gPad->PadtoY(gPad->AbsPixeltoY(pyold));
473  }
474  delete [] x; x = 0;
475  delete [] y; y = 0;
476  }
477  gPad->Modified(kTRUE);
478  gVirtualX->SetLineColor(-1);
479  }
480 }
481 
482 ////////////////////////////////////////////////////////////////////////////////
483 /// List this polyline with its attributes.
484 /// The option string is ignored.
485 
486 void TPolyLine::ls(Option_t *) const
487 {
489  printf("TPolyLine N=%d\n",fN);
490 }
491 
492 ////////////////////////////////////////////////////////////////////////////////
493 /// Merge polylines in the collection in this polyline
494 
496 {
497  if (!li) return 0;
498  TIter next(li);
499 
500  //first loop to count the number of entries
501  TPolyLine *pl;
502  Int_t npoints = 0;
503  while ((pl = (TPolyLine*)next())) {
504  if (!pl->InheritsFrom(TPolyLine::Class())) {
505  Error("Add","Attempt to add object of class: %s to a %s",pl->ClassName(),this->ClassName());
506  return -1;
507  }
508  npoints += pl->Size();
509  }
510 
511  //extend this polyline to hold npoints
512  if (npoints > 1) SetPoint(npoints-1,0,0);
513 
514  //merge all polylines
515  next.Reset();
516  while ((pl = (TPolyLine*)next())) {
517  Int_t np = pl->Size();
518  Double_t *x = pl->GetX();
519  Double_t *y = pl->GetY();
520  for (Int_t i=0;i<np;i++) {
521  SetPoint(i,x[i],y[i]);
522  }
523  }
524 
525  return npoints;
526 }
527 
528 ////////////////////////////////////////////////////////////////////////////////
529 /// Paint this polyline with its current attributes.
530 
532 {
533  if (TestBit(kPolyLineNDC)) {
534  if (strlen(option) > 0) PaintPolyLineNDC(fLastPoint+1, fX, fY, option);
536  } else {
537  if (strlen(option) > 0) PaintPolyLine(fLastPoint+1, fX, fY, option);
538  else PaintPolyLine(fLastPoint+1, fX, fY, fOption.Data());
539  }
540 }
541 
542 ////////////////////////////////////////////////////////////////////////////////
543 /// Draw this polyline with new coordinates.
544 ///
545 /// If option = 'f' or 'F' the fill area is drawn.
546 /// The default is to draw the lines only.
547 
549 {
550  if (n <= 0) return;
551  TAttLine::Modify(); //Change line attributes only if necessary
552  TAttFill::Modify(); //Change fill area attributes only if necessary
553  Double_t *xx = x;
554  Double_t *yy = y;
555  if (gPad->GetLogx()) {
556  xx = new Double_t[n];
557  for (Int_t ix=0;ix<n;ix++) xx[ix] = gPad->XtoPad(x[ix]);
558  }
559  if (gPad->GetLogy()) {
560  yy = new Double_t[n];
561  for (Int_t iy=0;iy<n;iy++) yy[iy] = gPad->YtoPad(y[iy]);
562  }
563  if (*option == 'f' || *option == 'F') gPad->PaintFillArea(n,xx,yy,option);
564  else gPad->PaintPolyLine(n,xx,yy,option);
565  if (x != xx) delete [] xx;
566  if (y != yy) delete [] yy;
567 }
568 
569 ////////////////////////////////////////////////////////////////////////////////
570 /// Draw this polyline with new coordinates in NDC.
571 
573 {
574  TAttLine::Modify(); //Change line attributes only if necessary
575  TAttFill::Modify(); //Change fill area attributes only if necessary
576  gPad->PaintPolyLineNDC(n,x,y,option);
577 }
578 
579 ////////////////////////////////////////////////////////////////////////////////
580 /// Dump this polyline with its attributes.
581 /// The option string is ignored.
582 
584 {
585  printf("PolyLine N=%d\n",fN);
586 }
587 
588 ////////////////////////////////////////////////////////////////////////////////
589 /// Save primitive as a C++ statement(s) on output stream out
590 
591 void TPolyLine::SavePrimitive(std::ostream &out, Option_t *option /*= ""*/)
592 {
593  char quote = '"';
594  out<<" "<<std::endl;
595  if (gROOT->ClassSaved(TPolyLine::Class())) {
596  out<<" ";
597  } else {
598  out<<" Double_t *dum = 0;"<<std::endl;
599  out<<" TPolyLine *";
600  }
601  out<<"pline = new TPolyLine("<<fN<<",dum,dum,"<<quote<<fOption<<quote<<");"<<std::endl;
602 
603  SaveFillAttributes(out,"pline",0,1001);
604  SaveLineAttributes(out,"pline",1,1,1);
605 
606  for (Int_t i=0;i<Size();i++) {
607  out<<" pline->SetPoint("<<i<<","<<fX[i]<<","<<fY[i]<<");"<<std::endl;
608  }
609  out<<" pline->Draw("
610  <<quote<<option<<quote<<");"<<std::endl;
611 }
612 
613 ////////////////////////////////////////////////////////////////////////////////
614 /// Set NDC mode on if isNDC = kTRUE, off otherwise
615 
617 {
619  if (isNDC) SetBit(kPolyLineNDC);
620 }
621 
622 ////////////////////////////////////////////////////////////////////////////////
623 /// Set point following LastPoint to x, y.
624 /// Returns index of the point (new last point).
625 
627 {
628  fLastPoint++;
629  SetPoint(fLastPoint, x, y);
630  return fLastPoint;
631 }
632 
633 ////////////////////////////////////////////////////////////////////////////////
634 /// Set point number n to (x, y)
635 /// If n is greater than the current size, the arrays are automatically
636 /// extended.
637 
639 {
640  if (n < 0) return;
641  if (!fX || !fY || n >= fN) {
642  // re-allocate the object
643  Int_t newN = TMath::Max(2*fN,n+1);
644  Double_t *savex = new Double_t [newN];
645  Double_t *savey = new Double_t [newN];
646  if (fX && fN){
647  memcpy(savex,fX,fN*sizeof(Double_t));
648  memset(&savex[fN],0,(newN-fN)*sizeof(Double_t));
649  delete [] fX;
650  }
651  if (fY && fN){
652  memcpy(savey,fY,fN*sizeof(Double_t));
653  memset(&savey[fN],0,(newN-fN)*sizeof(Double_t));
654  delete [] fY;
655  }
656  fX = savex;
657  fY = savey;
658  fN = newN;
659  }
660  fX[n] = x;
661  fY[n] = y;
663 }
664 
665 ////////////////////////////////////////////////////////////////////////////////
666 /// Resize this polyline to size n.
667 /// If n <= 0 the current arrays of points are deleted.
668 /// If n is greater than the current size, the new points are set to (0, 0)
669 
671 {
672  if (n <= 0) {
673  fN = 0;
674  fLastPoint = -1;
675  delete [] fX;
676  delete [] fY;
677  fX = fY = 0;
678  return;
679  }
680  if (n < fN) {
681  fN = n;
682  fLastPoint = n - 1;
683  } else {
684  SetPoint(n-1,0,0);
685  }
686 }
687 
688 ////////////////////////////////////////////////////////////////////////////////
689 /// Set new values for this polyline (single precision).
690 ///
691 /// If n <= 0 the current arrays of points are deleted.
692 
694 {
695  if (n <= 0) {
696  fN = 0;
697  fLastPoint = -1;
698  delete [] fX;
699  delete [] fY;
700  fX = fY = 0;
701  return;
702  }
703  fN =n;
704  if (fX) delete [] fX;
705  if (fY) delete [] fY;
706  fX = new Double_t[fN];
707  fY = new Double_t[fN];
708  for (Int_t i=0; i<fN;i++) {
709  if (x) fX[i] = (Double_t)x[i];
710  if (y) fY[i] = (Double_t)y[i];
711  }
712  fOption = option;
713  fLastPoint = fN-1;
714 }
715 
716 ////////////////////////////////////////////////////////////////////////////////
717 /// Set new values for this polyline (double precision).
718 ///
719 /// If n <= 0 the current arrays of points are deleted.
720 
722 {
723  if (n <= 0) {
724  fN = 0;
725  fLastPoint = -1;
726  delete [] fX;
727  delete [] fY;
728  fX = fY = 0;
729  return;
730  }
731  fN =n;
732  if (fX) delete [] fX;
733  if (fY) delete [] fY;
734  fX = new Double_t[fN];
735  fY = new Double_t[fN];
736  for (Int_t i=0; i<fN;i++) {
737  if (x) fX[i] = x[i];
738  if (y) fY[i] = y[i];
739  }
740  fOption = option;
741  fLastPoint = fN-1;
742 }
743 
744 ////////////////////////////////////////////////////////////////////////////////
745 /// Stream a class object.
746 
747 void TPolyLine::Streamer(TBuffer &b)
748 {
749  if (b.IsReading()) {
750  UInt_t R__s, R__c;
751  Version_t R__v = b.ReadVersion(&R__s, &R__c);
752  if (R__v > 1) {
753  b.ReadClassBuffer(TPolyLine::Class(), this, R__v, R__s, R__c);
754  return;
755  }
756  //====process old versions before automatic schema evolution
757  TObject::Streamer(b);
758  TAttLine::Streamer(b);
759  TAttFill::Streamer(b);
760  b >> fN;
761  fX = new Double_t[fN];
762  fY = new Double_t[fN];
763  Float_t *x = new Float_t[fN];
764  Float_t *y = new Float_t[fN];
765  b.ReadFastArray(x,fN);
766  b.ReadFastArray(y,fN);
767  for (Int_t i=0;i<fN;i++) {
768  fX[i] = x[i];
769  fY[i] = y[i];
770  }
771  fOption.Streamer(b);
772  b.CheckByteCount(R__s, R__c, TPolyLine::IsA());
773  //====end of old versions
774 
775  } else {
777  }
778 }
TPolyLine()
PolyLine default constructor.
Definition: TPolyLine.cxx:45
Bool_t IsReading() const
Definition: TBuffer.h:81
float xmin
Definition: THbookFile.cxx:93
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
Int_t fN
Number of points.
Definition: TPolyLine.h:26
Double_t * GetX() const
Definition: TPolyLine.h:54
short Version_t
Definition: RtypesCore.h:61
float Float_t
Definition: RtypesCore.h:53
const char Option_t
Definition: RtypesCore.h:62
float ymin
Definition: THbookFile.cxx:93
Bool_t TestBit(UInt_t f) const
Definition: TObject.h:159
Buffer base class used for serializing objects.
Definition: TBuffer.h:40
Double_t * fY
[fN] Array of Y coordinates
Definition: TPolyLine.h:29
TString fOption
options
Definition: TPolyLine.h:30
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)=0
#define gROOT
Definition: TROOT.h:375
TPolyLine & operator=(const TPolyLine &)
assignment operator
Definition: TPolyLine.cxx:124
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:168
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
virtual void SetPoint(Int_t point, Double_t x, Double_t y)
Set point number n to (x, y) If n is greater than the current size, the arrays are automatically exte...
Definition: TPolyLine.cxx:638
virtual Int_t Merge(TCollection *list)
Merge polylines in the collection in this polyline.
Definition: TPolyLine.cxx:495
virtual void Modify()
Change current line attributes if necessary.
Definition: TAttLine.cxx:232
Int_t fLastPoint
The index of the last filled point.
Definition: TPolyLine.h:27
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save primitive as a C++ statement(s) on output stream out.
Definition: TPolyLine.cxx:591
Short_t Abs(Short_t d)
Definition: TMathBase.h:108
void Reset()
Definition: TCollection.h:156
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:687
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:112
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition: TObject.cxx:135
Fill Area Attributes class.
Definition: TAttFill.h:19
Double_t x[n]
Definition: legend1.C:17
void Class()
Definition: Class.C:29
void Copy(TAttLine &attline) const
Copy this line attributes to a new TAttLine.
Definition: TAttLine.cxx:162
virtual void Copy(TObject &object) const
Copy this to obj.
Definition: TObject.cxx:60
virtual void SaveLineAttributes(std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1, Int_t widdef=1)
Save line attributes as C++ statement(s) on output stream out.
Definition: TAttLine.cxx:260
virtual void Modify()
Change current fill area attributes if necessary.
Definition: TAttFill.cxx:209
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute action corresponding to one event.
Definition: TPolyLine.cxx:261
TObject & operator=(const TObject &rhs)
TObject assignment operator.
Definition: TObject.h:250
Bool_t IsInside(T xp, T yp, Int_t np, T *x, T *y)
Definition: TMath.h:1205
float ymax
Definition: THbookFile.cxx:93
virtual void Draw(Option_t *option="")
Draw this polyline with its current attributes.
Definition: TPolyLine.cxx:232
virtual void ls(Option_t *option="") const
List this polyline with its attributes.
Definition: TPolyLine.cxx:486
Double_t * GetY() const
Definition: TPolyLine.h:55
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:436
Collection abstract base class.
Definition: TCollection.h:42
unsigned int UInt_t
Definition: RtypesCore.h:42
virtual void SaveFillAttributes(std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1001)
Save fill attributes as C++ statement(s) on output stream out.
Definition: TAttFill.cxx:232
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:873
static void IndentLevel()
Functions used by ls() to indent an object hierarchy.
Definition: TROOT.cxx:2632
virtual void PaintPolyLine(Int_t n, Double_t *x, Double_t *y, Option_t *option="")
Draw this polyline with new coordinates.
Definition: TPolyLine.cxx:548
float xmax
Definition: THbookFile.cxx:93
virtual void SetPolyLine(Int_t n)
Resize this polyline to size n.
Definition: TPolyLine.cxx:670
virtual void ReadFastArray(Bool_t *b, Int_t n)=0
#define gVirtualX
Definition: TVirtualX.h:350
const Bool_t kFALSE
Definition: RtypesCore.h:92
virtual ~TPolyLine()
PolyLine default destructor.
Definition: TPolyLine.cxx:142
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
virtual Int_t Size() const
Definition: TPolyLine.h:71
#define ClassImp(name)
Definition: Rtypes.h:336
Polyline coordinates are in NDC space.
Definition: TPolyLine.h:37
double Double_t
Definition: RtypesCore.h:55
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Returns closest distance in pixels from point (px, py) to a polyline.
Definition: TPolyLine.cxx:198
Double_t y[n]
Definition: legend1.C:17
Int_t DistancetoLine(Int_t px, Int_t py, Double_t xp1, Double_t yp1, Double_t xp2, Double_t yp2)
Compute distance from point px,py to a line.
Definition: TAttLine.cxx:196
virtual void Print(Option_t *option="") const
Dump this polyline with its attributes.
Definition: TPolyLine.cxx:583
Binding & operator=(OUT(*fun)(void))
Mother of all ROOT objects.
Definition: TObject.h:37
virtual void SetNDC(Bool_t isNDC=kTRUE)
Set NDC mode on if isNDC = kTRUE, off otherwise.
Definition: TPolyLine.cxx:616
Color_t fFillColor
Fill area color.
Definition: TAttFill.h:22
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:200
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
Defined by an array on N points in a 2-D space.
Definition: TPolyLine.h:23
#define gPad
Definition: TVirtualPad.h:284
virtual Int_t SetNextPoint(Double_t x, Double_t y)
Set point following LastPoint to x, y.
Definition: TPolyLine.cxx:626
virtual void Paint(Option_t *option="")
Paint this polyline with its current attributes.
Definition: TPolyLine.cxx:531
void ResetBit(UInt_t f)
Definition: TObject.h:158
Double_t * fX
[fN] Array of X coordinates
Definition: TPolyLine.h:28
virtual void Copy(TObject &polyline) const
Copy this polyline to polyline.
Definition: TPolyLine.cxx:163
virtual void PaintPolyLineNDC(Int_t n, Double_t *x, Double_t *y, Option_t *option="")
Draw this polyline with new coordinates in NDC.
Definition: TPolyLine.cxx:572
const Bool_t kTRUE
Definition: RtypesCore.h:91
const Int_t n
Definition: legend1.C:16
Line Attributes class.
Definition: TAttLine.h:18
virtual void DrawPolyLine(Int_t n, Double_t *x, Double_t *y, Option_t *option="")
Draw this polyline with new coordinates.
Definition: TPolyLine.cxx:240
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0
Style_t fFillStyle
Fill area style.
Definition: TAttFill.h:23
void Copy(TAttFill &attfill) const
Copy this fill attributes to a new TAttFill.
Definition: TAttFill.cxx:200
const char * Data() const
Definition: TString.h:347