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