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