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 "TMath.h"
14#include "TCrown.h"
15#include "TVirtualPad.h"
16
18
19/** \class TCrown
20\ingroup BasicGraphics
21
22To draw a Crown.
23
24A crown is specified with the position of its centre, its inner/outer radius
25a minimum and maximum angle. The attributes of the outline line are given via
26TAttLine. The attributes of the fill area are given via TAttFill.
27
28Example:
29
30Begin_Macro(source)
31{
32 TCanvas *c1 = new TCanvas("c1","c1",400,400);
33 auto cr1 = new TCrown(.5,.5,.3,.4);
34 cr1->SetLineStyle(2);
35 cr1->SetLineWidth(4);
36 cr1->Draw();
37 auto cr2 = new TCrown(.5,.5,.2,.3,45,315);
38 cr2->SetFillColor(38);
39 cr2->SetFillStyle(3010);
40 cr2->Draw();
41 auto cr3 = new TCrown(.5,.5,.2,.3,-45,45);
42 cr3->SetFillColor(50);
43 cr3->SetFillStyle(3025);
44 cr3->Draw();
45 auto cr4 = new TCrown(.5,.5,.0,.2);
46 cr4->SetFillColor(4);
47 cr4->SetFillStyle(3008);
48 cr4->Draw();
49 return c1;
50}
51End_Macro
52*/
53
54////////////////////////////////////////////////////////////////////////////////
55/// Crown default constructor.
56
60
61////////////////////////////////////////////////////////////////////////////////
62/// Crown normal constructor.
63///
64/// \param[in] x1,y1 coordinates of centre of crown
65/// \param[in] radin inner crown radius
66/// \param[in] radout outer crown radius
67/// \param[in] phimin min angle in degrees (default is 0)
68/// \param[in] phimax max angle in degrees (default is 360)
69///
70/// When a crown sector only is drawn, the lines connecting the center
71/// of the crown to the edges are drawn by default. One can specify
72/// the drawing option "only" to not draw these lines.
73
78
79////////////////////////////////////////////////////////////////////////////////
80/// Crown copy constructor.
81
83{
84 crown.TCrown::Copy(*this);
85}
86
87////////////////////////////////////////////////////////////////////////////////
88/// Crown default destructor.
89
93
94////////////////////////////////////////////////////////////////////////////////
95/// Copy this crown to crown.
96
98{
100 static_cast<TCrown &>(crown).fYXRatio = fYXRatio;
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) / fYXRatio;
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;
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
165
166////////////////////////////////////////////////////////////////////////////////
167/// Execute action corresponding to one event
168///
169/// For the time being TEllipse::ExecuteEvent is used.
170
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 const Double_t kPI = TMath::Pi();
183 const Int_t np = 40;
184 static Double_t xc[2 * np + 3], yc[2 * np + 3];
185
187 Double_t dphi = (fPhimax - fPhimin) * kPI / (180 * np);
188 Double_t ct = TMath::Cos(kPI * fTheta / 180);
189 Double_t st = TMath::Sin(kPI * fTheta / 180);
190 Int_t i;
191
192 // compute outer points
193 for (i = 0; i <= np; i++) {
194 angle = fPhimin * kPI / 180 + Double_t(i) * dphi;
195 dx = fR2 * TMath::Cos(angle);
197 xc[i] = fX1 + dx * ct - dy * st;
198 yc[i] = fY1 + dx * st + dy * ct;
199 }
200 // compute inner points
201 for (i = 0; i <= np; i++) {
202 angle = fPhimin * kPI / 180 + Double_t(i) * dphi;
203 dx = fR1 * TMath::Cos(angle);
205 xc[2 * np - i + 1] = fX1 + dx * ct - dy * st;
206 yc[2 * np - i + 1] = fY1 + dx * st + dy * ct;
207 }
208 xc[2 * np + 2] = xc[0];
209 yc[2 * np + 2] = yc[0];
210
211 return (Int_t)TMath::IsInside(x, y, 2 * np + 3, xc, yc);
212}
213
214////////////////////////////////////////////////////////////////////////////////
215/// Paint this crown with its current attributes.
216
218{
219 if (!gPad) return;
220
221 const Double_t kPI = TMath::Pi();
222 const Int_t np = 40;
223 static Double_t x[2 * np + 3], y[2 * np + 3];
226
228 Double_t dphi = (fPhimax - fPhimin) * kPI / (180 * np);
229 Double_t ct = TMath::Cos(kPI * fTheta / 180);
230 Double_t st = TMath::Sin(kPI * fTheta / 180);
231 Int_t i;
232 // compute outer points
233 for (i = 0; i <= np; i++) {
234 angle = fPhimin * kPI / 180 + Double_t(i) * dphi;
235 dx = fR2 * TMath::Cos(angle);
237 x[i] = fX1 + dx * ct - dy * st;
238 y[i] = fY1 + dx * st + dy * ct;
239 }
240 // compute inner points
241 for (i = 0; i <= np; i++) {
242 angle = fPhimin * kPI / 180 + Double_t(i) * dphi;
243 dx = fR1 * TMath::Cos(angle);
245 x[2 * np - i + 1] = fX1 + dx * ct - dy * st;
246 y[2 * np - i + 1] = fY1 + dx * st + dy * ct;
247 }
248 x[2 * np + 2] = x[0];
249 y[2 * np + 2] = y[0];
250 if (fPhimax - fPhimin >= 360 ) {
251 // a complete filled crown
252 if (GetFillColor() && GetFillStyle()) {
253 gPad->PaintFillArea(2 * np + 2,x,y);
254 }
255 // a complete empty crown
256 if (GetLineStyle()) {
257 gPad->PaintPolyLine(np + 1,x,y);
258 gPad->PaintPolyLine(np + 1,&x[np + 1],&y[np + 1]);
259 }
260 } else {
261 //crown segment
262 if (GetFillColor() && GetFillStyle()) gPad->PaintFillArea(2 * np + 2,x,y);
263 if (GetLineStyle()) gPad->PaintPolyLine(2 * np + 3,x,y);
264 }
265}
266
267////////////////////////////////////////////////////////////////////////////////
268/// Save primitive as a C++ statement(s) on output stream out.
269
270void TCrown::SavePrimitive(std::ostream &out, Option_t * /*= ""*/)
271{
272 SavePrimitiveConstructor(out, Class(), "crown",
273 TString::Format("%g, %g, %g, %g, %g, %g", fX1, fY1, fR1, fR2, fPhimin, fPhimax));
274
275 SaveFillAttributes(out, "crown", 0, 1001);
276 SaveLineAttributes(out, "crown", 1, 1, 1);
277
278 if (GetNoEdges())
279 out << " crown->SetNoEdges();" << std::endl;
280
281 if (fYXRatio != 1)
282 out << " crown->SetYXRatio(" << fYXRatio << ");" << std::endl;
283
284 out << " crown->Draw();" << std::endl;
285}
double Double_t
Definition RtypesCore.h:59
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:374
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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 gPad
virtual Color_t GetFillColor() const
Return the fill area color.
Definition TAttFill.h:31
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:32
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:36
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:90
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
Double_t GetYXRatio() const
Definition TCrown.h:30
void ExecuteEvent(Int_t event, Int_t px, Int_t py) override
Execute action corresponding to one event.
Definition TCrown.cxx:171
static TClass * Class()
TCrown()
Crown default constructor.
Definition TCrown.cxx:57
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save primitive as a C++ statement(s) on output stream out.
Definition TCrown.cxx:270
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:217
void Copy(TObject &crown) const override
Copy this crown to crown.
Definition TCrown.cxx:97
Double_t fYXRatio
Definition TCrown.h:20
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:642
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
static void SavePrimitiveConstructor(std::ostream &out, TClass *cl, const char *variable_name, const char *constructor_agrs="", Bool_t empty_line=kTRUE)
Save object constructor in the output stream "out".
Definition TObject.cxx:771
@ kCanDelete
if object in a list can be deleted
Definition TObject.h:66
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2378
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:636
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:1313
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:666
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Definition TMath.h:598
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:592
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123