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