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