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