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