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