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