Logo ROOT   6.08/07
Reference Guide
TGaxis.cxx
Go to the documentation of this file.
1 // @(#)root/graf:$Id$
2 // Author: Rene Brun, Olivier Couet 12/12/94
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 <string.h>
14 #include <time.h>
15 #include <math.h>
16 
17 #include "Riostream.h"
18 #include "TROOT.h"
19 #include "TGaxis.h"
20 #include "TAxisModLab.h"
21 #include "TVirtualPad.h"
22 #include "TVirtualX.h"
23 #include "TLine.h"
24 #include "TLatex.h"
25 #include "TStyle.h"
26 #include "TF1.h"
27 #include "TAxis.h"
28 #include "THashList.h"
29 #include "TObjString.h"
30 #include "TObject.h"
31 #include "TMath.h"
32 #include "THLimitsFinder.h"
33 #include "TColor.h"
34 #include "TClass.h"
35 #include "TTimeStamp.h"
36 #include "TSystem.h"
37 #include "TTimeStamp.h"
38 
40 Float_t TGaxis::fXAxisExpXOffset = 0.; //Exponent X offset for the X axis
41 Float_t TGaxis::fXAxisExpYOffset = 0.; //Exponent Y offset for the X axis
42 Float_t TGaxis::fYAxisExpXOffset = 0.; //Exponent X offset for the Y axis
43 Float_t TGaxis::fYAxisExpYOffset = 0.; //Exponent Y offset for the Y axis
44 const Int_t kHori = BIT(9); //defined in TPad
45 
47 
48 /** \class TGaxis
49 \ingroup BasicGraphics
50 
51 The axis painter class.
52 
53 Instances of this class are generated by the histograms and graphs painting
54 classes when `TAxis` are drawn. `TGaxis` is the "painter class" of
55 `TAxis`. Therefore it is mainly used via `TAxis`, even if is some
56 occasion it can be used directly to draw an axis which is not part of a graph
57 or an instance. For instance to draw an extra scale on a plot.
58 
59 - [Basic definition](#GA00)
60 - [Definition with a function](#GA01)
61 - [Logarithmic axis](#GA02)
62 - [Blank axis](#GA03)
63 - [Tick marks' orientation](#GA04)
64 - [Tick marks' size](#GA05)
65 - [Labels' positionning](#GA06)
66 - [Labels' orientation](#GA07)
67 - [Labels' position on tick marks](#GA08)
68 - [Labels' format](#GA09)
69 - [Alphanumeric labels](#GA10)
70 - [Changing axis labels](#GA10a)
71 - [Number of divisions optimisation](#GA11)
72 - [Maximum Number of Digits for the axis labels](#GA12)
73 - [Optional grid](#GA13)
74 - [Time axis](#GA14)
75 
76 ## <a name="GA00"></a> Basic definition
77 A `TGaxis` is defined the following way:
78 ~~~ {.cpp}
79  TGaxis::TGaxis(Double_t xmin, Double_t ymin, Double_t xmax, Double_t ymax,
80  Double_t wmin, Double_t wmax, Int_t ndiv, Option_t *chopt,
81  Double_t gridlength)
82 ~~~
83 Where:
84 
85 - xmin : X origin coordinate in user's coordinates space.
86 - xmax : X end axis coordinate in user's coordinates space.
87 - ymin : Y origin coordinate in user's coordinates space.
88 - ymax : Y end axis coordinate in user's coordinates space.
89 - wmin : Lowest value for the tick mark labels written on the axis.
90 - wmax : Highest value for the tick mark labels written on the axis.
91 - ndiv : Number of divisions.
92  - ndiv=N1 + 100*N2 + 10000*N3
93  - N1=number of 1st divisions.
94  - N2=number of 2nd divisions.
95  - N3=number of 3rd divisions. e.g.:
96  - ndiv=0 --> no tick marks.
97  - ndiv=2 --> 2 divisions, one tick mark in the middle of the axis.
98 - chopt : Drawing options (see below).
99 - gridlength: grid length on main tick marks.
100 
101 The example below generates various kind of axis.
102 
103 Begin_Macro(source)
104 {
105  TCanvas *c1 = new TCanvas("c1","Examples of TGaxis",10,10,700,500);
106 
107  c1->Range(-10,-1,10,1);
108 
109  TGaxis *axis1 = new TGaxis(-4.5,-0.2,5.5,-0.2,-6,8,510,"");
110  axis1->SetName("axis1");
111  axis1->Draw();
112 
113  TGaxis *axis2 = new TGaxis(-4.5,0.2,5.5,0.2,0.001,10000,510,"G");
114  axis2->SetName("axis2");
115  axis2->Draw();
116 
117  TGaxis *axis3 = new TGaxis(-9,-0.8,-9,0.8,-8,8,50510,"");
118  axis3->SetName("axis3");
119  axis3->Draw();
120 
121  TGaxis *axis4 = new TGaxis(-7,-0.8,-7,0.8,1,10000,50510,"G");
122  axis4->SetName("axis4");
123  axis4->Draw();
124 
125  TGaxis *axis5 = new TGaxis(-4.5,-0.6,5.5,-0.6,1.2,1.32,80506,"-+");
126  axis5->SetName("axis5");
127  axis5->SetLabelSize(0.03);
128  axis5->SetTextFont(72);
129  axis5->SetLabelOffset(0.025);
130 
131  axis5->Draw();
132 
133  TGaxis *axis6 = new TGaxis(-4.5,0.6,5.5,0.6,100,900,50510,"-");
134  axis6->SetName("axis6");
135  axis6->Draw();
136 
137  TGaxis *axis7 = new TGaxis(8,-0.8,8,0.8,0,9000,50510,"+L");
138  axis7->SetName("axis7");
139  axis7->SetLabelOffset(0.01);
140  axis7->Draw();
141 
142  //one can make axis going top->bottom. However because of a long standing
143  //problem, the two x values should not be equal
144  TGaxis *axis8 = new TGaxis(6.5,0.8,6.499,-0.8,0,90,50510,"-");
145  axis8->SetName("axis8");
146  axis8->Draw();
147  return c1;
148 }
149 End_Macro
150 
151 ## <a name="GA01"></a> Definition with a function
152 
153 Instead of the wmin,wmax arguments of the normal definition, the
154 name of a `TF1` function can be specified. This function will be used to
155 map the user coordinates to the axis values and ticks.
156 
157 A `TGaxis` is defined the following way:
158 ~~~ {.cpp}
159  TGaxis::TGaxis(Double_t xmin, Double_t ymin, Double_t xmax, Double_t ymax,
160  const char *func, Int_t ndiv, Option_t *chopt,
161  Double_t gridlength)
162 ~~~
163 Where:
164 
165 - xmin : X origin coordinate in user's coordinates space.
166 - xmax : X end axis coordinate in user's coordinates space.
167 - ymin : Y origin coordinate in user's coordinates space.
168 - ymax : Y end axis coordinate in user's coordinates space.
169 - func : function defining axis labels and tick marks.
170 - ndiv : Number of divisions.
171  - ndiv=N1 + 100*N2 + 10000*N3
172  - N1=number of 1st divisions.
173  - N2=number of 2nd divisions.
174  - N3=number of 3rd divisions. e.g.:
175  - ndiv=0 --> no tick marks.
176  - ndiv=2 --> 2 divisions, one tick mark in the middle of the axis.
177 - chopt : Drawing options (see below).
178 - gridlength: grid length on main tick marks.
179 
180 Examples:
181 
182 Begin_Macro(source)
183 {
184  TCanvas *c2 = new TCanvas("c2","c2",10,10,700,500);
185 
186  gStyle->SetOptStat(0);
187 
188  TH2F *h2 = new TH2F("h","Axes",100,0,10,100,-2,2);
189  h2->Draw();
190 
191  TF1 *f1=new TF1("f1","-x",-10,10);
192  TGaxis *A1 = new TGaxis(0,2,10,2,"f1",510,"-");
193  A1->SetTitle("axis with decreasing values");
194  A1->Draw();
195 
196  TF1 *f2=new TF1("f2","exp(x)",0,2);
197  TGaxis *A2 = new TGaxis(1,1,9,1,"f2");
198  A2->SetTitle("exponential axis");
199  A2->SetLabelSize(0.03);
200  A2->SetTitleSize(0.03);
201  A2->SetTitleOffset(1.2);
202  A2->Draw();
203 
204  TF1 *f3=new TF1("f3","log10(x)",1,1000);
205  TGaxis *A3 = new TGaxis(2,-2,2,0,"f3",505,"G");
206  A3->SetTitle("logarithmic axis");
207  A3->SetLabelSize(0.03);
208  A3->SetTitleSize(0.03);
209  A3->SetTitleOffset(1.2);
210  A3->Draw();
211  return c2;
212 }
213 End_Macro
214 
215 
216 ## <a name="GA02"></a> Logarithmic axis
217 
218 By default axis are linear. To define a `TGaxis` as logarithmic, it is
219 enough to create it with the option `"G"`.
220 
221 When plotting an histogram or a graph the logarithmic scale can be set using:
222 
223  - `gPad->SetLogx(1);` set the logarithmic scale on the X axis
224  - `gPad->SetLogy(1);` set the logarithmic scale on the Y axis
225 
226 When the `SetMoreLogLabels()` method is called more labels are drawn
227 when in logarithmic scale and there is a small number of decades (less than 3).
228 
229 ## <a name="GA03"></a> Blank axis
230 To draw only the axis tick marks without the axis body, it is enough to specify
231 the option `"B"`. It useful to superpose axis.
232 
233 ## <a name="GA04"></a> Tick marks' orientation
234 
235 By default tick marks are drawn on the positive side of the axis, except for
236 vertical axis for which the default is negative. The `chop` parameter
237 allows to control the tick marks orientation:
238 
239  - `chopt = "+"`: tick marks are drawn on Positive side. (default)
240  - `chopt ="-"`: tick mark are drawn on the negative side.
241  - `chopt = "+-"`: tick marks are drawn on both sides of the axis.
242  - `chopt = "U"`: Unlabelled axis, default is labeled.
243 
244 ## <a name="GA05"></a> Tick marks' size
245 
246 By default, tick marks have a length equal to 3 per cent of the axis length.
247 When the option "S" is specified, the length of the tick marks is equal to
248 `fTickSize*axis_length`, where `fTickSize` may be set via
249 `TGaxis::SetTickSize`.
250 
251 When plotting an histogram `h` the tick marks size can be changed using:
252 
253  - `h->GetXaxis()->SetTickLength(0.02);` set the tick length for the X axis
254  - `gStyle->SetTickLength(0.02,"x");` set the tick length for the X axis
255  of all histograms drawn after this instruction.
256 
257 A good way to remove tick marks on an axis is to set the tick length to 0:
258 `h->GetXaxis()->SetTickLength(0.);`
259 
260 ## <a name="GA06"></a> Labels' positionning
261 
262 Labels are normally drawn on side opposite to tick marks. However the option
263 `"="` allows to draw them on the same side.
264 
265 ## <a name="GA07"></a> Labels' orientation
266 
267 By default axis labels are drawn parallel to the axis. However if the axis is vertical
268 then are drawn perpendicular to the axis.
269 
270 ## <a name="GA08"></a> Labels' position on tick marks
271 
272 By default axis labels are centered on tick marks. However, for vertical axis,
273 they are right adjusted. The `chop` parameter allows to control the labels'
274 position on tick marks:
275 
276  - `chopt = "R"`: labels are Right adjusted on tick mark.(default is centered)
277  - `chopt = "L"`: labels are Left adjusted on tick mark.
278  - `chopt = "C"`: labels are Centered on tick mark.
279  - `chopt = "M"`: In the Middle of the divisions.
280 
281 ## <a name="GA09"></a> Labels' format
282 
283 Blank characters are stripped, and then the label is correctly aligned. the dot,
284 if last character of the string, is also stripped, unless the option `"."`
285 (a dot, or period) is specified. if `SetDecimals(kTRUE)` has been called
286 all labels have the same number of decimals after the `"."`
287 The same is true if `gStyle->SetStripDecimals(kFALSE)` has been called.
288 
289 In the following, we have some parameters, like tick marks length and characters
290 height (in percentage of the length of the axis (user's coordinates))
291 The default values are as follows:
292 
293  - Primary tick marks: 3.0 %
294  - Secondary tick marks: 1.5 %
295  - Third order tick marks: .75 %
296  - Characters height for labels: 4%
297  - Labels offset: 1.0 %
298 
299 By default, an exponent of the form 10^N is used when the label values are either
300 all very small or very large. One can disable the exponent by calling
301 `axis.SetNoExponent(kTRUE)`.
302 
303 `TGaxis::SetExponentOffset(Float_t xoff, Float_t yoff, Option_t *axis)` is
304 static function to set X and Y offset of the axis 10^n notation. It is in % of
305 the pad size. It can be negative. `axis` specifies which axis
306 (`"x"` or/and `"y"`), default is `"x"` if `axis = "xz"`
307 set the two axes
308 
309 ## <a name="GA10"></a> Alphanumeric labels
310 
311 Axis labels can be any alphanumeric character strings. Such axis can be produced
312 only with histograms because the labels'definition is stored in `TAxis`.
313 The following example demonstrates how to create such labels.
314 
315 Begin_Macro(source)
316 ../../../tutorials/hist/hlabels2.C
317 End_Macro
318 
319 Because the alphanumeric labels are usually longer that the numeric labels, their
320 size is by default equal to `0.66666 * the_numeric_labels_size`.
321 
322 ## <a name="GA10a"></a> Changing axis labels
323 \since **ROOT version 6.07/07:**
324 
325 After an axis has been created, TGaxis::ChangeLabel allows to define new text
326 attributes for a given label. A fine tuning of the labels can be done. All the
327 attributes can be changed as well as the text label itself.
328 
329 When plotting an histogram or a graph the labels can be changed like in the
330 following example which shows a way to produce \f$\pi\f$-axis :
331 
332 Begin_Macro(source)
333 {
334  Double_t pi = TMath::Pi();
335  TF1* f = new TF1("f","TMath::Cos(x/TMath::Pi())", -pi, pi);
336  TH1* h = f->GetHistogram();
337  TAxis* a = h->GetXaxis();
338  a->SetNdivisions(-502);
339  a->ChangeLabel(1,-1,-1,-1,-1,-1,"-#pi");
340  a->ChangeLabel(-1,-1,-1,-1,-1,-1,"#pi");
341  f->Draw();
342 }
343 End_Macro
344 
345 ## <a name="GA11"></a> Number of divisions optimisation
346 
347 By default the number of divisions on axis is optimised to show a coherent
348 labelling of the main tick marks. The number of division (`ndiv`) is a
349 composite integer given by:
350 
351 ` ndiv = N1 + 100*N2 + 10000*N3`
352 
353  - `N1` = number of 1st divisions.
354  - `N2` = number of 2nd divisions.
355  - `N3` = number of 3rd divisions.
356 
357 by default the value of `N1`, `N2` and `N3` are maximum
358 values. After optimisation the real number of divisions will be smaller or
359 equal to these value. If one wants to bypass the optimisation, the option `"N"`
360 should be given when the `TGaxis` is created. The option `"I"`
361 also act on the number of division as it will force an integer labelling of
362 the axis.
363 
364 On an histogram pointer `h` the number of divisions can be set in different ways:.
365 
366 - Directly on the histogram. The following will set the number of division
367  to 510 on the X axis of `h`. To avoid optimization the number of divisions
368  should be negative (ie: -510);
369 ~~~ {.cpp}
370  h->SetNdivisions(510, "X");
371 ~~~
372 - On the axis itself:
373 ~~~ {.cpp}
374  h->GetXaxis()->SetNdivisions(510, kTRUE);
375 ~~~
376 
377 The first parameter is the number of division. If it is negative of if the
378 second parameter is kFALSE then the number of divisions is not optimised.
379 And other signature is also allowed:
380 ~~~ {.cpp}
381  h->GetXaxis()->SetNdivisions(10, 5, 0, kTRUE);
382 ~~~
383 ## <a name="GA12"></a> Maximum Number of Digits for the axis labels
384 
385 The static function `TGaxis::SetMaxDigits` sets the maximum number of
386 digits permitted for the axis labels above which the notation with 10^N is used.
387 For example, to accept 6 digits number like 900000 on an axis call
388 `TGaxis::SetMaxDigits(6)`. The default value is 5.
389 `fgMaxDigits` must be greater than 0.
390 
391 ## <a name="GA13"></a> Optional grid
392 
393 The option `"W"` allows to draw a grid on the primary tick marks. In case
394 of a log axis, the grid is only drawn for the primary tick marks if the number
395 of secondary and tertiary divisions is 0. `SetGridLength()` allows to define
396 the length of the grid.
397 
398 When plotting an histogram or a graph the grid can be set ON or OFF using:
399 
400  - `gPad->SetGridy(1);` set the grid on the X axis
401  - `gPad->SetGridx(1);` set the grid on the Y axis
402  - `gPad->SetGrid(1,1);` set the grid on both axis.
403 
404 ## <a name="GA14"></a> Time axis
405 
406 Axis labels may be considered as times, plotted in a defined time format.
407 The format is set with `SetTimeFormat()`.
408 The `TGaxis` minimum (`wmin`) and maximum (`wmax`) values
409 are considered as two time values in seconds.
410 The time axis will be spread around the time offset value (set with
411 `SetTimeOffset()`). Actually it will go from `TimeOffset+wmin` to
412 `TimeOffset+wmax`
413 
414 Usually time axis are created automatically via histograms, but one may also
415 want to draw a time axis outside an "histogram context". This can be done
416 thanks to the option `"T"` of `TGaxis`.
417 
418 Begin_Macro(source)
419 {
420  c1 = new TCanvas("c1","Examples of TGaxis",10,10,700,100);
421  c1->Range(-10,-1,10,1);
422 
423  TGaxis *axis = new TGaxis(-8,0.,8,0.,-100000,150000,2405,"tS");
424  axis->SetLabelSize(0.3);
425  axis->SetTickSize(0.2);
426 
427  TDatime da(2003,02,28,12,00,00);
428  axis->SetTimeOffset(da.Convert());
429  axis->SetTimeFormat("%d-%m-%Y");
430  axis->Draw();
431  return c1;
432 }
433 End_Macro
434 
435 The following example compares what the system time function `gmtime`
436 and `localtime` give with what gives `TGaxis`. It can be used
437 as referenced test to check if the time option of `TGaxis` is working properly.
438 
439 Begin_Macro(source)
440 ../../../tutorials/graphs/timeonaxis3.C
441 End_Macro
442 
443 
444 The following macro illustrates the use, with histograms axis, of the time mode on the axis
445 with different time intervals and time formats.
446 
447 Begin_Macro(source)
448 ../../../tutorials/graphs/timeonaxis.C
449 End_Macro
450 
451 An other example showing how to define the time offset as 2003, January 1st
452 using histograms axis.
453 
454 Begin_Macro(source)
455 ../../../tutorials/graphs/timeonaxis2.C
456 End_Macro
457 */
458 
459 ////////////////////////////////////////////////////////////////////////////////
460 /// TGaxis default constructor.
461 
462 TGaxis::TGaxis(): TLine(), TAttText(11,0,1,62,0.040)
463 {
464 
465  fGridLength = 0.;
466  fLabelOffset = 0.005;
467  fLabelSize = 0.040;
468  fLabelFont = 62;
469  fLabelColor = 1;
470  fTickSize = 0.030;
471  fTitleOffset = 1;
472  fTitleSize = fLabelSize;
473  fChopt = "";
474  fName = "";
475  fTitle = "";
476  fTimeFormat = "";
477  fFunctionName= "";
478  fFunction = 0;
479  fAxis = 0;
480  fNdiv = 0;
481  fNModLabs = 0;
482  fModLabs = 0;
483  fWmin = 0.;
484  fWmax = 0.;
485 }
486 
487 ////////////////////////////////////////////////////////////////////////////////
488 /// TGaxis normal constructor.
489 
491  Double_t wmin, Double_t wmax, Int_t ndiv, Option_t *chopt,
492  Double_t gridlength)
493  : TLine(xmin,ymin,xmax,ymax), TAttText(11,0,1,62,0.040)
494 {
495 
496  fWmin = wmin;
497  fWmax = wmax;
498  fNdiv = ndiv;
499  fNModLabs = 0;
500  fModLabs = 0;
501  fGridLength = gridlength;
502  fLabelOffset = 0.005;
503  fLabelSize = 0.040;
504  fLabelFont = 62;
505  fLabelColor = 1;
506  fTickSize = 0.030;
507  fTitleOffset = 1;
508  fTitleSize = fLabelSize;
509  fChopt = chopt;
510  fName = "";
511  fTitle = "";
512  fTimeFormat = "";
513  fFunctionName= "";
514  fFunction = 0;
515  fAxis = 0;
516 }
517 
518 ////////////////////////////////////////////////////////////////////////////////
519 /// Constructor with a `TF1` to map axis values.
520 
522  const char *funcname, Int_t ndiv, Option_t *chopt,
523  Double_t gridlength)
524  : TLine(xmin,ymin,xmax,ymax), TAttText(11,0,1,62,0.040)
525 {
526 
527  fFunction = (TF1*)gROOT->GetFunction(funcname);
528  if (!fFunction) {
529  Error("TGaxis", "calling constructor with an unknown function: %s", funcname);
530  fWmin = 0;
531  fWmax = 1;
532  } else {
533  fWmin = fFunction->GetXmin();
534  fWmax = fFunction->GetXmax();
535  }
536  fFunctionName= funcname;
537  fNdiv = ndiv;
538  fNModLabs = 0;
539  fModLabs = 0;
540  fGridLength = gridlength;
541  fLabelOffset = 0.005;
542  fLabelSize = 0.040;
543  fLabelFont = 62;
544  fLabelColor = 1;
545  fTickSize = 0.030;
546  fTitleOffset = 1;
547  fTitleSize = fLabelSize;
548  fChopt = chopt;
549  fName = "";
550  fTitle = "";
551  fTimeFormat = "";
552  fAxis = 0;
553 }
554 
555 ////////////////////////////////////////////////////////////////////////////////
556 /// Copy constructor.
557 
558 TGaxis::TGaxis(const TGaxis& ax) :
559  TLine(ax),
560  TAttText(ax),
561  fWmin(ax.fWmin),
562  fWmax(ax.fWmax),
563  fGridLength(ax.fGridLength),
564  fTickSize(ax.fTickSize),
565  fLabelOffset(ax.fLabelOffset),
566  fLabelSize(ax.fLabelSize),
567  fTitleOffset(ax.fTitleOffset),
568  fTitleSize(ax.fTitleSize),
569  fNdiv(ax.fNdiv),
570  fLabelColor(ax.fLabelColor),
571  fLabelFont(ax.fLabelFont),
572  fNModLabs(ax.fNModLabs),
573  fChopt(ax.fChopt),
574  fName(ax.fName),
575  fTitle(ax.fTitle),
576  fTimeFormat(ax.fTimeFormat),
577  fFunctionName(ax.fFunctionName),
578  fFunction(ax.fFunction),
579  fAxis(ax.fAxis),
580  fModLabs(ax.fModLabs)
581 {
582 }
583 
584 ////////////////////////////////////////////////////////////////////////////////
585 /// Assignment operator.
586 
587 TGaxis& TGaxis::operator=(const TGaxis& ax)
588 {
589 
590  if(this!=&ax) {
592  TAttText::operator=(ax);
593  fWmin=ax.fWmin;
594  fWmax=ax.fWmax;
596  fTickSize=ax.fTickSize;
601  fNdiv=ax.fNdiv;
602  fModLabs=ax.fModLabs;
605  fChopt=ax.fChopt;
606  fName=ax.fName;
607  fTitle=ax.fTitle;
610  fFunction=ax.fFunction;
611  fAxis=ax.fAxis;
612  fNModLabs=ax.fNModLabs;
613  }
614  return *this;
615 }
616 
617 ////////////////////////////////////////////////////////////////////////////////
618 /// TGaxis default destructor.
619 
621 {
622 }
623 
624 ////////////////////////////////////////////////////////////////////////////////
625 /// If center = kTRUE axis labels are centered in the center of the bin.
626 /// The default is to center on the primary tick marks.
627 /// This option does not make sense if there are more bins than tick marks.
628 
629 void TGaxis::CenterLabels(Bool_t center)
630 {
631 
632  if (center) SetBit(TAxis::kCenterLabels);
634 }
635 
636 ////////////////////////////////////////////////////////////////////////////////
637 /// If center = kTRUE axis title will be centered. The default is right adjusted.
638 
639 void TGaxis::CenterTitle(Bool_t center)
640 {
641 
642  if (center) SetBit(TAxis::kCenterTitle);
644 }
645 
646 ////////////////////////////////////////////////////////////////////////////////
647 /// Draw this axis with new attributes.
648 
650  Double_t wmin, Double_t wmax, Int_t ndiv, Option_t *chopt,
651  Double_t gridlength)
652 {
654  TGaxis *newaxis = new TGaxis(xmin,ymin,xmax,ymax,wmin,wmax,ndiv,chopt,gridlength);
655  newaxis->SetLineColor(fLineColor);
656  newaxis->SetLineWidth(fLineWidth);
657  newaxis->SetLineStyle(fLineStyle);
658  newaxis->SetTextAlign(fTextAlign);
659  newaxis->SetTextAngle(fTextAngle);
660  newaxis->SetTextColor(fTextColor);
661  newaxis->SetTextFont(fTextFont);
662  newaxis->SetTextSize(fTextSize);
663  newaxis->SetTitleSize(fTitleSize);
664  newaxis->SetTitleOffset(fTitleOffset);
665  newaxis->SetLabelFont(fLabelFont);
666  newaxis->SetLabelColor(fLabelColor);
667  newaxis->SetLabelSize(fLabelSize);
668  newaxis->SetLabelOffset(fLabelOffset);
669  newaxis->SetTickSize(fTickSize);
670  newaxis->SetBit(kCanDelete);
671  newaxis->SetTitle(GetTitle());
673  newaxis->AppendPad();
674 }
675 
676 ////////////////////////////////////////////////////////////////////////////////
677 /// Static function returning fgMaxDigits (See SetMaxDigits).
678 
680 {
681 
682  return fgMaxDigits;
683 }
684 
685 ////////////////////////////////////////////////////////////////////////////////
686 /// Internal method to import TAxis attributes to this TGaxis.
687 
689 {
690 
691  fAxis = axis;
693  SetTextColor(axis->GetTitleColor());
694  SetTextFont(axis->GetTitleFont());
695  SetLabelColor(axis->GetLabelColor());
696  SetLabelFont(axis->GetLabelFont());
697  SetLabelSize(axis->GetLabelSize());
699  SetTickSize(axis->GetTickLength());
700  SetTitle(axis->GetTitle());
702  SetTitleSize(axis->GetTitleSize());
710  if (axis->GetDecimals()) SetBit(TAxis::kDecimals); //the bit is in TAxis::fAxis2
711  SetTimeFormat(axis->GetTimeFormat());
712 }
713 
714 ////////////////////////////////////////////////////////////////////////////////
715 /// Draw this axis with its current attributes.
716 
717 void TGaxis::Paint(Option_t *)
718 {
719 
720  Double_t wmin = fWmin;
721  Double_t wmax = fWmax;
722  Int_t ndiv = fNdiv;
723 
724  // following code required to support toggle of lin/log scales
725  Double_t x1 = gPad->XtoPad(fX1);
726  Double_t y1 = gPad->YtoPad(fY1);
727  Double_t x2 = gPad->XtoPad(fX2);
728  Double_t y2 = gPad->YtoPad(fY2);
729 
730  PaintAxis(x1,y1,x2,y2,wmin,wmax,ndiv,fChopt.Data(),fGridLength);
731 }
732 
733 ////////////////////////////////////////////////////////////////////////////////
734 /// Control function to draw an axis.
735 /// Original authors: O.Couet C.E.Vandoni N.Cremel-Somon.
736 /// Modified and converted to C++ class by Rene Brun.
737 
738 void TGaxis::PaintAxis(Double_t xmin, Double_t ymin, Double_t xmax, Double_t ymax,
739  Double_t &wmin, Double_t &wmax, Int_t &ndiv, Option_t *chopt,
740  Double_t gridlength, Bool_t drawGridOnly)
741 {
743  const char *where = "PaintAxis";
744 
745  Double_t alfa, beta, ratio1, ratio2, grid_side;
746  Double_t axis_lengthN = 0;
747  Double_t axis_length0 = 0;
748  Double_t axis_length1 = 0;
749  Double_t axis_length;
750  Double_t atick[3];
751  Double_t tick_side;
752  Double_t charheight;
753  Double_t phil, phi, sinphi, cosphi, asinphi, acosphi;
754  Double_t binLow = 0., binLow2 = 0., binLow3 = 0.;
755  Double_t binHigh = 0., binHigh2 = 0., binHigh3 = 0.;
756  Double_t binWidth = 0., binWidth2 = 0., binWidth3 = 0.;
757  Double_t xpl1, xpl2, ypl1, ypl2;
758  Double_t xtick = 0;
759  Double_t xtick0, xtick1, dxtick=0;
760  Double_t ytick, ytick0, ytick1;
761  Double_t wlabel, dwlabel;
762  Double_t xfactor, yfactor;
763  Double_t xlabel, ylabel, dxlabel;
764  Double_t xone, xtwo;
765  Double_t rlab;
766  Double_t x0, x1, y0, y1, xx0, xx1, yy0, yy1;
767  xx0 = xx1 = yy0 = yy1 = 0;
768  Double_t xxmin, xxmax, yymin, yymax;
769  xxmin = xxmax = yymin = yymax = 0;
770  Double_t xlside, xmside;
771  Double_t ww, af, rne;
772  Double_t xx, yy;
773  Double_t xmnlog, x00, x11, h2, h2sav, axmul, y;
774  Float_t chupxvsav, chupyvsav;
775  Double_t rtxw, rtyw;
776  Int_t nlabels, nticks, nticks0, nticks1;
777  Int_t i, j, k, l, decade, ltick;
778  Int_t mside, lside;
779  Int_t nexe = 0;
780  Int_t lnlen = 0;
781  Int_t iexe, if1, if2, na, nf, ih1, ih2, nbinin, nch, kmod;
782  Int_t optionLog,optionBlank,optionVert,optionPlus,optionMinus,optionUnlab,optionPara;
783  Int_t optionDown,optionRight,optionLeft,optionCent,optionEqual,optionDecimals=0,optionDot;
784  Int_t optionY,optionText,optionGrid,optionSize,optionNoopt,optionInt,optionM,optionUp,optionX;
785  Int_t optionTime;
786  Int_t first=0,last=0,labelnumber;
787  Int_t xalign, yalign;
788  Int_t nn1, nn2, nn3, n1a, n2a, n3a, nb2, nb3;
789  Int_t nbins=10, n1aold, nn1old;
790  n1aold = nn1old = 0;
791  Int_t ndyn;
792  Int_t nhilab = 0;
793  Int_t idn;
794  Bool_t flexe = 0;
795  Bool_t flexpo,flexne;
796  char *label;
797  char *chtemp;
798  char *coded;
799  char chlabel[256];
800  char kchtemp[256];
801  char chcoded[8];
802  TLine *linegrid;
803  TString timeformat;
804  TString typolabel;
805  time_t timelabel;
806  Double_t timed, wTimeIni;
807  struct tm* utctis;
808  Double_t rangeOffset = 0;
809 
810  Double_t epsilon = 1e-5;
811  const Double_t kPI = TMath::Pi();
812 
813  Double_t rwmi = wmin;
814  Double_t rwma = wmax;
815  chtemp = &kchtemp[0];
816  label = &chlabel[0];
817  linegrid = 0;
818 
819  fFunction = (TF1*)gROOT->GetFunction(fFunctionName.Data());
820 
821  Bool_t noExponent = TestBit(TAxis::kNoExponent);
822 
823 // If moreLogLabels = kTRUE more Log Intermediate Labels are drawn.
824  Bool_t moreLogLabels = TestBit(TAxis::kMoreLogLabels);
825 
826 // the following parameters correspond to the pad range in NDC
827 // and the user's coordinates in the pad
828 
829  Double_t padh = gPad->GetWh()*gPad->GetAbsHNDC();
830  Double_t rwxmin = gPad->GetX1();
831  Double_t rwxmax = gPad->GetX2();
832  Double_t rwymin = gPad->GetY1();
833  Double_t rwymax = gPad->GetY2();
834 
835  if(strchr(chopt,'G')) optionLog = 1; else optionLog = 0;
836  if(strchr(chopt,'B')) optionBlank= 1; else optionBlank= 0;
837  if(strchr(chopt,'V')) optionVert = 1; else optionVert = 0;
838  if(strchr(chopt,'+')) optionPlus = 1; else optionPlus = 0;
839  if(strchr(chopt,'-')) optionMinus= 1; else optionMinus= 0;
840  if(strchr(chopt,'U')) optionUnlab= 1; else optionUnlab= 0;
841  if(strchr(chopt,'P')) optionPara = 1; else optionPara = 0;
842  if(strchr(chopt,'O')) optionDown = 1; else optionDown = 0;
843  if(strchr(chopt,'R')) optionRight= 1; else optionRight= 0;
844  if(strchr(chopt,'L')) optionLeft = 1; else optionLeft = 0;
845  if(strchr(chopt,'C')) optionCent = 1; else optionCent = 0;
846  if(strchr(chopt,'=')) optionEqual= 1; else optionEqual= 0;
847  if(strchr(chopt,'Y')) optionY = 1; else optionY = 0;
848  if(strchr(chopt,'T')) optionText = 1; else optionText = 0;
849  if(strchr(chopt,'W')) optionGrid = 1; else optionGrid = 0;
850  if(strchr(chopt,'S')) optionSize = 1; else optionSize = 0;
851  if(strchr(chopt,'N')) optionNoopt= 1; else optionNoopt= 0;
852  if(strchr(chopt,'I')) optionInt = 1; else optionInt = 0;
853  if(strchr(chopt,'M')) optionM = 1; else optionM = 0;
854  if(strchr(chopt,'0')) optionUp = 1; else optionUp = 0;
855  if(strchr(chopt,'X')) optionX = 1; else optionX = 0;
856  if(strchr(chopt,'t')) optionTime = 1; else optionTime = 0;
857  if(strchr(chopt,'.')) optionDot = 1; else optionDot = 0;
858  if (TestBit(TAxis::kTickPlus)) optionPlus = 2;
859  if (TestBit(TAxis::kTickMinus)) optionMinus = 2;
860  if (TestBit(TAxis::kCenterLabels)) optionM = 1;
861  if (TestBit(TAxis::kDecimals)) optionDecimals = 1;
862  if (!gStyle->GetStripDecimals()) optionDecimals = 1;
863  if (fAxis) {
864  if (fAxis->GetLabels()) {
865  optionM = 1;
866  optionText = 1;
867  ndiv = fAxis->GetLast()-fAxis->GetFirst()+1;
868  }
869  TList *ml = fAxis->GetModifiedLabels();
870  if (ml) {
871  fModLabs = ml;
873  } else {
874  fModLabs = 0;
875  fNModLabs = 0;
876  }
877  }
878  if (ndiv < 0) {
879  Error(where, "Invalid number of divisions: %d",ndiv);
880  return;
881  }
882 
883 // Set the grid length
884 
885  if (optionGrid) {
886  if (gridlength == 0) gridlength = 0.8;
887  linegrid = new TLine();
888  linegrid->SetLineColor(gStyle->GetGridColor());
889  if (linegrid->GetLineColor() == 0) linegrid->SetLineColor(GetLineColor());
890  linegrid->SetLineStyle(gStyle->GetGridStyle());
891  linegrid->SetLineWidth(gStyle->GetGridWidth());
892  }
893 
894 // No labels if the axis label offset is big.
895 // In that case the labels are not visible anyway.
896 
897  if (GetLabelOffset() > 1.1 ) optionUnlab = 1;
898 
899 // Determine time format
900 
901  Int_t idF = fTimeFormat.Index("%F");
902  if (idF>=0) {
903  timeformat = fTimeFormat(0,idF);
904  } else {
905  timeformat = fTimeFormat;
906  }
907 
908  //GMT option
909  if (fTimeFormat.Index("GMT")>=0) optionTime =2;
910 
911  // Determine the time offset and correct for time offset not being integer.
912  Double_t timeoffset =0;
913  if (optionTime) {
914  if (idF>=0) {
915  Int_t lnF = fTimeFormat.Length();
916  TString stringtimeoffset = fTimeFormat(idF+2,lnF);
917  Int_t year, mm, dd, hh, mi, ss;
918  if (sscanf(stringtimeoffset.Data(), "%d-%d-%d %d:%d:%d", &year, &mm, &dd, &hh, &mi, &ss) == 6) {
919  //Get time offset in seconds since EPOCH:
920  struct tm tp;
921  tp.tm_year = year-1900;
922  tp.tm_mon = mm-1;
923  tp.tm_mday = dd;
924  tp.tm_hour = hh;
925  tp.tm_min = mi;
926  tp.tm_sec = ss;
927  tp.tm_isdst = 0; //no DST for UTC (and forced to 0 in MktimeFromUTC function)
928  timeoffset = TTimeStamp::MktimeFromUTC(&tp);
929 
930  // Add the time offset's decimal part if it is there
931  Int_t ids = stringtimeoffset.Index("s");
932  if (ids >= 0) {
933  Float_t dp;
934  Int_t lns = stringtimeoffset.Length();
935  TString sdp = stringtimeoffset(ids+1,lns);
936  sscanf(sdp.Data(),"%g",&dp);
937  timeoffset += dp;
938  }
939  } else {
940  Error(where, "Time offset has not the right format");
941  }
942  } else {
943  timeoffset = gStyle->GetTimeOffset();
944  }
945  wmin += timeoffset - (int)(timeoffset);
946  wmax += timeoffset - (int)(timeoffset);
947 
948  // correct for time offset at a good limit (min, hour, day, month, year)
949  struct tm* tp0;
950  time_t timetp = (time_t)((Long_t)(timeoffset));
951  Double_t range = wmax - wmin;
952  Long_t rangeBase = 60;
953  if (range>60) rangeBase = 60*20; // minutes
954  if (range>3600) rangeBase = 3600*20; // hours
955  if (range>86400) rangeBase = 86400*20; // days
956  if (range>2419200) rangeBase = 31556736; // months (average # days)
957  rangeOffset = (Double_t) ((Long_t)(timeoffset)%rangeBase);
958  if (range>31536000) {
959  tp0 = gmtime(&timetp);
960  tp0->tm_mon = 0;
961  tp0->tm_mday = 1;
962  tp0->tm_hour = 0;
963  tp0->tm_min = 0;
964  tp0->tm_sec = 0;
965  tp0->tm_isdst = 1; // daylight saving time is on.
966  rangeBase = (timetp-mktime(tp0)); // years
967  rangeOffset = (Double_t) (rangeBase);
968  }
969  wmax += rangeOffset;
970  wmin += rangeOffset;
971  }
972 
973 // Determine number of divisions 1, 2 and 3
974  n1a = ndiv%100;
975  n2a = (ndiv%10000 - n1a)/100;
976  n3a = ndiv/10000;
977  nn3 = TMath::Max(n3a,1);
978  nn2 = TMath::Max(n2a,1)*nn3;
979  nn1 = TMath::Max(n1a,1)*nn2+1;
980  nticks= nn1;
981 
982 // Axis bining optimisation is ignored if:
983 // - the first and the last label are equal
984 // - the number of divisions is 0
985 // - less than 1 primary division is requested
986 // - logarithmic scale is requested
987 
988  if (wmin == wmax || ndiv == 0 || n1a <= 1 || optionLog) {
989  optionNoopt = 1;
990  optionInt = 0;
991  }
992 
993 // Axis bining optimisation
994  if ( (wmax-wmin) < 1 && optionInt) {
995  Error(where, "option I not available");
996  optionInt = 0;
997  }
998  if (!optionNoopt || optionInt ) {
999 
1000 // Primary divisions optimisation
1001 // When integer labelling is required, Optimize is invoked first
1002 // and only if the result is not an integer labelling, AdjustBinSize is invoked.
1003 
1004  THLimitsFinder::Optimize(wmin,wmax,n1a,binLow,binHigh,nbins,binWidth,fChopt.Data());
1005  if (optionInt) {
1006  if (binLow != Double_t(int(binLow)) || binWidth != Double_t(int(binWidth))) {
1007  AdjustBinSize(wmin,wmax,n1a,binLow,binHigh,nbins,binWidth);
1008  }
1009  }
1010  if ((wmin-binLow) > epsilon) { binLow += binWidth; nbins--; }
1011  if ((binHigh-wmax) > epsilon) { binHigh -= binWidth; nbins--; }
1012  if (xmax == xmin) {
1013  rtyw = (ymax-ymin)/(wmax-wmin);
1014  xxmin = xmin;
1015  xxmax = xmax;
1016  yymin = rtyw*(binLow-wmin) + ymin;
1017  yymax = rtyw*(binHigh-wmin) + ymin;
1018  }
1019  else {
1020  rtxw = (xmax-xmin)/(wmax-wmin);
1021  xxmin = rtxw*(binLow-wmin) + xmin;
1022  xxmax = rtxw*(binHigh-wmin) + xmin;
1023  if (ymax == ymin) {
1024  yymin = ymin;
1025  yymax = ymax;
1026  }
1027  else {
1028  alfa = (ymax-ymin)/(xmax-xmin);
1029  beta = (ymin*xmax-ymax*xmin)/(xmax-xmin);
1030  yymin = alfa*xxmin + beta;
1031  yymax = alfa*xxmax + beta;
1032  }
1033  }
1034  if (fFunction) {
1035  yymin = ymin;
1036  yymax = ymax;
1037  xxmin = xmin;
1038  xxmax = xmax;
1039  } else {
1040  wmin = binLow;
1041  wmax = binHigh;
1042  }
1043 
1044 // Secondary divisions optimisation
1045  nb2 = n2a;
1046  if (!optionNoopt && n2a > 1 && binWidth > 0) {
1047  THLimitsFinder::Optimize(wmin,wmin+binWidth,n2a,binLow2,binHigh2,nb2,binWidth2,fChopt.Data());
1048  }
1049 
1050 // Tertiary divisions optimisation
1051  nb3 = n3a;
1052  if (!optionNoopt && n3a > 1 && binWidth2 > 0) {
1053  THLimitsFinder::Optimize(binLow2,binLow2+binWidth2,n3a,binLow3,binHigh3,nb3,binWidth3,fChopt.Data());
1054  }
1055  n1aold = n1a;
1056  nn1old = nn1;
1057  n1a = nbins;
1058  nn3 = TMath::Max(nb3,1);
1059  nn2 = TMath::Max(nb2,1)*nn3;
1060  nn1 = TMath::Max(n1a,1)*nn2+1;
1061  nticks = nn1;
1062  }
1063 
1064 // Coordinates are normalized
1065 
1066  ratio1 = 1/(rwxmax-rwxmin);
1067  ratio2 = 1/(rwymax-rwymin);
1068  x0 = ratio1*(xmin-rwxmin);
1069  x1 = ratio1*(xmax-rwxmin);
1070  y0 = ratio2*(ymin-rwymin);
1071  y1 = ratio2*(ymax-rwymin);
1072  if (!optionNoopt || optionInt ) {
1073  xx0 = ratio1*(xxmin-rwxmin);
1074  xx1 = ratio1*(xxmax-rwxmin);
1075  yy0 = ratio2*(yymin-rwymin);
1076  yy1 = ratio2*(yymax-rwymin);
1077  }
1078 
1079  if ((x0 == x1) && (y0 == y1)) {
1080  Error(where, "length of axis is 0");
1081  return;
1082  }
1083 
1084 // Return wmin, wmax and the number of primary divisions
1085  if (optionX) {
1086  ndiv = n1a;
1087  return;
1088  }
1089 
1090  Int_t maxDigits = fgMaxDigits;
1091 
1092  TLatex *textaxis = new TLatex();
1093  SetLineStyle(1); // axis line style
1094  textaxis->SetTextColor(GetTextColor());
1095  textaxis->SetTextFont(GetTextFont());
1096 
1097  if (!gPad->IsBatch()) {
1098  gVirtualX->GetCharacterUp(chupxvsav, chupyvsav);
1099  gVirtualX->SetClipOFF(gPad->GetCanvasID());
1100  }
1101 
1102 // Compute length of axis
1103  axis_length = TMath::Sqrt((x1-x0)*(x1-x0)+(y1-y0)*(y1-y0));
1104  if (axis_length == 0) {
1105  Error(where, "length of axis is 0");
1106  goto L210;
1107  }
1108  if (!optionNoopt || optionInt) {
1109  axis_lengthN = TMath::Sqrt((xx1-xx0)*(xx1-xx0)+(yy1-yy0)*(yy1-yy0));
1110  axis_length0 = TMath::Sqrt((xx0-x0)*(xx0-x0)+(yy0-y0)*(yy0-y0));
1111  axis_length1 = TMath::Sqrt((x1-xx1)*(x1-xx1)+(y1-yy1)*(y1-yy1));
1112  if (axis_lengthN < epsilon) {
1113  optionNoopt = 1;
1114  optionInt = 0;
1115  wmin = rwmi;
1116  wmax = rwma;
1117  n1a = n1aold;
1118  nn1 = nn1old;
1119  nticks = nn1;
1120  if (optionTime) {
1121  wmin += timeoffset - (int)(timeoffset) + rangeOffset;
1122  wmax += timeoffset - (int)(timeoffset) + rangeOffset;
1123  }
1124  }
1125  }
1126 
1127  if (x0 == x1) {
1128  phi = 0.5*kPI;
1129  phil = phi;
1130  } else {
1131  phi = TMath::ATan2((y1-y0),(x1-x0));
1132  Int_t px0 = gPad->UtoPixel(x0);
1133  Int_t py0 = gPad->VtoPixel(y0);
1134  Int_t px1 = gPad->UtoPixel(x1);
1135  Int_t py1 = gPad->VtoPixel(y1);
1136  if (x0 < x1) phil = TMath::ATan2(Double_t(py0-py1), Double_t(px1-px0));
1137  else phil = TMath::ATan2(Double_t(py1-py0), Double_t(px0-px1));
1138  }
1139  cosphi = TMath::Cos(phi);
1140  sinphi = TMath::Sin(phi);
1141  acosphi = TMath::Abs(cosphi);
1142  asinphi = TMath::Abs(sinphi);
1143  if (acosphi <= epsilon) { acosphi = 0; cosphi = 0; }
1144  if (asinphi <= epsilon) { asinphi = 0; sinphi = 0; }
1145 
1146 // mside positive, tick marks on positive side
1147 // mside negative, tick marks on negative side
1148 // mside zero, tick marks on both sides
1149 // Default is positive except for vertical axis
1150 
1151  mside=1;
1152  if (x0 == x1 && y1 > y0) mside = -1;
1153  if (optionPlus) mside = 1;
1154  if (optionMinus) mside = -1;
1155  if (optionPlus && optionMinus) mside = 0;
1156  xmside = mside;
1157  lside = -mside;
1158  if (optionEqual) lside = mside;
1159  if (optionPlus && optionMinus) {
1160  lside = -1;
1161  if (optionEqual) lside=1;
1162  }
1163  xlside = lside;
1164 
1165 // Tick marks size
1166  if(xmside >= 0) tick_side = 1;
1167  else tick_side = -1;
1168  if (optionSize) atick[0] = tick_side*axis_length*fTickSize;
1169  else atick[0] = tick_side*axis_length*0.03;
1170 
1171  atick[1] = 0.5*atick[0];
1172  atick[2] = 0.5*atick[1];
1173 
1174 // Set the side of the grid
1175  if ((x0 == x1) && (y1 > y0)) grid_side =-1;
1176  else grid_side = 1;
1177 
1178 // Compute Values if Function is given
1179  if(fFunction) {
1180  rwmi = fFunction->Eval(wmin);
1181  rwma = fFunction->Eval(wmax);
1182  if(rwmi > rwma) {
1183  Double_t t = rwma;
1184  rwma = rwmi;
1185  rwmi = t;
1186  }
1187  }
1188 
1189 // Draw the axis if needed...
1190  if (!optionBlank) {
1191  xpl1 = x0;
1192  xpl2 = x1;
1193  ypl1 = y0;
1194  ypl2 = y1;
1195  PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1196  }
1197 
1198 // Draw axis title if it exists
1199  if (!drawGridOnly && strlen(GetTitle())) {
1200  textaxis->SetTextSize (GetTitleSize());
1201  charheight = GetTitleSize();
1202  if ((GetTextFont() % 10) > 2) {
1203  charheight = charheight/gPad->GetWh();
1204  }
1205  Double_t toffset = GetTitleOffset();
1206 //////if (toffset < 0.1) toffset = 1; // Negative offset should be allowed
1207  if (x1 == x0) ylabel = xlside*1.6*charheight*toffset;
1208  else ylabel = xlside*1.3*charheight*toffset;
1209  if (y1 == y0) ylabel = xlside*1.6*charheight*toffset;
1210  Double_t axispos;
1211  if (TestBit(TAxis::kCenterTitle)) axispos = 0.5*axis_length;
1212  else axispos = axis_length;
1214  if (x1 >= x0) {
1215  if (TestBit(TAxis::kCenterTitle)) textaxis->SetTextAlign(22);
1216  else textaxis->SetTextAlign(12);
1217  Rotate(axispos,ylabel,cosphi,sinphi,x0,y0,xpl1,ypl1);
1218  } else {
1219  if (TestBit(TAxis::kCenterTitle)) textaxis->SetTextAlign(22);
1220  else textaxis->SetTextAlign(32);
1221  Rotate(axispos,ylabel,cosphi,sinphi,x0,y0,xpl1,ypl1);
1222  }
1223  textaxis->PaintLatex(gPad->GetX1() + xpl1*(gPad->GetX2() - gPad->GetX1()),
1224  gPad->GetY1() + ypl1*(gPad->GetY2() - gPad->GetY1()),
1225  phil=(kPI+phil)*180/kPI,
1226  GetTitleSize(),
1227  GetTitle());
1228  } else {
1229  if (x1 >= x0) {
1230  if (TestBit(TAxis::kCenterTitle)) textaxis->SetTextAlign(22);
1231  else textaxis->SetTextAlign(32);
1232  Rotate(axispos,ylabel,cosphi,sinphi,x0,y0,xpl1,ypl1);
1233  } else {
1234  if (TestBit(TAxis::kCenterTitle)) textaxis->SetTextAlign(22);
1235  else textaxis->SetTextAlign(12);
1236  Rotate(axispos,ylabel,cosphi,sinphi,x0,y0,xpl1,ypl1);
1237  }
1238  textaxis->PaintLatex(gPad->GetX1() + xpl1*(gPad->GetX2() - gPad->GetX1()),
1239  gPad->GetY1() + ypl1*(gPad->GetY2() - gPad->GetY1()),
1240  phil*180/kPI,
1241  GetTitleSize(),
1242  GetTitle());
1243  }
1244  }
1245 
1246 // No bining
1247 
1248  if (ndiv == 0)goto L210;
1249  if (wmin == wmax) {
1250  Error(where, "wmin (%f) == wmax (%f)", wmin, wmax);
1251  goto L210;
1252  }
1253 
1254 // Labels preparation:
1255 // Get character height
1256 // Compute the labels orientation in case of overlaps
1257 // (with alphanumeric labels for horizontal axis).
1258 
1259  charheight = GetLabelSize();
1260  if (optionText && GetLabelFont()%10 != 3) charheight *= 0.66666;
1261  textaxis->SetTextFont(GetLabelFont());
1262  if ((GetLabelFont()%10 < 2) && optionLog) // force TLatex mode in PaintLatex
1263  textaxis->SetTextFont((Int_t)(GetLabelFont()/10)*10+2);
1264  textaxis->SetTextColor(GetLabelColor());
1265  textaxis->SetTextSize (charheight);
1266  textaxis->SetTextAngle(GetTextAngle());
1267  if (GetLabelFont()%10 > 2) {
1268  charheight /= padh;
1269  }
1270  if (!optionUp && !optionDown && !optionY && !optionUnlab) {
1271  if (!drawGridOnly && optionText && ((ymin == ymax) || (xmin == xmax))) {
1272  textaxis->SetTextAlign(32);
1273  optionText = 2;
1274  Int_t nl = fAxis->GetLast()-fAxis->GetFirst()+1;
1275  Double_t angle = 0;
1276  for (i=fAxis->GetFirst(); i<=fAxis->GetLast(); i++) {
1277  textaxis->SetText(0,0,fAxis->GetBinLabel(i));
1278  if (textaxis->GetXsize() < (xmax-xmin)/nl) continue;
1279  angle = -20;
1280  break;
1281  }
1282  for (i=fAxis->GetFirst(); i<=fAxis->GetLast(); i++) {
1283  if ((!strcmp(fAxis->GetName(),"xaxis") && !gPad->TestBit(kHori))
1284  ||(!strcmp(fAxis->GetName(),"yaxis") && gPad->TestBit(kHori))) {
1285  if (nl > 50) angle = 90;
1286  if (fAxis->TestBit(TAxis::kLabelsHori)) angle = 0;
1287  if (fAxis->TestBit(TAxis::kLabelsVert)) angle = 90;
1288  if (fAxis->TestBit(TAxis::kLabelsUp)) angle = 20;
1289  if (fAxis->TestBit(TAxis::kLabelsDown)) angle =-20;
1290  if (angle == 0) textaxis->SetTextAlign(23);
1291  if (angle == -20) textaxis->SetTextAlign(12);
1292  Double_t s = -3;
1293  if (ymin == gPad->GetUymax()) {
1294  if (angle == 0) textaxis->SetTextAlign(21);
1295  s = 3;
1296  }
1297  textaxis->PaintLatex(fAxis->GetBinCenter(i),
1298  ymin + s*fAxis->GetLabelOffset()*(gPad->GetUymax()-gPad->GetUymin()),
1299  angle,
1300  textaxis->GetTextSize(),
1301  fAxis->GetBinLabel(i));
1302  } else if ((!strcmp(fAxis->GetName(),"yaxis") && !gPad->TestBit(kHori))
1303  || (!strcmp(fAxis->GetName(),"xaxis") && gPad->TestBit(kHori))) {
1304  Double_t s = -3;
1305  if (xmin == gPad->GetUxmax()) {
1306  textaxis->SetTextAlign(12);
1307  s = 3;
1308  }
1309  textaxis->PaintLatex(xmin + s*fAxis->GetLabelOffset()*(gPad->GetUxmax()-gPad->GetUxmin()),
1310  fAxis->GetBinCenter(i),
1311  0,
1312  textaxis->GetTextSize(),
1313  fAxis->GetBinLabel(i));
1314  } else {
1315  textaxis->PaintLatex(xmin - 3*fAxis->GetLabelOffset()*(gPad->GetUxmax()-gPad->GetUxmin()),
1316  ymin +(i-0.5)*(ymax-ymin)/nl,
1317  0,
1318  textaxis->GetTextSize(),
1319  fAxis->GetBinLabel(i));
1320  }
1321  }
1322  }
1323  }
1324 
1325 // Now determine orientation of labels on axis
1326  if (!gPad->IsBatch()) {
1327  if (cosphi > 0) gVirtualX->SetCharacterUp(-sinphi,cosphi);
1328  else gVirtualX->SetCharacterUp(sinphi,-cosphi);
1329  if (x0 == x1) gVirtualX->SetCharacterUp(0,1);
1330  if (optionVert) gVirtualX->SetCharacterUp(0,1);
1331  if (optionPara) gVirtualX->SetCharacterUp(-sinphi,cosphi);
1332  if (optionDown) gVirtualX->SetCharacterUp(cosphi,sinphi);
1333  }
1334 
1335 // Now determine text alignment
1336  xalign = 2;
1337  yalign = 1;
1338  if (x0 == x1) xalign = 3;
1339  if (y0 != y1) yalign = 2;
1340  if (optionCent) xalign = 2;
1341  if (optionRight) xalign = 3;
1342  if (optionLeft) xalign = 1;
1343  if (TMath::Abs(cosphi) > 0.9) {
1344  xalign = 2;
1345  } else {
1346  if (cosphi*sinphi > 0) xalign = 1;
1347  if (cosphi*sinphi < 0) xalign = 3;
1348  }
1349  textaxis->SetTextAlign(10*xalign+yalign);
1350 
1351 // Position of labels in Y
1352  if (x0 == x1) {
1353  if (optionPlus && !optionMinus) {
1354  if (optionEqual) ylabel = fLabelOffset/2 + atick[0];
1355  else ylabel = -fLabelOffset;
1356  } else {
1357  ylabel = fLabelOffset;
1358  if (lside < 0) ylabel += atick[0];
1359  }
1360  } else if (y0 == y1) {
1361  if (optionMinus && !optionPlus) {
1362  if ((GetLabelFont() % 10) == 3 ) {
1363  ylabel = fLabelOffset+0.5*
1364  ((gPad->AbsPixeltoY(0)-gPad->AbsPixeltoY((Int_t)fLabelSize))/
1365  (gPad->GetY2() - gPad->GetY1()));
1366  } else {
1367  ylabel = fLabelOffset+0.5*fLabelSize;
1368  }
1369  ylabel += TMath::Abs(atick[0]);
1370  } else {
1371  ylabel = -fLabelOffset;
1372  if (mside <= 0) ylabel -= TMath::Abs(atick[0]);
1373  }
1374  if (optionLog) ylabel -= 0.5*charheight;
1375  } else {
1376  if (mside+lside >= 0) ylabel = fLabelOffset;
1377  else ylabel = -fLabelOffset;
1378  }
1379  if (optionText) ylabel /= 2;
1380 
1381 // Draw the linear tick marks if needed...
1382  if (!optionLog) {
1383  if (ndiv) {
1384  if (fFunction) {
1385  dxtick=(binHigh-binLow)/Double_t(nticks-1);
1386  } else {
1387  if (optionNoopt && !optionInt) dxtick=axis_length/Double_t(nticks-1);
1388  else dxtick=axis_lengthN/Double_t(nticks-1);
1389  }
1390  for (k=0;k<nticks; k++) {
1391  ltick = 2;
1392  if (k%nn3 == 0) ltick = 1;
1393  if (k%nn2 == 0) ltick = 0;
1394  if (fFunction) {
1395  Double_t xf = binLow+Double_t(k)*dxtick;
1396  Double_t zz = fFunction->Eval(xf)-rwmi;
1397  xtick = zz* axis_length / TMath::Abs(rwma-rwmi);
1398  } else {
1399  xtick = Double_t(k)*dxtick;
1400  }
1401  ytick = 0;
1402  if (!mside) ytick -= atick[ltick];
1403  if ( optionNoopt && !optionInt) {
1404  Rotate(xtick,ytick,cosphi,sinphi,x0,y0,xpl2,ypl2);
1405  Rotate(xtick,atick[ltick],cosphi,sinphi,x0,y0,xpl1,ypl1);
1406  }
1407  else {
1408  Rotate(xtick,ytick,cosphi,sinphi,xx0,yy0,xpl2,ypl2);
1409  Rotate(xtick,atick[ltick],cosphi,sinphi,xx0,yy0,xpl1,ypl1);
1410  }
1411  if (optionVert) {
1412  if ((x0 != x1) && (y0 != y1)) {
1413  if (mside) {
1414  xpl1 = xpl2;
1415  if (cosphi > 0) ypl1 = ypl2 + atick[ltick];
1416  else ypl1 = ypl2 - atick[ltick];
1417  }
1418  else {
1419  xpl1 = 0.5*(xpl1 + xpl2);
1420  xpl2 = xpl1;
1421  ypl1 = 0.5*(ypl1 + ypl2) + atick[ltick];
1422  ypl2 = 0.5*(ypl1 + ypl2) - atick[ltick];
1423  }
1424  }
1425  }
1426  if (!drawGridOnly) PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1427 
1428  if (optionGrid) {
1429  if (ltick == 0) {
1430  if (optionNoopt && !optionInt) {
1431  Rotate(xtick,0,cosphi,sinphi,x0,y0 ,xpl2,ypl2);
1432  Rotate(xtick,grid_side*gridlength ,cosphi,sinphi,x0,y0 ,xpl1,ypl1);
1433  }
1434  else {
1435  Rotate(xtick,0,cosphi ,sinphi,xx0,yy0 ,xpl2,ypl2);
1436  Rotate(xtick,grid_side*gridlength ,cosphi,sinphi,xx0,yy0 ,xpl1,ypl1);
1437  }
1438  linegrid->PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1439  }
1440  }
1441  }
1442  xtick0 = 0;
1443  xtick1 = xtick;
1444 
1445  if (fFunction) axis_length0 = binLow-wmin;
1446  if ((!optionNoopt || optionInt) && axis_length0) {
1447  nticks0 = Int_t(axis_length0/dxtick);
1448  if (nticks0 > 1000) nticks0 = 1000;
1449  for (k=0; k<=nticks0; k++) {
1450  ltick = 2;
1451  if (k%nn3 == 0) ltick = 1;
1452  if (k%nn2 == 0) ltick = 0;
1453  ytick0 = 0;
1454  if (!mside) ytick0 -= atick[ltick];
1455  if (fFunction) {
1456  xtick0 = (fFunction->Eval(binLow - Double_t(k)*dxtick)-rwmi)
1457  * axis_length / TMath::Abs(rwma-rwmi);
1458  }
1459  Rotate(xtick0,ytick0,cosphi,sinphi,xx0,yy0 ,xpl2,ypl2);
1460  Rotate(xtick0,atick[ltick],cosphi,sinphi,xx0,yy0 ,xpl1,ypl1);
1461  if (optionVert) {
1462  if ((x0 != x1) && (y0 != y1)) {
1463  if (mside) {
1464  xpl1 = xpl2;
1465  if (cosphi > 0) ypl1 = ypl2 + atick[ltick];
1466  else ypl1 = ypl2 - atick[ltick];
1467  }
1468  else {
1469  xpl1 = 0.5*(xpl1 + xpl2);
1470  xpl2 = xpl1;
1471  ypl1 = 0.5*(ypl1 + ypl2) + atick[ltick];
1472  ypl2 = 0.5*(ypl1 + ypl2) - atick[ltick];
1473  }
1474  }
1475  }
1476  if (!drawGridOnly) PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1477 
1478  if (optionGrid) {
1479  if (ltick == 0) {
1480  Rotate(xtick0,0,cosphi,sinphi,xx0,yy0,xpl2,ypl2);
1481  Rotate(xtick0,grid_side*gridlength ,cosphi,sinphi,xx0,yy0 ,xpl1,ypl1);
1482  linegrid->PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1483  }
1484  }
1485  xtick0 -= dxtick;
1486  }
1487  }
1488 
1489  if (fFunction) axis_length1 = wmax-binHigh;
1490  if ((!optionNoopt || optionInt) && axis_length1) {
1491  nticks1 = int(axis_length1/dxtick);
1492  if (nticks1 > 1000) nticks1 = 1000;
1493  for (k=0; k<=nticks1; k++) {
1494  ltick = 2;
1495  if (k%nn3 == 0) ltick = 1;
1496  if (k%nn2 == 0) ltick = 0;
1497  ytick1 = 0;
1498  if (!mside) ytick1 -= atick[ltick];
1499  if (fFunction) {
1500  xtick1 = (fFunction->Eval(binHigh + Double_t(k)*dxtick)-rwmi)
1501  * axis_length / TMath::Abs(rwma-rwmi);
1502  }
1503  Rotate(xtick1,ytick1,cosphi,sinphi,xx0,yy0 ,xpl2,ypl2);
1504  Rotate(xtick1,atick[ltick],cosphi,sinphi,xx0,yy0 ,xpl1,ypl1);
1505  if (optionVert) {
1506  if ((x0 != x1) && (y0 != y1)) {
1507  if (mside) {
1508  xpl1 = xpl2;
1509  if (cosphi > 0) ypl1 = ypl2 + atick[ltick];
1510  else ypl1 = ypl2 - atick[ltick];
1511  }
1512  else {
1513  xpl1 = 0.5*(xpl1 + xpl2);
1514  xpl2 = xpl1;
1515  ypl1 = 0.5*(ypl1 + ypl2) + atick[ltick];
1516  ypl2 = 0.5*(ypl1 + ypl2) - atick[ltick];
1517  }
1518  }
1519  }
1520  if (!drawGridOnly) PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1521  if (optionGrid) {
1522  if (ltick == 0) {
1523  Rotate(xtick1,0,cosphi,sinphi,xx0,yy0 ,xpl2,ypl2);
1524  Rotate(xtick1,grid_side*gridlength,cosphi,sinphi,xx0,yy0,xpl1,ypl1);
1525  linegrid->PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1526  }
1527  }
1528  xtick1 += dxtick;
1529  }
1530  }
1531  }
1532  }
1533 
1534 // Draw the numeric labels if needed...
1535  if (!drawGridOnly && !optionUnlab) {
1536  if (!optionLog) {
1537  if (n1a) {
1538 // Spacing of labels
1539  if ((wmin == wmax) || (ndiv == 0)) {
1540  Error(where, "wmin (%f) == wmax (%f), or ndiv == 0", wmin, wmax);
1541  goto L210;
1542  }
1543  wlabel = wmin;
1544  dwlabel = (wmax-wmin)/Double_t(n1a);
1545  if (optionNoopt && !optionInt) dxlabel = axis_length/Double_t(n1a);
1546  else dxlabel = axis_lengthN/Double_t(n1a);
1547 
1548  if (!optionText && !optionTime) {
1549 
1550 // We have to decide what format to generate
1551 // (for numeric labels only)
1552 // Test the magnitude, decide format
1553  flexe = kFALSE;
1554  nexe = 0;
1555  flexpo = kFALSE;
1556  flexne = kFALSE;
1557  ww = TMath::Max(TMath::Abs(wmin),TMath::Abs(wmax));
1558 
1559 // First case : (wmax-wmin)/n1a less than 0.001
1560 // (0.001 fgMaxDigits of 5 (fgMaxDigits) characters). Then we use x 10 n
1561 // format. If af >=0 x10 n cannot be used
1562  Double_t xmicros = 0.00099;
1563  if (maxDigits) xmicros = TMath::Power(10,-maxDigits);
1564  if (!noExponent && (TMath::Abs(wmax-wmin)/Double_t(n1a)) < xmicros) {
1565  af = TMath::Log10(ww) + epsilon;
1566  if (af < 0) {
1567  flexe = kTRUE;
1568  nexe = int(af);
1569  iexe = TMath::Abs(nexe);
1570  if (iexe%3 == 1) iexe += 2;
1571  else if(iexe%3 == 2) iexe += 1;
1572  if (nexe < 0) nexe = -iexe;
1573  else nexe = iexe;
1574  wlabel = wlabel*TMath::Power(10,iexe);
1575  dwlabel = dwlabel*TMath::Power(10,iexe);
1576  if1 = maxDigits;
1577  if2 = maxDigits-2;
1578  goto L110;
1579  }
1580  }
1581  if (ww >= 1) af = TMath::Log10(ww);
1582  else af = TMath::Log10(ww*0.0001);
1583  af += epsilon;
1584  nf = Int_t(af)+1;
1585  if (!noExponent && nf > maxDigits) flexpo = kTRUE;
1586  if (!noExponent && nf < -maxDigits) flexne = kTRUE;
1587 
1588 // Use x 10 n format. (only powers of 3 allowed)
1589 
1590  if (flexpo) {
1591  flexe = kTRUE;
1592  while (1) {
1593  nexe++;
1594  ww /= 10;
1595  wlabel /= 10;
1596  dwlabel /= 10;
1597  if (nexe%3 == 0 && ww <= TMath::Power(10,maxDigits-1)) break;
1598  }
1599  }
1600 
1601  if (flexne) {
1602  flexe = kTRUE;
1603  rne = 1/TMath::Power(10,maxDigits-2);
1604  while (1) {
1605  nexe--;
1606  ww *= 10;
1607  wlabel *= 10;
1608  dwlabel *= 10;
1609  if (nexe%3 == 0 && ww >= rne) break;
1610  }
1611  }
1612 
1613  na = 0;
1614  for (i=maxDigits-1; i>0; i--) {
1615  if (TMath::Abs(ww) < TMath::Power(10,i)) na = maxDigits-i;
1616  }
1617  ndyn = n1a;
1618  while (ndyn) {
1619  Double_t wdyn = TMath::Abs((wmax-wmin)/ndyn);
1620  if (wdyn <= 0.999 && na < maxDigits-2) {
1621  na++;
1622  ndyn /= 10;
1623  }
1624  else break;
1625  }
1626 
1627  if2 = na;
1628  if1 = TMath::Max(nf+na,maxDigits)+1;
1629 L110:
1630  if (TMath::Min(wmin,wmax) < 0)if1 = if1+1;
1631  if1 = TMath::Min(if1,32);
1632 
1633 // In some cases, if1 and if2 are too small....
1634  while (dwlabel < TMath::Power(10,-if2)) {
1635  if1++;
1636  if2++;
1637  }
1638  coded = &chcoded[0];
1639  if (if1 > 14) if1=14;
1640  if (if2 > 14) if2=14;
1641  if (if2>0) snprintf(coded,8,"%%%d.%df",if1,if2);
1642  else snprintf(coded,8,"%%%d.%df",if1+1,1);
1643  }
1644 
1645 // We draw labels
1646 
1647  snprintf(chtemp,256,"%g",dwlabel);
1648  Int_t ndecimals = 0;
1649  if (optionDecimals) {
1650  char *dot = strchr(chtemp,'.');
1651  if (dot) {
1652  ndecimals = chtemp + strlen(chtemp) -dot;
1653  } else {
1654  char *exp;
1655  exp = strstr(chtemp,"e-");
1656  if (exp) {
1657  sscanf(&exp[2],"%d",&ndecimals);
1658  ndecimals++;
1659  }
1660  }
1661  }
1662  if (optionM) nlabels = n1a-1;
1663  else nlabels = n1a;
1664  wTimeIni = wlabel;
1665  for ( k=0; k<=nlabels; k++) {
1666  if (fFunction) {
1667  Double_t xf = binLow+Double_t(k*nn2)*dxtick;
1668  Double_t zz = fFunction->Eval(xf)-rwmi;
1669  wlabel = xf;
1670  xlabel = zz* axis_length / TMath::Abs(rwma-rwmi);
1671  } else {
1672  xlabel = dxlabel*k;
1673  }
1674  if (optionM) xlabel += 0.5*dxlabel;
1675 
1676  if (!optionText && !optionTime) {
1677  snprintf(label,256,&chcoded[0],wlabel);
1678  label[28] = 0;
1679  wlabel += dwlabel;
1680 
1681  LabelsLimits(label,first,last); //Eliminate blanks
1682 
1683  if (label[first] == '.') { //check if '.' is preceded by a digit
1684  strncpy(chtemp, "0",256);
1685  strlcat(chtemp, &label[first],256);
1686  strncpy(label, chtemp,256);
1687  first = 1; last = strlen(label);
1688  }
1689  if (label[first] == '-' && label[first+1] == '.') {
1690  strncpy(chtemp, "-0",256);
1691  strlcat(chtemp, &label[first+1],256);
1692  strncpy(label, chtemp, 256);
1693  first = 1; last = strlen(label);
1694  }
1695 
1696 // We eliminate the non significant 0 after '.'
1697  if (ndecimals) {
1698  char *adot = strchr(label,'.');
1699  if (adot) adot[ndecimals] = 0;
1700  } else {
1701  while (label[last] == '0') { label[last] = 0; last--;}
1702  }
1703 
1704 // We eliminate the dot, unless dot is forced.
1705  if (label[last] == '.') {
1706  if (!optionDot) { label[last] = 0; last--;}
1707  }
1708 
1709 // Make sure the label is not "-0"
1710  if (last-first == 1 && label[first] == '-'
1711  && label[last] == '0') {
1712  strncpy(label, "0", 256);
1713  label[last] = 0;
1714  }
1715  }
1716 
1717 // Generate the time labels
1718 
1719  if (optionTime) {
1720  timed = wlabel + (int)(timeoffset) - rangeOffset;
1721  timelabel = (time_t)((Long_t)(timed));
1722  if (optionTime == 1) {
1723  utctis = localtime(&timelabel);
1724  } else {
1725  utctis = gmtime(&timelabel);
1726  }
1727  TString timeformattmp;
1728  if (timeformat.Length() < 220) timeformattmp = timeformat;
1729  else timeformattmp = "#splitline{Format}{too long}";
1730 
1731 // Appends fractional part if seconds displayed
1732  if (dwlabel<0.9) {
1733  double tmpdb;
1734  int tmplast;
1735  snprintf(label, 256, "%%S%7.5f", modf(timed,&tmpdb));
1736  tmplast = strlen(label)-1;
1737 
1738 // We eliminate the non significant 0 after '.'
1739  while (label[tmplast] == '0') {
1740  label[tmplast] = 0; tmplast--;
1741  }
1742 
1743  timeformattmp.ReplaceAll("%S",label);
1744 // replace the "0." at the beginning by "s"
1745  timeformattmp.ReplaceAll("%S0.","%Ss");
1746 
1747  }
1748 
1749  if (utctis != nullptr) {
1750  strftime(label, 256, timeformattmp.Data(), utctis);
1751  } else {
1752  strncpy(label, "invalid", 256);
1753  }
1754  strncpy(chtemp, &label[0], 256);
1755  first = 0; last=strlen(label)-1;
1756  wlabel = wTimeIni + (k+1)*dwlabel;
1757  }
1758 
1759 // We generate labels (numeric or alphanumeric).
1760 
1761  if (optionNoopt && !optionInt)
1762  Rotate (xlabel,ylabel,cosphi,sinphi,x0,y0,xx,yy);
1763  else Rotate (xlabel,ylabel,cosphi,sinphi,xx0,yy0,xx,yy);
1764  if (y0 == y1 && !optionDown && !optionUp) {
1765  yy -= 0.80*charheight;
1766  }
1767  if (optionVert) {
1768  if (x0 != x1 && y0 != y1) {
1769  if (optionNoopt && !optionInt)
1770  Rotate (xlabel,0,cosphi,sinphi,x0,y0,xx,yy);
1771  else Rotate (xlabel,0,cosphi,sinphi,xx0,yy0,xx,yy);
1772  if (cosphi > 0 ) yy += ylabel;
1773  if (cosphi < 0 ) yy -= ylabel;
1774  }
1775  }
1776  if (!optionY || (x0 == x1)) {
1777  if (!optionText) {
1778  if (first > last) strncpy(chtemp, " ", 256);
1779  else strncpy(chtemp, &label[first], 256);
1780  if (fNModLabs) ChangeLabelAttributes(k+1, nlabels, textaxis, chtemp);
1781  typolabel = chtemp;
1782  if (!optionTime) typolabel.ReplaceAll("-", "#minus");
1783  textaxis->PaintLatex(gPad->GetX1() + xx*(gPad->GetX2() - gPad->GetX1()),
1784  gPad->GetY1() + yy*(gPad->GetY2() - gPad->GetY1()),
1785  textaxis->GetTextAngle(),
1786  textaxis->GetTextSize(),
1787  typolabel.Data());
1788  if (fNModLabs) ResetLabelAttributes(textaxis);
1789  }
1790  else {
1791  if (optionText == 1) textaxis->PaintLatex(gPad->GetX1() + xx*(gPad->GetX2() - gPad->GetX1()),
1792  gPad->GetY1() + yy*(gPad->GetY2() - gPad->GetY1()),
1793  0,
1794  textaxis->GetTextSize(),
1795  fAxis->GetBinLabel(k+fAxis->GetFirst()));
1796  }
1797  }
1798  else {
1799 
1800 // Text alignment is down
1801  if (!optionText) lnlen = last-first+1;
1802  else {
1803  if (k+1 > nhilab) lnlen = 0;
1804  }
1805  for ( l=1; l<=lnlen; l++) {
1806  if (!optionText) *chtemp = label[first+l-2];
1807  else {
1808  if (lnlen == 0) strncpy(chtemp, " ", 256);
1809  else strncpy(chtemp, "1", 256);
1810  }
1811  typolabel = chtemp;
1812  typolabel.ReplaceAll("-", "#minus");
1813  textaxis->PaintLatex(gPad->GetX1() + xx*(gPad->GetX2() - gPad->GetX1()),
1814  gPad->GetY1() + yy*(gPad->GetY2() - gPad->GetY1()),
1815  0,
1816  textaxis->GetTextSize(),
1817  typolabel.Data());
1818  yy -= charheight*1.3;
1819  }
1820  }
1821  }
1822 
1823 // We use the format x 10 ** n
1824 
1825  if (flexe && !optionText && nexe) {
1826  snprintf(label,256,"#times10^{%d}", nexe);
1827  if (x0 != x1) { xfactor = axis_length+0.1*charheight; yfactor = 0; }
1828  else { xfactor = y1-y0+0.1*charheight; yfactor = 0; }
1829  Rotate (xfactor,yfactor,cosphi,sinphi,x0,y0,xx,yy);
1830  textaxis->SetTextAlign(11);
1831  if (GetLabelFont()%10 < 2) // force TLatex mode in PaintLatex
1832  textaxis->SetTextFont((Int_t)(GetLabelFont()/10)*10+2);
1833  if (fAxis && !strcmp(fAxis->GetName(),"xaxis")) {
1834  xx = xx + fXAxisExpXOffset;
1835  yy = yy + fXAxisExpYOffset;
1836  }
1837  if (fAxis && !strcmp(fAxis->GetName(),"yaxis")) {
1838  xx = xx + fYAxisExpXOffset;
1839  yy = yy + fYAxisExpYOffset;
1840  }
1841  typolabel = label;
1842  typolabel.ReplaceAll("-", "#minus");
1843  textaxis->PaintLatex(gPad->GetX1() + xx*(gPad->GetX2() - gPad->GetX1()),
1844  gPad->GetY1() + yy*(gPad->GetY2() - gPad->GetY1()),
1845  0,
1846  textaxis->GetTextSize(),
1847  typolabel.Data());
1848  }
1849  }
1850  }
1851  }
1852 
1853 // Log axis
1854 
1855  if (optionLog && ndiv) {
1856  UInt_t xi1=0,xi2=0,wi=0,yi1=0,yi2=0,hi=0,xl=0,xh=0;
1857  Bool_t firstintlab = kTRUE, overlap = kFALSE;
1858  if ((wmin == wmax) || (ndiv == 0)) {
1859  Error(where, "wmin (%f) == wmax (%f), or ndiv == 0", wmin, wmax);
1860  goto L210;
1861  }
1862  if (wmin <= 0) {
1863  Error(where, "negative logarithmic axis");
1864  goto L210;
1865  }
1866  if (wmax <= 0) {
1867  Error(where, "negative logarithmic axis");
1868  goto L210;
1869  }
1870  xmnlog = TMath::Log10(wmin);
1871  if (xmnlog > 0) xmnlog += 1.E-6;
1872  else xmnlog -= 1.E-6;
1873  x00 = 0;
1874  x11 = axis_length;
1875  h2 = TMath::Log10(wmax);
1876  h2sav = h2;
1877  if (h2 > 0) h2 += 1.E-6;
1878  else h2 -= 1.E-6;
1879  ih1 = int(xmnlog);
1880  ih2 = 1+int(h2);
1881  nbinin = ih2-ih1+1;
1882  axmul = (x11-x00)/(h2sav-xmnlog);
1883 
1884 // Plot decade and intermediate tick marks
1885  decade = ih1-2;
1886  labelnumber = ih1;
1887  if ( xmnlog > 0 && (xmnlog-Double_t(ih1) > 0) ) labelnumber++;
1888  for (j=1; j<=nbinin; j++) {
1889 
1890 // Plot decade
1891  firstintlab = kTRUE, overlap = kFALSE;
1892  decade++;
1893  if (x0 == x1 && j == 1) ylabel += charheight*0.33;
1894  if (y0 == y1 && j == 1) ylabel -= charheight*0.65;
1895  xone = x00+axmul*(Double_t(decade)-xmnlog);
1896  //the following statement is a trick to circumvent a gcc bug
1897  if (j < 0) printf("j=%d\n",j);
1898  if (x00 > xone) goto L160;
1899  if ((xone-x11)>epsilon) break;
1900  xtwo = xone;
1901  y = 0;
1902  if (!mside) y -= atick[0];
1903  Rotate(xone,y,cosphi,sinphi,x0,y0,xpl2,ypl2);
1904  Rotate(xtwo,atick[0],cosphi,sinphi,x0,y0,xpl1,ypl1);
1905  if (optionVert) {
1906  if ((x0 != x1) && (y0 != y1)) {
1907  if (mside) {
1908  xpl1=xpl2;
1909  if (cosphi > 0) ypl1 = ypl2 + atick[0];
1910  else ypl1 = ypl2 - atick[0];
1911  }
1912  else {
1913  xpl1 = 0.5*(xpl1 + xpl2);
1914  xpl2 = xpl1;
1915  ypl1 = 0.5*(ypl1 + ypl2) + atick[0];
1916  ypl2 = 0.5*(ypl1 + ypl2) - atick[0];
1917  }
1918  }
1919  }
1920  if (!drawGridOnly) PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1921 
1922  if (optionGrid) {
1923  Rotate(xone,0,cosphi,sinphi,x0,y0,xpl2,ypl2);
1924  Rotate(xone,grid_side*gridlength,cosphi,sinphi,x0,y0,xpl1,ypl1);
1925  linegrid->PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1926  }
1927 
1928  if (!drawGridOnly && !optionUnlab) {
1929 
1930 // We generate labels (numeric only).
1931  if (noExponent) {
1932  rlab = TMath::Power(10,labelnumber);
1933  snprintf(label,256, "%f", rlab);
1934  LabelsLimits(label,first,last);
1935  while (last > first) {
1936  if (label[last] != '0') break;
1937  label[last] = 0;
1938  last--;
1939  }
1940  if (label[last] == '.') {label[last] = 0; last--;}
1941  } else {
1942  snprintf(label,256, "%d", labelnumber);
1943  LabelsLimits(label,first,last);
1944  }
1945  Rotate (xone,ylabel,cosphi,sinphi,x0,y0,xx,yy);
1946  if ((x0 == x1) && !optionPara) {
1947  if (lside < 0) {
1948  if (mside < 0) {
1949  if (labelnumber == 0) nch=1;
1950  else nch=2;
1951  xx += nch*charheight;
1952  } else {
1953  xx += 0.25*charheight;
1954  }
1955  }
1956  xx += 0.25*charheight;
1957  }
1958  if ((y0 == y1) && !optionDown && !optionUp) {
1959  if (noExponent) yy += 0.33*charheight;
1960  }
1961  if (n1a == 0)goto L210;
1962  kmod = nbinin/n1a;
1963  if (kmod == 0) kmod=1000000;
1964  if ((nbinin <= n1a) || (j == 1) || (j == nbinin) || ((nbinin > n1a)
1965  && (j%kmod == 0))) {
1966  if (labelnumber == 0) {
1967  textaxis->PaintTextNDC(xx,yy,"1");
1968  } else if (labelnumber == 1) {
1969  textaxis->PaintTextNDC(xx,yy,"10");
1970  } else {
1971  if (noExponent) {
1972  textaxis->PaintTextNDC(xx,yy,&label[first]);
1973  } else {
1974  snprintf(chtemp,256, "10^{%d}", labelnumber);
1975  typolabel = chtemp;
1976  typolabel.ReplaceAll("-", "#minus");
1977  textaxis->PaintLatex(gPad->GetX1() + xx*(gPad->GetX2() - gPad->GetX1()),
1978  gPad->GetY1() + yy*(gPad->GetY2() - gPad->GetY1()),
1979  0, textaxis->GetTextSize(), typolabel.Data());
1980 
1981  }
1982  }
1983  }
1984  labelnumber++;
1985  }
1986 L160:
1987  for (k=2;k<10;k++) {
1988 
1989 // Plot intermediate tick marks
1990  xone = x00+axmul*(TMath::Log10(Double_t(k))+Double_t(decade)-xmnlog);
1991  if (x00 > xone) continue;
1992  if (xone > x11) goto L200;
1993  y = 0;
1994  if (!mside) y -= atick[1];
1995  xtwo = xone;
1996  Rotate(xone,y,cosphi,sinphi,x0,y0,xpl2,ypl2);
1997  Rotate(xtwo,atick[1],cosphi,sinphi,x0,y0,xpl1,ypl1);
1998  if (optionVert) {
1999  if ((x0 != x1) && (y0 != y1)) {
2000  if (mside) {
2001  xpl1 = xpl2;
2002  if (cosphi > 0) ypl1 = ypl2 + atick[1];
2003  else ypl1 = ypl2 - atick[1];
2004  }
2005  else {
2006  xpl1 = 0.5*(xpl1+xpl2);
2007  xpl2 = xpl1;
2008  ypl1 = 0.5*(ypl1+ypl2) + atick[1];
2009  ypl2 = 0.5*(ypl1+ypl2) - atick[1];
2010  }
2011  }
2012  }
2013  idn = n1a*2;
2014  if ((nbinin <= idn) || ((nbinin > idn) && (k == 5))) {
2015  if (!drawGridOnly) PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
2016 
2017 // Draw the intermediate LOG labels if requested
2018 
2019  if (moreLogLabels && !optionUnlab && !drawGridOnly && !overlap) {
2020  if (noExponent) {
2021  rlab = Double_t(k)*TMath::Power(10,labelnumber-1);
2022  snprintf(chtemp,256, "%g", rlab);
2023  } else {
2024  if (labelnumber-1 == 0) {
2025  snprintf(chtemp,256, "%d", k);
2026  } else if (labelnumber-1 == 1) {
2027  snprintf(chtemp,256, "%d", 10*k);
2028  } else {
2029  snprintf(chtemp,256, "%d#times10^{%d}", k, labelnumber-1);
2030  }
2031  }
2032  Rotate (xone,ylabel,cosphi,sinphi,x0,y0,xx,yy);
2033  if ((x0 == x1) && !optionPara) {
2034  if (lside < 0) {
2035  if (mside < 0) {
2036  if (labelnumber == 0) nch=1;
2037  else nch=2;
2038  xx += nch*charheight;
2039  } else {
2040  if (labelnumber >= 0) xx += 0.25*charheight;
2041  else xx += 0.50*charheight;
2042  }
2043  }
2044  xx += 0.25*charheight;
2045  }
2046  if ((y0 == y1) && !optionDown && !optionUp) {
2047  if (noExponent) yy += 0.33*charheight;
2048  }
2049  if (optionVert) {
2050  if ((x0 != x1) && (y0 != y1)) {
2051  Rotate(xone,ylabel,cosphi,sinphi,x0,y0,xx,yy);
2052  if (cosphi > 0) yy += ylabel;
2053  else yy -= ylabel;
2054  }
2055  }
2056  textaxis->SetTitle(chtemp);
2057  Double_t u = gPad->GetX1() + xx*(gPad->GetX2() - gPad->GetX1());
2058  Double_t v = gPad->GetY1() + yy*(gPad->GetY2() - gPad->GetY1());
2059  if (firstintlab) {
2060  textaxis->GetBoundingBox(wi, hi); wi=(UInt_t)(wi*1.3); hi=(UInt_t)(hi*1.3);
2061  xi1 = gPad->XtoAbsPixel(u);
2062  yi1 = gPad->YtoAbsPixel(v);
2063  firstintlab = kFALSE;
2064  typolabel = chtemp;
2065  typolabel.ReplaceAll("-", "#minus");
2066  textaxis->PaintLatex(u,v,0,textaxis->GetTextSize(),typolabel.Data());
2067  } else {
2068  xi2 = gPad->XtoAbsPixel(u);
2069  yi2 = gPad->YtoAbsPixel(v);
2070  xl = TMath::Min(xi1,xi2);
2071  xh = TMath::Max(xi1,xi2);
2072  if ((x0 == x1 && yi1-hi <= yi2) || (y0 == y1 && xl+wi >= xh)){
2073  overlap = kTRUE;
2074  } else {
2075  xi1 = xi2;
2076  yi1 = yi2;
2077  textaxis->GetBoundingBox(wi, hi); wi=(UInt_t)(wi*1.3); hi=(UInt_t)(hi*1.3);
2078  typolabel = chtemp;
2079  typolabel.ReplaceAll("-", "#minus");
2080  textaxis->PaintLatex(u,v,0,textaxis->GetTextSize(),typolabel.Data());
2081  }
2082  }
2083  }
2084 
2085 // Draw the intermediate LOG grid if only three decades are requested
2086  if (optionGrid && nbinin <= 5 && ndiv > 100) {
2087  Rotate(xone,0,cosphi,sinphi,x0,y0,xpl2, ypl2);
2088  Rotate(xone,grid_side*gridlength,cosphi,sinphi,x0,y0, xpl1,ypl1);
2089  linegrid->PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
2090  }
2091  } //endif ((nbinin <= idn) ||
2092  } //endfor (k=2;k<10;k++)
2093  } //endfor (j=1; j<=nbinin; j++)
2094 L200:
2095  Int_t dummy = 0; if (dummy) { }
2096  } //endif (optionLog && ndiv)
2097 
2098 
2099 L210:
2100  if (optionGrid) delete linegrid;
2101  delete textaxis;
2102 }
2103 
2104 ////////////////////////////////////////////////////////////////////////////////
2105 /// Internal method for axis labels optimisation. This method adjusts the bining
2106 /// of the axis in order to have integer values for the labels.
2107 ///
2108 /// \param[in] A1,A2 Old WMIN,WMAX
2109 /// \param[out] binLow,binHigh New WMIN,WMAX
2110 /// \param[in] nold Old NDIV (primary divisions)
2111 /// \param[out] nbins New NDIV
2112 /// \param[out] binWidth Bin width
2113 
2114 void TGaxis::AdjustBinSize(Double_t A1, Double_t A2, Int_t nold
2115  ,Double_t &binLow, Double_t &binHigh, Int_t &nbins, Double_t &binWidth)
2116 {
2117 
2118  binWidth = TMath::Abs(A2-A1)/Double_t(nold);
2119  if (binWidth <= 1) { binWidth = 1; binLow = int(A1); }
2120  else {
2121  Int_t width = int(binWidth/5) + 1;
2122  binWidth = 5*width;
2123  binLow = int(A1/binWidth)*binWidth;
2124 
2125 // We determine binLow to have one tick mark at 0
2126 // if there are negative labels.
2127 
2128  if (A1 < 0) {
2129  for (Int_t ic=0; ic<1000; ic++) {
2130  Double_t rbl = binLow/binWidth;
2131  Int_t ibl = int(binLow/binWidth);
2132  if ( (rbl-ibl) == 0 || ic > width) { binLow -= 5; break;}
2133  }
2134  }
2135  }
2136  binHigh = int(A2);
2137  nbins = 0;
2138  Double_t xb = binLow;
2139  while (xb <= binHigh) {
2140  xb += binWidth;
2141  nbins++;
2142  }
2143  binHigh = xb - binWidth;
2144 }
2145 
2146 ////////////////////////////////////////////////////////////////////////////////
2147 /// Internal method to find first and last character of a label.
2148 
2149 void TGaxis::LabelsLimits(const char *label, Int_t &first, Int_t &last)
2150 {
2151  last = strlen(label)-1;
2152  for (Int_t i=0; i<=last; i++) {
2153  if (strchr("1234567890-+.", label[i]) ) { first = i; return; }
2154  }
2155  Error("LabelsLimits", "attempt to draw a blank label");
2156 }
2157 
2158 ////////////////////////////////////////////////////////////////////////////////
2159 /// Internal method to rotate axis coordinates.
2160 
2162  ,Double_t XT, Double_t YT, Double_t &U, Double_t &V)
2163 {
2164  U = CFI*X-SFI*Y+XT;
2165  V = SFI*X+CFI*Y+YT;
2166 }
2167 
2168 ////////////////////////////////////////////////////////////////////////////////
2169 /// Save primitive as a C++ statement(s) on output stream out
2170 
2171 void TGaxis::SavePrimitive(std::ostream &out, Option_t * /*= ""*/)
2172 {
2173  char quote = '"';
2174  if (gROOT->ClassSaved(TGaxis::Class())) {
2175  out<<" ";
2176  } else {
2177  out<<" TGaxis *";
2178  }
2179  out<<"gaxis = new TGaxis("<<fX1<<","<<fY1<<","<<fX2<<","<<fY2
2180  <<","<<fWmin<<","<<fWmax<<","<<fNdiv<<","<<quote<<fChopt.Data()<<quote<<");"<<std::endl;
2181  out<<" gaxis->SetLabelOffset("<<GetLabelOffset()<<");"<<std::endl;
2182  out<<" gaxis->SetLabelSize("<<GetLabelSize()<<");"<<std::endl;
2183  out<<" gaxis->SetTickSize("<<GetTickSize()<<");"<<std::endl;
2184  out<<" gaxis->SetGridLength("<<GetGridLength()<<");"<<std::endl;
2185  out<<" gaxis->SetTitleOffset("<<GetTitleOffset()<<");"<<std::endl;
2186  out<<" gaxis->SetTitleSize("<<GetTitleSize()<<");"<<std::endl;
2187  out<<" gaxis->SetTitleColor("<<GetTextColor()<<");"<<std::endl;
2188  out<<" gaxis->SetTitleFont("<<GetTextFont()<<");"<<std::endl;
2189 
2190  if (strlen(GetName())) {
2191  out<<" gaxis->SetName("<<quote<<GetName()<<quote<<");"<<std::endl;
2192  }
2193  if (strlen(GetTitle())) {
2194  out<<" gaxis->SetTitle("<<quote<<GetTitle()<<quote<<");"<<std::endl;
2195  }
2196 
2197  if (fLabelColor != 1) {
2198  if (fLabelColor > 228) {
2200  out<<" gaxis->SetLabelColor(ci);" << std::endl;
2201  } else
2202  out<<" gaxis->SetLabelColor("<<GetLabelColor()<<");"<<std::endl;
2203  }
2204  if (fLineColor != 1) {
2205  if (fLineColor > 228) {
2207  out<<" gaxis->SetLineColor(ci);" << std::endl;
2208  } else
2209  out<<" gaxis->SetLineColor("<<GetLineColor()<<");"<<std::endl;
2210  }
2211  if (fLineStyle != 1) {
2212  out<<" gaxis->SetLineStyle("<<GetLineStyle()<<");"<<std::endl;
2213  }
2214  if (fLineWidth != 1) {
2215  out<<" gaxis->SetLineWidth("<<GetLineWidth()<<");"<<std::endl;
2216  }
2217  if (fLabelFont != 62) {
2218  out<<" gaxis->SetLabelFont("<<GetLabelFont()<<");"<<std::endl;
2219  }
2221  out<<" gaxis->SetMoreLogLabels();"<<std::endl;
2222  }
2223  if (TestBit(TAxis::kNoExponent)) {
2224  out<<" gaxis->SetNoExponent();"<<std::endl;
2225  }
2226 
2227  out<<" gaxis->Draw();"<<std::endl;
2228 }
2229 
2230 ////////////////////////////////////////////////////////////////////////////////
2231 /// Set the decimals flag. By default, blank characters are stripped, and then the
2232 /// label is correctly aligned. The dot, if last character of the string, is also
2233 /// stripped, unless this option is specified. One can disable the option by
2234 /// calling `axis.SetDecimals(kTRUE)`.
2235 /// Note the bit is set in fBits (as opposed to fBits2 in TAxis!)
2236 
2237 void TGaxis::SetDecimals(Bool_t dot)
2238 {
2239  if (dot) SetBit(TAxis::kDecimals);
2240  else ResetBit(TAxis::kDecimals);
2242 
2243 ////////////////////////////////////////////////////////////////////////////////
2244 /// Specify a function to map the axis values.
2245 
2246 void TGaxis::SetFunction(const char *funcname)
2247 {
2248  fFunctionName = funcname;
2249  if (!funcname[0]) {
2251  return;
2252  }
2253  fFunction = (TF1*)gROOT->GetFunction(funcname);
2254  if (!fFunction) {
2255  Error("SetFunction", "unknown function: %s", funcname);
2256  } else {
2257  fWmin = fFunction->GetXmin();
2258  fWmax = fFunction->GetXmax();
2259  }
2260 }
2261 
2262 ////////////////////////////////////////////////////////////////////////////////
2263 /// Define new text attributes for the label number "labNum". It allows to do a
2264 /// fine tuning of the labels. All the attributes can be changed and even the
2265 /// label text itself.
2266 ///
2267 /// \param[in] labNum Number of the label to be changed, negative numbers start from the end
2268 /// \param[in] labAngle New angle value
2269 /// \param[in] labSize New size (0 erase the label)
2270 /// \param[in] labAlign New alignment value
2271 /// \param[in] labColor New label color
2272 /// \param[in] labText New label text
2273 ///
2274 /// If an attribute should not be changed just give the value
2275 /// "-1".The following macro gives an example:
2276 ///
2277 /// Begin_Macro(source)
2278 /// {
2279 /// c1 = new TCanvas("c1","Examples of Gaxis",10,10,900,500);
2280 /// c1->Range(-6,-0.1,6,0.1);
2281 /// TGaxis *axis1 = new TGaxis(-5.5,0.,5.5,0.,0.0,100,510,"");
2282 /// axis1->SetName("axis1");
2283 /// axis1->SetTitle("Axis Title");
2284 /// axis1->SetTitleSize(0.05);
2285 /// axis1->SetTitleColor(kBlue);
2286 /// axis1->SetTitleFont(42);
2287 /// axis1->ChangeLabel(1,-1,-1,-1,2);
2288 /// axis1->ChangeLabel(3,-1,0.);
2289 /// axis1->ChangeLabel(5,30.,-1,0);
2290 /// axis1->ChangeLabel(6,-1,-1,-1,3,-1,"6th label");
2291 /// axis1->ChangeLabel(-2,-1,-1,-1,3,-1,"2nd to last label");
2292 /// axis1->Draw();
2293 /// }
2294 /// End_Macro
2295 ///
2296 /// If labnum=0 the list of modified labels is reset.
2297 
2298 void TGaxis::ChangeLabel(Int_t labNum, Double_t labAngle, Double_t labSize,
2299  Int_t labAlign, Int_t labColor, Int_t labFont,
2300  TString labText)
2301 {
2303  if (!fModLabs) fModLabs = new TList();
2304 
2305  // Reset the list of modified labels.
2306  if (labNum == 0) {
2307  delete fModLabs;
2308  fModLabs = 0;
2309  fNModLabs = 0;
2310  return;
2311  }
2312 
2313  TAxisModLab *ml = new TAxisModLab();
2314  ml->SetLabNum(labNum);
2315  ml->SetAngle(labAngle);
2316  ml->SetSize(labSize);
2317  ml->SetAlign(labAlign);
2318  ml->SetColor(labColor);
2319  ml->SetFont(labFont);
2320  ml->SetText(labText);
2321 
2322  fModLabs->Add((TObject*)ml);
2323 }
2324 
2325 ////////////////////////////////////////////////////////////////////////////////
2326 /// Change the label attributes of label number i. If needed.
2327 static Double_t SavedTextAngle;
2328 static Double_t SavedTextSize;
2329 static Int_t SavedTextAlign;
2330 static Int_t SavedTextColor;
2333 void TGaxis::ChangeLabelAttributes(Int_t i, Int_t nlabels, TLatex* t, char* c)
2335  if (!fModLabs) return;
2336 
2338  TAxisModLab *ml;
2339  Int_t labNum;
2340  while ( (ml = (TAxisModLab*)next()) ) {
2341  SavedTextAngle = t->GetTextAngle();
2342  SavedTextSize = t->GetTextSize();
2343  SavedTextAlign = t->GetTextAlign();
2344  SavedTextColor = t->GetTextColor();
2345  SavedTextFont = t->GetTextFont();
2346  labNum = ml->GetLabNum();
2347  if (labNum < 0) labNum = nlabels + labNum + 2;
2348  if (i == labNum) {
2349  if (ml->GetAngle()>=0.) t->SetTextAngle(ml->GetAngle());
2350  if (ml->GetSize()>=0.) t->SetTextSize(ml->GetSize());
2351  if (ml->GetAlign()>0) t->SetTextAlign(ml->GetAlign());
2352  if (ml->GetColor()>=0) t->SetTextColor(ml->GetColor());
2353  if (ml->GetFont()>0) t->SetTextFont(ml->GetFont());
2354  if (!(ml->GetText().IsNull())) strncpy(c, (ml->GetText()).Data(), 256);
2355  return;
2356  }
2357  }
2358 }
2359 
2360 ////////////////////////////////////////////////////////////////////////////////
2361 /// Reset the label attributes to the value they have before the last call to
2362 /// ChangeLabelAttributes.
2363 
2365 {
2366  t->SetTextAngle(SavedTextAngle);
2367  t->SetTextSize(SavedTextSize);
2368  t->SetTextAlign(SavedTextAlign);
2369  t->SetTextColor(SavedTextColor);
2370  t->SetTextFont(SavedTextFont);
2371 }
2372 
2373 ////////////////////////////////////////////////////////////////////////////////
2374 /// Static function to set `fgMaxDigits` for axis.`fgMaxDigits` is
2375 /// the maximum number of digits permitted for the axis labels above which the
2376 /// notation with 10^N is used.For example, to accept 6 digits number like 900000
2377 /// on an axis call `TGaxis::SetMaxDigits(6)`. The default value is 5.
2378 /// `fgMaxDigits` must be greater than 0.
2379 
2380 void TGaxis::SetMaxDigits(Int_t maxd)
2381 {
2382  fgMaxDigits = maxd;
2383  if (maxd < 1) fgMaxDigits = 1;
2385 
2386 ////////////////////////////////////////////////////////////////////////////////
2387 /// Change the name of the axis.
2388 
2389 void TGaxis::SetName(const char *name)
2390 {
2391  fName = name;
2392 }
2394 ////////////////////////////////////////////////////////////////////////////////
2395 /// Set the kMoreLogLabels bit flag. When this option is selected more labels are
2396 /// drawn when in logarithmic scale and there is a small number of decades (less than 3).
2397 /// Note that this option is automatically inherited from TAxis
2398 
2400 {
2401  if (more) SetBit(TAxis::kMoreLogLabels);
2404 
2405 ////////////////////////////////////////////////////////////////////////////////
2406 /// Set the NoExponent flag. By default, an exponent of the form 10^N is used
2407 /// when the label values are either all very small or very large. One can disable
2408 /// the exponent by calling axis.SetNoExponent(kTRUE).
2409 
2410 void TGaxis::SetNoExponent(Bool_t noExponent)
2411 {
2412  if (noExponent) SetBit(TAxis::kNoExponent);
2415 
2416 ////////////////////////////////////////////////////////////////////////////////
2417 /// To set axis options.
2418 
2419 void TGaxis::SetOption(Option_t *option)
2420 {
2421  fChopt = option;
2422 }
2424 ////////////////////////////////////////////////////////////////////////////////
2425 /// Change the title of the axis.
2426 
2427 void TGaxis::SetTitle(const char *title)
2428 {
2429  fTitle = title;
2430 }
2432 ////////////////////////////////////////////////////////////////////////////////
2433 /// Change the format used for time plotting.
2434 /// The format string for date and time use the same options as the one used
2435 /// in the standard strftime C function, i.e. :
2436 ///
2437 /// for date :
2438 ///
2439 /// - `%a` abbreviated weekday name
2440 /// - `%b` abbreviated month name
2441 /// - `%d` day of the month (01-31)
2442 /// - `%m` month (01-12)
2443 /// - `%y` year without century
2444 ///
2445 /// for time :
2446 ///
2447 /// - `%H` hour (24-hour clock)
2448 /// - `%I` hour (12-hour clock)
2449 /// - `%p` local equivalent of AM or PM
2450 /// - `%M` minute (00-59)
2451 /// - `%S` seconds (00-61)
2452 /// - `%%` %
2453 
2454 void TGaxis::SetTimeFormat(const char *tformat)
2455 {
2456  TString timeformat = tformat;
2457 
2458  if (timeformat.Index("%F")>=0 || timeformat.IsNull()) {
2459  fTimeFormat = timeformat;
2460  return;
2461  }
2462 
2463  Int_t idF = fTimeFormat.Index("%F");
2464  if (idF>=0) {
2465  Int_t lnF = fTimeFormat.Length();
2466  TString stringtimeoffset = fTimeFormat(idF,lnF);
2467  fTimeFormat = tformat;
2468  fTimeFormat.Append(stringtimeoffset);
2469  } else {
2470  fTimeFormat = tformat;
2472  }
2473 }
2474 
2475 ////////////////////////////////////////////////////////////////////////////////
2476 /// Change the time offset. If option = "gmt", set display mode to GMT.
2477 
2478 void TGaxis::SetTimeOffset(Double_t toffset, Option_t *option)
2479 {
2480  TString opt = option;
2481  opt.ToLower();
2483  char tmp[20];
2484  time_t timeoff;
2485  struct tm* utctis;
2486  Int_t idF = fTimeFormat.Index("%F");
2487  if (idF>=0) fTimeFormat.Remove(idF);
2488  fTimeFormat.Append("%F");
2489 
2490  timeoff = (time_t)((Long_t)(toffset));
2491 
2492  // offset is always saved in GMT to allow file transport
2493  // to different time zones
2494  utctis = gmtime(&timeoff);
2495 
2496  if (utctis != nullptr) {
2497  strftime(tmp, 20,"%Y-%m-%d %H:%M:%S",utctis);
2498  fTimeFormat.Append(tmp);
2499  } else {
2500  fTimeFormat.Append("1970-01-01 00:00:00");
2501  }
2502 
2503  // append the decimal part of the time offset
2504  Double_t ds = toffset-(Int_t)toffset;
2505  snprintf(tmp,20,"s%g",ds);
2506  fTimeFormat.Append(tmp);
2507 
2508  // add GMT/local option
2509  if (opt.Contains("gmt")) fTimeFormat.Append(" GMT");
2510 }
2511 
2512 ////////////////////////////////////////////////////////////////////////////////
2513 /// Static function to set X and Y offset of the axis 10^n notation.
2514 /// It is in % of the pad size. It can be negative.
2515 /// axis specifies which axis ("x","y"), default = "x"
2516 /// if axis="xz" set the two axes
2517 
2518 void TGaxis::SetExponentOffset(Float_t xoff, Float_t yoff, Option_t *axis)
2519 {
2520  TString opt = axis;
2521  opt.ToLower();
2523  if (opt.Contains("x")) {
2524  fXAxisExpXOffset = xoff;
2525  fXAxisExpYOffset = yoff;
2526  }
2527  if (opt.Contains("y")) {
2528  fYAxisExpXOffset = xoff;
2529  fYAxisExpYOffset = yoff;
2530  }
2531 }
2532 
2533 ////////////////////////////////////////////////////////////////////////////////
2534 /// Stream an object of class TGaxis.
2535 
2536 void TGaxis::Streamer(TBuffer &R__b)
2537 {
2538  if (R__b.IsReading()) {
2539  UInt_t R__s, R__c;
2540  Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
2541  if (R__v > 3) {
2542  R__b.ReadClassBuffer(TGaxis::Class(), this, R__v, R__s, R__c);
2543  return;
2544  }
2545  //====process old versions before automatic schema evolution
2546  TLine::Streamer(R__b);
2547  TAttText::Streamer(R__b);
2548  R__b >> fNdiv;
2549  R__b >> fWmin;
2550  R__b >> fWmax;
2551  R__b >> fGridLength;
2552  R__b >> fTickSize;
2553  R__b >> fLabelOffset;
2554  R__b >> fLabelSize;
2555  R__b >> fTitleOffset;
2556  R__b >> fTitleSize;
2557  R__b >> fLabelFont;
2558  if (R__v > 2) {
2559  R__b >> fLabelColor;
2560  }
2561  fChopt.Streamer(R__b);
2562  fName.Streamer(R__b);
2563  fTitle.Streamer(R__b);
2564  fTimeFormat.Streamer(R__b);
2565  if (R__v > 1) {
2566  fFunctionName.Streamer(R__b);
2567  fFunction = (TF1*)gROOT->GetFunction(fFunctionName.Data());
2568  }
2569  R__b.CheckByteCount(R__s, R__c, TGaxis::IsA());
2570  //====end of old versions
2571 
2572  } else {
2573  R__b.WriteClassBuffer(TGaxis::Class(),this);
2574  }
2575 }
static time_t MktimeFromUTC(tm_t *tmstruct)
Equivalent of standard routine "mktime" but using the assumption that tm struct is filled with UTC...
Definition: TTimeStamp.cxx:767
Int_t GetLabelFont() const
Definition: TGaxis.h:86
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition: TAttLine.h:49
Bool_t IsReading() const
Definition: TBuffer.h:83
virtual Float_t GetTickLength() const
Definition: TAttAxis.h:49
virtual void SetName(const char *name)
Change the name of the axis.
Definition: TGaxis.cxx:2393
virtual void SetMoreLogLabels(Bool_t more=kTRUE)
Set the kMoreLogLabels bit flag.
Definition: TGaxis.cxx:2403
float xmin
Definition: THbookFile.cxx:93
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
Int_t GetColor()
Definition: TAxisModLab.h:50
const char * GetBinLabel(Int_t bin) const
Return label for bin.
Definition: TAxis.cxx:426
void SetText(TString t="")
Set modified label text.
Definition: TAxisModLab.cxx:84
Int_t GetFirst() const
Return first bin on the axis i.e.
Definition: TAxis.cxx:444
Double_t fX1
X of 1st point.
Definition: TLine.h:36
short Version_t
Definition: RtypesCore.h:61
virtual Color_t GetTextColor() const
Return the text color.
Definition: TAttText.h:40
float Float_t
Definition: RtypesCore.h:53
Style_t GetGridStyle() const
Definition: TStyle.h:224
virtual Float_t GetLabelOffset() const
Definition: TAttAxis.h:45
virtual Short_t GetTextAlign() const
Return the text alignment.
Definition: TAttText.h:38
return c
const char Option_t
Definition: RtypesCore.h:62
float ymin
Definition: THbookFile.cxx:93
static Int_t fgMaxDigits
! Number of digits above which the 10>N notation is used
Definition: TGaxis.h:55
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:635
R__EXTERN TStyle * gStyle
Definition: TStyle.h:418
static void SetExponentOffset(Float_t xoff=0., Float_t yoff=0., Option_t *axis="xy")
Static function to set X and Y offset of the axis 10^n notation.
Definition: TGaxis.cxx:2522
Bool_t TestBit(UInt_t f) const
Definition: TObject.h:157
#define BIT(n)
Definition: Rtypes.h:120
virtual Color_t GetAxisColor() const
Definition: TAttAxis.h:42
static void SaveColor(std::ostream &out, Int_t ci)
Save a color with index > 228 as a C++ statement(s) on output stream out.
Definition: TColor.cxx:2021
Float_t GetLabelSize() const
Definition: TGaxis.h:88
virtual Float_t GetTextAngle() const
Return the text angle.
Definition: TAttText.h:39
static void Optimize(Double_t A1, Double_t A2, Int_t nold, Double_t &BinLow, Double_t &BinHigh, Int_t &nbins, Double_t &BWID, Option_t *option="")
static function to compute reasonable axis limits
virtual void CenterLabels(Bool_t center=kTRUE)
If center = kTRUE axis labels are centered in the center of the bin.
Definition: TGaxis.cxx:633
TString fTimeFormat
Time format, ex: 09/12/99 12:34:00.
Definition: TGaxis.h:49
Buffer base class used for serializing objects.
Definition: TBuffer.h:42
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)=0
#define gROOT
Definition: TROOT.h:364
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition: TString.h:582
TGaxis & operator=(const TGaxis &)
Assignment operator.
Definition: TGaxis.cxx:591
virtual void SetTitle(const char *title="")
Change the title of the axis.
Definition: TGaxis.cxx:2431
Basic string class.
Definition: TString.h:137
virtual void ImportAxisAttributes(TAxis *axis)
Internal method to import TAxis attributes to this TGaxis.
Definition: TGaxis.cxx:692
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:170
void ToLower()
Change string to lower-case.
Definition: TString.cxx:1089
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
Double_t fX2
X of 2nd point.
Definition: TLine.h:38
void SetFont(Int_t f=-1)
Set modified label font.
Definition: TAxisModLab.cxx:77
int nbins[3]
virtual void PaintTextNDC(Double_t u, Double_t v, const char *text)
Draw this text with new coordinates in NDC.
Definition: TText.cxx:762
virtual Float_t GetLabelSize() const
Definition: TAttAxis.h:46
TAxis * fAxis
! Pointer to original TAxis axis (if any)
Definition: TGaxis.h:52
const char * Class
Definition: TXMLSetup.cxx:64
void SetAlign(Int_t a=-1)
Set modified label alignment.
Definition: TAxisModLab.cxx:63
TString fChopt
Axis options.
Definition: TGaxis.h:46
Float_t fLabelOffset
Offset of label wrt axis.
Definition: TGaxis.h:38
Bool_t GetDecimals() const
Definition: TAxis.h:122
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
virtual void DrawAxis(Double_t xmin, Double_t ymin, Double_t xmax, Double_t ymax, Double_t wmin, Double_t wmax, Int_t ndiv=510, Option_t *chopt="", Double_t gridlength=0)
Draw this axis with new attributes.
Definition: TGaxis.cxx:653
virtual Width_t GetLineWidth() const
Return the line width.
Definition: TAttLine.h:41
void ChangeLabelAttributes(Int_t i, Int_t nlabels, TLatex *t, char *c)
Definition: TGaxis.cxx:2337
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:501
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:739
double beta(double x, double y)
Calculates the beta function.
Color_t GetGridColor() const
Definition: TStyle.h:223
if object in a list can be deleted
Definition: TObject.h:56
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition: TObject.cxx:165
static Double_t SavedTextSize
Definition: TGaxis.cxx:2332
virtual Style_t GetTitleFont() const
Definition: TAttAxis.h:51
void SetTitleSize(Float_t titlesize)
Definition: TGaxis.h:132
Double_t GetSize()
Definition: TAxisModLab.h:48
virtual void SetTextFont(Font_t tfont=62)
Set the text font.
Definition: TAttText.h:51
static Float_t fXAxisExpYOffset
! Exponent Y offset for the X axis
Definition: TGaxis.h:57
virtual Style_t GetLineStyle() const
Return the line style.
Definition: TAttLine.h:40
static const double x2[5]
virtual void PaintLineNDC(Double_t u1, Double_t v1, Double_t u2, Double_t v2)
Draw this line with new coordinates in NDC.
Definition: TLine.cxx:389
static Float_t fYAxisExpYOffset
! Exponent Y offset for the Y axis
Definition: TGaxis.h:59
virtual void SetText(Double_t x, Double_t y, const char *text)
Definition: TText.h:82
TString fTitle
Axis title.
Definition: TGaxis.h:48
Float_t fGridLength
Length of the grid in NDC.
Definition: TGaxis.h:36
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save primitive as a C++ statement(s) on output stream out.
Definition: TGaxis.cxx:2175
virtual Float_t GetTextSize() const
Return the text size.
Definition: TAttText.h:42
Int_t fNModLabs
Number of modified labels.
Definition: TGaxis.h:45
To draw Mathematical Formula.
Definition: TLatex.h:25
THashList * GetLabels() const
Definition: TAxis.h:123
Double_t Log10(Double_t x)
Definition: TMath.h:529
void SetOption(Option_t *option="")
To set axis options.
Definition: TGaxis.cxx:2423
virtual Double_t GetBinCenter(Int_t bin) const
Return center of bin.
Definition: TAxis.cxx:464
virtual const char * GetName() const
Returns name of object.
Definition: TGaxis.h:91
TString & Append(const char *cs)
Definition: TString.h:492
std::vector< std::vector< double > > Data
void SetTimeFormat(const char *tformat)
Change the format used for time plotting.
Definition: TGaxis.cxx:2458
void SetAngle(Double_t a=-1.)
Set modified label angle.
Definition: TAxisModLab.cxx:49
void SetColor(Int_t c=-1)
Set modified label color.
Definition: TAxisModLab.cxx:70
TObject & operator=(const TObject &rhs)
TObject assignment operator.
Definition: TObject.cxx:103
TString fName
Axis name.
Definition: TGaxis.h:47
Double_t ATan2(Double_t, Double_t)
Definition: TMath.h:454
TList * GetModifiedLabels() const
Definition: TAxis.h:124
virtual void PaintLatex(Double_t x, Double_t y, Double_t angle, Double_t size, const char *text)
Main drawing function.
Definition: TLatex.cxx:2038
virtual ~TGaxis()
TGaxis default destructor.
Definition: TGaxis.cxx:624
void Error(const char *location, const char *msgfmt,...)
virtual const char * GetTimeFormat() const
Definition: TAxis.h:133
const Int_t kHori
Definition: TGaxis.cxx:44
virtual Color_t GetLabelColor() const
Definition: TAttAxis.h:43
void SetLabelSize(Float_t labelsize)
Definition: TGaxis.h:114
virtual void SetTextAlign(Short_t align=11)
Set the text alignment.
Definition: TAttText.h:47
A doubly linked list.
Definition: TList.h:47
TList * fModLabs
List of modified labels.
Definition: TGaxis.h:53
Float_t GetGridLength() const
Definition: TGaxis.h:83
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:46
Float_t fTextAngle
Text angle.
Definition: TAttText.h:27
float ymax
Definition: THbookFile.cxx:93
Style_t fLineStyle
Line style.
Definition: TAttLine.h:28
virtual void PaintAxis(Double_t xmin, Double_t ymin, Double_t xmax, Double_t ymax, Double_t &wmin, Double_t &wmax, Int_t &ndiv, Option_t *chopt="", Double_t gridlength=0, Bool_t drawGridOnly=kFALSE)
Control function to draw an axis.
Definition: TGaxis.cxx:742
Int_t GetLast() const
Return last bin on the axis i.e.
Definition: TAxis.cxx:455
const char * GetTitle() const
Returns title of object.
Definition: TAxis.h:135
Double_t GetXsize()
Return size of the formula along X in pad coordinates.
Definition: TLatex.cxx:2507
Class to manage histogram axis.
Definition: TAxis.h:36
static void SetMaxDigits(Int_t maxd=5)
Static function to set fgMaxDigits for axis.
Definition: TGaxis.cxx:2384
SVector< double, 2 > v
Definition: Dict.h:5
void SetFunction(const char *funcname="")
Specify a function to map the axis values.
Definition: TGaxis.cxx:2250
Double_t fY1
Y of 1st point.
Definition: TLine.h:37
Double_t fWmin
Lowest value on the axis.
Definition: TGaxis.h:34
virtual Font_t GetTextFont() const
Return the text font.
Definition: TAttText.h:41
Color_t fLineColor
Line color.
Definition: TAttLine.h:27
Int_t GetAlign()
Definition: TAxisModLab.h:49
virtual void SetTextAngle(Float_t tangle=0)
Set the text angle.
Definition: TAttText.h:48
Text Attributes class.
Definition: TAttText.h:24
unsigned int UInt_t
Definition: RtypesCore.h:42
virtual Float_t GetTitleOffset() const
Definition: TAttAxis.h:47
Width_t fLineWidth
Line width.
Definition: TAttLine.h:29
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:925
Ssiz_t Length() const
Definition: TString.h:390
A simple line.
Definition: TLine.h:33
static Int_t SavedTextAlign
Definition: TGaxis.cxx:2333
TLine * l
Definition: textangle.C:4
void GetBoundingBox(UInt_t &w, UInt_t &h, Bool_t angle=kFALSE)
Return text size in pixels.
Definition: TLatex.cxx:2538
The axis painter class.
Definition: TGaxis.h:30
void SetSize(Double_t s=-1.)
Set modified label size.
Definition: TAxisModLab.cxx:56
float xmax
Definition: THbookFile.cxx:93
virtual Color_t GetTitleColor() const
Definition: TAttAxis.h:50
Font_t fTextFont
Text font.
Definition: TAttText.h:31
const Double_t kPI
Definition: TEllipse.cxx:22
virtual Double_t GetXmin() const
Definition: TF1.h:388
#define gVirtualX
Definition: TVirtualX.h:362
REAL epsilon
Definition: triangle.c:617
TF1 * fFunction
! Pointer to function computing axis values
Definition: TGaxis.h:51
Double_t Cos(Double_t)
Definition: TMath.h:424
void SetLabNum(Int_t n=0)
Set modified label number.
Definition: TAxisModLab.cxx:42
static Float_t fYAxisExpXOffset
! Exponent X offset for the Y axis
Definition: TGaxis.h:58
virtual Color_t GetLineColor() const
Return the line color.
Definition: TAttLine.h:39
Double_t Pi()
Definition: TMath.h:44
void ResetLabelAttributes(TLatex *t)
Reset the label attributes to the value they have before the last call to ChangeLabelAttributes.
Definition: TGaxis.cxx:2368
TString & Remove(Ssiz_t pos)
Definition: TString.h:616
long Long_t
Definition: RtypesCore.h:50
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
virtual Double_t Eval(Double_t x, Double_t y=0, Double_t z=0, Double_t t=0) const
Evaluate this function.
Definition: TF1.cxx:1196
Int_t GetLabNum()
Definition: TAxisModLab.h:46
Int_t fLabelFont
Font for labels.
Definition: TGaxis.h:44
virtual void AdjustBinSize(Double_t A1, Double_t A2, Int_t nold, Double_t &BinLow, Double_t &BinHigh, Int_t &nbins, Double_t &BinWidth)
Internal method for axis labels optimisation.
Definition: TGaxis.cxx:2118
void SetLabelOffset(Float_t labeloffset)
Definition: TGaxis.h:113
static const double x1[5]
#define ClassImp(name)
Definition: Rtypes.h:279
virtual void SetDecimals(Bool_t dot=kTRUE)
Set the decimals flag.
Definition: TGaxis.cxx:2241
TAxis helper class used to store the modified labels.
Definition: TAxisModLab.h:27
double Double_t
Definition: RtypesCore.h:55
Double_t GetAngle()
Definition: TAxisModLab.h:47
Float_t fTitleOffset
Offset of title wrt axis.
Definition: TGaxis.h:40
static RooMathCoreReg dummy
TString fFunctionName
Name of mapping function pointed by fFunction.
Definition: TGaxis.h:50
Double_t y[n]
Definition: legend1.C:17
virtual const char * GetTitle() const
Returns title of object.
Definition: TGaxis.h:93
virtual Double_t GetXmax() const
Definition: TF1.h:389
void SetLabelFont(Int_t labelfont)
Definition: TGaxis.h:112
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:567
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
Definition: TRolke.cxx:630
virtual Float_t GetTitleSize() const
Definition: TAttAxis.h:48
Float_t fTextSize
Text size.
Definition: TAttText.h:28
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition: TAttLine.h:48
Bool_t IsNull() const
Definition: TString.h:387
static Int_t GetMaxDigits()
Static function returning fgMaxDigits (See SetMaxDigits).
Definition: TGaxis.cxx:683
Mother of all ROOT objects.
Definition: TObject.h:37
TGaxis()
TGaxis default constructor.
Definition: TGaxis.cxx:466
virtual void Paint(Option_t *chopt="")
Draw this axis with its current attributes.
Definition: TGaxis.cxx:721
Float_t GetTitleOffset() const
Definition: TGaxis.h:89
TString GetText()
Definition: TAxisModLab.h:52
virtual void Rotate(Double_t X, Double_t Y, Double_t CFI, Double_t SFI, Double_t XT, Double_t YT, Double_t &U, Double_t &V)
Internal method to rotate axis coordinates.
Definition: TGaxis.cxx:2165
virtual void Add(TObject *obj)
Definition: TList.h:81
static Int_t SavedTextColor
Definition: TGaxis.cxx:2334
Float_t fLabelSize
Size of labels in NDC.
Definition: TGaxis.h:39
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:202
1-Dim function class
Definition: TF1.h:149
static Int_t SavedTextFont
Definition: TGaxis.cxx:2335
void SetTimeOffset(Double_t toffset, Option_t *option="local")
Change the time offset. If option = "gmt", set display mode to GMT.
Definition: TGaxis.cxx:2482
Double_t Sin(Double_t)
Definition: TMath.h:421
void SetTickSize(Float_t ticksize)
Definition: TGaxis.h:125
Int_t fLabelColor
Color for labels.
Definition: TGaxis.h:43
Float_t GetLabelOffset() const
Definition: TGaxis.h:87
#define snprintf
Definition: civetweb.c:822
void SetLabelColor(Int_t labelcolor)
Definition: TGaxis.h:111
#define gPad
Definition: TVirtualPad.h:289
Double_t fY2
Y of 2nd point.
Definition: TLine.h:39
virtual void SetTextColor(Color_t tcolor=1)
Set the text color.
Definition: TAttText.h:49
float type_of_call hi(const int &, const int &)
void ResetBit(UInt_t f)
Definition: TObject.h:156
static Double_t SavedTextAngle
Change the label attributes of label number i. If needed.
Definition: TGaxis.cxx:2331
Width_t GetGridWidth() const
Definition: TStyle.h:225
Definition: first.py:1
Double_t Sqrt(Double_t x)
Definition: TMath.h:464
virtual void SetTextSize(Float_t tsize=1)
Set the text size.
Definition: TAttText.h:52
Float_t fTitleSize
Size of title in NDC.
Definition: TGaxis.h:41
virtual Int_t GetSize() const
Definition: TCollection.h:95
double exp(double)
void SetTitleOffset(Float_t titleoffset=1)
Definition: TGaxis.h:131
const Bool_t kTRUE
Definition: Rtypes.h:91
Double_t fWmax
Highest value on the axis.
Definition: TGaxis.h:35
Float_t GetTickSize() const
Definition: TGaxis.h:98
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition: TNamed.cxx:155
Color_t fTextColor
Text color.
Definition: TAttText.h:30
Int_t GetStripDecimals() const
Definition: TStyle.h:268
Int_t GetFont()
Definition: TAxisModLab.h:51
virtual void SetNoExponent(Bool_t noExponent=kTRUE)
Set the NoExponent flag.
Definition: TGaxis.cxx:2414
virtual void CenterTitle(Bool_t center=kTRUE)
If center = kTRUE axis title will be centered. The default is right adjusted.
Definition: TGaxis.cxx:643
Float_t fTickSize
Size of primary tick mark in NDC.
Definition: TGaxis.h:37
char name[80]
Definition: TGX11.cxx:109
static Float_t fXAxisExpXOffset
! Exponent X offset for the X axis
Definition: TGaxis.h:56
virtual Style_t GetLabelFont() const
Definition: TAttAxis.h:44
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0
void LabelsLimits(const char *label, Int_t &first, Int_t &last)
Internal method to find first and last character of a label.
Definition: TGaxis.cxx:2153
TLine()
Line default constructor.
Definition: TLine.cxx:34
Int_t GetLabelColor() const
Definition: TGaxis.h:85
Int_t fNdiv
Number of divisions.
Definition: TGaxis.h:42
Short_t fTextAlign
Text alignment.
Definition: TAttText.h:29
void ChangeLabel(Int_t labNum=0, Double_t labAngle=-1., Double_t labSize=-1., Int_t labAlign=-1, Int_t labColor=-1, Int_t labFont=-1, TString labText="")
Define new text attributes for the label number "labNum".
Definition: TGaxis.cxx:2302
Float_t GetTitleSize() const
Definition: TGaxis.h:90
const char * Data() const
Definition: TString.h:349
Double_t GetTimeOffset() const
Definition: TStyle.h:269