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
17
18/** \class TCrown
19\ingroup BasicGraphics
20
21To draw a Crown.
22
23A crown is specified with the position of its centre, its inner/outer radius
24a minimum and maximum angle. The attributes of the outline line are given via
25TAttLine. The attributes of the fill area are given via TAttFill.
26
27Example:
28
29Begin_Macro(source)
30{
31 TCanvas *c1 = new TCanvas("c1","c1",400,400);
32 auto cr1 = new TCrown(.5,.5,.3,.4);
33 cr1->SetLineStyle(2);
34 cr1->SetLineWidth(4);
35 cr1->Draw();
36 auto cr2 = new TCrown(.5,.5,.2,.3,45,315);
37 cr2->SetFillColor(38);
38 cr2->SetFillStyle(3010);
39 cr2->Draw();
40 auto cr3 = new TCrown(.5,.5,.2,.3,-45,45);
41 cr3->SetFillColor(50);
42 cr3->SetFillStyle(3025);
43 cr3->Draw();
44 auto cr4 = new TCrown(.5,.5,.0,.2);
45 cr4->SetFillColor(4);
46 cr4->SetFillStyle(3008);
47 cr4->Draw();
48 return c1;
49}
50End_Macro
51*/
52
53////////////////////////////////////////////////////////////////////////////////
54/// Crown default constructor.
55
59
60////////////////////////////////////////////////////////////////////////////////
61/// Crown normal constructor.
62///
63/// \param[in] x1,y1 coordinates of centre of crown
64/// \param[in] radin inner crown radius
65/// \param[in] radout outer crown radius
66/// \param[in] phimin min angle in degrees (default is 0)
67/// \param[in] phimax max angle in degrees (default is 360)
68///
69/// When a crown sector only is drawn, the lines connecting the center
70/// of the crown to the edges are drawn by default. One can specify
71/// the drawing option "only" to not draw these lines.
72
77
78////////////////////////////////////////////////////////////////////////////////
79/// Crown copy constructor.
80
82{
83 crown.TCrown::Copy(*this);
84}
85
86////////////////////////////////////////////////////////////////////////////////
87/// Crown default destructor.
88
92
93////////////////////////////////////////////////////////////////////////////////
94/// Copy this crown to crown.
95
97{
99 static_cast<TCrown &>(crown).fYXRatio = fYXRatio;
100}
101
102////////////////////////////////////////////////////////////////////////////////
103/// Compute distance from point px,py to a crown.
104///
105/// If crown is filled, return OK if we are inside
106/// otherwise, crown is found if near the crown edges.
107
109{
110 if (!gPad) return 9999;
111 const Double_t kPI = TMath::Pi();
112 Double_t x = gPad->PadtoX(gPad->AbsPixeltoX(px)) - fX1;
113 Double_t y = (gPad->PadtoY(gPad->AbsPixeltoY(py)) - fY1) / fYXRatio;
114
115 Double_t r1 = fR1;
116 Double_t r2 = fR2;
118
119 if (r1>r2) {
120 r1 = fR2;
121 r2 = fR1;
122 }
123
124 Int_t dist = 9999;
125 if (r > r2) return dist;
126 if (r < r1) return dist;
127 if (fPhimax-fPhimin < 360) {
128 Double_t phi = 180*TMath::ACos(x/r)/kPI;
129 if (y<0) phi = 360-phi;
132 if (phi1<0) phi1=phi1+360;
133 if (phi2<0) phi2=phi2+360;
134 if (phi2<phi1) {
135 if (phi < phi1 && phi > phi2) return dist;
136 } else {
137 if (phi < phi1) return dist;
138 if (phi > phi2) return dist;
139 }
140 }
141
142 if (GetFillColor() && GetFillStyle()) {
143 return 0;
144 } else {
145 if (TMath::Abs(r2-r)/r2 < 0.02) return 0;
146 if (TMath::Abs(r1-r)/r1 < 0.02) return 0;
147 }
148 return dist;
149}
150
151////////////////////////////////////////////////////////////////////////////////
152/// Draw this crown with new coordinates.
153
164
165////////////////////////////////////////////////////////////////////////////////
166/// Execute action corresponding to one event
167///
168/// For the time being TEllipse::ExecuteEvent is used.
169
171{
172 TEllipse::ExecuteEvent(event, px, py);
173}
174
175////////////////////////////////////////////////////////////////////////////////
176/// Return 1 if the point (x,y) is inside the polygon defined by
177/// the crown 0 otherwise.
178
180{
181 const Double_t kPI = TMath::Pi();
182 const Int_t np = 40;
183 static Double_t xc[2 * np + 3], yc[2 * np + 3];
184
186 Double_t dphi = (fPhimax - fPhimin) * kPI / (180 * np);
187 Double_t ct = TMath::Cos(kPI * fTheta / 180);
188 Double_t st = TMath::Sin(kPI * fTheta / 180);
189 Int_t i;
190
191 // compute outer points
192 for (i = 0; i <= np; i++) {
193 angle = fPhimin * kPI / 180 + Double_t(i) * dphi;
194 dx = fR2 * TMath::Cos(angle);
196 xc[i] = fX1 + dx * ct - dy * st;
197 yc[i] = fY1 + dx * st + dy * ct;
198 }
199 // compute inner points
200 for (i = 0; i <= np; i++) {
201 angle = fPhimin * kPI / 180 + Double_t(i) * dphi;
202 dx = fR1 * TMath::Cos(angle);
204 xc[2 * np - i + 1] = fX1 + dx * ct - dy * st;
205 yc[2 * np - i + 1] = fY1 + dx * st + dy * ct;
206 }
207 xc[2 * np + 2] = xc[0];
208 yc[2 * np + 2] = yc[0];
209
210 return (Int_t)TMath::IsInside(x, y, 2 * np + 3, xc, yc);
211}
212
213////////////////////////////////////////////////////////////////////////////////
214/// Paint this crown with its current attributes.
215
217{
218 if (!gPad) return;
219
220 const Double_t kPI = TMath::Pi();
221 const Int_t np = 40;
222 static Double_t x[2 * np + 3], y[2 * np + 3];
225
227 Double_t dphi = (fPhimax - fPhimin) * kPI / (180 * np);
228 Double_t ct = TMath::Cos(kPI * fTheta / 180);
229 Double_t st = TMath::Sin(kPI * fTheta / 180);
230 Int_t i;
231 // compute outer points
232 for (i = 0; i <= np; i++) {
233 angle = fPhimin * kPI / 180 + Double_t(i) * dphi;
234 dx = fR2 * TMath::Cos(angle);
236 x[i] = fX1 + dx * ct - dy * st;
237 y[i] = fY1 + dx * st + dy * ct;
238 }
239 // compute inner points
240 for (i = 0; i <= np; i++) {
241 angle = fPhimin * kPI / 180 + Double_t(i) * dphi;
242 dx = fR1 * TMath::Cos(angle);
244 x[2 * np - i + 1] = fX1 + dx * ct - dy * st;
245 y[2 * np - i + 1] = fY1 + dx * st + dy * ct;
246 }
247 x[2 * np + 2] = x[0];
248 y[2 * np + 2] = y[0];
249 if (fPhimax - fPhimin >= 360 ) {
250 // a complete filled crown
251 if (GetFillColor() && GetFillStyle()) {
252 gPad->PaintFillArea(2 * np + 2,x,y);
253 }
254 // a complete empty crown
255 if (GetLineStyle()) {
256 gPad->PaintPolyLine(np + 1,x,y);
257 gPad->PaintPolyLine(np + 1,&x[np + 1],&y[np + 1]);
258 }
259 } else {
260 //crown segment
261 if (GetFillColor() && GetFillStyle()) gPad->PaintFillArea(2 * np + 2,x,y);
262 if (GetLineStyle()) gPad->PaintPolyLine(2 * np + 3,x,y);
263 }
264}
265
266////////////////////////////////////////////////////////////////////////////////
267/// Save primitive as a C++ statement(s) on output stream out.
268
269void TCrown::SavePrimitive(std::ostream &out, Option_t *option)
270{
271 SavePrimitiveConstructor(out, Class(), "crown",
272 TString::Format("%g, %g, %g, %g, %g, %g", fX1, fY1, fR1, fR2, fPhimin, fPhimax));
273
274 SaveFillAttributes(out, "crown", 0, 1001);
275 SaveLineAttributes(out, "crown", 1, 1, 1);
276
277 if (GetNoEdges())
278 out << " crown->SetNoEdges();\n";
279
280 if (fYXRatio != 1)
281 out << " crown->SetYXRatio(" << fYXRatio << ");\n";
282
283 SavePrimitiveDraw(out, "crown", option);
284}
double Double_t
Double 8 bytes.
Definition RtypesCore.h:73
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
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:206
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:215
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:238
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:246
void Copy(TAttLine &attline) const
Copy this line attributes to a new TAttLine.
Definition TAttLine.cxx:176
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:274
To draw a Crown.
Definition TCrown.h:19
~TCrown() override
Crown default destructor.
Definition TCrown.cxx:89
Int_t DistancetoPrimitive(Int_t px, Int_t py) override
Compute distance from point px,py to a crown.
Definition TCrown.cxx:108
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:154
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:170
static TClass * Class()
TCrown()
Crown default constructor.
Definition TCrown.cxx:56
void SavePrimitive(std::ostream &out, Option_t *option="") override
Save primitive as a C++ statement(s) on output stream out.
Definition TCrown.cxx:269
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:179
void Paint(Option_t *option="") override
Paint this crown with its current attributes.
Definition TCrown.cxx:216
void Copy(TObject &crown) const override
Copy this crown to crown.
Definition TCrown.cxx:96
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:200
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:641
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:109
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 SavePrimitiveDraw(std::ostream &out, const char *variable_name, Option_t *option=nullptr)
Save invocation of primitive Draw() method Skipped if option contains "nodraw" string.
Definition TObject.cxx:822
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:68
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:2384
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:643
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:1320
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:673
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Definition TMath.h:605
constexpr Double_t Pi()
Definition TMath.h:40
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:599
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:124