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 // Process options
103 TString opt = options;
104 opt.ToUpper();
105
106 // Ignore same
107 opt.ReplaceAll("SAME","");
108
109 // ReDraw polargram if required by options
110 if (opt.Contains("A")) fOptionAxis = kTRUE;
111 opt.ReplaceAll("A","");
112
113 AppendPad(opt);
114}
115
116////////////////////////////////////////////////////////////////////////////////
117/// Return points in polar coordinates.
118
120{
121 if (!fXpol) fXpol = new Double_t[fNpoints];
122 return fXpol;
123}
124
125////////////////////////////////////////////////////////////////////////////////
126/// Return points in polar coordinates.
127
129{
130 if (!fYpol) fYpol = new Double_t[fNpoints];
131 return fYpol;
132}
133
134////////////////////////////////////////////////////////////////////////////////
135/// Set maximum Polar.
136
138{
140}
141
142////////////////////////////////////////////////////////////////////////////////
143/// Set maximum radial at the intersection of the positive X axis part and
144/// the circle.
145
147{
149}
150
151////////////////////////////////////////////////////////////////////////////////
152/// Set minimum Polar.
153
155{
157}
158
159////////////////////////////////////////////////////////////////////////////////
160/// Set minimum radial in the center of the circle.
161
163{
165}
constexpr Bool_t kFALSE
Definition RtypesCore.h:101
constexpr Bool_t kTRUE
Definition RtypesCore.h:100
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:377
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.
To draw a polar graph.
Definition TGraphPolar.h:23
void SetMinRadial(Double_t minimum=0)
Set minimum radial in the center of the circle.
Bool_t fOptionAxis
Force drawing of new coord system.
Definition TGraphPolar.h:26
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.
void SetMinPolar(Double_t minimum=0)
Set minimum Polar.
Double_t * fYpol
[fNpoints] points in polar coordinates
Definition TGraphPolar.h:31
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
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:2274
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition TObject.cxx:184
Basic string class.
Definition TString.h:139
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition TString.h:704
void ToUpper()
Change string to upper case.
Definition TString.cxx:1183
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:636
const Int_t n
Definition legend1.C:16