Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
TCandle.cxx
Go to the documentation of this file.
1 // @(#)root/graf:$Id$
2 // Author: Georg Troska 19/05/16
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 <stdlib.h>
13 #include <iostream>
14 #include "TROOT.h"
15 #include "TCandle.h"
16 #include "TClass.h"
17 #include "TPad.h"
18 #include "TRandom2.h"
19 
20 Double_t TCandle::fWhiskerRange = 1.0;
21 Double_t TCandle::fBoxRange = 0.5;
22 Bool_t TCandle::fScaledCandle = false;
23 Bool_t TCandle::fScaledViolin = true;
24 
25 ClassImp(TCandle);
26 
27 /** \class TCandle
28 \ingroup BasicGraphics
29 
30 The candle plot painter class.
31 
32 Instances of this class are generated by the histograms painting
33 classes (THistPainter and THStack) when an candle plot (box plot) is drawn.
34 TCandle is the "painter class" of the box plots. Therefore it is never used
35 directly to draw a candle.
36 */
37 
38 ////////////////////////////////////////////////////////////////////////////////
39 /// TCandle default constructor.
40 
41 TCandle::TCandle()
42 {
43  fIsCalculated = 0;
44  fIsRaw = 0;
45  fPosCandleAxis = 0.;
46  fCandleWidth = 1.0;
47  fHistoWidth = 1.0;
48  fMean = 0.;
49  fMedian = 0.;
50  fMedianErr = 0;
51  fBoxUp = 0.;
52  fBoxDown = 0.;
53  fWhiskerUp = 0.;
54  fWhiskerDown = 0.;
55  fNDatapoints = 0;
56  fDismiss = 0;
57  fLogX = 0;
58  fLogY = 0;
59  fLogZ = 0;
60  fNDrawPoints = 0;
61  fNHistoPoints = 0;
62  fAxisMin = 0.;
63  fAxisMax = 0.;
64  fOption = kNoOption;
65  fProj = NULL;
66  fDatapoints = 0;
67 
68 }
69 
70 ////////////////////////////////////////////////////////////////////////////////
71 /// TCandle constructor passing a draw-option.
72 
73 TCandle::TCandle(const char *opt)
74 {
75  fIsCalculated = 0;
76  fIsRaw = 0;
77  fPosCandleAxis = 0.;
78  fCandleWidth = 1.0;
79  fHistoWidth = 1.0;
80  fMean = 0.;
81  fMedian = 0.;
82  fMedianErr = 0;
83  fBoxUp = 0.;
84  fBoxDown = 0.;
85  fWhiskerUp = 0.;
86  fWhiskerDown = 0.;
87  fNDatapoints = 0;
88  fDismiss = 0;
89  fLogX = 0;
90  fLogY = 0;
91  fLogZ = 0;
92  fNDrawPoints = 0;
93  fNHistoPoints = 0;
94  fAxisMin = 0.;
95  fAxisMax = 0.;
96  fOption = kNoOption;
97  fProj = NULL;
98  fDatapoints = 0;
99 
100 
101  // Conversion necessary in order to cast from const char* to char*
102  char myopt[128];
103  strlcpy(myopt,opt,128);
104 
105 
106  ParseOption(myopt);
107 }
108 
109 
110 ////////////////////////////////////////////////////////////////////////////////
111 /// TCandle constructor for raw-data candles.
112 
113 TCandle::TCandle(const Double_t candlePos, const Double_t candleWidth, Long64_t n, Double_t * points)
114  : TAttLine(), TAttFill(), TAttMarker()
115 {
116  //Preliminary values only, need to be calculated before paint
117  fMean = 0;
118  fMedian = 0;
119  fMedianErr = 0;
120  fBoxUp = 0;
121  fBoxDown = 0;
122  fWhiskerUp = 0;
123  fWhiskerDown = 0;
124  fNDatapoints = n;
125  fIsCalculated = 0;
126  fIsRaw = true;
127  fPosCandleAxis = candlePos;
128  fCandleWidth = candleWidth;
129  fHistoWidth = candleWidth;
130  fDatapoints = points;
131  fProj = NULL;
132  fDismiss = 0;
133  fOption = kNoOption;
134  fLogX = 0;
135  fLogY = 0;
136  fLogZ = 0;
137  fNDrawPoints = 0;
138  fNHistoPoints = 0;
139  fAxisMin = 0.;
140  fAxisMax = 0.;
141  sprintf(fOptionStr," ");
142 }
143 
144 ////////////////////////////////////////////////////////////////////////////////
145 /// TCandle TH1 data constructor.
146 
147 TCandle::TCandle(const Double_t candlePos, const Double_t candleWidth, TH1D *proj)
148  : TAttLine(), TAttFill(), TAttMarker()
149 {
150  //Preliminary values only, need to be calculated before paint
151  fMean = 0;
152  fMedian = 0;
153  fMedianErr = 0;
154  fBoxUp = 0;
155  fBoxDown = 0;
156  fWhiskerUp = 0;
157  fWhiskerDown = 0;
158  fNDatapoints = 0;
159  fIsCalculated = 0;
160  fIsRaw = 0;
161  fPosCandleAxis = candlePos;
162  fCandleWidth = candleWidth;
163  fHistoWidth = candleWidth;
164  fDatapoints = 0;
165  fProj = proj;
166  fDismiss = 0;
167  fOption = kNoOption;
168  fLogX = 0;
169  fLogY = 0;
170  fLogZ = 0;
171  fNDrawPoints = 0;
172  fNHistoPoints = 0;
173  fAxisMin = 0.;
174  fAxisMax = 0.;
175  sprintf(fOptionStr," ");
176 }
177 
178 ////////////////////////////////////////////////////////////////////////////////
179 /// TCandle default destructor.
180 
181 TCandle::~TCandle() {
182  if (fIsRaw && fProj) delete fProj;
183 }
184 
185 Bool_t TCandle::IsCandleScaled()
186 {
187  return fScaledCandle;
188 }
189 
190 Bool_t TCandle::IsViolinScaled()
191 {
192  return fScaledViolin;
193 }
194 
195 ////////////////////////////////////////////////////////////////////////////////
196 /// Static function to set fWhiskerRange, by setting whisker-range, one can force
197 /// the whiskers to cover the fraction of the distribution.
198 /// Set wRange between 0 and 1. Default is 1
199 /// TCandle::SetWhiskerRange(0.95) will set all candle-charts to cover 95% of
200 /// the distribution with the whiskers.
201 /// Can only be used with the standard-whisker definition
202 
203 void TCandle::SetWhiskerRange(const Double_t wRange) {
204  if (wRange < 0) fWhiskerRange = 0;
205  else if (wRange > 1) fWhiskerRange = 1;
206  else fWhiskerRange = wRange;
207 
208 }
209 
210 ////////////////////////////////////////////////////////////////////////////////
211 /// Static function to set fBoxRange, by setting whisker-range, one can force the
212 /// box of the candle-chart to cover that given fraction of the distribution.
213 /// Set bRange between 0 and 1. Default is 0.5
214 /// TCandle::SetBoxRange(0.68) will set all candle-charts to cover 68% of the
215 /// distribution by the box
216 
217 void TCandle::SetBoxRange(const Double_t bRange) {
218  if (bRange < 0) fBoxRange = 0;
219  else if (bRange > 1) fBoxRange = 1;
220  else fBoxRange = bRange;
221 }
222 
223 ////////////////////////////////////////////////////////////////////////////////
224 /// Static function to set scaling between candles-withs. A candle containing
225 /// 100 entries with be two times wider than a candle containing 50 entries
226 
227 void TCandle::SetScaledCandle(const Bool_t cScale) {
228  fScaledCandle = cScale;
229 }
230 
231 ////////////////////////////////////////////////////////////////////////////////
232 /// Static function to set scaling between violin-withs. A violin or histo chart
233 /// with a maximum bin content to 100 will be two times as high as a violin with
234 /// a maximum bin content of 50
235 
236 void TCandle::SetScaledViolin(const Bool_t vScale) {
237  fScaledViolin = vScale;
238 }
239 
240 ////////////////////////////////////////////////////////////////////////////////
241 /// Parsing of the option-string.
242 /// The option-string will be empty at the end (by-reference).
243 
244 int TCandle::ParseOption(char * opt) {
245  fOption = kNoOption;
246  char *l;
247 
248  l = strstr(opt,"CANDLE");
249  if (l) {
250  const CandleOption fallbackCandle = (CandleOption)(kBox + kMedianLine + kMeanCircle + kWhiskerAll + kAnchor);
251 
252  char direction = ' ';
253  char preset = ' ';
254 
255 
256 
257  if (l[6] >= 'A' && l[6] <= 'Z') direction = l[6];
258  if (l[6] >= '1' && l[6] <= '9') preset = l[6];
259  if (l[7] >= 'A' && l[7] <= 'Z' && preset != ' ') direction = l[7];
260  if (l[7] >= '1' && l[7] <= '9' && direction != ' ') preset = l[7];
261 
262  if (direction == 'X' || direction == 'V') { /* nothing */ }
263  if (direction == 'Y' || direction == 'H') { fOption = (CandleOption)(fOption + kHorizontal); }
264  if (preset == '1') //Standard candle using old candle-definition
265  fOption = (CandleOption)(fOption + fallbackCandle);
266  else if (preset == '2') //New standard candle with better whisker definition + outlier
267  fOption = (CandleOption)(fOption + kBox + kMeanLine + kMedianLine + kWhisker15 + kAnchor + kPointsOutliers);
268  else if (preset == '3') //Like candle2 but with a fMean as a circle
269  fOption = (CandleOption)(fOption + kBox + kMeanCircle + kMedianLine + kWhisker15 + kAnchor + kPointsOutliers);
270  else if (preset == '4') //Like candle3 but showing the uncertainty of the fMedian as well
271  fOption = (CandleOption)(fOption + kBox + kMeanCircle + kMedianNotched + kWhisker15 + kAnchor + kPointsOutliers);
272  else if (preset == '5') //Like candle2 but showing all datapoints
273  fOption = (CandleOption)(fOption + kBox + kMeanLine + kMedianLine + kWhisker15 + kAnchor + kPointsAll);
274  else if (preset == '6') //Like candle2 but showing all datapoints scattered
275  fOption = (CandleOption)(fOption + kBox + kMeanCircle + kMedianLine + kWhisker15 + kAnchor + kPointsAllScat);
276  else if (preset != ' ') //For all other presets not implemented yet used fallback candle
277  fOption = (CandleOption)(fOption + fallbackCandle);
278 
279  if (preset != ' ' && direction != ' ')
280  memcpy(l," ",8);
281  else if (preset != ' ' || direction != ' ')
282  memcpy(l," ",7);
283  else
284  memcpy(l," ",6);
285 
286  Bool_t useIndivOption = false;
287 
288  if (direction == ' ') direction = 'X';
289  if (preset == ' ') { // Check if the user wants to set the properties individually
290  char *brOpen = strstr(opt,"(");
291  char *brClose = strstr(opt,")");
292  char indivOption[32];
293  if (brOpen && brClose) {
294  useIndivOption = true;
295  bool isHorizontal = IsHorizontal();
296  strlcpy(indivOption, brOpen, brClose-brOpen+2); //Now the string "(....)" including brackets is in this array
297  sscanf(indivOption,"(%d)", (int*) &fOption);
298  if (isHorizontal) {fOption = (CandleOption)(fOption + kHorizontal);}
299  memcpy(brOpen," ",brClose-brOpen+1); //Cleanup
300 
301  sprintf(fOptionStr,"CANDLE%c(%ld)",direction,(long)fOption);
302  } else {
303  preset = 1;
304  fOption = (CandleOption)(fOption + fallbackCandle);
305  }
306  } else {
307  sprintf(fOptionStr,"CANDLE%c%c",direction,preset);
308  }
309  //Handle option "CANDLE" ,"CANDLEX" or "CANDLEY" to behave like "CANDLEX1" or "CANDLEY1"
310  if (!useIndivOption && !fOption ) {
311  fOption = fallbackCandle;
312  sprintf(fOptionStr,"CANDLE%c2",direction);
313  }
314  }
315 
316  l = strstr(opt,"VIOLIN");
317  if (l) {
318  const CandleOption fallbackCandle = (CandleOption)(kMeanCircle + kWhiskerAll + kHistoViolin + kHistoZeroIndicator);
319 
320  char direction = ' ';
321  char preset = ' ';
322 
323  if (l[6] >= 'A' && l[6] <= 'Z') direction = l[6];
324  if (l[6] >= '1' && l[6] <= '9') preset = l[6];
325  if (l[7] >= 'A' && l[7] <= 'Z' && preset != ' ') direction = l[7];
326  if (l[7] >= '1' && l[7] <= '9' && direction != ' ') preset = l[7];
327 
328  if (direction == 'X' || direction == 'V') { /* nothing */ }
329  if (direction == 'Y' || direction == 'H') { fOption = (CandleOption)(fOption + kHorizontal); }
330  if (preset == '1') //Standard candle using old candle-definition
331  fOption = (CandleOption)(fOption + fallbackCandle);
332  else if (preset == '2') //New standard candle with better whisker definition + outlier
333  fOption = (CandleOption)(fOption + kMeanCircle + kWhisker15 + kHistoViolin + kHistoZeroIndicator + kPointsOutliers);
334  else if (preset != ' ') //For all other presets not implemented yet used fallback candle
335  fOption = (CandleOption)(fOption + fallbackCandle);
336 
337  if (preset != ' ' && direction != ' ')
338  memcpy(l," ",8);
339  else if (preset != ' ' || direction != ' ')
340  memcpy(l," ",7);
341  else
342  memcpy(l," ",6);
343 
344  Bool_t useIndivOption = false;
345 
346  if (direction == ' ') direction = 'X';
347  if (preset == ' ') { // Check if the user wants to set the properties individually
348  char *brOpen = strstr(opt,"(");
349  char *brClose = strstr(opt,")");
350  char indivOption[32];
351  if (brOpen && brClose) {
352  useIndivOption = true;
353  bool isHorizontal = IsHorizontal();
354  strlcpy(indivOption, brOpen, brClose-brOpen +2); //Now the string "(....)" including brackets is in this array
355  sscanf(indivOption,"(%d)", (int*) &fOption);
356  if (isHorizontal) {fOption = (CandleOption)(fOption + kHorizontal);}
357  memcpy(brOpen," ",brClose-brOpen+1); //Cleanup
358 
359  sprintf(fOptionStr,"VIOLIN%c(%ld)",direction,(long)fOption);
360 
361  } else {
362  preset = 1;
363  fOption = (CandleOption)(fOption + fallbackCandle);
364  }
365  } else {
366  sprintf(fOptionStr,"VIOLIN%c%c",direction,preset);
367  }
368  //Handle option "VIOLIN" ,"VIOLINX" or "VIOLINY" to behave like "VIOLINX1" or "VIOLINY1"
369  if (!useIndivOption && !fOption ) {
370  fOption = fallbackCandle;
371  sprintf(fOptionStr,"VIOLIN%c1",direction);
372  }
373  }
374 
375  fIsCalculated = false;
376 
377  return fOption;
378 
379 }
380 
381 ////////////////////////////////////////////////////////////////////////////////
382 /// Calculates all values needed by the candle definition depending on the
383 /// candle options.
384 
385 void TCandle::Calculate() {
386  //Reset everything
387  fNDrawPoints = 0;
388  fNHistoPoints = 0;
389 
390  Bool_t swapXY = IsOption(kHorizontal);
391  Bool_t doLogY = (!(swapXY) && fLogY) || (swapXY && fLogX);
392  Bool_t doLogX = (!(swapXY) && fLogX) || (swapXY && fLogY);
393  Bool_t doLogZ = fLogZ;
394 
395  //Will be min and max values of raw-data
396  Double_t min = 1e15;
397  Double_t max = -1e15;
398 
399  // Determining the quantiles
400  Double_t *prob = new Double_t[5];
401 
402  if (fWhiskerRange >= 1) {
403  prob[0] = 1e-15;
404  prob[4] = 1-1e-15;
405  } else {
406  prob[0] = 0.5 - fWhiskerRange/2.;
407  prob[4] = 0.5 + fWhiskerRange/2.;
408  }
409 
410 
411  if (fBoxRange >= 1) {
412  prob[1] = 1E-14;
413  prob[3] = 1-1E-14;
414  } else {
415  prob[1] = 0.5 - fBoxRange/2.;
416  prob[3] = 0.5 + fBoxRange/2.;
417  }
418 
419  prob[2]=0.5;
420  Double_t *quantiles = new Double_t[5];
421  quantiles[0]=0.; quantiles[1]=0.; quantiles[2] = 0.; quantiles[3] = 0.; quantiles[4] = 0.;
422  if (!fIsRaw && fProj) { //Need a calculation for a projected histo
423  if (((IsOption(kHistoLeft)) || (IsOption(kHistoRight)) || (IsOption(kHistoViolin))) && fProj->GetNbinsX() > 500) {
424  // When using the histooption the number of bins of the projection is
425  // limited because of the array space defined by kNMAXPOINTS.
426  // So the histo is rebinned, that it can be displayed at any time.
427  // Finer granularity is not useful anyhow
428  int divideBy = ((fProj->GetNbinsX() - 1)/((kNMAXPOINTS-10)/4))+1;
429  fProj->RebinX(divideBy);
430  }
431  fProj->GetQuantiles(5, quantiles, prob);
432  } else { //Need a calculation for a raw-data candle
433  TMath::Quantiles(fNDatapoints,5,fDatapoints,quantiles,prob,kFALSE);
434  }
435 
436  // Check if the quantiles are valid, seems the under- and overflow is taken
437  // into account as well, we need to ignore this!
438  if (quantiles[0] >= quantiles[4] ||
439  quantiles[1] >= quantiles[3]) {
440  delete [] prob;
441  delete [] quantiles;
442  return;
443  }
444 
445  // Definition of the candle in the standard case
446  fBoxUp = quantiles[3];
447  fBoxDown = quantiles[1];
448  fWhiskerUp = quantiles[4]; //Standard case
449  fWhiskerDown = quantiles[0]; //Standard case
450  fMedian = quantiles[2];
451  Double_t iqr = fBoxUp-fBoxDown;
452  Int_t nOutliers = 0;
453 
454  if (IsOption(kWhisker15)) { // Improved whisker definition, with 1.5*iqr
455  if (!fIsRaw && fProj) { //Need a calculation for a projected histo
456  int bin = fProj->FindBin(fBoxDown-1.5*iqr);
457  // extending only to the lowest data value within this range
458  while (fProj->GetBinContent(bin) == 0 && bin <= fProj->GetNbinsX()) bin++;
459  fWhiskerDown = fProj->GetBinCenter(bin);
460 
461  bin = fProj->FindBin(fBoxUp+1.5*iqr);
462  while (fProj->GetBinContent(bin) == 0 && bin >= 1) bin--;
463  fWhiskerUp = fProj->GetBinCenter(bin);
464  } else { //Need a calculation for a raw-data candle
465  fWhiskerUp = fBoxDown;
466  fWhiskerDown = fBoxUp;
467 
468  //Need to find highest value up to 1.5*iqr from the BoxUp-pos, and the lowest value up to -1.5*iqr from the boxLow-pos
469  for (Long64_t i = 0; i < fNDatapoints; ++i) {
470  Double_t myData = fDatapoints[i];
471  if (myData > fWhiskerUp && myData <= fBoxUp + 1.5*iqr) fWhiskerUp = myData;
472  if (myData < fWhiskerDown && myData >= fBoxDown - 1.5*iqr) fWhiskerDown = myData;
473  }
474  }
475  }
476 
477  if (!fIsRaw && fProj) { //Need a calculation for a projected histo
478  fMean = fProj->GetMean();
479  fMedianErr = 1.57*iqr/sqrt(fProj->GetEntries());
480  fAxisMin = fProj->GetXaxis()->GetXmin();
481  fAxisMax = fProj->GetXaxis()->GetXmax();
482  } else { //Need a calculation for a raw-data candle
483  //Calculate the Mean
484  fMean = 0;
485  for (Long64_t i = 0; i < fNDatapoints; ++i) {
486  fMean += fDatapoints[i];
487  if (fDatapoints[i] < min) min = fDatapoints[i];
488  if (fDatapoints[i] > max) max = fDatapoints[i];
489  if (fDatapoints[i] < fWhiskerDown || fDatapoints[i] > fWhiskerUp) nOutliers++;
490  }
491  fMean /= fNDatapoints;
492  fMedianErr = 1.57*iqr/sqrt(fNDatapoints);
493  }
494 
495  delete [] prob;
496  delete [] quantiles;
497 
498  //Doing the outliers and other single points to show
499  if (GetCandleOption(5) > 0) { //Draw outliers
500  TRandom2 random;
501  const int maxOutliers = kNMAXPOINTS;
502  Double_t myScale = 1.;
503  if (!fIsRaw && fProj) { //Need a calculation for a projected histo
504  if (fProj->GetEntries() > maxOutliers/2) myScale = fProj->GetEntries()/(maxOutliers/2.);
505  fNDrawPoints = 0;
506  for (int bin = 0; bin < fProj->GetNbinsX(); bin++) {
507  // Either show them only outside the whiskers, or all of them
508  if (fProj->GetBinContent(bin) > 0 && (fProj->GetBinCenter(bin) < fWhiskerDown || fProj->GetBinCenter(bin) > fWhiskerUp || (GetCandleOption(5) > 1)) ) {
509  Double_t scaledBinContent = fProj->GetBinContent(bin)/myScale;
510  if (scaledBinContent >0 && scaledBinContent < 1) scaledBinContent = 1; //Outliers have a typical bin content between 0 and 1, when scaling they would disappear
511  for (int j=0; j < (int)scaledBinContent; j++) {
512  if (fNDrawPoints > maxOutliers) break;
513  if (IsOption(kPointsAllScat)) { //Draw outliers and "all" values scattered
514  fDrawPointsX[fNDrawPoints] = fPosCandleAxis - fCandleWidth/2. + fCandleWidth*random.Rndm();
515  fDrawPointsY[fNDrawPoints] = fProj->GetBinLowEdge(bin) + fProj->GetBinWidth(bin)*random.Rndm();
516  } else { //Draw them in the "candle line"
517  fDrawPointsX[fNDrawPoints] = fPosCandleAxis;
518  if ((int)scaledBinContent == 1) //If there is only one datapoint available put it in the middle of the bin
519  fDrawPointsY[fNDrawPoints] = fProj->GetBinCenter(bin);
520  else //If there is more than one datapoint scatter it along the bin, otherwise all marker would be (invisibly) stacked on top of each other
521  fDrawPointsY[fNDrawPoints] = fProj->GetBinLowEdge(bin) + fProj->GetBinWidth(bin)*random.Rndm();
522  }
523  if (swapXY) {
524  //Swap X and Y
525  Double_t keepCurrently;
526  keepCurrently = fDrawPointsX[fNDrawPoints];
527  fDrawPointsX[fNDrawPoints] = fDrawPointsY[fNDrawPoints];
528  fDrawPointsY[fNDrawPoints] = keepCurrently;
529  }
530  // Continue fMeans, that fNDrawPoints is not increased, so that value will not be shown
531  if (doLogX) {
532  if (fDrawPointsX[fNDrawPoints] > 0) fDrawPointsX[fNDrawPoints] = TMath::Log10(fDrawPointsX[fNDrawPoints]); else continue;
533  }
534  if (doLogY) {
535  if (fDrawPointsY[fNDrawPoints] > 0) fDrawPointsY[fNDrawPoints] = TMath::Log10(fDrawPointsY[fNDrawPoints]); else continue;
536  }
537  fNDrawPoints++;
538  }
539  }
540  if (fNDrawPoints > maxOutliers) { //Should never happen, due to myScale!!!
541  Error ("PaintCandlePlot","Not possible to draw all outliers.");
542  break;
543  }
544  }
545  } else { //Raw data candle
546  //If only outliers are shown, calculate myScale only based on nOutliers, use fNDatapoints (all) instead
547  if (IsOption(kPointsOutliers) && nOutliers > maxOutliers/2) {
548  myScale = nOutliers/(maxOutliers/2.);
549  } else {
550  if (fNDatapoints > maxOutliers/2) myScale = fNDatapoints/(maxOutliers/2.);
551  }
552  fNDrawPoints = 0;
553  for (int i = 0; i < fNDatapoints; i++ ) {
554  Double_t myData = fDatapoints[i];
555  Double_t maxScatter = (fWhiskerUp-fWhiskerDown)/100;
556  if (!(i % (int) myScale == 0 )) continue; //If the amount of data is too large take only every 2nd or 3rd to reduce the amount
557  // Either show them only outside the whiskers, or all of them
558  if (myData < fWhiskerDown || myData > fWhiskerUp || (GetCandleOption(5) > 1)) {
559  if (IsOption(kPointsAllScat)) { //Draw outliers and "all" values scattered
560  fDrawPointsX[fNDrawPoints] = fPosCandleAxis - fCandleWidth/2. + fCandleWidth*random.Rndm();
561  fDrawPointsY[fNDrawPoints] = myData + (random.Rndm() - 0.5)*maxScatter; //random +- 0.5 of candle-height
562  } else { //Draw them in the "candle line"
563  fDrawPointsX[fNDrawPoints] = fPosCandleAxis;
564  fDrawPointsY[fNDrawPoints] = myData + (random.Rndm() - 0.5)*maxScatter; //random +- 0.5 of candle-height
565  }
566  if (swapXY) {
567  //Swap X and Y
568  Double_t keepCurrently;
569  keepCurrently = fDrawPointsX[fNDrawPoints];
570  fDrawPointsX[fNDrawPoints] = fDrawPointsY[fNDrawPoints];
571  fDrawPointsY[fNDrawPoints] = keepCurrently;
572  }
573  // Continue fMeans, that fNDrawPoints is not increased, so that value will not be shown
574  if (doLogX) {
575  if (fDrawPointsX[fNDrawPoints] > 0) fDrawPointsX[fNDrawPoints] = TMath::Log10(fDrawPointsX[fNDrawPoints]);
576  else continue;
577  }
578  if (doLogY) {
579  if (fDrawPointsY[fNDrawPoints] > 0) fDrawPointsY[fNDrawPoints] = TMath::Log10(fDrawPointsY[fNDrawPoints]);
580  else continue;
581  }
582  fNDrawPoints++;
583  if (fNDrawPoints > maxOutliers) { //Should never happen, due to myScale!!!
584  Error ("PaintCandlePlotRaw","Not possible to draw all outliers.");
585  break;
586  }
587  }
588  }
589  }
590  }
591  if (IsOption(kHistoRight) || IsOption(kHistoLeft) || IsOption(kHistoViolin)) {
592  //We are starting with kHistoRight, left will be modified from right later
593  if (fIsRaw) { //This is a raw-data candle
594  if (!fProj) {
595  fProj = new TH1D("hpa","hpa",100,min,max+0.0001*(max-min));
596  for (Long64_t i = 0; i < fNDatapoints; ++i) {
597  fProj->Fill(fDatapoints[i]);
598  }
599  }
600  }
601 
602  fNHistoPoints = 0;
603  Double_t maxContent = fProj->GetMaximum();
604  Double_t maxHistoHeight = fHistoWidth;
605  if (IsOption(kHistoViolin)) maxHistoHeight *= 0.5;
606 
607  bool isFirst = true;
608  int lastNonZero = 0;
609  for (int bin = 1; bin <= fProj->GetNbinsX(); bin++) {
610  if (isFirst) {
611  if (fProj->GetBinContent(bin) > 0) {
612  fHistoPointsX[fNHistoPoints] = fPosCandleAxis;
613  fHistoPointsY[fNHistoPoints] = fProj->GetBinLowEdge(bin);
614  if (doLogX) {
615  if (fHistoPointsX[fNHistoPoints] > 0) fHistoPointsX[fNHistoPoints] = TMath::Log10(fHistoPointsX[fNHistoPoints]); else continue;
616  }
617  if (doLogY) {
618  if (fHistoPointsY[fNHistoPoints] > 0) fHistoPointsY[fNHistoPoints] = TMath::Log10(fHistoPointsY[fNHistoPoints]); else continue;
619  }
620  fNHistoPoints++;
621  isFirst = false;
622  } else {
623  continue;
624  }
625  }
626 
627  Double_t myBinValue = fProj->GetBinContent(bin);
628  if (doLogZ) {
629  if (myBinValue > 0) myBinValue = TMath::Log10(myBinValue); else myBinValue = 0;
630  }
631  fHistoPointsX[fNHistoPoints] = fPosCandleAxis + myBinValue/maxContent*maxHistoHeight;
632  fHistoPointsY[fNHistoPoints] = fProj->GetBinLowEdge(bin);
633  fNHistoPoints++;
634  fHistoPointsX[fNHistoPoints] = fPosCandleAxis + myBinValue/maxContent*maxHistoHeight;
635  fHistoPointsY[fNHistoPoints] = fProj->GetBinLowEdge(bin)+fProj->GetBinWidth(bin);
636  if (doLogX) {
637  if (fHistoPointsX[fNHistoPoints -1] > 0) fHistoPointsX[fNHistoPoints - 1] = TMath::Log10(fHistoPointsX[fNHistoPoints - 1]); else continue;
638  if (fHistoPointsX[fNHistoPoints] > 0) fHistoPointsX[fNHistoPoints] = TMath::Log10(fHistoPointsX[fNHistoPoints]); else continue;
639  }
640  if (doLogY) {
641  if (fHistoPointsY[fNHistoPoints -1] > 0) fHistoPointsY[fNHistoPoints - 1] = TMath::Log10(fHistoPointsY[fNHistoPoints - 1]); else continue;
642  if (fHistoPointsY[fNHistoPoints] > 0) fHistoPointsY[fNHistoPoints] = TMath::Log10(fHistoPointsY[fNHistoPoints]); else continue;
643  }
644 
645  fNHistoPoints++;
646  if (fProj->GetBinContent(bin) > 0) lastNonZero = fNHistoPoints;
647  }
648 
649  fHistoPointsX[fNHistoPoints] = fPosCandleAxis;
650  fHistoPointsY[fNHistoPoints] = fHistoPointsY[fNHistoPoints-1];
651  fNHistoPoints = lastNonZero+1; //+1 so that the line down to 0 is added as well
652 
653  if (IsOption(kHistoLeft)) {
654  for (int i = 0; i < fNHistoPoints; i++) {
655  fHistoPointsX[i] = 2*fPosCandleAxis - fHistoPointsX[i];
656  }
657  }
658  if (IsOption(kHistoViolin)) {
659  for (int i = 0; i < fNHistoPoints; i++) {
660  fHistoPointsX[fNHistoPoints + i] = 2*fPosCandleAxis - fHistoPointsX[fNHistoPoints -i-1];
661  fHistoPointsY[fNHistoPoints + i] = fHistoPointsY[fNHistoPoints -i-1];
662  }
663  fNHistoPoints *= 2;
664  }
665  }
666 
667 
668  fIsCalculated = true;
669 }
670 
671 ////////////////////////////////////////////////////////////////////////////////
672 /// Paint one candle with its current attributes.
673 
674 void TCandle::Paint(Option_t *)
675 {
676  //If something was changed before, we need to recalculate some values
677  if (!fIsCalculated) Calculate();
678 
679  // Save the attributes as they were set originally
680  Style_t saveLine = GetLineStyle();
681  Style_t saveMarker = GetMarkerStyle();
682  Style_t saveFillStyle = GetFillStyle();
683  Style_t saveFillColor = GetFillColor();
684  Style_t saveLineColor = GetLineColor();
685 
686  Double_t dimLeft = fPosCandleAxis-0.5*fCandleWidth;
687  Double_t dimRight = fPosCandleAxis+0.5*fCandleWidth;
688 
689  TAttLine::Modify();
690  TAttFill::Modify();
691  TAttMarker::Modify();
692 
693  Bool_t swapXY = IsOption(kHorizontal);
694  Bool_t doLogY = (!(swapXY) && fLogY) || (swapXY && fLogX);
695  Bool_t doLogX = (!(swapXY) && fLogX) || (swapXY && fLogY);
696 
697  // From now on this is real painting only, no calculations anymore
698 
699  if (IsOption(kHistoZeroIndicator)) {
700  SetLineColor(saveFillColor);
701  TAttLine::Modify();
702  PaintLine(fPosCandleAxis, fAxisMin, fPosCandleAxis, fAxisMax, swapXY);
703  SetLineColor(saveLineColor);
704  TAttLine::Modify();
705  }
706 
707 
708  if (IsOption(kHistoRight) || IsOption(kHistoLeft) || IsOption(kHistoViolin)) {
709  if (IsOption(kHistoZeroIndicator) && (saveFillStyle != 0)) {
710  SetLineColor(saveFillColor);
711  TAttLine::Modify();
712  }
713  if (!swapXY) {
714  gPad->PaintFillArea(fNHistoPoints, fHistoPointsX, fHistoPointsY);
715  gPad->PaintPolyLine(fNHistoPoints, fHistoPointsX, fHistoPointsY);
716  } else {
717  gPad->PaintFillArea(fNHistoPoints, fHistoPointsY, fHistoPointsX);
718  gPad->PaintPolyLine(fNHistoPoints, fHistoPointsY, fHistoPointsX);
719  }
720  if (IsOption(kHistoZeroIndicator) && (saveFillStyle != 0)) {
721  SetLineColor(saveLineColor);
722  TAttLine::Modify();
723  }
724  }
725 
726  if (IsOption(kBox)) { // Draw a simple box
727  if (IsOption(kMedianNotched)) { // Check if we have to draw a box with notches
728  Double_t x[] = {dimLeft, dimLeft, dimLeft+fCandleWidth/3., dimLeft, dimLeft, dimRight,
729  dimRight, dimRight-fCandleWidth/3., dimRight, dimRight, dimLeft};
730  Double_t y[] = {fBoxDown, fMedian-fMedianErr, fMedian, fMedian+fMedianErr, fBoxUp, fBoxUp,
731  fMedian+fMedianErr, fMedian, fMedian-fMedianErr, fBoxDown, fBoxDown};
732  PaintBox(11, x, y, swapXY);
733  } else { // draw a simple box
734  Double_t x[] = {dimLeft, dimLeft, dimRight, dimRight, dimLeft};
735  Double_t y[] = {fBoxDown, fBoxUp, fBoxUp, fBoxDown, fBoxDown};
736  PaintBox(5, x, y, swapXY);
737  }
738  }
739 
740  if (IsOption(kAnchor)) { // Draw the anchor line
741  PaintLine(dimLeft, fWhiskerUp, dimRight, fWhiskerUp, swapXY);
742  PaintLine(dimLeft, fWhiskerDown, dimRight, fWhiskerDown, swapXY);
743  }
744 
745  if (IsOption(kWhiskerAll) && !IsOption(kHistoZeroIndicator)) { // Whiskers are dashed
746  SetLineStyle(2);
747  TAttLine::Modify();
748  PaintLine(fPosCandleAxis, fWhiskerUp, fPosCandleAxis, fBoxUp, swapXY);
749  PaintLine(fPosCandleAxis, fBoxDown, fPosCandleAxis, fWhiskerDown, swapXY);
750  SetLineStyle(saveLine);
751  TAttLine::Modify();
752  } else if ((IsOption(kWhiskerAll) && IsOption(kHistoZeroIndicator)) || IsOption(kWhisker15) ) { // Whiskers without dashing, better whisker definition, or forced when using zero line
753  PaintLine(fPosCandleAxis, fWhiskerUp, fPosCandleAxis, fBoxUp, swapXY);
754  PaintLine(fPosCandleAxis, fBoxDown, fPosCandleAxis, fWhiskerDown, swapXY);
755  }
756 
757  if (IsOption(kMedianLine)) { // Paint fMedian as a line
758  PaintLine(dimLeft, fMedian, dimRight, fMedian, swapXY);
759  } else if (IsOption(kMedianNotched)) { // Paint fMedian as a line (using notches, fMedian line is shorter)
760  PaintLine(dimLeft+fCandleWidth/3, fMedian, dimRight-fCandleWidth/3., fMedian, swapXY);
761  } else if (IsOption(kMedianCircle)) { // Paint fMedian circle
762  Double_t myMedianX[1], myMedianY[1];
763  if (!swapXY) {
764  myMedianX[0] = fPosCandleAxis;
765  myMedianY[0] = fMedian;
766  } else {
767  myMedianX[0] = fMedian;
768  myMedianY[0] = fPosCandleAxis;
769  }
770 
771  Bool_t isValid = true;
772  if (doLogX) {
773  if (myMedianX[0] > 0) myMedianX[0] = TMath::Log10(myMedianX[0]); else isValid = false;
774  }
775  if (doLogY) {
776  if (myMedianY[0] > 0) myMedianY[0] = TMath::Log10(myMedianY[0]); else isValid = false;
777  }
778 
779  SetMarkerStyle(24);
780  TAttMarker::Modify();
781 
782  if (isValid) gPad->PaintPolyMarker(1,myMedianX,myMedianY); // A circle for the fMedian
783 
784  SetMarkerStyle(saveMarker);
785  TAttMarker::Modify();
786 
787  }
788 
789  if (IsOption(kMeanCircle)) { // Paint fMean as a circle
790  Double_t myMeanX[1], myMeanY[1];
791  if (!swapXY) {
792  myMeanX[0] = fPosCandleAxis;
793  myMeanY[0] = fMean;
794  } else {
795  myMeanX[0] = fMean;
796  myMeanY[0] = fPosCandleAxis;
797  }
798 
799  Bool_t isValid = true;
800  if (doLogX) {
801  if (myMeanX[0] > 0) myMeanX[0] = TMath::Log10(myMeanX[0]); else isValid = false;
802  }
803  if (doLogY) {
804  if (myMeanY[0] > 0) myMeanY[0] = TMath::Log10(myMeanY[0]); else isValid = false;
805  }
806 
807  SetMarkerStyle(24);
808  TAttMarker::Modify();
809 
810  if (isValid) gPad->PaintPolyMarker(1,myMeanX,myMeanY); // A circle for the fMean
811 
812  SetMarkerStyle(saveMarker);
813  TAttMarker::Modify();
814 
815  } else if (IsOption(kMeanLine)) { // Paint fMean as a dashed line
816  SetLineStyle(2);
817  TAttLine::Modify();
818 
819  PaintLine(dimLeft, fMean, dimRight, fMean, swapXY);
820  SetLineStyle(saveLine);
821  TAttLine::Modify();
822 
823  }
824 
825  if (IsOption(kAnchor)) { //Draw standard anchor
826  PaintLine(dimLeft, fWhiskerDown, dimRight, fWhiskerDown, swapXY); // the lower anchor line
827  PaintLine(dimLeft, fWhiskerUp, dimRight, fWhiskerUp, swapXY); // the upper anchor line
828  }
829 
830  // This is a bit complex. All values here are handled as outliers. Usually
831  // only the datapoints outside the whiskers are shown.
832  // One can show them in one row as crosses, or scattered randomly. If activated
833  // all datapoint are shown in the same way
834 
835  if (GetCandleOption(5) > 0) { //Draw outliers
836  if (IsOption(kPointsAllScat)) { //Draw outliers and "all" values scattered
837  SetMarkerStyle(0);
838  } else {
839  SetMarkerStyle(5);
840  }
841  TAttMarker::Modify();
842  gPad->PaintPolyMarker(fNDrawPoints,fDrawPointsX, fDrawPointsY);
843  }
844 }
845 
846 ////////////////////////////////////////////////////////////////////////////////
847 /// Return true is this option is activated in fOption
848 
849 bool TCandle::IsOption(CandleOption opt) {
850  long myOpt = 9;
851  int pos = 0;
852  for (pos = 0; pos < 16; pos++) {
853  if (myOpt > opt) break;
854  else myOpt *=10;
855  }
856  myOpt /= 9;
857  int thisOpt = GetCandleOption(pos);
858 
859  return ((thisOpt * myOpt) == opt);
860 }
861 
862 ////////////////////////////////////////////////////////////////////////////////
863 /// Paint a box for candle.
864 
865 void TCandle::PaintBox(Int_t nPoints, Double_t *x, Double_t *y, Bool_t swapXY)
866 {
867  Bool_t doLogY = (!(swapXY) && fLogY) || (swapXY && fLogX);
868  Bool_t doLogX = (!(swapXY) && fLogX) || (swapXY && fLogY);
869  if (doLogY) {
870  for (int i=0; i<nPoints; i++) {
871  if (y[i] > 0) y[i] = TMath::Log10(y[i]);
872  else return;
873  }
874  }
875  if (doLogX) {
876  for (int i=0; i<nPoints; i++) {
877  if (x[i] > 0) x[i] = TMath::Log10(x[i]);
878  else return;
879  }
880  }
881  if (!swapXY) {
882  gPad->PaintFillArea(nPoints, x, y);
883  gPad->PaintPolyLine(nPoints, x, y);
884  } else {
885  gPad->PaintFillArea(nPoints, y, x);
886  gPad->PaintPolyLine(nPoints, y, x);
887  }
888 }
889 
890 ////////////////////////////////////////////////////////////////////////////////
891 /// Paint a line for candle.
892 
893 void TCandle::PaintLine(Double_t x1, Double_t y1, Double_t x2, Double_t y2, Bool_t swapXY)
894 {
895  Bool_t doLogY = (!(swapXY) && fLogY) || (swapXY && fLogX);
896  Bool_t doLogX = (!(swapXY) && fLogX) || (swapXY && fLogY);
897  if (doLogY) {
898  if (y1 > 0) y1 = TMath::Log10(y1); else return;
899  if (y2 > 0) y2 = TMath::Log10(y2); else return;
900  }
901  if (doLogX) {
902  if (x1 > 0) x1 = TMath::Log10(x1); else return;
903  if (x2 > 0) x2 = TMath::Log10(x2); else return;
904  }
905  if (!swapXY) {
906  gPad->PaintLine(x1, y1, x2, y2);
907  } else {
908  gPad->PaintLine(y1, x1, y2, x2);
909  }
910 }
911 
912 ////////////////////////////////////////////////////////////////////////////////
913 /// Stream an object of class TCandle.
914 
915 void TCandle::Streamer(TBuffer &R__b)
916 {
917  if (R__b.IsReading()) {
918  UInt_t R__s, R__c;
919  Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
920  if (R__v > 3) {
921  R__b.ReadClassBuffer(TCandle::Class(), this, R__v, R__s, R__c);
922  return;
923  }
924  } else {
925  R__b.WriteClassBuffer(TCandle::Class(),this);
926  }
927 }
928 
929 ////////////////////////////////////////////////////////////////////////////////
930 /// The coordinates in the TParallelCoordVar-class are in Pad-Coordinates, so we need to convert them
931 
932 void TCandle::ConvertToPadCoords(Double_t minAxis, Double_t maxAxis, Double_t axisMinCoord, Double_t axisMaxCoord)
933 {
934  if (!fIsCalculated) Calculate();
935  Double_t a,b;
936  if (fLogY) {
937  a = TMath::Log10(minAxis);
938  b = TMath::Log10(maxAxis/minAxis);
939  } else {
940  a = minAxis;
941  b = maxAxis-minAxis;
942  }
943 
944  fMean = axisMinCoord + ((fMean-a)/b)*(axisMaxCoord-axisMinCoord);
945  fMedian = axisMinCoord + ((fMedian-a)/b)*(axisMaxCoord-axisMinCoord);
946  fMedianErr = axisMinCoord + ((fMedianErr-a)/b)*(axisMaxCoord-axisMinCoord);
947  fBoxUp = axisMinCoord + ((fBoxUp-a)/b)*(axisMaxCoord-axisMinCoord);
948  fBoxDown = axisMinCoord + ((fBoxDown-a)/b)*(axisMaxCoord-axisMinCoord);
949  fWhiskerUp = axisMinCoord + ((fWhiskerUp-a)/b)*(axisMaxCoord-axisMinCoord);
950  fWhiskerDown = axisMinCoord + ((fWhiskerDown-a)/b)*(axisMaxCoord-axisMinCoord);
951 
952  for (int i = 0; i < fNDrawPoints; i++) {
953  fDrawPointsY[i] = axisMinCoord + ((fDrawPointsY[i]-a)/b)*(axisMaxCoord-axisMinCoord);
954  }
955  for (int i = 0; i < fNHistoPoints; i++) {
956  fHistoPointsY[i] = axisMinCoord + ((fHistoPointsY[i]-a)/b)*(axisMaxCoord-axisMinCoord);
957  }
958 }