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