Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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 <iostream>
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
23To draw a Crown.
24
25A crown is specified with the position of its centre, its inner/outer radius
26a minimum and maximum angle. The attributes of the outline line are given via
27TAttLine. The attributes of the fill area are given via TAttFill.
28
29Example:
30
31Begin_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}
52End_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
83TCrown::TCrown(const TCrown &crown) : TEllipse(crown)
84{
85 crown.TCrown::Copy(*this);
86}
87
88////////////////////////////////////////////////////////////////////////////////
89/// Crown default destructor.
90
92{
93}
94
95////////////////////////////////////////////////////////////////////////////////
96/// Copy this crown to crown.
97
98void TCrown::Copy(TObject &crown) const
99{
100 TEllipse::Copy(crown);
101}
102
103////////////////////////////////////////////////////////////////////////////////
104/// Compute distance from point px,py to a crown.
105///
106/// If crown is filled, return OK if we are inside
107/// otherwise, crown is found if near the crown edges.
108
110{
111 if (!gPad) return 9999;
112 const Double_t kPI = TMath::Pi();
113 Double_t x = gPad->PadtoX(gPad->AbsPixeltoX(px)) - fX1;
114 Double_t y = gPad->PadtoY(gPad->AbsPixeltoY(py)) - fY1;
115
116 Double_t r1 = fR1;
117 Double_t r2 = fR2;
119
120 if (r1>r2) {
121 r1 = fR2;
122 r2 = fR1;
123 }
124
125 Int_t dist = 9999;
126 if (r > r2) return dist;
127 if (r < r1) return dist;
128 if (fPhimax-fPhimin < 360) {
129 Double_t phi = 180*TMath::ACos(x/r)/kPI;
130 if (y<0) phi = 360-phi;
131 Double_t phi1 = fPhimin;
132 Double_t phi2 = fPhimax;
133 if (phi1<0) phi1=phi1+360;
134 if (phi2<0) phi2=phi2+360;
135 if (phi2<phi1) {
136 if (phi < phi1 && phi > phi2) return dist;
137 } else {
138 if (phi < phi1) return dist;
139 if (phi > phi2) return dist;
140 }
141 }
142
143 if (GetFillColor() && GetFillStyle()) {
144 return 0;
145 } else {
146 if (TMath::Abs(r2-r)/r2 < 0.02) return 0;
147 if (TMath::Abs(r1-r)/r1 < 0.02) return 0;
148 }
149 return dist;
150}
151
152////////////////////////////////////////////////////////////////////////////////
153/// Draw this crown with new coordinates.
154
156{
157 TCrown *newcrown = new TCrown(x1, y1, radin, radout, phimin, phimax);
158 TAttLine::Copy(*newcrown);
159 TAttFill::Copy(*newcrown);
160 newcrown->SetBit(kCanDelete);
161 newcrown->AppendPad(option);
162 return newcrown;
163}
164
165////////////////////////////////////////////////////////////////////////////////
166/// Execute action corresponding to one event
167///
168/// For the time being TEllipse::ExecuteEvent is used.
169
171{
172
173 TEllipse::ExecuteEvent(event,px,py);
174}
175
176////////////////////////////////////////////////////////////////////////////////
177/// Return 1 if the point (x,y) is inside the polygon defined by
178/// the crown 0 otherwise.
179
181{
182
183 const Double_t kPI = TMath::Pi();
184 const Int_t np = 40;
185 static Double_t xc[2 * np + 3], yc[2 * np +3];
186
187 Double_t angle, dx, dy;
188 Double_t dphi = (fPhimax - fPhimin) * kPI / (180 * np);
189 Double_t ct = TMath::Cos(kPI * fTheta / 180);
190 Double_t st = TMath::Sin(kPI * fTheta / 180);
191 Int_t i;
192
193 // compute outer points
194 for (i = 0; i <= np; i++) {
195 angle = fPhimin * kPI / 180 + Double_t(i) * dphi;
196 dx = fR2 * TMath::Cos(angle);
197 dy = fR2 * TMath::Sin(angle);
198 xc[i] = fX1 + dx * ct - dy * st;
199 yc[i] = fY1 + dx * st + dy * ct;
200 }
201 // compute inner points
202 for (i = 0; i <= np; i++) {
203 angle = fPhimin * kPI / 180 + Double_t(i) * dphi;
204 dx = fR1 * TMath::Cos(angle);
205 dy = fR1 * TMath::Sin(angle);
206 xc[2 * np - i + 1] = fX1 + dx * ct - dy * st;
207 yc[2 * np - i + 1] = fY1 + dx * st + dy * ct;
208 }
209 xc[2 * np + 2] = xc[0];
210 yc[2 * np + 2] = yc[0];
211
212 return (Int_t)TMath::IsInside(x, y, 2 * np + 3, xc, yc);
213}
214
215////////////////////////////////////////////////////////////////////////////////
216/// Paint this crown with its current attributes.
217
219{
220 if (!gPad) return;
221
222 const Double_t kPI = TMath::Pi();
223 const Int_t np = 40;
224 static Double_t x[2 * np + 3], y[2 * np + 3];
227
228 Double_t angle,dx,dy;
229 Double_t dphi = (fPhimax - fPhimin) * kPI / (180 * np);
230 Double_t ct = TMath::Cos(kPI * fTheta / 180);
231 Double_t st = TMath::Sin(kPI * fTheta / 180);
232 Int_t i;
233 // compute outer points
234 for (i = 0; i <= np; i++) {
235 angle = fPhimin * kPI / 180 + Double_t(i) * dphi;
236 dx = fR2 * TMath::Cos(angle);
237 dy = fR2 * TMath::Sin(angle);
238 x[i] = fX1 + dx * ct - dy * st;
239 y[i] = fY1 + dx * st + dy * ct;
240 }
241 // compute inner points
242 for (i = 0; i <= np; i++) {
243 angle = fPhimin * kPI / 180 + Double_t(i) * dphi;
244 dx = fR1 * TMath::Cos(angle);
245 dy = fR1 * TMath::Sin(angle);
246 x[2 * np - i + 1] = fX1 + dx * ct - dy * st;
247 y[2 * np - i + 1] = fY1 + dx * st + dy * ct;
248 }
249 x[2 * np + 2] = x[0];
250 y[2 * np + 2] = y[0];
251 if (fPhimax - fPhimin >= 360 ) {
252 // a complete filled crown
253 if (GetFillColor() && GetFillStyle()) {
254 gPad->PaintFillArea(2 * np + 2,x,y);
255 }
256 // a complete empty crown
257 if (GetLineStyle()) {
258 gPad->PaintPolyLine(np + 1,x,y);
259 gPad->PaintPolyLine(np + 1,&x[np + 1],&y[np + 1]);
260 }
261 } else {
262 //crown segment
263 if (GetFillColor() && GetFillStyle()) gPad->PaintFillArea(2 * np + 2,x,y);
264 if (GetLineStyle()) gPad->PaintPolyLine(2 * np + 3,x,y);
265 }
266}
267
268////////////////////////////////////////////////////////////////////////////////
269/// Save primitive as a C++ statement(s) on output stream out.
270
271void TCrown::SavePrimitive(std::ostream &out, Option_t * /*= ""*/)
272{
273
274 out<<" "<<std::endl;
275 if (gROOT->ClassSaved(TCrown::Class())) {
276 out<<" ";
277 } else {
278 out<<" TCrown *";
279 }
280 out<<"crown = new TCrown("<<fX1<<","<<fY1<<","<<fR1<<","<<fR2
281 <<","<<fPhimin<<","<<fPhimax<<");"<<std::endl;
282
283 SaveFillAttributes(out,"crown",0,1001);
284 SaveLineAttributes(out,"crown",1,1,1);
285
286 if (GetNoEdges()) out<<" crown->SetNoEdges();"<<std::endl;
287
288 out<<" crown->Draw();"<<std::endl;
289}
double Double_t
Definition RtypesCore.h:59
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:377
constexpr Double_t kPI
Definition TEllipse.cxx:24
Option_t Option_t option
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t np
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint angle
Option_t Option_t TPoint TPoint const char y1
#define gROOT
Definition TROOT.h:406
#define gPad
virtual Color_t GetFillColor() const
Return the fill area color.
Definition TAttFill.h:30
void Copy(TAttFill &attfill) const
Copy this fill attributes to a new TAttFill.
Definition TAttFill.cxx:207
virtual Style_t GetFillStyle() const
Return the fill area style.
Definition TAttFill.h:31
virtual void Modify()
Change current fill area attributes if necessary.
Definition TAttFill.cxx:216
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:239
virtual Style_t GetLineStyle() const
Return the line style.
Definition TAttLine.h:34
virtual void Modify()
Change current line attributes if necessary.
Definition TAttLine.cxx:247
void Copy(TAttLine &attline) const
Copy this line attributes to a new TAttLine.
Definition TAttLine.cxx:177
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:275
To draw a Crown.
Definition TCrown.h:19
~TCrown() override
Crown default destructor.
Definition TCrown.cxx:91
Int_t DistancetoPrimitive(Int_t px, Int_t py) override
Compute distance from point px,py to a crown.
Definition TCrown.cxx:109
virtual TCrown * 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:155
void ExecuteEvent(Int_t event, Int_t px, Int_t py) override
Execute action corresponding to one event.
Definition TCrown.cxx:170
static TClass * Class()
TCrown()
Crown default constructor.
Definition TCrown.cxx:58
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save primitive as a C++ statement(s) on output stream out.
Definition TCrown.cxx:271
Int_t IsInside(Double_t x, Double_t y) const
Return 1 if the point (x,y) is inside the polygon defined by the crown 0 otherwise.
Definition TCrown.cxx:180
void Paint(Option_t *option="") override
Paint this crown with its current attributes.
Definition TCrown.cxx:218
void Copy(TObject &crown) const override
Copy this crown to crown.
Definition TCrown.cxx:98
Draw Ellipses.
Definition TEllipse.h:23
void ExecuteEvent(Int_t event, Int_t px, Int_t py) override
Execute action corresponding to one event.
Definition TEllipse.cxx:201
Double_t fPhimax
Maximum angle (degrees)
Definition TEllipse.h:31
Double_t fX1
X coordinate of centre.
Definition TEllipse.h:26
Bool_t GetNoEdges() const
Return kTRUE if kNoEdges bit is set, kFALSE otherwise.
Definition TEllipse.cxx:646
Double_t fY1
Y coordinate of centre.
Definition TEllipse.h:27
Double_t fTheta
Rotation angle (degrees)
Definition TEllipse.h:32
Double_t fR2
second radius
Definition TEllipse.h:29
void Copy(TObject &ellipse) const override
Copy this ellipse to ellipse.
Definition TEllipse.cxx:110
Double_t fPhimin
Minimum angle (degrees)
Definition TEllipse.h:30
Double_t fR1
first radius
Definition TEllipse.h:28
Mother of all ROOT objects.
Definition TObject.h:41
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition TObject.cxx:184
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:780
@ kCanDelete
if object in a list can be deleted
Definition TObject.h:62
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
Double_t ACos(Double_t)
Returns the principal value of the arc cosine of x, expressed in radians.
Definition TMath.h:632
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:1233
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:662
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Definition TMath.h:594
constexpr Double_t Pi()
Definition TMath.h:37
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:588
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123