Logo ROOT   6.08/07
Reference Guide
TCrown.cxx
Go to the documentation of this file.
1 // @(#)root/graf:$Id$
2 // Author: Rene Brun 108/08/2002
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 "TCrown.h"
16 #include "TVirtualPad.h"
17 
19 
20 /** \class TCrown
21 \ingroup BasicGraphics
22 
23 To draw a Crown.
24 
25 A crown is specified with the position of its centre, its inner/outer radius
26 a minimum and maximum angle. The attributes of the outline line are given via
27 TAttLine. The attributes of the fill area are given via TAttFill.
28 
29 Example:
30 
31 Begin_Macro(source)
32 {
33  TCanvas *c1 = new TCanvas("c1","c1",400,400);
34  TCrown cr1(.5,.5,.3,.4);
35  cr1.SetLineStyle(2);
36  cr1.SetLineWidth(4);
37  cr1.Draw();
38  TCrown cr2(.5,.5,.2,.3,45,315);
39  cr2.SetFillColor(38);
40  cr2.SetFillStyle(3010);
41  cr2.Draw();
42  TCrown cr3(.5,.5,.2,.3,-45,45);
43  cr3.SetFillColor(50);
44  cr3.SetFillStyle(3025);
45  cr3.Draw();
46  TCrown cr4(.5,.5,.0,.2);
47  cr4.SetFillColor(4);
48  cr4.SetFillStyle(3008);
49  cr4.Draw();
50  return c1;
51 }
52 End_Macro
53 */
54 
55 ////////////////////////////////////////////////////////////////////////////////
56 /// Crown default constructor.
57 
59 {
60 }
61 
62 ////////////////////////////////////////////////////////////////////////////////
63 /// Crown normal constructor.
64 ///
65 /// \param[in] x1,y1 coordinates of centre of crown
66 /// \param[in] radin inner crown radius
67 /// \param[in] radout outer crown radius
68 /// \param[in] phimin min angle in degrees (default is 0)
69 /// \param[in] phimax max angle in degrees (default is 360)
70 ///
71 /// When a crown sector only is drawn, the lines connecting the center
72 /// of the crown to the edges are drawn by default. One can specify
73 /// the drawing option "only" to not draw these lines.
74 
76  :TEllipse(x1,y1,radin,radout,phimin,phimax,0)
77 {
78 }
79 
80 ////////////////////////////////////////////////////////////////////////////////
81 /// Crown copy constructor.
82 
83 TCrown::TCrown(const TCrown &crown) : TEllipse(crown)
84 {
85 
86  ((TCrown&)crown).Copy(*this);
87 }
88 
89 ////////////////////////////////////////////////////////////////////////////////
90 /// Crown default destructor.
91 
93 {
94 }
95 
96 ////////////////////////////////////////////////////////////////////////////////
97 /// Copy this crown to crown.
98 
99 void TCrown::Copy(TObject &crown) const
100 {
101 
102  TEllipse::Copy(crown);
103 }
104 
105 ////////////////////////////////////////////////////////////////////////////////
106 /// Compute distance from point px,py to a crown.
107 ///
108 /// If crown is filled, return OK if we are inside
109 /// otherwise, crown is found if near the crown edges.
110 
112 {
113 
114  const Double_t kPI = TMath::Pi();
115  Double_t x = gPad->PadtoX(gPad->AbsPixeltoX(px)) - fX1;
116  Double_t y = gPad->PadtoY(gPad->AbsPixeltoY(py)) - fY1;
117 
118  Double_t r1 = fR1;
119  Double_t r2 = fR2;
120  Double_t r = TMath::Sqrt(x*x+y*y);
121 
122  if (r1>r2) {
123  r1 = fR2;
124  r2 = fR1;
125  }
126 
127  Int_t dist = 9999;
128  if (r > r2) return dist;
129  if (r < r1) return dist;
130  if (fPhimax-fPhimin < 360) {
131  Double_t phi = 180*TMath::ACos(x/r)/kPI;
132  if (y<0) phi = 360-phi;
133  Double_t phi1 = fPhimin;
134  Double_t phi2 = fPhimax;
135  if (phi1<0) phi1=phi1+360;
136  if (phi2<0) phi2=phi2+360;
137  if (phi2<phi1) {
138  if (phi < phi1 && phi > phi2) return dist;
139  } else {
140  if (phi < phi1) return dist;
141  if (phi > phi2) return dist;
142  }
143  }
144 
145  if (GetFillColor() && GetFillStyle()) {
146  return 0;
147  } else {
148  if (TMath::Abs(r2-r)/r2 < 0.02) return 0;
149  if (TMath::Abs(r1-r)/r1 < 0.02) return 0;
150  }
151  return dist;
152 }
153 
154 ////////////////////////////////////////////////////////////////////////////////
155 /// Draw this crown with new coordinates.
156 
157 void TCrown::DrawCrown(Double_t x1, Double_t y1,Double_t radin,Double_t radout,Double_t phimin,Double_t phimax,Option_t *option)
158 {
159 
160  TCrown *newcrown = new TCrown(x1, y1, radin, radout, phimin, phimax);
161  TAttLine::Copy(*newcrown);
162  TAttFill::Copy(*newcrown);
163  newcrown->SetBit(kCanDelete);
164  newcrown->AppendPad(option);
165 }
166 
167 ////////////////////////////////////////////////////////////////////////////////
168 /// Execute action corresponding to one event
169 ///
170 /// For the time being TEllipse::ExecuteEvent is used.
171 
173 {
174 
175  TEllipse::ExecuteEvent(event,px,py);
176 }
177 
178 ////////////////////////////////////////////////////////////////////////////////
179 /// Paint this crown with its current attributes.
180 
182 {
183 
184  const Double_t kPI = TMath::Pi();
185  const Int_t np = 40;
186  static Double_t x[2*np+3], y[2*np+3];
189 
190  Double_t angle,dx,dy;
191  Double_t dphi = (fPhimax-fPhimin)*kPI/(180*np);
192  Double_t ct = TMath::Cos(kPI*fTheta/180);
193  Double_t st = TMath::Sin(kPI*fTheta/180);
194  Int_t i;
195  //compute outer points
196  for (i=0;i<=np;i++) {
197  angle = fPhimin*kPI/180 + Double_t(i)*dphi;
198  dx = fR2*TMath::Cos(angle);
199  dy = fR2*TMath::Sin(angle);
200  x[i] = fX1 + dx*ct - dy*st;
201  y[i] = fY1 + dx*st + dy*ct;
202  }
203  //compute inner points
204  for (i=0;i<=np;i++) {
205  angle = fPhimin*kPI/180 + Double_t(i)*dphi;
206  dx = fR1*TMath::Cos(angle);
207  dy = fR1*TMath::Sin(angle);
208  x[2*np-i+1] = fX1 + dx*ct - dy*st;
209  y[2*np-i+1] = fY1 + dx*st + dy*ct;
210  }
211  x[2*np+2] = x[0];
212  y[2*np+2] = y[0];
213  if (fPhimax-fPhimin >= 360 ) {
214  // a complete filled crown
215  if (GetFillColor() && GetFillStyle()) {
216  gPad->PaintFillArea(2*np+2,x,y);
217  }
218  // a complete empty crown
219  if (GetLineStyle()) {
220  gPad->PaintPolyLine(np+1,x,y);
221  gPad->PaintPolyLine(np+1,&x[np+1],&y[np+1]);
222  }
223  } else {
224  //crown segment
225  if (GetFillColor() && GetFillStyle()) gPad->PaintFillArea(2*np+2,x,y);
226  if (GetLineStyle()) gPad->PaintPolyLine(2*np+3,x,y);
227  }
228 }
229 
230 ////////////////////////////////////////////////////////////////////////////////
231 /// Save primitive as a C++ statement(s) on output stream out.
232 
233 void TCrown::SavePrimitive(std::ostream &out, Option_t * /*= ""*/)
234 {
235 
236  out<<" "<<std::endl;
237  if (gROOT->ClassSaved(TCrown::Class())) {
238  out<<" ";
239  } else {
240  out<<" TCrown *";
241  }
242  out<<"crown = new TCrown("<<fX1<<","<<fY1<<","<<fR1<<","<<fR2
243  <<","<<fPhimin<<","<<fPhimax<<");"<<std::endl;
244 
245  SaveFillAttributes(out,"crown",0,1001);
246  SaveLineAttributes(out,"crown",1,1,1);
247 
248  if (GetNoEdges()) out<<" crown->SetNoEdges();"<<std::endl;
249 
250  out<<" crown->Draw();"<<std::endl;
251 }
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute action corresponding to one event.
Definition: TEllipse.cxx:198
double dist(Rotation3D const &r1, Rotation3D const &r2)
Definition: 3DDistances.cxx:48
To draw a Crown.
Definition: TCrown.h:21
Double_t fX1
X coordinate of centre.
Definition: TEllipse.h:39
const char Option_t
Definition: RtypesCore.h:62
Double_t fY1
Y coordinate of centre.
Definition: TEllipse.h:40
virtual void Paint(Option_t *option="")
Paint this crown with its current attributes.
Definition: TCrown.cxx:181
#define gROOT
Definition: TROOT.h:364
int Int_t
Definition: RtypesCore.h:41
virtual void Modify()
Change current line attributes if necessary.
Definition: TAttLine.cxx:232
const char * Class
Definition: TXMLSetup.cxx:64
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:739
virtual Int_t DistancetoPrimitive(Int_t px, Int_t py)
Compute distance from point px,py to a crown.
Definition: TCrown.cxx:111
if object in a list can be deleted
Definition: TObject.h:56
TCrown()
Crown default constructor.
Definition: TCrown.cxx:58
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition: TObject.cxx:165
virtual Style_t GetLineStyle() const
Return the line style.
Definition: TAttLine.h:40
Double_t x[n]
Definition: legend1.C:17
void Copy(TAttLine &attline) const
Copy this line attributes to a new TAttLine.
Definition: TAttLine.cxx:162
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
void Copy(TObject &crown) const
Copy this crown to crown.
Definition: TCrown.cxx:99
virtual void DrawCrown(Double_t x1, Double_t y1, Double_t radin, Double_t radout, Double_t phimin=0, Double_t phimax=360, Option_t *option="")
Draw this crown with new coordinates.
Definition: TCrown.cxx:157
Double_t fTheta
Rotation angle (degrees)
Definition: TEllipse.h:45
Double_t fR1
first radius
Definition: TEllipse.h:41
Double_t fPhimax
Maximum angle (degrees)
Definition: TEllipse.h:44
TRandom2 r(17)
unsigned int r1[N_CITIES]
Definition: simanTSP.cxx:321
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save primitive as a C++ statement(s) on output stream out.
Definition: TCrown.cxx:233
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
Double_t ACos(Double_t)
Definition: TMath.h:445
const Double_t kPI
Definition: TEllipse.cxx:22
Double_t fR2
second radius
Definition: TEllipse.h:42
Double_t Cos(Double_t)
Definition: TMath.h:424
Double_t Pi()
Definition: TMath.h:44
static const double x1[5]
#define ClassImp(name)
Definition: Rtypes.h:279
double Double_t
Definition: RtypesCore.h:55
Double_t y[n]
Definition: legend1.C:17
virtual Color_t GetFillColor() const
Return the fill area color.
Definition: TAttFill.h:35
Draw Ellipses.
Definition: TEllipse.h:36
virtual ~TCrown()
Crown default destructor.
Definition: TCrown.cxx:92
Mother of all ROOT objects.
Definition: TObject.h:37
void Copy(TObject &ellipse) const
Copy this ellipse to ellipse.
Definition: TEllipse.cxx:109
Double_t Sin(Double_t)
Definition: TMath.h:421
virtual void ExecuteEvent(Int_t event, Int_t px, Int_t py)
Execute action corresponding to one event.
Definition: TCrown.cxx:172
#define gPad
Definition: TVirtualPad.h:289
Bool_t GetNoEdges() const
Return kTRUE if kNoEdges bit is set, kFALSE otherwise.
Definition: TEllipse.cxx:610
virtual Style_t GetFillStyle() const
Return the fill area style.
Definition: TAttFill.h:36
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
Double_t fPhimin
Minimum angle (degrees)
Definition: TEllipse.h:43
unsigned int r2[N_CITIES]
Definition: simanTSP.cxx:322
void Copy(TAttFill &attfill) const
Copy this fill attributes to a new TAttFill.
Definition: TAttFill.cxx:200