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