Logo ROOT   6.08/07
Reference Guide
TCandle.h
Go to the documentation of this file.
1 // @(#)root/graf:$Id$
2 // Author: Georg Troska 2016/04/14
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 #ifndef ROOT_TCandle
13 #define ROOT_TCandle
14 
15 #ifndef ROOT_TObject
16 #include "TObject.h"
17 #endif
18 #ifndef ROOT_TAttLine
19 #include "TAttLine.h"
20 #endif
21 #ifndef ROOT_TAttFill
22 #include "TAttFill.h"
23 #endif
24 #ifndef ROOT_TAttMarket
25 #include "TAttMarker.h"
26 #endif
27 
28 #include "TH1D.h"
29 #include "TMath.h"
30 
31 class TCandle : public TAttLine, public TAttFill, public TAttMarker {
32 public:
33  //Candle Option
34  enum CandleOption : int {
35  kNoOption = 0,
36  kBox = 1,
41  kMeanLine = 100,
42  kMeanCircle = 300,
43  kWhiskerAll = 1000,
44  kWhisker15 = 2000,
45  kAnchor = 10000,
46  kPointsOutliers = 100000,
47  kPointsAll = 200000,
48  kPointsAllScat = 300000,
49  kHorizontal = 1000000 // if this bit is not set it is vertical!
50  };
51 
52 protected:
53 
54  bool fIsRaw; ///< 0: for TH1 projection, 1: using raw data
57  bool fDismiss; ///< True if the candle cannot be painted
58 
59  Double_t fPosCandleAxis; ///< x-pos for a vertical candle
60  Double_t fCandleWidth; ///< The candle width
61 
62  Double_t fMean; ///< Position of the mean
63  Double_t fMedian; ///< Position of the median
64  Double_t fBoxUp; ///< Position of the upper box end
65  Double_t fBoxDown; ///< Position of the lower box end
66  Double_t fWhiskerUp; ///< Position of the upper whisker end
67  Double_t fWhiskerDown; ///< Position of the lower whisker end
68 
69  Double_t * fDatapoints; ///< position of all Datapoints within this candle
70  Double_t fNDatapoints; ///< Number of Datapoints within this candle
71 
72  CandleOption fOption; ///< Setting the style of the candle
73  int fLogX;
74  int fLogY;
75 
76  void Calculate();
77 
78  int GetCandleOption(const int pos) {return (fOption/(int)TMath::Power(10,pos))%10;}
79  bool IsOption(CandleOption opt);
80  void PaintBox(Int_t nPoints, Double_t *x, Double_t *y, Bool_t swapXY, Bool_t fill);
81  void PaintLine(Double_t x1, Double_t y1, Double_t x2, Double_t y2, Bool_t swapXY);
82 
83 public:
84 
85  TCandle();
86  TCandle(const Double_t candlePos, const Double_t candleWidth, const Int_t n, const Double_t * points);
87  TCandle(const Double_t candlePos, const Double_t candleWidth, TH1D *proj);
88  TCandle(const TCandle &candle);
89  virtual ~TCandle();
90 
91  Double_t GetMean() const {return fMean;}
92  Double_t GetMedian() const {return fMedian;}
93  Double_t GetQ1() const {return fBoxUp;}
94  Double_t GetQ2() const {return fMedian;}
95  Double_t GetQ3() const {return fBoxDown;}
98 
99  void SetOption(CandleOption opt) { fOption = opt; }
100  void SetLog(int x, int y) { fLogX = x; fLogY = y; }
101  void SetAxisPosition(const Double_t candlePos) { fPosCandleAxis = candlePos; }
102 
103  void SetWidth(const Double_t width) { fCandleWidth = width; }
104  void SetHistogram(TH1D *proj) { fProj = proj; fIsCalculated = false;}
105 
106  virtual void Paint(Option_t *option="");
107 
108  virtual void SetMean(Double_t mean) { fMean = mean; }
109  virtual void SetMedian(Double_t median) { fMedian = median; }
110  virtual void SetQ1(Double_t q1) { fBoxUp = q1; }
111  virtual void SetQ2(Double_t q2) { fMedian = q2; }
112  virtual void SetQ3(Double_t q3) { fBoxDown = q3; }
113 
114  int ParseOption(char *optin);
115 
116  ClassDef(TCandle,1) //A Candle
117 };
118 #endif
Double_t GetMean() const
Definition: TCandle.h:91
virtual void SetMean(Double_t mean)
Definition: TCandle.h:108
TCandle()
TCandle default constructor.
Definition: TCandle.cxx:38
Double_t fMedian
Position of the median.
Definition: TCandle.h:63
Bool_t IsHorizontal()
Definition: TCandle.h:96
Double_t fWhiskerUp
Position of the upper whisker end.
Definition: TCandle.h:66
bool IsOption(CandleOption opt)
Return true is this option is activated in fOption.
Definition: TCandle.cxx:412
const char Option_t
Definition: RtypesCore.h:62
Double_t fNDatapoints
Number of Datapoints within this candle.
Definition: TCandle.h:70
fill
Definition: fit1_py.py:6
int ParseOption(char *optin)
Parsing of the option-string.
Definition: TCandle.cxx:95
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
void SetHistogram(TH1D *proj)
Definition: TCandle.h:104
Double_t * fDatapoints
position of all Datapoints within this candle
Definition: TCandle.h:69
void SetLog(int x, int y)
Definition: TCandle.h:100
void PaintLine(Double_t x1, Double_t y1, Double_t x2, Double_t y2, Bool_t swapXY)
Paint a line for candle.
Definition: TCandle.cxx:457
void SetOption(CandleOption opt)
Definition: TCandle.h:99
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:501
Double_t fMean
Position of the mean.
Definition: TCandle.h:62
Marker Attributes class.
Definition: TAttMarker.h:24
virtual void SetMedian(Double_t median)
Definition: TCandle.h:109
static const double x2[5]
Fill Area Attributes class.
Definition: TAttFill.h:24
bool fIsCalculated
Definition: TCandle.h:55
Double_t x[n]
Definition: legend1.C:17
#define ClassDef(name, id)
Definition: Rtypes.h:254
Double_t fCandleWidth
The candle width.
Definition: TCandle.h:60
int fLogX
Definition: TCandle.h:73
Double_t GetQ3() const
Definition: TCandle.h:95
Double_t GetMedian() const
Definition: TCandle.h:92
Bool_t IsVertical()
Definition: TCandle.h:97
Double_t fBoxDown
Position of the lower box end.
Definition: TCandle.h:65
point * points
Definition: X3DBuffer.c:20
Double_t fBoxUp
Position of the upper box end.
Definition: TCandle.h:64
Double_t GetQ1() const
Definition: TCandle.h:93
int fLogY
Definition: TCandle.h:74
The candle plot painter class.
Definition: TCandle.h:31
void Calculate()
Calculates all values needed by the candle definition depending on the candle options.
Definition: TCandle.cxx:161
void SetWidth(const Double_t width)
Definition: TCandle.h:103
tomato 1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:618
virtual void Paint(Option_t *option="")
Paint one candle with its current attributes.
Definition: TCandle.cxx:207
void PaintBox(Int_t nPoints, Double_t *x, Double_t *y, Bool_t swapXY, Bool_t fill)
Paint a box for candle.
Definition: TCandle.cxx:428
bool fDismiss
True if the candle cannot be painted.
Definition: TCandle.h:57
static const double x1[5]
TH1D * fProj
Definition: TCandle.h:56
double Double_t
Definition: RtypesCore.h:55
CandleOption fOption
Setting the style of the candle.
Definition: TCandle.h:72
Double_t y[n]
Definition: legend1.C:17
int GetCandleOption(const int pos)
Definition: TCandle.h:78
CandleOption
Definition: TCandle.h:34
virtual void SetQ2(Double_t q2)
Definition: TCandle.h:111
Double_t fWhiskerDown
Position of the lower whisker end.
Definition: TCandle.h:67
Double_t fPosCandleAxis
x-pos for a vertical candle
Definition: TCandle.h:59
virtual void SetQ1(Double_t q1)
Definition: TCandle.h:110
virtual ~TCandle()
TCandle default destructor.
Definition: TCandle.cxx:88
virtual void SetQ3(Double_t q3)
Definition: TCandle.h:112
bool fIsRaw
0: for TH1 projection, 1: using raw data
Definition: TCandle.h:54
const Int_t n
Definition: legend1.C:16
Line Attributes class.
Definition: TAttLine.h:24
void SetAxisPosition(const Double_t candlePos)
Definition: TCandle.h:101
Double_t GetQ2() const
Definition: TCandle.h:94