Logo ROOT   6.10/09
Reference Guide
TEllipse.cxx
Go to the documentation of this file.
1 // @(#)root/graf:$Id$
2 // Author: Rene Brun 16/10/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 <stdlib.h>
13 
14 #include "Riostream.h"
15 #include "TROOT.h"
16 #include "TEllipse.h"
17 #include "TVirtualPad.h"
18 #include "TMath.h"
19 #include "TClass.h"
20 #include "TPoint.h"
21 
22 const Double_t kPI = 3.14159265358979323846;
23 
25 
26 /** \class TEllipse
27 \ingroup BasicGraphics
28 
29 Draw Ellipses.
30 
31 The ellipse can be truncated and rotated. It is defined by its center `(x1,y1)`
32 and two radius`r1` and `r2`.
33 
34 A minimum and maximum angle may be specified `(phimin, phimax)`.
35 The ellipse may be rotated with an angle `theta`. All these
36 angles are in degrees.
37 The attributes of the outline line are given via `TAttLine`.
38 The attributes of the fill area are given via `TAttFill`.
39 The picture below illustrates different types of ellipses.
40 
41 When an ellipse sector only is drawn, the lines connecting the center
42 of the ellipse to the edges are drawn by default. One can specify
43 the drawing option "only" to not draw these lines or alternatively
44 call the function `SetNoEdges()`. To remove completely the ellipse
45 outline it is enough to specify 0 as line style.
46 
47 Begin_Macro(source)
48 ../../../tutorials/graphics/ellipse.C
49 End_Macro
50 */
51 
52 ////////////////////////////////////////////////////////////////////////////////
53 /// Ellipse default constructor.
54 
56 {
57  fX1 = 0;
58  fY1 = 0;
59  fR1 = 1;
60  fR2 = 1;
61  fPhimin = 0;
62  fPhimax = 360;
63  fTheta = 0;
64 }
65 
66 ////////////////////////////////////////////////////////////////////////////////
67 /// Ellipse normal constructor.
68 
70  :TObject(), TAttLine(), TAttFill(0,1001)
71 {
72  fX1 = x1;
73  fY1 = y1;
74  fR1 = r1;
75  fR2 = r2;
76  fPhimin = phimin;
77  fPhimax = phimax;
78  fTheta = theta;
79  if (r2 <= 0) fR2 = fR1;
80 }
81 
82 ////////////////////////////////////////////////////////////////////////////////
83 /// Ellipse default destructor.
84 
86 {
87 }
88 
89 ////////////////////////////////////////////////////////////////////////////////
90 /// Copy constructor.
91 
92 TEllipse::TEllipse(const TEllipse &ellipse) : TObject(ellipse), TAttLine(ellipse), TAttFill(ellipse), TAttBBox2D(ellipse)
93 {
94  fX1 = 0;
95  fY1 = 0;
96  fR1 = 1;
97  fR2 = 1;
98  fPhimin = 0;
99  fPhimax = 360;
100  fTheta = 0;
101 
102  ((TEllipse&)ellipse).Copy(*this);
103 }
104 
105 ////////////////////////////////////////////////////////////////////////////////
106 /// Copy this ellipse to ellipse.
107 
108 void TEllipse::Copy(TObject &obj) const
109 {
110  TObject::Copy(obj);
111  TAttLine::Copy(((TEllipse&)obj));
112  TAttFill::Copy(((TEllipse&)obj));
113  ((TEllipse&)obj).fX1 = fX1;
114  ((TEllipse&)obj).fY1 = fY1;
115  ((TEllipse&)obj).fR1 = fR1;
116  ((TEllipse&)obj).fR2 = fR2;
117  ((TEllipse&)obj).fPhimin = fPhimin;
118  ((TEllipse&)obj).fPhimax = fPhimax;
119  ((TEllipse&)obj).fTheta = fTheta;
120 }
121 
122 ////////////////////////////////////////////////////////////////////////////////
123 /// Compute distance from point px,py to an ellipse.
124 ///
125 /// Compute the closest distance of approach from point px,py to this
126 /// ellipse. The distance is computed in pixels units.
127 ///
128 /// In case of a filled ellipse the distance returned is 0 if the point
129 /// (px,py) is inside the ellipse, and is huge if the point is outside.
130 
132 {
133  Double_t x = gPad->PadtoX(gPad->AbsPixeltoX(px));
134  Double_t y = gPad->PadtoY(gPad->AbsPixeltoY(py));
135 
136  Double_t dxnr = x - fX1;
137  Double_t dynr = y - fY1;
138 
139  Double_t ct = TMath::Cos(kPI*GetTheta()/180.0);
140  Double_t st = TMath::Sin(kPI*GetTheta()/180.0);
141 
142  Double_t dx = dxnr*ct + dynr*st;
143  Double_t dy = -dxnr*st + dynr*ct;
144 
145  Double_t r1 = fR1;
146  Double_t r2 = fR2;
147 
148  if (dx == 0 || r1 == 0 || r2 == 0) return 9999;
149  Double_t distp = TMath::Sqrt(dx*dx + dy*dy);
150 
151  Double_t tana = dy/dx;
152  tana *= tana;
153  Double_t distr = TMath::Sqrt((1+tana)/(1.0/(r1*r1) + tana/(r2*r2)));
154  Int_t dist = 9999;
155  if (GetFillColor() && GetFillStyle()) {
156  if (distr > distp) dist = 0;
157  } else {
158  if (TMath::Abs(distr-distp)/(r1+r2) < 0.01) dist = 0;
159  }
160  return dist;
161 }
162 
163 ////////////////////////////////////////////////////////////////////////////////
164 /// Draw this ellipse with its current attributes.
165 
166 void TEllipse::Draw(Option_t *option)
167 {
168  AppendPad(option);
169 }
170 
171 ////////////////////////////////////////////////////////////////////////////////
172 /// Draw this ellipse with new coordinates.
173 
175 {
176  TEllipse *newellipse = new TEllipse(x1, y1, r1, r2, phimin, phimax,theta);
177  TAttLine::Copy(*newellipse);
178  TAttFill::Copy(*newellipse);
179  newellipse->SetBit(kCanDelete);
180  newellipse->AppendPad(option);
181  if (TestBit(kNoEdges)) newellipse->SetBit(kNoEdges);
182 }
183 
184 ////////////////////////////////////////////////////////////////////////////////
185 /// Execute action corresponding to one event.
186 ///
187 /// This member function is called when a line is clicked with the locator
188 ///
189 /// If Left button clicked on one of the line end points, this point
190 /// follows the cursor until button is released.
191 ///
192 /// if Middle button clicked, the line is moved parallel to itself
193 /// until the button is released.
194 ///
195 /// NOTE that support for log scale is not implemented
196 
197 void TEllipse::ExecuteEvent(Int_t event, Int_t px, Int_t py)
198 {
199  if (!gPad) return;
200 
201  Int_t kMaxDiff = 10;
202 
203  Int_t i, dpx, dpy;
204  Double_t angle,dx,dy,dphi,ct,st,fTy,fBy,fLx,fRx;
205  static Int_t px1,py1,npe,r1,r2,sav1,sav2;
206  const Int_t kMinSize = 25;
207  const Int_t np = 40;
208  static Bool_t pTop, pL, pR, pBot, pINSIDE;
209  static Int_t pTx,pTy,pLx,pLy,pRx,pRy,pBx,pBy;
210  static Int_t x[np+2], y[np+2];
211  static Int_t pxold, pyold;
212  static Int_t sig,impair;
213  static Double_t sdx, sdy;
214  static Double_t oldX1, oldY1, oldR1, oldR2;
215 
216  Bool_t opaque = gPad->OpaqueMoving();
217 
218  if (!gPad->IsEditable()) return;
219 
220  switch (event) {
221 
222  case kArrowKeyPress:
223  case kButton1Down:
224  oldX1 = fX1;
225  oldY1 = fY1;
226  oldR1 = fR1;
227  oldR2 = fR2;
228  dphi = (fPhimax-fPhimin)*kPI/(180*np);
229  ct = TMath::Cos(kPI*fTheta/180);
230  st = TMath::Sin(kPI*fTheta/180);
231  for (i=0;i<np;i++) {
232  angle = fPhimin*kPI/180 + Double_t(i)*dphi;
233  dx = fR1*TMath::Cos(angle);
234  dy = fR2*TMath::Sin(angle);
235  x[i] = gPad->XtoAbsPixel(fX1 + dx*ct - dy*st);
236  y[i] = gPad->YtoAbsPixel(fY1 + dx*st + dy*ct);
237  }
238  if (fPhimax-fPhimin >= 360 ) {
239  x[np] = x[0];
240  y[np] = y[0];
241  npe = np;
242  } else {
243  x[np] = gPad->XtoAbsPixel(fX1);
244  y[np] = gPad->YtoAbsPixel(fY1);
245  x[np+1] = x[0];
246  y[np+1] = y[0];
247  npe = np + 1;
248  }
249  impair = 0;
250  px1 = gPad->XtoAbsPixel(fX1);
251  py1 = gPad->YtoAbsPixel(fY1);
252  pTx = pBx = px1;
253  pLy = pRy = py1;
254  pTy = gPad->YtoAbsPixel(fR2+fY1);
255  pBy = gPad->YtoAbsPixel(-fR2+fY1);
256  pLx = gPad->XtoAbsPixel(-fR1+fX1);
257  pRx = gPad->XtoAbsPixel(fR1+fX1);
258  r2 = (pBy-pTy)/2;
259  r1 = (pRx-pLx)/2;
260  if (!opaque) {
261  gVirtualX->SetLineColor(-1);
263  gVirtualX->DrawLine(pRx+4, py1+4, pRx-4, py1+4);
264  gVirtualX->DrawLine(pRx-4, py1+4, pRx-4, py1-4);
265  gVirtualX->DrawLine(pRx-4, py1-4, pRx+4, py1-4);
266  gVirtualX->DrawLine(pRx+4, py1-4, pRx+4, py1+4);
267  gVirtualX->DrawLine(pLx+4, py1+4, pLx-4, py1+4);
268  gVirtualX->DrawLine(pLx-4, py1+4, pLx-4, py1-4);
269  gVirtualX->DrawLine(pLx-4, py1-4, pLx+4, py1-4);
270  gVirtualX->DrawLine(pLx+4, py1-4, pLx+4, py1+4);
271  gVirtualX->DrawLine(px1+4, pBy+4, px1-4, pBy+4);
272  gVirtualX->DrawLine(px1-4, pBy+4, px1-4, pBy-4);
273  gVirtualX->DrawLine(px1-4, pBy-4, px1+4, pBy-4);
274  gVirtualX->DrawLine(px1+4, pBy-4, px1+4, pBy+4);
275  gVirtualX->DrawLine(px1+4, pTy+4, px1-4, pTy+4);
276  gVirtualX->DrawLine(px1-4, pTy+4, px1-4, pTy-4);
277  gVirtualX->DrawLine(px1-4, pTy-4, px1+4, pTy-4);
278  gVirtualX->DrawLine(px1+4, pTy-4, px1+4, pTy+4);
279  }
280  else {
281  sdx = this->GetX1()-gPad->AbsPixeltoX(px);
282  sdy = this->GetY1()-gPad->AbsPixeltoY(py);
283  }
284  // No break !!!
285 
286  case kMouseMotion:
287  px1 = gPad->XtoAbsPixel(fX1);
288  py1 = gPad->YtoAbsPixel(fY1);
289  pTx = pBx = px1;
290  pLy = pRy = py1;
291  pTy = gPad->YtoAbsPixel(fR2+fY1);
292  pBy = gPad->YtoAbsPixel(-fR2+fY1);
293  pLx = gPad->XtoAbsPixel(-fR1+fX1);
294  pRx = gPad->XtoAbsPixel(fR1+fX1);
295  pTop = pL = pR = pBot = pINSIDE = kFALSE;
296  if ((TMath::Abs(px - pTx) < kMaxDiff) &&
297  (TMath::Abs(py - pTy) < kMaxDiff)) { // top edge
298  pTop = kTRUE;
299  gPad->SetCursor(kTopSide);
300  }
301  else
302  if ((TMath::Abs(px - pBx) < kMaxDiff) &&
303  (TMath::Abs(py - pBy) < kMaxDiff)) { // bottom edge
304  pBot = kTRUE;
305  gPad->SetCursor(kBottomSide);
306  }
307  else
308  if ((TMath::Abs(py - pLy) < kMaxDiff) &&
309  (TMath::Abs(px - pLx) < kMaxDiff)) { // left edge
310  pL = kTRUE;
311  gPad->SetCursor(kLeftSide);
312  }
313  else
314  if ((TMath::Abs(py - pRy) < kMaxDiff) &&
315  (TMath::Abs(px - pRx) < kMaxDiff)) { // right edge
316  pR = kTRUE;
317  gPad->SetCursor(kRightSide);
318  }
319  else {pINSIDE= kTRUE; gPad->SetCursor(kMove); }
320  pxold = px; pyold = py;
321 
322  break;
323 
324  case kArrowKeyRelease:
325  case kButton1Motion:
326  if (!opaque)
327  {
328  gVirtualX->DrawLine(pRx+4, py1+4, pRx-4, py1+4);
329  gVirtualX->DrawLine(pRx-4, py1+4, pRx-4, py1-4);
330  gVirtualX->DrawLine(pRx-4, py1-4, pRx+4, py1-4);
331  gVirtualX->DrawLine(pRx+4, py1-4, pRx+4, py1+4);
332  gVirtualX->DrawLine(pLx+4, py1+4, pLx-4, py1+4);
333  gVirtualX->DrawLine(pLx-4, py1+4, pLx-4, py1-4);
334  gVirtualX->DrawLine(pLx-4, py1-4, pLx+4, py1-4);
335  gVirtualX->DrawLine(pLx+4, py1-4, pLx+4, py1+4);
336  gVirtualX->DrawLine(px1+4, pBy+4, px1-4, pBy+4);
337  gVirtualX->DrawLine(px1-4, pBy+4, px1-4, pBy-4);
338  gVirtualX->DrawLine(px1-4, pBy-4, px1+4, pBy-4);
339  gVirtualX->DrawLine(px1+4, pBy-4, px1+4, pBy+4);
340  gVirtualX->DrawLine(px1+4, pTy+4, px1-4, pTy+4);
341  gVirtualX->DrawLine(px1-4, pTy+4, px1-4, pTy-4);
342  gVirtualX->DrawLine(px1-4, pTy-4, px1+4, pTy-4);
343  gVirtualX->DrawLine(px1+4, pTy-4, px1+4, pTy+4);
344  for (i=0;i<npe;i++) gVirtualX->DrawLine(x[i], y[i], x[i+1], y[i+1]);
345  }
346  if (pTop) {
347  sav1 = py1;
348  sav2 = r2;
349  py1 += (py - pyold)/2;
350  r2 -= (py - pyold)/2;
351  if (TMath::Abs(pyold-py)%2==1) impair++;
352  if (py-pyold>0) sig=+1;
353  else sig=-1;
354  if (impair==2) { impair = 0; py1 += sig; r2 -= sig;}
355  if (py1 > pBy-kMinSize) {py1 = sav1; r2 = sav2; py = pyold;}
356  }
357  if (pBot) {
358  sav1 = py1;
359  sav2 = r2;
360  py1 += (py - pyold)/2;
361  r2 += (py - pyold)/2;
362  if (TMath::Abs(pyold-py)%2==1) impair++;
363  if (py-pyold>0) sig=+1;
364  else sig=-1;
365  if (impair==2) { impair = 0; py1 += sig; r2 += sig;}
366  if (py1 < pTy+kMinSize) {py1 = sav1; r2 = sav2; py = pyold;}
367  }
368  if (pL) {
369  sav1 = px1;
370  sav2 = r1;
371  px1 += (px - pxold)/2;
372  r1 -= (px - pxold)/2;
373  if (TMath::Abs(pxold-px)%2==1) impair++;
374  if (px-pxold>0) sig=+1;
375  else sig=-1;
376  if (impair==2) { impair = 0; px1 += sig; r1 -= sig;}
377  if (px1 > pRx-kMinSize) {px1 = sav1; r1 = sav2; px = pxold;}
378  }
379  if (pR) {
380  sav1 = px1;
381  sav2 = r1;
382  px1 += (px - pxold)/2;
383  r1 += (px - pxold)/2;
384  if (TMath::Abs(pxold-px)%2==1) impair++;
385  if (px-pxold>0) sig=+1;
386  else sig=-1;
387  if (impair==2) { impair = 0; px1 += sig; r1 += sig;}
388  if (px1 < pLx+kMinSize) {px1 = sav1; r1 = sav2; px = pxold;}
389  }
390  if (pTop || pBot || pL || pR) {
391  if (!opaque) {
392  dphi = (fPhimax-fPhimin)*kPI/(180*np);
393  ct = TMath::Cos(kPI*fTheta/180);
394  st = TMath::Sin(kPI*fTheta/180);
395  for (i=0;i<np;i++) {
396  angle = fPhimin*kPI/180 + Double_t(i)*dphi;
397  dx = r1*TMath::Cos(angle);
398  dy = r2*TMath::Sin(angle);
399  x[i] = px1 + Int_t(dx*ct - dy*st);
400  y[i] = py1 + Int_t(dx*st + dy*ct);
401  }
402  if (fPhimax-fPhimin >= 360 ) {
403  x[np] = x[0];
404  y[np] = y[0];
405  npe = np;
406  } else {
407  x[np] = px1;
408  y[np] = py1;
409  x[np+1] = x[0];
410  y[np+1] = y[0];
411  npe = np + 1;
412  }
413  gVirtualX->SetLineColor(-1);
415  for (i=0;i<npe;i++)
416  gVirtualX->DrawLine(x[i], y[i], x[i+1], y[i+1]);
417  }
418  else
419  {
420  this->SetX1(gPad->AbsPixeltoX(px1));
421  this->SetY1(gPad->AbsPixeltoY(py1));
422  this->SetR1(TMath::Abs(gPad->AbsPixeltoX(px1-r1)-gPad->AbsPixeltoX(px1+r1))/2);
423  this->SetR2(TMath::Abs(gPad->AbsPixeltoY(py1-r2)-gPad->AbsPixeltoY(py1+r2))/2);
424  if (pTop) gPad->ShowGuidelines(this, event, 't', true);
425  if (pBot) gPad->ShowGuidelines(this, event, 'b', true);
426  if (pL) gPad->ShowGuidelines(this, event, 'l', true);
427  if (pR) gPad->ShowGuidelines(this, event, 'r', true);
428  gPad->Modified(kTRUE);
429  gPad->Update();
430  }
431  }
432  if (pINSIDE) {
433  if (!opaque){
434  dpx = px-pxold; dpy = py-pyold;
435  px1 += dpx; py1 += dpy;
436  for (i=0;i<=npe;i++) { x[i] += dpx; y[i] += dpy;}
437  for (i=0;i<npe;i++) gVirtualX->DrawLine(x[i], y[i], x[i+1], y[i+1]);
438  }
439  else {
440  this->SetX1(gPad->AbsPixeltoX(px)+sdx);
441  this->SetY1(gPad->AbsPixeltoY(py)+sdy);
442  gPad->ShowGuidelines(this, event, 'i', true);
443  gPad->Modified(kTRUE);
444  gPad->Update();
445  }
446  }
447  if (!opaque){
448  pTx = pBx = px1;
449  pRx = px1+r1;
450  pLx = px1-r1;
451  pRy = pLy = py1;
452  pTy = py1-r2;
453  pBy = py1+r2;
454  gVirtualX->DrawLine(pRx+4, py1+4, pRx-4, py1+4);
455  gVirtualX->DrawLine(pRx-4, py1+4, pRx-4, py1-4);
456  gVirtualX->DrawLine(pRx-4, py1-4, pRx+4, py1-4);
457  gVirtualX->DrawLine(pRx+4, py1-4, pRx+4, py1+4);
458  gVirtualX->DrawLine(pLx+4, py1+4, pLx-4, py1+4);
459  gVirtualX->DrawLine(pLx-4, py1+4, pLx-4, py1-4);
460  gVirtualX->DrawLine(pLx-4, py1-4, pLx+4, py1-4);
461  gVirtualX->DrawLine(pLx+4, py1-4, pLx+4, py1+4);
462  gVirtualX->DrawLine(px1+4, pBy+4, px1-4, pBy+4);
463  gVirtualX->DrawLine(px1-4, pBy+4, px1-4, pBy-4);
464  gVirtualX->DrawLine(px1-4, pBy-4, px1+4, pBy-4);
465  gVirtualX->DrawLine(px1+4, pBy-4, px1+4, pBy+4);
466  gVirtualX->DrawLine(px1+4, pTy+4, px1-4, pTy+4);
467  gVirtualX->DrawLine(px1-4, pTy+4, px1-4, pTy-4);
468  gVirtualX->DrawLine(px1-4, pTy-4, px1+4, pTy-4);
469  gVirtualX->DrawLine(px1+4, pTy-4, px1+4, pTy+4);
470  }
471  pxold = px;
472  pyold = py;
473  break;
474 
475  case kButton1Up:
476  if (gROOT->IsEscaped()) {
477  gROOT->SetEscape(kFALSE);
478  if (opaque) {
479  this->SetX1(oldX1);
480  this->SetY1(oldY1);
481  this->SetR1(oldR1);
482  this->SetR2(oldR2);
483  gPad->Modified(kTRUE);
484  gPad->Update();
485  }
486  break;
487  }
488 
489  if (opaque) {
490  gPad->ShowGuidelines(this, event);
491  } else {
492  fX1 = gPad->AbsPixeltoX(px1);
493  fY1 = gPad->AbsPixeltoY(py1);
494  fBy = gPad->AbsPixeltoY(py1+r2);
495  fTy = gPad->AbsPixeltoY(py1-r2);
496  fLx = gPad->AbsPixeltoX(px1+r1);
497  fRx = gPad->AbsPixeltoX(px1-r1);
498  fR1 = TMath::Abs(fRx-fLx)/2;
499  fR2 = TMath::Abs(fTy-fBy)/2;
500  gPad->Modified(kTRUE);
501  gVirtualX->SetLineColor(-1);
502  }
503  }
504 }
505 
506 ////////////////////////////////////////////////////////////////////////////////
507 /// List this ellipse with its attributes.
508 
509 void TEllipse::ls(Option_t *) const
510 {
512  printf("%s: X1= %f Y1=%f R1=%f R2=%f\n",GetName(),fX1,fY1,fR1,fR2);
513 }
514 
515 ////////////////////////////////////////////////////////////////////////////////
516 /// Paint this ellipse with its current attributes.
517 
518 void TEllipse::Paint(Option_t *option)
519 {
521 }
522 
523 ////////////////////////////////////////////////////////////////////////////////
524 /// Draw this ellipse with new coordinates.
525 
527  Double_t phimin, Double_t phimax, Double_t theta,
528  Option_t *option)
529 {
530  const Int_t np = 200;
531  static Double_t x[np+3], y[np+3];
532  TAttLine::Modify(); //Change line attributes only if necessary
533  TAttFill::Modify(); //Change fill attributes only if necessary
534 
535  Double_t phi1 = TMath::Min(phimin,phimax);
536  Double_t phi2 = TMath::Max(phimin,phimax);
537 
538  //set number of points approximatively proportional to the ellipse circumference
539  Double_t circ = kPI*(r1+r2)*(phi2-phi1)/360;
540  Int_t n = (Int_t)(np*circ/((gPad->GetX2()-gPad->GetX1())+(gPad->GetY2()-gPad->GetY1())));
541  if (n < 8) n= 8;
542  if (n > np) n = np;
543  Double_t angle,dx,dy;
544  Double_t dphi = (phi2-phi1)*kPI/(180*n);
545  Double_t ct = TMath::Cos(kPI*theta/180);
546  Double_t st = TMath::Sin(kPI*theta/180);
547  for (Int_t i=0;i<=n;i++) {
548  angle = phi1*kPI/180 + Double_t(i)*dphi;
549  dx = r1*TMath::Cos(angle);
550  dy = r2*TMath::Sin(angle);
551  x[i] = gPad->XtoPad(x1 + dx*ct - dy*st);
552  y[i] = gPad->YtoPad(y1 + dx*st + dy*ct);
553  }
554  TString opt = option;
555  opt.ToLower();
556  if (phi2-phi1 >= 360 ) {
557  if (GetFillStyle()) gPad->PaintFillArea(n,x,y);
558  if (GetLineStyle()) gPad->PaintPolyLine(n+1,x,y);
559  } else {
560  x[n+1] = gPad->XtoPad(x1);
561  y[n+1] = gPad->YtoPad(y1);
562  x[n+2] = x[0];
563  y[n+2] = y[0];
564  if (GetFillStyle()) gPad->PaintFillArea(n+2,x,y);
565  if (GetLineStyle()) {
566  if (TestBit(kNoEdges) || opt.Contains("only")) gPad->PaintPolyLine(n+1,x,y);
567  else gPad->PaintPolyLine(n+3,x,y);
568  }
569  }
570 }
571 
572 ////////////////////////////////////////////////////////////////////////////////
573 /// Dump this ellipse with its attributes.
574 
575 void TEllipse::Print(Option_t *) const
576 {
577  printf("Ellipse: X1=%f Y1=%f R1=%f R2=%f",fX1,fY1,fR1,fR2);
578  if (GetLineColor() != 1) printf(" Color=%d",GetLineColor());
579  if (GetLineStyle() != 1) printf(" Style=%d",GetLineStyle());
580  if (GetLineWidth() != 1) printf(" Width=%d",GetLineWidth());
581  printf("\n");
582 }
583 
584 ////////////////////////////////////////////////////////////////////////////////
585 /// Save primitive as a C++ statement(s) on output stream out
586 
587 void TEllipse::SavePrimitive(std::ostream &out, Option_t * /*= ""*/)
588 {
589  out<<" "<<std::endl;
590  if (gROOT->ClassSaved(TEllipse::Class())) {
591  out<<" ";
592  } else {
593  out<<" TEllipse *";
594  }
595  out<<"ellipse = new TEllipse("<<fX1<<","<<fY1<<","<<fR1<<","<<fR2
596  <<","<<fPhimin<<","<<fPhimax<<","<<fTheta<<");"<<std::endl;
597 
598  SaveFillAttributes(out,"ellipse",0,1001);
599  SaveLineAttributes(out,"ellipse",1,1,1);
600 
601  if (GetNoEdges()) out<<" ellipse->SetNoEdges();"<<std::endl;
602 
603  out<<" ellipse->Draw();"<<std::endl;
604 }
605 
606 ////////////////////////////////////////////////////////////////////////////////
607 /// Return kTRUE if kNoEdges bit is set, kFALSE otherwise.
608 
610 {
611  return TestBit(kNoEdges) ? kTRUE : kFALSE;
612 }
613 
614 ////////////////////////////////////////////////////////////////////////////////
615 /// if noEdges = kTRUE the lines connecting the center to the edges
616 /// will not be drawn.
617 /// default is to draw the edges.
618 
619 void TEllipse::SetNoEdges(Bool_t noEdges)
620 {
621  if (noEdges) SetBit(kNoEdges);
622  else ResetBit(kNoEdges);
623 }
624 
625 ////////////////////////////////////////////////////////////////////////////////
626 /// Stream an object of class TEllipse.
627 
628 void TEllipse::Streamer(TBuffer &R__b)
629 {
630  if (R__b.IsReading()) {
631  UInt_t R__s, R__c;
632  Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
633  if (R__v > 1) {
634  R__b.ReadClassBuffer(TEllipse::Class(), this, R__v, R__s, R__c);
635  return;
636  }
637  //====process old versions before automatic schema evolution
638  TObject::Streamer(R__b);
639  TAttLine::Streamer(R__b);
640  TAttFill::Streamer(R__b);
641  Float_t x1,y1,r1,r2,phimin,phimax,theta;
642  R__b >> x1; fX1 = x1;
643  R__b >> y1; fY1 = y1;
644  R__b >> r1; fR1 = r1;
645  R__b >> r2; fR2 = r2;
646  R__b >> phimin; fPhimin = phimin;
647  R__b >> phimax; fPhimax = phimax;
648  R__b >> theta; fTheta = theta;
649  R__b.CheckByteCount(R__s, R__c, TEllipse::IsA());
650  //====end of old versions
651 
652  } else {
653  R__b.WriteClassBuffer(TEllipse::Class(),this);
654  }
655 }
656 
657 ////////////////////////////////////////////////////////////////////////////////
658 /// Return the bounding Box of the Ellipse, currently not taking into
659 /// account the rotating angle.
660 
662 {
663  Rectangle_t BBox;
664  BBox.fX = gPad->XtoPixel(fX1-fR1);
665  BBox.fY = gPad->YtoPixel(fY1+fR2);
666  BBox.fWidth = gPad->XtoPixel(fX1+fR1)-gPad->XtoPixel(fX1-fR1);
667  BBox.fHeight = gPad->YtoPixel(fY1-fR2)-gPad->YtoPixel(fY1+fR2);
668  return (BBox);
669 }
670 
671 ////////////////////////////////////////////////////////////////////////////////
672 /// Return the center of the Ellipse as TPoint in pixels
673 
675 {
676  TPoint p;
677  p.SetX(gPad->XtoPixel(fX1));
678  p.SetY(gPad->YtoPixel(fY1));
679  return(p);
680 }
681 
682 ////////////////////////////////////////////////////////////////////////////////
683 /// Set center of the Ellipse
684 
685 void TEllipse::SetBBoxCenter(const TPoint &p)
686 {
687  fX1 = gPad->PixeltoX(p.GetX());
688  fY1 = gPad->PixeltoY(p.GetY()-gPad->VtoPixel(0));
689 }
690 
691 ////////////////////////////////////////////////////////////////////////////////
692 /// Set X coordinate of the center of the Ellipse
693 
694 void TEllipse::SetBBoxCenterX(const Int_t x)
695 {
696  fX1 = gPad->PixeltoX(x);
697 }
698 
699 ////////////////////////////////////////////////////////////////////////////////
700 /// Set Y coordinate of the center of the Ellipse
701 
702 void TEllipse::SetBBoxCenterY(const Int_t y)
703 {
704  fY1 = gPad->PixeltoY(y-gPad->VtoPixel(0));
705 }
706 
707 ////////////////////////////////////////////////////////////////////////////////
708 /// Set left hand side of BoundingBox to a value
709 /// (resize in x direction on left)
710 
711 void TEllipse::SetBBoxX1(const Int_t x)
712 {
713  Double_t x1 = gPad->PixeltoX(x);
714  if (x1>fX1+fR1) return;
715 
716  fR1 = (fX1+fR1-x1)*0.5;
717  fX1 = x1 + fR1;
718 }
719 
720 ////////////////////////////////////////////////////////////////////////////////
721 /// Set right hand side of BoundingBox to a value
722 /// (resize in x direction on right)
723 
724 void TEllipse::SetBBoxX2(const Int_t x)
725 {
726  Double_t x2 = gPad->PixeltoX(x);
727  if (x2<fX1-fR1) return;
728 
729  fR1 = (x2-fX1+fR1)*0.5;
730  fX1 = x2-fR1;
731 }
732 
733 ////////////////////////////////////////////////////////////////////////////////
734 /// Set top of BoundingBox to a value (resize in y direction on top)
735 
736 void TEllipse::SetBBoxY1(const Int_t y)
737 {
738  Double_t y1 = gPad->PixeltoY(y-gPad->VtoPixel(0));
739  if (y1<fY1-fR2) return;
740 
741  fR2 = (y1-fY1+fR2)*0.5;
742  fY1 = y1-fR2;
743 }
744 
745 ////////////////////////////////////////////////////////////////////////////////
746 /// Set bottom of BoundingBox to a value
747 /// (resize in y direction on bottom)
748 
749 void TEllipse::SetBBoxY2(const Int_t y)
750 {
751  Double_t y2 = gPad->PixeltoY(y-gPad->VtoPixel(0));
752 
753  if (y2>fY1+fR2) return;
754 
755  fR2 = (fY1+fR2-y2)*0.5;
756  fY1 = y2+fR2;
757 }
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute action corresponding to one event.
Definition: TEllipse.cxx:198
UShort_t fWidth
Definition: GuiTypes.h:362
Bool_t IsReading() const
Definition: TBuffer.h:81
double dist(Rotation3D const &r1, Rotation3D const &r2)
Definition: 3DDistances.cxx:48
void SetX(SCoord_t x)
Definition: TPoint.h:49
Short_t fY
Definition: GuiTypes.h:361
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
virtual void PaintEllipse(Double_t x1, Double_t y1, Double_t r1, Double_t r2, Double_t phimin, Double_t phimax, Double_t theta, Option_t *option="")
Draw this ellipse with new coordinates.
Definition: TEllipse.cxx:527
virtual void DrawEllipse(Double_t x1, Double_t y1, Double_t r1, Double_t r2, Double_t phimin, Double_t phimax, Double_t theta, Option_t *option="")
Draw this ellipse with new coordinates.
Definition: TEllipse.cxx:175
Double_t GetTheta() const
Definition: TEllipse.h:55
UShort_t fHeight
Definition: GuiTypes.h:362
TAttLine()
AttLine default constructor.
Definition: TAttLine.cxx:128
Double_t fX1
X coordinate of centre.
Definition: TEllipse.h:27
short Version_t
Definition: RtypesCore.h:61
float Float_t
Definition: RtypesCore.h:53
const char Option_t
Definition: RtypesCore.h:62
Double_t fY1
Y coordinate of centre.
Definition: TEllipse.h:28
virtual void ls(Option_t *option="") const
List this ellipse with its attributes.
Definition: TEllipse.cxx:510
virtual void SetBBoxCenterX(const Int_t x)
Set X coordinate of the center of the Ellipse.
Definition: TEllipse.cxx:695
Bool_t TestBit(UInt_t f) const
Definition: TObject.h:159
Double_t GetY1() const
Definition: TEllipse.h:50
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save primitive as a C++ statement(s) on output stream out.
Definition: TEllipse.cxx:588
Buffer base class used for serializing objects.
Definition: TBuffer.h:40
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)=0
#define gROOT
Definition: TROOT.h:375
TEllipse()
Ellipse default constructor.
Definition: TEllipse.cxx:56
Basic string class.
Definition: TString.h:129
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:168
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1099
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute distance from point px,py to an ellipse.
Definition: TEllipse.cxx:132
SCoord_t GetY() const
Definition: TPoint.h:48
SCoord_t GetX() const
Definition: TPoint.h:47
virtual void SetBBoxX1(const Int_t x)
Set left hand side of BoundingBox to a value (resize in x direction on left)
Definition: TEllipse.cxx:712
virtual void Modify()
Change current line attributes if necessary.
Definition: TAttLine.cxx:232
virtual void SetR1(Double_t r1)
Definition: TEllipse.h:65
void SetY(SCoord_t y)
Definition: TPoint.h:50
Short_t Abs(Short_t d)
Definition: TMathBase.h:108
virtual Width_t GetLineWidth() const
Return the line width.
Definition: TAttLine.h:35
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 Paint(Option_t *option="")
Paint this ellipse with its current attributes.
Definition: TEllipse.cxx:519
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition: TObject.cxx:112
virtual Style_t GetLineStyle() const
Return the line style.
Definition: TAttLine.h:34
static const double x2[5]
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 SetBBoxY1(const Int_t y)
Set top of BoundingBox to a value (resize in y direction on top)
Definition: TEllipse.cxx:737
TObject()
TObject constructor.
Definition: TObject.h:213
virtual void SetNoEdges(Bool_t noEdges=kTRUE)
if noEdges = kTRUE the lines connecting the center to the edges will not be drawn.
Definition: TEllipse.cxx:620
Short_t fX
Definition: GuiTypes.h:361
virtual void Draw(Option_t *option="")
Draw this ellipse with its current attributes.
Definition: TEllipse.cxx:167
Definition: TPoint.h:31
Double_t fTheta
Rotation angle (degrees)
Definition: TEllipse.h:33
Double_t fR1
first radius
Definition: TEllipse.h:29
Double_t fPhimax
Maximum angle (degrees)
Definition: TEllipse.h:32
virtual void SetY1(Double_t y1)
Definition: TEllipse.h:69
virtual void SetBBoxX2(const Int_t x)
Set right hand side of BoundingBox to a value (resize in x direction on right)
Definition: TEllipse.cxx:725
virtual void SetBBoxCenterY(const Int_t y)
Set Y coordinate of the center of the Ellipse.
Definition: TEllipse.cxx:703
virtual void Print(Option_t *option="") const
Dump this ellipse with its attributes.
Definition: TEllipse.cxx:576
virtual void SetX1(Double_t x1)
Definition: TEllipse.h:68
unsigned int r1[N_CITIES]
Definition: simanTSP.cxx:321
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 SetR2(Double_t r2)
Definition: TEllipse.h:66
static void IndentLevel()
Functions used by ls() to indent an object hierarchy.
Definition: TROOT.cxx:2632
const Double_t kPI
Definition: TEllipse.cxx:22
Double_t fR2
second radius
Definition: TEllipse.h:30
#define gVirtualX
Definition: TVirtualX.h:350
Double_t Cos(Double_t)
Definition: TMath.h:551
const Bool_t kFALSE
Definition: RtypesCore.h:92
virtual Color_t GetLineColor() const
Return the line color.
Definition: TAttLine.h:33
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
virtual ~TEllipse()
Ellipse default destructor.
Definition: TEllipse.cxx:86
static const double x1[5]
#define ClassImp(name)
Definition: Rtypes.h:336
double Double_t
Definition: RtypesCore.h:55
virtual TPoint GetBBoxCenter()
Return the center of the Ellipse as TPoint in pixels.
Definition: TEllipse.cxx:675
virtual void SetBBoxCenter(const TPoint &p)
Set center of the Ellipse.
Definition: TEllipse.cxx:686
Double_t y[n]
Definition: legend1.C:17
virtual Color_t GetFillColor() const
Return the fill area color.
Definition: TAttFill.h:30
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:572
Double_t GetX1() const
Definition: TEllipse.h:49
Draw Ellipses.
Definition: TEllipse.h:24
Mother of all ROOT objects.
Definition: TObject.h:37
void Copy(TObject &ellipse) const
Copy this ellipse to ellipse.
Definition: TEllipse.cxx:109
TAttFill()
AttFill default constructor.
Definition: TAttFill.cxx:172
Abstract base class for elements drawn in the editor.
Definition: TAttBBox2D.h:19
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:200
virtual void SetBBoxY2(const Int_t y)
Set bottom of BoundingBox to a value (resize in y direction on bottom)
Definition: TEllipse.cxx:750
Double_t Sin(Double_t)
Definition: TMath.h:548
virtual Rectangle_t GetBBox()
Return the bounding Box of the Ellipse, currently not taking into account the rotating angle...
Definition: TEllipse.cxx:662
#define gPad
Definition: TVirtualPad.h:284
Bool_t GetNoEdges() const
Return kTRUE if kNoEdges bit is set, kFALSE otherwise.
Definition: TEllipse.cxx:610
void ResetBit(UInt_t f)
Definition: TObject.h:158
virtual Style_t GetFillStyle() const
Return the fill area style.
Definition: TAttFill.h:31
Double_t Sqrt(Double_t x)
Definition: TMath.h:591
virtual const char * GetName() const
Returns name of object.
Definition: TObject.cxx:364
Double_t fPhimin
Minimum angle (degrees)
Definition: TEllipse.h:31
const Bool_t kTRUE
Definition: RtypesCore.h:91
const Int_t n
Definition: legend1.C:16
Line Attributes class.
Definition: TAttLine.h:18
unsigned int r2[N_CITIES]
Definition: simanTSP.cxx:322
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0
void Copy(TAttFill &attfill) const
Copy this fill attributes to a new TAttFill.
Definition: TAttFill.cxx:200