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