ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 
15 To draw a polar graph.
16 
17 TGraphPolar creates a polar graph (including error bars). A TGraphPolar is
18 a TGraphErrors represented in polar coordinates.
19 It uses the class TGraphPolargram to draw the polar axis.
20 
21 Example:
22 
23 Begin_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 }
55 End_Macro
56 */
57 
58 #include "TROOT.h"
59 #include "TGraphPolar.h"
60 #include "TGraphPolargram.h"
61 #include "TVirtualPad.h"
62 
64 
65 ////////////////////////////////////////////////////////////////////////////////
66 /// TGraphPolar default constructor.
67 
69  fOptionAxis(kFALSE),fPolargram(0),fXpol(0),fYpol(0)
70 {
71 }
72 
73 ////////////////////////////////////////////////////////////////////////////////
74 /// TGraphPolar constructor.
75 ///
76 /// \param[in] n number of points.
77 /// \param[in] theta angular values.
78 /// \param[in] r radial values.
79 /// \param[in] etheta errors on angular values.
80 /// \param[in] er errors on radial values.
81 
83  const Double_t *etheta, const Double_t* er)
84  : TGraphErrors(n,theta,r,etheta,er),
85  fOptionAxis(kFALSE),fPolargram(0),fXpol(0),fYpol(0)
86 {
88 }
89 
90 ////////////////////////////////////////////////////////////////////////////////
91 /// TGraphPolar destructor.
92 
94 {
95  delete []fXpol;
96  delete []fYpol;
97 }
98 
99 ////////////////////////////////////////////////////////////////////////////////
100 /// Draw TGraphPolar.
101 
103 {
104  // Process options
105  TString opt = options;
106  opt.ToUpper();
107 
108  // Ignore same
109  opt.ReplaceAll("SAME","");
110 
111  // ReDraw polargram if required by options
112  if (opt.Contains("A")) fOptionAxis = kTRUE;
113  opt.ReplaceAll("A","");
114 
115  AppendPad(opt);
116 }
117 
118 ////////////////////////////////////////////////////////////////////////////////
119 /// Return points in polar coordinates.
120 
122 {
123  if (!fXpol) fXpol = new Double_t[fNpoints];
124  return fXpol;
125 }
126 
127 ////////////////////////////////////////////////////////////////////////////////
128 /// Return points in polar coordinates.
129 
131 {
132  if (!fYpol) fYpol = new Double_t[fNpoints];
133  return fYpol;
134 }
135 
136 ////////////////////////////////////////////////////////////////////////////////
137 /// Set maximum Polar.
138 
140 {
142 }
143 
144 ////////////////////////////////////////////////////////////////////////////////
145 /// Set maximum radial at the intersection of the positive X axis part and
146 /// the circle.
147 
149 {
151 }
152 
153 ////////////////////////////////////////////////////////////////////////////////
154 /// Set minimum Polar.
155 
157 {
159 }
160 
161 ////////////////////////////////////////////////////////////////////////////////
162 /// Set minimum radial in the center of the circle.
163 
165 {
167 }
void SetMaxRadial(Double_t maximum=1)
Set maximum radial at the intersection of the positive X axis part and the circle.
Int_t fNpoints
Current dimension of arrays fX and fY.
Definition: TGraph.h:58
void SetMinRadial(Double_t minimum=0)
Set minimum radial in the center of the circle.
TGraphPolar()
TGraphPolar default constructor.
Definition: TGraphPolar.cxx:68
Double_t GetRMin()
Float_t theta
Definition: shapesAnim.C:5
const char Option_t
Definition: RtypesCore.h:62
Bool_t fOptionAxis
Definition: TGraphPolar.h:41
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:635
void SetMinPolar(Double_t minimum=0)
Set minimum Polar.
Double_t * fXpol
Definition: TGraphPolar.h:45
void ToUpper()
Change string to upper case.
Definition: TString.cxx:1088
Basic string class.
Definition: TString.h:137
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition: TObject.cxx:164
void ChangeRangePolar(Double_t tmin, Double_t tmax)
Set the Polar range.
Double_t GetTMin()
To draw a polar graph.
Definition: TGraphPolar.h:38
ClassImp(TGraphPolar)
ROOT::R::TRInterface & r
Definition: Object.C:4
void SetRangeRadial(Double_t rmin, Double_t rmax)
Set the radial range.
Double_t GetTMax()
void SetMaxPolar(Double_t maximum=6.28318530717958623)
Set maximum Polar.
double Double_t
Definition: RtypesCore.h:55
virtual ~TGraphPolar()
TGraphPolar destructor.
Definition: TGraphPolar.cxx:93
Double_t * fYpol
Definition: TGraphPolar.h:46
Double_t GetRMax()
TGraphPolargram * fPolargram
Definition: TGraphPolar.h:44
void Draw(Option_t *options="")
Draw TGraphPolar.
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:567
A TGraphErrors is a TGraph with error bars.
Definition: TGraphErrors.h:28
Double_t * GetXpol()
Return points in polar coordinates.
const Bool_t kTRUE
Definition: Rtypes.h:91
Double_t * GetYpol()
Return points in polar coordinates.
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:2100
const Int_t n
Definition: legend1.C:16