Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGraphPolar.cxx
Go to the documentation of this file.
1// @(#)root/graf:$Id$
2// Author: Sebastian Boser, Mathieu Demaret 02/02/06
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/** \class TGraphPolar
13\ingroup BasicGraphics
14
15To draw a polar graph.
16
17TGraphPolar creates a polar graph (including error bars). A TGraphPolar is
18a TGraphErrors represented in polar coordinates.
19It uses the class TGraphPolargram to draw the polar axis.
20
21Example:
22
23Begin_Macro(source)
24{
25 TCanvas * CPol = new TCanvas("CPol","TGraphPolar Example",500,500);
26
27 Double_t theta[8];
28 Double_t radius[8];
29 Double_t etheta[8];
30 Double_t eradius[8];
31
32 for (int i=0; i<8; i++) {
33 theta[i] = (i+1)*(TMath::Pi()/4.);
34 radius[i] = (i+1)*0.05;
35 etheta[i] = TMath::Pi()/8.;
36 eradius[i] = 0.05;
37 }
38
39 TGraphPolar * grP1 = new TGraphPolar(8, theta, radius, etheta, eradius);
40 grP1->SetTitle("TGraphPolar Example");
41
42 grP1->SetMarkerStyle(20);
43 grP1->SetMarkerSize(2.);
44 grP1->SetMarkerColor(4);
45 grP1->SetLineColor(2);
46 grP1->SetLineWidth(3);
47 grP1->Draw("PE");
48
49 // Update, otherwise GetPolargram returns 0
50 CPol->Update();
51 grP1->GetPolargram()->SetToRadian();
52
53 return CPol;
54}
55End_Macro
56*/
57
58#include "TGraphPolar.h"
59#include "TGraphPolargram.h"
60
62
63////////////////////////////////////////////////////////////////////////////////
64/// TGraphPolar default constructor.
65
67 fOptionAxis(kFALSE),fPolargram(nullptr),fXpol(nullptr),fYpol(nullptr)
68{
69}
70
71////////////////////////////////////////////////////////////////////////////////
72/// TGraphPolar constructor.
73///
74/// \param[in] n number of points.
75/// \param[in] theta angular values.
76/// \param[in] r radial values.
77/// \param[in] etheta errors on angular values.
78/// \param[in] er errors on radial values.
79
81 const Double_t *etheta, const Double_t* er)
82 : TGraphErrors(n,theta,r,etheta,er),
83 fOptionAxis(kFALSE),fPolargram(nullptr),fXpol(nullptr),fYpol(nullptr)
84{
86}
87
88////////////////////////////////////////////////////////////////////////////////
89/// TGraphPolar destructor.
90
92{
93 delete [] fXpol;
94 delete [] fYpol;
95}
96
97////////////////////////////////////////////////////////////////////////////////
98/// Draw TGraphPolar.
99
101{
102 AppendPad(options);
103}
104
105////////////////////////////////////////////////////////////////////////////////
106/// Return points in polar coordinates.
107
109{
110 if (!fXpol) fXpol = new Double_t[fNpoints];
111 return fXpol;
112}
113
114////////////////////////////////////////////////////////////////////////////////
115/// Return points in polar coordinates.
116
118{
119 if (!fYpol) fYpol = new Double_t[fNpoints];
120 return fYpol;
121}
122
123////////////////////////////////////////////////////////////////////////////////
124/// Set maximum Polar.
125
127{
129}
130
131////////////////////////////////////////////////////////////////////////////////
132/// Set maximum radial at the intersection of the positive X axis part and
133/// the circle.
134
136{
138}
139
140////////////////////////////////////////////////////////////////////////////////
141/// Set minimum Polar.
142
144{
146}
147
148////////////////////////////////////////////////////////////////////////////////
149/// Set minimum radial in the center of the circle.
150
152{
154}
155
156////////////////////////////////////////////////////////////////////////////////
157/// Create polargram object for given draw options
158
160{
161 Int_t theNpoints = GetN();
162
163 if (theNpoints < 1)
164 return nullptr;
165
166 Double_t *theX = GetX();
167 Double_t *theY = GetY();
168 Double_t *theEX = GetEX();
169 Double_t *theEY = GetEY();
170
171 // Get range, initialize with first/last value
172 Double_t rwrmin = theY[0], rwrmax = theY[theNpoints - 1], rwtmin = theX[0], rwtmax = theX[theNpoints - 1];
173
174 for (Int_t ipt = 0; ipt < theNpoints; ipt++) {
175 // Check for errors if available
176 if (theEX) {
177 if (theX[ipt] - theEX[ipt] < rwtmin)
178 rwtmin = theX[ipt] - theEX[ipt];
179 if (theX[ipt] + theEX[ipt] > rwtmax)
180 rwtmax = theX[ipt] + theEX[ipt];
181 } else {
182 if (theX[ipt] < rwtmin)
183 rwtmin = theX[ipt];
184 if (theX[ipt] > rwtmax)
185 rwtmax = theX[ipt];
186 }
187 if (theEY) {
188 if (theY[ipt] - theEY[ipt] < rwrmin)
189 rwrmin = theY[ipt] - theEY[ipt];
190 if (theY[ipt] + theEY[ipt] > rwrmax)
191 rwrmax = theY[ipt] + theEY[ipt];
192 } else {
193 if (theY[ipt] < rwrmin)
194 rwrmin = theY[ipt];
195 if (theY[ipt] > rwrmax)
196 rwrmax = theY[ipt];
197 }
198 }
199
200 // Add radial and Polar margins.
201 if (rwrmin == rwrmax)
202 rwrmax += 1.;
203 if (rwtmin == rwtmax)
204 rwtmax += 1.;
205
206 Double_t dr = rwrmax - rwrmin, dt = rwtmax - rwtmin;
207
208 rwrmax += 0.1 * dr;
209 rwrmin -= 0.1 * dr;
210
211 // Assume equally spaced points for full 2*Pi.
212 rwtmax += dt / theNpoints;
213
214 return new TGraphPolargram("Polargram", rwrmin, rwrmax, rwtmin, rwtmax, opt);
215}
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:382
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
A TGraphErrors is a TGraph with error bars.
Double_t * GetEX() const override
Double_t * GetEY() const override
To draw a polar graph.
Definition TGraphPolar.h:23
void SetMinRadial(Double_t minimum=0)
Set minimum radial in the center of the circle.
Double_t * fXpol
[fNpoints] points in polar coordinates
Definition TGraphPolar.h:30
void SetMaxPolar(Double_t maximum=6.28318530717958623)
Set maximum Polar.
Double_t * GetYpol()
Return points in polar coordinates.
TGraphPolargram * fPolargram
The polar coordinates system.
Definition TGraphPolar.h:29
void SetMaxRadial(Double_t maximum=1)
Set maximum radial at the intersection of the positive X axis part and the circle.
TGraphPolar()
TGraphPolar default constructor.
void Draw(Option_t *options="") override
Draw TGraphPolar.
~TGraphPolar() override
TGraphPolar destructor.
Double_t * GetXpol()
Return points in polar coordinates.
TGraphPolargram * CreatePolargram(const char *opt)
Create polargram object for given draw options.
void SetMinPolar(Double_t minimum=0)
Set minimum Polar.
Double_t * fYpol
[fNpoints] points in polar coordinates
Definition TGraphPolar.h:31
To draw polar axis.
void SetRangeRadial(Double_t rmin, Double_t rmax)
Set the radial range.
void ChangeRangePolar(Double_t tmin, Double_t tmax)
Set the Polar range.
Int_t fNpoints
Number of points <= fMaxSize.
Definition TGraph.h:46
Double_t * GetY() const
Definition TGraph.h:140
Int_t GetN() const
Definition TGraph.h:132
Double_t * GetX() const
Definition TGraph.h:139
virtual void SetEditable(Bool_t editable=kTRUE)
if editable=kFALSE, the graph cannot be modified with the mouse by default a TGraph is editable
Definition TGraph.cxx:2306
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition TObject.cxx:202
const Int_t n
Definition legend1.C:16