Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGaxis.cxx
Go to the documentation of this file.
1// @(#)root/graf:$Id$
2// Author: Rene Brun, Olivier Couet 12/12/94
3
4/*************************************************************************
5 * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12#include <cstdlib>
13#include <cstring>
14#include <ctime>
15#include <cmath>
16#include <iostream>
17
18#include "TROOT.h"
19#include "TBuffer.h"
20#include "TGaxis.h"
21#include "TAxisModLab.h"
22#include "TVirtualPad.h"
23#include "TVirtualX.h"
24#include "TLine.h"
25#include "TLatex.h"
26#include "TStyle.h"
27#include "TF1.h"
28#include "TAxis.h"
29#include "THashList.h"
30#include "TObject.h"
31#include "TMath.h"
32#include "THLimitsFinder.h"
33#include "TColor.h"
34#include "TTimeStamp.h"
35#include "strlcpy.h"
36#include "snprintf.h"
37
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) {
808 TAttText::operator=(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 optionNoopt = 1;
1084 ndiv = fAxis->GetLast()-fAxis->GetFirst()+1;
1085 }
1086 TList *ml = fAxis->GetModifiedLabels();
1087 if (ml) {
1088 fModLabs = ml;
1090 } else {
1091 fModLabs = 0;
1092 fNModLabs = 0;
1093 }
1094 }
1095 if (ndiv < 0) {
1096 Error(where, "Invalid number of divisions: %d",ndiv);
1097 return;
1098 }
1099
1100// Set the grid length
1101
1102 if (optionGrid) {
1103 if (gridlength == 0) gridlength = 0.8;
1104 linegrid = new TLine();
1105 linegrid->SetLineColor(gStyle->GetGridColor());
1106 if (linegrid->GetLineColor() == 0) linegrid->SetLineColor(GetLineColor());
1107 linegrid->SetLineStyle(gStyle->GetGridStyle());
1108 linegrid->SetLineWidth(gStyle->GetGridWidth());
1109 }
1110
1111// No labels if the axis label offset is big.
1112// In that case the labels are not visible anyway.
1113
1114 if (GetLabelOffset() > 1.1 ) optionUnlab = 1;
1115
1116// Determine time format
1117
1118 Int_t idF = fTimeFormat.Index("%F");
1119 if (idF>=0) {
1120 timeformat = fTimeFormat(0,idF);
1121 } else {
1122 timeformat = fTimeFormat;
1123 }
1124
1125 //GMT option
1126 if (fTimeFormat.Index("GMT")>=0) optionTime =2;
1127
1128 // Determine the time offset and correct for time offset not being integer.
1129 Double_t timeoffset =0;
1130 if (optionTime) {
1131 if (idF>=0) {
1132 Int_t lnF = fTimeFormat.Length();
1133 TString stringtimeoffset = fTimeFormat(idF+2,lnF);
1134 Int_t year, mm, dd, hh, mi, ss;
1135 if (sscanf(stringtimeoffset.Data(), "%d-%d-%d %d:%d:%d", &year, &mm, &dd, &hh, &mi, &ss) == 6) {
1136 //Get time offset in seconds since EPOCH:
1137 struct tm tp;
1138 tp.tm_year = year-1900;
1139 tp.tm_mon = mm-1;
1140 tp.tm_mday = dd;
1141 tp.tm_hour = hh;
1142 tp.tm_min = mi;
1143 tp.tm_sec = ss;
1144 tp.tm_isdst = 0; //no DST for UTC (and forced to 0 in MktimeFromUTC function)
1145 timeoffset = TTimeStamp::MktimeFromUTC(&tp);
1146
1147 // Add the time offset's decimal part if it is there
1148 Int_t ids = stringtimeoffset.Index("s");
1149 if (ids >= 0) {
1150 Float_t dp;
1151 Int_t lns = stringtimeoffset.Length();
1152 TString sdp = stringtimeoffset(ids+1,lns);
1153 sscanf(sdp.Data(),"%g",&dp);
1154 timeoffset += dp;
1155 }
1156 } else {
1157 Error(where, "Time offset has not the right format");
1158 }
1159 } else {
1160 timeoffset = gStyle->GetTimeOffset();
1161 }
1162 wmin += timeoffset - (int)(timeoffset);
1163 wmax += timeoffset - (int)(timeoffset);
1164
1165 // correct for time offset at a good limit (min, hour, day, month, year)
1166 struct tm* tp0;
1167 time_t timetp = (time_t)((Long_t)(timeoffset));
1168 Double_t range = wmax - wmin;
1169 Long_t rangeBase = 60;
1170 if (range>60) rangeBase = 60*20; // minutes
1171 if (range>3600) rangeBase = 3600*20; // hours
1172 if (range>86400) rangeBase = 86400*20; // days
1173 if (range>2419200) rangeBase = 31556736; // months (average # days)
1174 rangeOffset = (Double_t) ((Long_t)(timeoffset)%rangeBase);
1175 if (range>31536000) {
1176 tp0 = gmtime(&timetp);
1177 tp0->tm_mon = 0;
1178 tp0->tm_mday = 1;
1179 tp0->tm_hour = 0;
1180 tp0->tm_min = 0;
1181 tp0->tm_sec = 0;
1182 tp0->tm_isdst = 1; // daylight saving time is on.
1183 rangeBase = (timetp-mktime(tp0)); // years
1184 rangeOffset = (Double_t) (rangeBase);
1185 }
1186 wmax += rangeOffset;
1187 wmin += rangeOffset;
1188 }
1189
1190// Determine number of divisions 1, 2 and 3 and the maximum digits for this axis
1191 n1a = (ndiv%100);
1192 n2a = (ndiv%10000 - n1a)/100;
1193 n3a = (ndiv%1000000 - n2a -n1a)/10000;
1194 nn3 = TMath::Max(n3a,1);
1195 nn2 = TMath::Max(n2a,1)*nn3;
1196 nn1 = TMath::Max(n1a,1)*nn2+1;
1197 nticks = nn1;
1198 maxDigits = (ndiv/1000000);
1199 if (maxDigits==0) maxDigits = fgMaxDigits;
1200
1201// Axis bining optimisation is ignored if:
1202// - the first and the last label are equal
1203// - the number of divisions is 0
1204// - less than 1 primary division is requested
1205// - logarithmic scale is requested
1206
1207 if (wmin == wmax || ndiv == 0 || n1a <= 1 || optionLog) {
1208 optionNoopt = 1;
1209 optionInt = 0;
1210 }
1211
1212// Axis bining optimisation
1213 if ( (wmax-wmin) < 1 && optionInt) {
1214 Error(where, "option I not available");
1215 optionInt = 0;
1216 }
1217 if (!optionNoopt || optionInt ) {
1218
1219// Primary divisions optimisation
1220// When integer labelling is required, Optimize is invoked first
1221// and only if the result is not an integer labelling, AdjustBinSize is invoked.
1222
1223 THLimitsFinder::Optimize(wmin,wmax,n1a,binLow,binHigh,nbins,binWidth,fChopt.Data());
1224 if (optionInt) {
1225 if (binLow != Double_t(int(binLow)) || binWidth != Double_t(int(binWidth))) {
1226 AdjustBinSize(wmin,wmax,n1a,binLow,binHigh,nbins,binWidth);
1227 }
1228 }
1229 if ((wmin-binLow) > epsilon) { binLow += binWidth; nbins--; }
1230 if ((binHigh-wmax) > epsilon) { binHigh -= binWidth; nbins--; }
1231 if (xmax == xmin) {
1232 rtyw = (ymax-ymin)/(wmax-wmin);
1233 xxmin = xmin;
1234 xxmax = xmax;
1235 yymin = rtyw*(binLow-wmin) + ymin;
1236 yymax = rtyw*(binHigh-wmin) + ymin;
1237 }
1238 else {
1239 rtxw = (xmax-xmin)/(wmax-wmin);
1240 xxmin = rtxw*(binLow-wmin) + xmin;
1241 xxmax = rtxw*(binHigh-wmin) + xmin;
1242 if (ymax == ymin) {
1243 yymin = ymin;
1244 yymax = ymax;
1245 }
1246 else {
1247 alfa = (ymax-ymin)/(xmax-xmin);
1248 beta = (ymin*xmax-ymax*xmin)/(xmax-xmin);
1249 yymin = alfa*xxmin + beta;
1250 yymax = alfa*xxmax + beta;
1251 }
1252 }
1253 if (fFunction) {
1254 yymin = ymin;
1255 yymax = ymax;
1256 xxmin = xmin;
1257 xxmax = xmax;
1258 } else {
1259 wmin = binLow;
1260 wmax = binHigh;
1261 }
1262
1263// Secondary divisions optimisation
1264 nb2 = n2a;
1265 if (!optionNoopt && n2a > 1 && binWidth > 0) {
1266 THLimitsFinder::Optimize(wmin,wmin+binWidth,n2a,binLow2,binHigh2,nb2,binWidth2,fChopt.Data());
1267 }
1268
1269// Tertiary divisions optimisation
1270 nb3 = n3a;
1271 if (!optionNoopt && n3a > 1 && binWidth2 > 0) {
1272 THLimitsFinder::Optimize(binLow2,binLow2+binWidth2,n3a,binLow3,binHigh3,nb3,binWidth3,fChopt.Data());
1273 }
1274 n1aold = n1a;
1275 nn1old = nn1;
1276 n1a = nbins;
1277 nn3 = TMath::Max(nb3,1);
1278 nn2 = TMath::Max(nb2,1)*nn3;
1279 nn1 = TMath::Max(n1a,1)*nn2+1;
1280 nticks = nn1;
1281 }
1282
1283// Coordinates are normalized
1284
1285 ratio1 = 1/(rwxmax-rwxmin);
1286 ratio2 = 1/(rwymax-rwymin);
1287 x0 = ratio1*(xmin-rwxmin);
1288 x1 = ratio1*(xmax-rwxmin);
1289 y0 = ratio2*(ymin-rwymin);
1290 y1 = ratio2*(ymax-rwymin);
1291 if (!optionNoopt || optionInt ) {
1292 xx0 = ratio1*(xxmin-rwxmin);
1293 xx1 = ratio1*(xxmax-rwxmin);
1294 yy0 = ratio2*(yymin-rwymin);
1295 yy1 = ratio2*(yymax-rwymin);
1296 }
1297
1298 if ((x0 == x1) && (y0 == y1)) {
1299 Error(where, "length of axis is 0");
1300 return;
1301 }
1302
1303// Title offset. If 0 it is automatically computed
1304 Double_t toffset = GetTitleOffset();
1305 Bool_t autotoff = kFALSE;
1306 if (toffset==0 && x1 == x0) autotoff = kTRUE;
1307
1308// Return wmin, wmax and the number of primary divisions
1309 if (optionX) {
1310 ndiv = n1a;
1311 return;
1312 }
1313
1314 TLatex *textaxis = new TLatex();
1315 SetLineStyle(1); // axis line style
1316 Int_t TitleColor = GetTextColor();
1317 Int_t TitleFont = GetTextFont();
1318
1319 if (!gPad->IsBatch()) {
1320 gVirtualX->GetCharacterUp(chupxvsav, chupyvsav);
1321 gVirtualX->SetClipOFF(gPad->GetCanvasID());
1322 }
1323
1324// Compute length of axis
1325 axis_length = TMath::Sqrt((x1-x0)*(x1-x0)+(y1-y0)*(y1-y0));
1326 if (axis_length == 0) {
1327 Error(where, "length of axis is 0");
1328 goto L210;
1329 }
1330 if (!optionNoopt || optionInt) {
1331 axis_lengthN = TMath::Sqrt((xx1-xx0)*(xx1-xx0)+(yy1-yy0)*(yy1-yy0));
1332 axis_length0 = TMath::Sqrt((xx0-x0)*(xx0-x0)+(yy0-y0)*(yy0-y0));
1333 axis_length1 = TMath::Sqrt((x1-xx1)*(x1-xx1)+(y1-yy1)*(y1-yy1));
1334 if (axis_lengthN < epsilon) {
1335 optionNoopt = 1;
1336 optionInt = 0;
1337 wmin = rwmi;
1338 wmax = rwma;
1339 n1a = n1aold;
1340 nn1 = nn1old;
1341 nticks = nn1;
1342 if (optionTime) {
1343 wmin += timeoffset - (int)(timeoffset) + rangeOffset;
1344 wmax += timeoffset - (int)(timeoffset) + rangeOffset;
1345 }
1346 }
1347 }
1348
1349 if (x0 == x1) {
1350 if (y1>=y0) phi = 0.5*kPI;
1351 else phi = 1.5*kPI;
1352 phil = phi;
1353 } else {
1354 phi = TMath::ATan2((y1-y0),(x1-x0));
1355 Int_t px0 = gPad->UtoPixel(x0);
1356 Int_t py0 = gPad->VtoPixel(y0);
1357 Int_t px1 = gPad->UtoPixel(x1);
1358 Int_t py1 = gPad->VtoPixel(y1);
1359 if (x0 < x1) phil = TMath::ATan2(Double_t(py0-py1), Double_t(px1-px0));
1360 else phil = TMath::ATan2(Double_t(py1-py0), Double_t(px0-px1));
1361 }
1362 cosphi = TMath::Cos(phi);
1363 sinphi = TMath::Sin(phi);
1364 acosphi = TMath::Abs(cosphi);
1365 asinphi = TMath::Abs(sinphi);
1366 if (acosphi <= epsilon) { acosphi = 0; cosphi = 0; }
1367 if (asinphi <= epsilon) { asinphi = 0; sinphi = 0; }
1368
1369// mside positive, tick marks on positive side
1370// mside negative, tick marks on negative side
1371// mside zero, tick marks on both sides
1372// Default is positive except for vertical axis
1373
1374 mside=1;
1375 if (x0 == x1 && y1 > y0) mside = -1;
1376 if (optionPlus) mside = 1;
1377 if (optionMinus) mside = -1;
1378 if (optionPlus && optionMinus) mside = 0;
1379 xmside = mside;
1380 lside = -mside;
1381 if (optionEqual) lside = mside;
1382 if (optionPlus && optionMinus) {
1383 lside = -1;
1384 if (optionEqual) lside=1;
1385 }
1386 xlside = lside;
1387
1388// Tick marks size
1389 if(xmside >= 0) tick_side = 1;
1390 else tick_side = -1;
1391 if (optionSize) atick[0] = tick_side*axis_length*fTickSize;
1392 else atick[0] = tick_side*axis_length*0.03;
1393
1394 atick[1] = 0.5*atick[0];
1395 atick[2] = 0.5*atick[1];
1396
1397// Set the side of the grid
1398 if ((x0 == x1) && (y1 > y0)) grid_side =-1;
1399 else grid_side = 1;
1400
1401// Compute Values if Function is given
1402 if(fFunction) {
1403 rwmi = fFunction->Eval(wmin);
1404 rwma = fFunction->Eval(wmax);
1405 if(rwmi > rwma) {
1406 Double_t t = rwma;
1407 rwma = rwmi;
1408 rwmi = t;
1409 }
1410 }
1411
1412// Draw the axis if needed...
1413 if (!optionBlank) {
1414 xpl1 = x0;
1415 xpl2 = x1;
1416 ypl1 = y0;
1417 ypl2 = y1;
1418 PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1419 }
1420
1421// No bining
1422
1423 if (ndiv == 0)goto L210;
1424 if (wmin == wmax) {
1425 Error(where, "wmin (%f) == wmax (%f)", wmin, wmax);
1426 goto L210;
1427 }
1428
1429// Labels preparation:
1430// Get character height
1431// Compute the labels orientation in case of overlaps
1432// (with alphanumeric labels for horizontal axis).
1433
1434 charheight = GetLabelSize();
1435 if (optionText && GetLabelFont()%10 != 3) charheight *= 0.66666;
1436 textaxis->SetTextFont(GetLabelFont());
1437 if ((GetLabelFont()%10 < 2) && optionLog) // force TLatex mode in PaintLatex
1438 textaxis->SetTextFont((Int_t)(GetLabelFont()/10)*10+2);
1439 textaxis->SetTextColor(GetLabelColor());
1440 textaxis->SetTextSize (charheight);
1441 textaxis->SetTextAngle(GetTextAngle());
1442 if (GetLabelFont()%10 > 2) {
1443 charheight /= padh;
1444 }
1445 if (!optionUp && !optionDown && !optionY && !optionUnlab) {
1446 if (!drawGridOnly && optionText && ((ymin == ymax) || (xmin == xmax))) {
1447 textaxis->SetTextAlign(32);
1448 optionText = 2;
1449 Int_t nl = fAxis->GetLast()-fAxis->GetFirst()+1;
1450 Double_t angle = 0;
1451 for (i=fAxis->GetFirst(); i<=fAxis->GetLast(); i++) {
1452 textaxis->SetText(0,0,fAxis->GetBinLabel(i));
1453 if (textaxis->GetXsize() < (xmax-xmin)/nl) continue;
1454 angle = -20;
1455 break;
1456 }
1457 for (i=fAxis->GetFirst(); i<=fAxis->GetLast(); i++) {
1458 if ((!strcmp(fAxis->GetName(),"xaxis") && !gPad->TestBit(kHori))
1459 ||(!strcmp(fAxis->GetName(),"yaxis") && gPad->TestBit(kHori))) {
1460 if (nl > 50) angle = 90;
1461 if (fAxis->TestBit(TAxis::kLabelsHori)) angle = 0;
1462 if (fAxis->TestBit(TAxis::kLabelsVert)) angle = 90;
1463 if (fAxis->TestBit(TAxis::kLabelsUp)) angle = 20;
1464 if (fAxis->TestBit(TAxis::kLabelsDown)) angle =-20;
1465 if (angle == 0) textaxis->SetTextAlign(23);
1466 if (angle == -20) textaxis->SetTextAlign(12);
1467 textaxis->SetTextAngle(angle);
1468 Double_t s = -3;
1469 if (ymin == gPad->GetUymax()) {
1470 if (angle == 0) textaxis->SetTextAlign(21);
1471 s = 3;
1472 }
1473 strncpy(chtemp, fAxis->GetBinLabel(i), 255);
1474 if (fNModLabs) ChangeLabelAttributes(i, fAxis->GetLabels()->GetSize()-1, textaxis, chtemp);
1475 textaxis->PaintLatex(fAxis->GetBinCenter(i),
1476 ymin + s*fAxis->GetLabelOffset()*(gPad->GetUymax()-gPad->GetUymin()),
1477 textaxis->GetTextAngle(),
1478 textaxis->GetTextSize(),
1479 chtemp);
1480 if (fNModLabs) ResetLabelAttributes(textaxis);
1481 } else if ((!strcmp(fAxis->GetName(),"yaxis") && !gPad->TestBit(kHori))
1482 || (!strcmp(fAxis->GetName(),"xaxis") && gPad->TestBit(kHori))) {
1483 Double_t s = -3;
1484 if (xmin == gPad->GetUxmax()) {
1485 textaxis->SetTextAlign(12);
1486 s = 3;
1487 }
1488 if (autotoff) {
1489 UInt_t w,h;
1490 textaxis->SetText(0.,0., fAxis->GetBinLabel(i));
1491 textaxis->GetBoundingBox(w,h);
1492 double scale=gPad->GetWw()*gPad->GetWNDC();
1493 if (scale>0.0) toffset = TMath::Max(toffset,(double)w/scale);
1494 }
1495 strncpy(chtemp, fAxis->GetBinLabel(i), 255);
1496 if (fNModLabs) ChangeLabelAttributes(i, fAxis->GetLabels()->GetSize()-1, textaxis, chtemp);
1497 textaxis->PaintLatex(xmin + s*fAxis->GetLabelOffset()*(gPad->GetUxmax()-gPad->GetUxmin()),
1498 fAxis->GetBinCenter(i),
1499 0,
1500 textaxis->GetTextSize(),
1501 chtemp);
1502 if (fNModLabs) ResetLabelAttributes(textaxis);
1503 } else {
1504 strncpy(chtemp, fAxis->GetBinLabel(i), 255);
1505 if (fNModLabs) ChangeLabelAttributes(i, fAxis->GetLabels()->GetSize()-1, textaxis, chtemp);
1506 textaxis->PaintLatex(xmin - 3*fAxis->GetLabelOffset()*(gPad->GetUxmax()-gPad->GetUxmin()),
1507 ymin +(i-0.5)*(ymax-ymin)/nl,
1508 0,
1509 textaxis->GetTextSize(),
1510 chtemp);
1511 if (fNModLabs) ResetLabelAttributes(textaxis);
1512 }
1513 }
1514 }
1515 }
1516
1517// Now determine orientation of labels on axis
1518 if (!gPad->IsBatch()) {
1519 if (cosphi > 0) gVirtualX->SetCharacterUp(-sinphi,cosphi);
1520 else gVirtualX->SetCharacterUp(sinphi,-cosphi);
1521 if (x0 == x1) gVirtualX->SetCharacterUp(0,1);
1522 if (optionVert) gVirtualX->SetCharacterUp(0,1);
1523 if (optionPara) gVirtualX->SetCharacterUp(-sinphi,cosphi);
1524 if (optionDown) gVirtualX->SetCharacterUp(cosphi,sinphi);
1525 }
1526
1527// Now determine text alignment
1528 xalign = 2;
1529 yalign = 1;
1530 if (x0 == x1) xalign = 3;
1531 if (y0 != y1) yalign = 2;
1532 if (optionCent) xalign = 2;
1533 if (optionRight) xalign = 3;
1534 if (optionLeft) xalign = 1;
1535 if (TMath::Abs(cosphi) > 0.9) {
1536 xalign = 2;
1537 } else {
1538 if (cosphi*sinphi > 0) xalign = 1;
1539 if (cosphi*sinphi < 0) xalign = 3;
1540 }
1541 textaxis->SetTextAlign(10*xalign+yalign);
1542
1543// Position of labels in Y
1544 if (x0 == x1) {
1545 if (optionPlus && !optionMinus) {
1546 if (optionEqual) ylabel = fLabelOffset/2 + atick[0];
1547 else ylabel = -fLabelOffset;
1548 } else {
1549 ylabel = fLabelOffset;
1550 if (lside < 0) ylabel += atick[0];
1551 }
1552 } else if (y0 == y1) {
1553 if (optionMinus && !optionPlus) {
1554 if ((GetLabelFont() % 10) == 3 ) {
1555 ylabel = fLabelOffset+0.5*
1556 ((gPad->AbsPixeltoY(0)-gPad->AbsPixeltoY((Int_t)fLabelSize))/
1557 (gPad->GetY2() - gPad->GetY1()));
1558 } else {
1559 ylabel = fLabelOffset+0.5*fLabelSize;
1560 }
1561 ylabel += TMath::Abs(atick[0]);
1562 } else {
1563 ylabel = -fLabelOffset;
1564 if (mside <= 0) ylabel -= TMath::Abs(atick[0]);
1565 }
1566 if (optionLog) ylabel -= 0.5*charheight;
1567 } else {
1568 if (mside+lside >= 0) ylabel = fLabelOffset;
1569 else ylabel = -fLabelOffset;
1570 }
1571 if (optionText) ylabel /= 2;
1572
1573// Draw the linear tick marks if needed...
1574 if (!optionLog) {
1575 if (ndiv) {
1576 if (fFunction) {
1577 dxtick=(binHigh-binLow)/Double_t(nticks-1);
1578 } else {
1579 if (optionNoopt && !optionInt) dxtick=axis_length/Double_t(nticks-1);
1580 else dxtick=axis_lengthN/Double_t(nticks-1);
1581 }
1582 for (k=0;k<nticks; k++) {
1583 ltick = 2;
1584 if (k%nn3 == 0) ltick = 1;
1585 if (k%nn2 == 0) ltick = 0;
1586 if (fFunction) {
1587 Double_t xf = binLow+Double_t(k)*dxtick;
1588 Double_t zz = fFunction->Eval(xf)-rwmi;
1589 xtick = zz* axis_length / TMath::Abs(rwma-rwmi);
1590 } else {
1591 xtick = Double_t(k)*dxtick;
1592 }
1593 ytick = 0;
1594 if (!mside) ytick -= atick[ltick];
1595 if ( optionNoopt && !optionInt) {
1596 Rotate(xtick,ytick,cosphi,sinphi,x0,y0,xpl2,ypl2);
1597 Rotate(xtick,atick[ltick],cosphi,sinphi,x0,y0,xpl1,ypl1);
1598 }
1599 else {
1600 Rotate(xtick,ytick,cosphi,sinphi,xx0,yy0,xpl2,ypl2);
1601 Rotate(xtick,atick[ltick],cosphi,sinphi,xx0,yy0,xpl1,ypl1);
1602 }
1603 if (optionVert) {
1604 if ((x0 != x1) && (y0 != y1)) {
1605 if (mside) {
1606 xpl1 = xpl2;
1607 if (cosphi > 0) ypl1 = ypl2 + atick[ltick];
1608 else ypl1 = ypl2 - atick[ltick];
1609 }
1610 else {
1611 xpl1 = 0.5*(xpl1 + xpl2);
1612 xpl2 = xpl1;
1613 ypl1 = 0.5*(ypl1 + ypl2) + atick[ltick];
1614 ypl2 = 0.5*(ypl1 + ypl2) - atick[ltick];
1615 }
1616 }
1617 }
1618 if (!drawGridOnly) PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1619
1620 if (optionGrid) {
1621 if (ltick == 0) {
1622 if (optionNoopt && !optionInt) {
1623 Rotate(xtick,0,cosphi,sinphi,x0,y0 ,xpl2,ypl2);
1624 Rotate(xtick,grid_side*gridlength ,cosphi,sinphi,x0,y0 ,xpl1,ypl1);
1625 }
1626 else {
1627 Rotate(xtick,0,cosphi ,sinphi,xx0,yy0 ,xpl2,ypl2);
1628 Rotate(xtick,grid_side*gridlength ,cosphi,sinphi,xx0,yy0 ,xpl1,ypl1);
1629 }
1630 linegrid->PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1631 }
1632 }
1633 }
1634 xtick0 = 0;
1635 xtick1 = xtick;
1636
1637 if (fFunction) axis_length0 = binLow-wmin;
1638 if ((!optionNoopt || optionInt) && axis_length0) {
1639 nticks0 = Int_t(axis_length0/dxtick);
1640 if (nticks0 > 1000) nticks0 = 1000;
1641 for (k=0; k<=nticks0; k++) {
1642 ltick = 2;
1643 if (k%nn3 == 0) ltick = 1;
1644 if (k%nn2 == 0) ltick = 0;
1645 ytick0 = 0;
1646 if (!mside) ytick0 -= atick[ltick];
1647 if (fFunction) {
1648 xtick0 = (fFunction->Eval(binLow - Double_t(k)*dxtick)-rwmi)
1649 * axis_length / TMath::Abs(rwma-rwmi);
1650 }
1651 Rotate(xtick0,ytick0,cosphi,sinphi,xx0,yy0 ,xpl2,ypl2);
1652 Rotate(xtick0,atick[ltick],cosphi,sinphi,xx0,yy0 ,xpl1,ypl1);
1653 if (optionVert) {
1654 if ((x0 != x1) && (y0 != y1)) {
1655 if (mside) {
1656 xpl1 = xpl2;
1657 if (cosphi > 0) ypl1 = ypl2 + atick[ltick];
1658 else ypl1 = ypl2 - atick[ltick];
1659 }
1660 else {
1661 xpl1 = 0.5*(xpl1 + xpl2);
1662 xpl2 = xpl1;
1663 ypl1 = 0.5*(ypl1 + ypl2) + atick[ltick];
1664 ypl2 = 0.5*(ypl1 + ypl2) - atick[ltick];
1665 }
1666 }
1667 }
1668 if (!drawGridOnly) PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1669
1670 if (optionGrid) {
1671 if (ltick == 0) {
1672 Rotate(xtick0,0,cosphi,sinphi,xx0,yy0,xpl2,ypl2);
1673 Rotate(xtick0,grid_side*gridlength ,cosphi,sinphi,xx0,yy0 ,xpl1,ypl1);
1674 linegrid->PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1675 }
1676 }
1677 xtick0 -= dxtick;
1678 }
1679 }
1680
1681 if (fFunction) axis_length1 = wmax-binHigh;
1682 if ((!optionNoopt || optionInt) && axis_length1) {
1683 nticks1 = int(axis_length1/dxtick);
1684 if (nticks1 > 1000) nticks1 = 1000;
1685 for (k=0; k<=nticks1; k++) {
1686 ltick = 2;
1687 if (k%nn3 == 0) ltick = 1;
1688 if (k%nn2 == 0) ltick = 0;
1689 ytick1 = 0;
1690 if (!mside) ytick1 -= atick[ltick];
1691 if (fFunction) {
1692 xtick1 = (fFunction->Eval(binHigh + Double_t(k)*dxtick)-rwmi)
1693 * axis_length / TMath::Abs(rwma-rwmi);
1694 }
1695 Rotate(xtick1,ytick1,cosphi,sinphi,xx0,yy0 ,xpl2,ypl2);
1696 Rotate(xtick1,atick[ltick],cosphi,sinphi,xx0,yy0 ,xpl1,ypl1);
1697 if (optionVert) {
1698 if ((x0 != x1) && (y0 != y1)) {
1699 if (mside) {
1700 xpl1 = xpl2;
1701 if (cosphi > 0) ypl1 = ypl2 + atick[ltick];
1702 else ypl1 = ypl2 - atick[ltick];
1703 }
1704 else {
1705 xpl1 = 0.5*(xpl1 + xpl2);
1706 xpl2 = xpl1;
1707 ypl1 = 0.5*(ypl1 + ypl2) + atick[ltick];
1708 ypl2 = 0.5*(ypl1 + ypl2) - atick[ltick];
1709 }
1710 }
1711 }
1712 if (!drawGridOnly) PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1713 if (optionGrid) {
1714 if (ltick == 0) {
1715 Rotate(xtick1,0,cosphi,sinphi,xx0,yy0 ,xpl2,ypl2);
1716 Rotate(xtick1,grid_side*gridlength,cosphi,sinphi,xx0,yy0,xpl1,ypl1);
1717 linegrid->PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1718 }
1719 }
1720 xtick1 += dxtick;
1721 }
1722 }
1723 }
1724 }
1725
1726// Draw the numeric labels if needed...
1727 if (!drawGridOnly && !optionUnlab) {
1728 if (!optionLog) {
1729 if (n1a) {
1730// Spacing of labels
1731 if ((wmin == wmax) || (ndiv == 0)) {
1732 Error(where, "wmin (%f) == wmax (%f), or ndiv == 0", wmin, wmax);
1733 goto L210;
1734 }
1735 wlabel = wmin;
1736 dwlabel = (wmax-wmin)/Double_t(n1a);
1737 if (optionNoopt && !optionInt) dxlabel = axis_length/Double_t(n1a);
1738 else dxlabel = axis_lengthN/Double_t(n1a);
1739
1740 if (!optionText && !optionTime) {
1741
1742// We have to decide what format to generate
1743// (for numeric labels only)
1744// Test the magnitude, decide format
1745 flexe = kFALSE;
1746 nexe = 0;
1747 flexpo = kFALSE;
1748 flexne = kFALSE;
1749 ww = TMath::Max(TMath::Abs(wmin),TMath::Abs(wmax));
1750
1751// First case : (wmax-wmin)/n1a less than 0.001
1752// (0.001 fgMaxDigits of 5 (fgMaxDigits) characters). Then we use x 10 n
1753// format. If af >=0 x10 n cannot be used
1754 Double_t xmicros = 0.00099;
1755 if (maxDigits) xmicros = TMath::Power(10,-maxDigits);
1756 if (!noExponent && (TMath::Abs(wmax-wmin)/Double_t(n1a)) < xmicros) {
1757 af = TMath::Log10(ww) + epsilon;
1758 if (af < 0) {
1759 flexe = kTRUE;
1760 nexe = int(af);
1761 iexe = TMath::Abs(nexe);
1762 if (iexe%3 == 1) iexe += 2;
1763 else if(iexe%3 == 2) iexe += 1;
1764 if (nexe < 0) nexe = -iexe;
1765 else nexe = iexe;
1766 wlabel = wlabel*TMath::Power(10,iexe);
1767 dwlabel = dwlabel*TMath::Power(10,iexe);
1768 if1 = maxDigits;
1769 if2 = maxDigits-2;
1770 goto L110;
1771 }
1772 }
1773 if (ww >= 1) af = TMath::Log10(ww);
1774 else af = TMath::Log10(ww*0.0001);
1775 af += epsilon;
1776 nf = Int_t(af)+1;
1777 if (!noExponent && nf > maxDigits) flexpo = kTRUE;
1778 if (!noExponent && nf < -maxDigits) flexne = kTRUE;
1779
1780// Use x 10 n format. (only powers of 3 allowed)
1781
1782 if (flexpo) {
1783 flexe = kTRUE;
1784 while (1) {
1785 nexe++;
1786 ww /= 10;
1787 wlabel /= 10;
1788 dwlabel /= 10;
1789 if (nexe%3 == 0 && ww <= TMath::Power(10,maxDigits-1)) break;
1790 }
1791 }
1792
1793 if (flexne) {
1794 flexe = kTRUE;
1795 rne = 1/TMath::Power(10,maxDigits-2);
1796 while (1) {
1797 nexe--;
1798 ww *= 10;
1799 wlabel *= 10;
1800 dwlabel *= 10;
1801 if (nexe%3 == 0 && ww >= rne) break;
1802 }
1803 }
1804
1805 na = 0;
1806 for (i=maxDigits-1; i>0; i--) {
1807 if (TMath::Abs(ww) < TMath::Power(10,i)) na = maxDigits-i;
1808 }
1809 ndyn = n1a;
1810 while (ndyn) {
1811 Double_t wdyn = TMath::Abs((wmax-wmin)/ndyn);
1812 if (wdyn <= 0.999 && na < maxDigits-2) {
1813 na++;
1814 ndyn /= 10;
1815 }
1816 else break;
1817 }
1818// if1 and if2 are the two digits defining the format used to produce the
1819// labels. The format used will be %[if1].[if2]f .
1820// if1 and if2 are positive (small) integers.
1821 if2 = na;
1822 if1 = TMath::Max(nf+na,maxDigits)+1;
1823L110:
1824 if (TMath::Min(wmin,wmax) < 0)if1 = if1+1;
1825 if1 = TMath::Min(if1,32);
1826
1827// In some cases, if1 and if2 are too small....
1828 while (dwlabel < TMath::Power(10,-if2)) {
1829 if1++;
1830 if2++;
1831 }
1832 if (if1 > 14) if1 = 14;
1833 if (if2 > 14) if2 = 14;
1834 if (if1 < 0) if1 = 0;
1835 int len = 0;
1836 if (if2 > 0) {
1837 len = snprintf(chcoded,sizeof(chcoded),"%%%d.%df",if1,if2);
1838 } else {
1839 len = snprintf(chcoded,sizeof(chcoded),"%%%d.%df",if1+1,1);
1840 }
1841 // check improbable error condition, suppress gcc9 warnings
1842 if ((len < 0) || (len >= (int) sizeof(chcoded)))
1843 strcpy(chcoded,"%7.3f");
1844 }
1845
1846// We draw labels
1847
1848 snprintf(chtemp,256,"%g",dwlabel);
1849 Int_t ndecimals = 0;
1850 if (optionDecimals) {
1851 char *dot = strchr(chtemp,'.');
1852 if (dot) {
1853 ndecimals = chtemp + strlen(chtemp) -dot;
1854 } else {
1855 char *exp;
1856 exp = strstr(chtemp,"e-");
1857 if (exp) {
1858 sscanf(&exp[2],"%d",&ndecimals);
1859 ndecimals++;
1860 }
1861 }
1862 }
1863 if (optionM) nlabels = n1a-1;
1864 else nlabels = n1a;
1865 wTimeIni = wlabel;
1866 for ( k=0; k<=nlabels; k++) {
1867 if (fFunction) {
1868 Double_t xf = binLow+Double_t(k*nn2)*dxtick;
1869 Double_t zz = fFunction->Eval(xf)-rwmi;
1870 wlabel = xf;
1871 xlabel = zz* axis_length / TMath::Abs(rwma-rwmi);
1872 } else {
1873 xlabel = dxlabel*k;
1874 }
1875 if (optionM) xlabel += 0.5*dxlabel;
1876
1877 if (!optionText && !optionTime) {
1878 snprintf(label,256,chcoded,wlabel);
1879 label[28] = 0;
1880 wlabel += dwlabel;
1881
1882 LabelsLimits(label,first,last); //Eliminate blanks
1883
1884 if (label[first] == '.') { //check if '.' is preceded by a digit
1885 strncpy(chtemp, "0",256);
1886 strlcat(chtemp, &label[first],256);
1887 strncpy(label, chtemp,256);
1888 first = 1; last = strlen(label);
1889 }
1890 if (label[first] == '-' && label[first+1] == '.') {
1891 strncpy(chtemp, "-0",256);
1892 strlcat(chtemp, &label[first+1],256);
1893 strncpy(label, chtemp, 256);
1894 first = 1; last = strlen(label);
1895 }
1896
1897// We eliminate the non significant 0 after '.'
1898 if (ndecimals) {
1899 char *adot = strchr(label,'.');
1900 if (adot) adot[ndecimals] = 0;
1901 } else {
1902 while (label[last] == '0') { label[last] = 0; last--;}
1903 }
1904
1905// We eliminate the dot, unless dot is forced.
1906 if (label[last] == '.') {
1907 if (!optionDot) { label[last] = 0; last--;}
1908 }
1909
1910// Make sure the label is not "-0"
1911 if (last-first == 1 && label[first] == '-'
1912 && label[last] == '0') {
1913 strncpy(label, "0", 256);
1914 label[last] = 0;
1915 }
1916 }
1917
1918// Generate the time labels
1919
1920 if (optionTime) {
1921 timed = wlabel + (int)(timeoffset) - rangeOffset;
1922 timelabel = (time_t)((Long_t)(timed));
1923 if (optionTime == 1) {
1924 utctis = localtime(&timelabel);
1925 } else {
1926 utctis = gmtime(&timelabel);
1927 }
1928 TString timeformattmp;
1929 if (timeformat.Length() < 220) timeformattmp = timeformat;
1930 else timeformattmp = "#splitline{Format}{too long}";
1931
1932// Appends fractional part if seconds displayed
1933 if (dwlabel<0.9) {
1934 double tmpdb;
1935 int tmplast;
1936 snprintf(label, 256, "%%S%7.5f", modf(timed,&tmpdb));
1937 tmplast = strlen(label)-1;
1938
1939// We eliminate the non significant 0 after '.'
1940 while (label[tmplast] == '0') {
1941 label[tmplast] = 0; tmplast--;
1942 }
1943
1944 timeformattmp.ReplaceAll("%S",label);
1945// replace the "0." at the beginning by "s"
1946 timeformattmp.ReplaceAll("%S0.","%Ss");
1947
1948 }
1949
1950 if (utctis != nullptr) {
1951 strftime(label, 256, timeformattmp.Data(), utctis);
1952 } else {
1953 strncpy(label, "invalid", 256);
1954 }
1955 strncpy(chtemp, &label[0], 256);
1956 first = 0; last=strlen(label)-1;
1957 wlabel = wTimeIni + (k+1)*dwlabel;
1958 }
1959
1960// We generate labels (numeric or alphanumeric).
1961
1962 if (optionNoopt && !optionInt)
1963 Rotate (xlabel,ylabel,cosphi,sinphi,x0,y0,xx,yy);
1964 else Rotate (xlabel,ylabel,cosphi,sinphi,xx0,yy0,xx,yy);
1965 if (y0 == y1 && !optionDown && !optionUp) {
1966 yy -= 0.80*charheight;
1967 }
1968 if (optionVert) {
1969 if (x0 != x1 && y0 != y1) {
1970 if (optionNoopt && !optionInt)
1971 Rotate (xlabel,0,cosphi,sinphi,x0,y0,xx,yy);
1972 else Rotate (xlabel,0,cosphi,sinphi,xx0,yy0,xx,yy);
1973 if (cosphi > 0 ) yy += ylabel;
1974 if (cosphi < 0 ) yy -= ylabel;
1975 }
1976 }
1977 if (!optionY || (x0 == x1)) {
1978 if (!optionText) {
1979 if (first > last) strncpy(chtemp, " ", 256);
1980 else strncpy(chtemp, &label[first], 255);
1981 if (fNModLabs) ChangeLabelAttributes(k+1, nlabels, textaxis, chtemp);
1982 typolabel = chtemp;
1983 if (!optionTime) typolabel.ReplaceAll("-", "#minus");
1984 if (autotoff) {
1985 UInt_t w,h;
1986 textaxis->SetText(0.,0., typolabel.Data());
1987 textaxis->GetBoundingBox(w,h);
1988 double scale=gPad->GetWw()*gPad->GetWNDC();
1989 if (scale>0.0) toffset = TMath::Max(toffset,(double)w/scale);
1990 }
1991 textaxis->PaintLatex(gPad->GetX1() + xx*(gPad->GetX2() - gPad->GetX1()),
1992 gPad->GetY1() + yy*(gPad->GetY2() - gPad->GetY1()),
1993 textaxis->GetTextAngle(),
1994 textaxis->GetTextSize(),
1995 typolabel.Data());
1996 if (fNModLabs) ResetLabelAttributes(textaxis);
1997 } else {
1998 strncpy(chtemp, fAxis->GetBinLabel(k+fAxis->GetFirst()), 255);
1999 if (fNModLabs) ChangeLabelAttributes(k+fAxis->GetFirst(), fAxis->GetLabels()->GetSize()-1, textaxis, chtemp);
2000 if (optionText == 1) textaxis->PaintLatex(gPad->GetX1() + xx*(gPad->GetX2() - gPad->GetX1()),
2001 gPad->GetY1() + yy*(gPad->GetY2() - gPad->GetY1()),
2002 0,
2003 textaxis->GetTextSize(),
2004 chtemp);
2005 if (fNModLabs) ResetLabelAttributes(textaxis);
2006 }
2007 } else {
2008
2009// Text alignment is down
2010 if (!optionText) lnlen = last-first+1;
2011 else {
2012 if (k+1 > nhilab) lnlen = 0;
2013 }
2014 for ( l=1; l<=lnlen; l++) {
2015 if (!optionText) *chtemp = label[first+l-2];
2016 else {
2017 if (lnlen == 0) strncpy(chtemp, " ", 256);
2018 else strncpy(chtemp, "1", 256);
2019 }
2020 typolabel = chtemp;
2021 typolabel.ReplaceAll("-", "#minus");
2022 textaxis->PaintLatex(gPad->GetX1() + xx*(gPad->GetX2() - gPad->GetX1()),
2023 gPad->GetY1() + yy*(gPad->GetY2() - gPad->GetY1()),
2024 0,
2025 textaxis->GetTextSize(),
2026 typolabel.Data());
2027 yy -= charheight*1.3;
2028 }
2029 }
2030 }
2031
2032// We use the format x 10 ** n
2033
2034 if (flexe && !optionText && nexe) {
2035 snprintf(label,256,"#times10^{%d}", nexe);
2036 if (x0 != x1) { xfactor = axis_length+0.1*charheight; yfactor = 0; }
2037 else { xfactor = y1-y0+0.1*charheight; yfactor = 0; }
2038 Rotate (xfactor,yfactor,cosphi,sinphi,x0,y0,xx,yy);
2039 textaxis->SetTextAlign(11);
2040 if (GetLabelFont()%10 < 2) // force TLatex mode in PaintLatex
2041 textaxis->SetTextFont((Int_t)(GetLabelFont()/10)*10+2);
2042 if (fAxis && !strcmp(fAxis->GetName(),"xaxis")) {
2043 xx = xx + fXAxisExpXOffset;
2044 yy = yy + fXAxisExpYOffset;
2045 }
2046 if (fAxis && !strcmp(fAxis->GetName(),"yaxis")) {
2047 xx = xx + fYAxisExpXOffset;
2048 yy = yy + fYAxisExpYOffset;
2049 }
2050 typolabel = label;
2051 typolabel.ReplaceAll("-", "#minus");
2052 textaxis->PaintLatex(gPad->GetX1() + xx*(gPad->GetX2() - gPad->GetX1()),
2053 gPad->GetY1() + yy*(gPad->GetY2() - gPad->GetY1()),
2054 0,
2055 textaxis->GetTextSize(),
2056 typolabel.Data());
2057 }
2058 }
2059 }
2060 }
2061
2062// Log axis
2063
2064 if (optionLog && ndiv) {
2065 UInt_t xi1=0,xi2=0,wi=0,yi1=0,yi2=0,hi=0,xl=0,xh=0;
2066 Bool_t firstintlab = kTRUE, overlap = kFALSE;
2067 if ((wmin == wmax) || (ndiv == 0)) {
2068 Error(where, "wmin (%f) == wmax (%f), or ndiv == 0", wmin, wmax);
2069 goto L210;
2070 }
2071 if (wmin <= 0) {
2072 Error(where, "negative logarithmic axis");
2073 goto L210;
2074 }
2075 if (wmax <= 0) {
2076 Error(where, "negative logarithmic axis");
2077 goto L210;
2078 }
2079 xmnlog = TMath::Log10(wmin);
2080 if (xmnlog > 0) xmnlog += 1.E-6;
2081 else xmnlog -= 1.E-6;
2082 x00 = 0;
2083 x11 = axis_length;
2084 h2 = TMath::Log10(wmax);
2085 h2sav = h2;
2086 if (h2 > 0) h2 += 1.E-6;
2087 else h2 -= 1.E-6;
2088 ih1 = int(xmnlog);
2089 ih2 = 1+int(h2);
2090 nbinin = ih2-ih1+1;
2091 axmul = (x11-x00)/(h2sav-xmnlog);
2092
2093// Plot decade and intermediate tick marks
2094 decade = ih1-2;
2095 labelnumber = ih1;
2096 if ( xmnlog > 0 && (xmnlog-Double_t(ih1) > 0) ) labelnumber++;
2097 Int_t changelablogid = 0;
2098 Int_t changelablognum = 0;
2099 for (j=1; j<=nbinin; j++) {
2100
2101// Plot decade
2102 firstintlab = kTRUE, overlap = kFALSE;
2103 decade++;
2104 if (x0 == x1 && j == 1) ylabel += charheight*0.33;
2105 if (y0 == y1 && j == 1) ylabel -= charheight*0.65;
2106 xone = x00+axmul*(Double_t(decade)-xmnlog);
2107 //the following statement is a trick to circumvent a gcc bug
2108 if (j < 0) printf("j=%d\n",j);
2109 if (x00 > xone) goto L160;
2110 if ((xone-x11)>epsilon) break;
2111 xtwo = xone;
2112 y = 0;
2113 if (!mside) y -= atick[0];
2114 Rotate(xone,y,cosphi,sinphi,x0,y0,xpl2,ypl2);
2115 Rotate(xtwo,atick[0],cosphi,sinphi,x0,y0,xpl1,ypl1);
2116 if (optionVert) {
2117 if ((x0 != x1) && (y0 != y1)) {
2118 if (mside) {
2119 xpl1=xpl2;
2120 if (cosphi > 0) ypl1 = ypl2 + atick[0];
2121 else ypl1 = ypl2 - atick[0];
2122 }
2123 else {
2124 xpl1 = 0.5*(xpl1 + xpl2);
2125 xpl2 = xpl1;
2126 ypl1 = 0.5*(ypl1 + ypl2) + atick[0];
2127 ypl2 = 0.5*(ypl1 + ypl2) - atick[0];
2128 }
2129 }
2130 }
2131 if (!drawGridOnly) PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
2132
2133 if (optionGrid) {
2134 Rotate(xone,0,cosphi,sinphi,x0,y0,xpl2,ypl2);
2135 Rotate(xone,grid_side*gridlength,cosphi,sinphi,x0,y0,xpl1,ypl1);
2136 linegrid->PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
2137 }
2138
2139 if (!drawGridOnly && !optionUnlab) {
2140
2141// We generate labels (numeric only).
2142 if (noExponent) {
2143 rlab = TMath::Power(10,labelnumber);
2144 snprintf(label,256, "%f", rlab);
2145 LabelsLimits(label,first,last);
2146 while (last > first) {
2147 if (label[last] != '0') break;
2148 label[last] = 0;
2149 last--;
2150 }
2151 if (label[last] == '.') {label[last] = 0; last--;}
2152 } else {
2153 snprintf(label,256, "%d", labelnumber);
2154 LabelsLimits(label,first,last);
2155 }
2156 Rotate (xone,ylabel,cosphi,sinphi,x0,y0,xx,yy);
2157 if ((x0 == x1) && !optionPara) {
2158 if (lside < 0) {
2159 if (mside < 0) {
2160 if (labelnumber == 0) nch=1;
2161 else nch=2;
2162 xx += nch*charheight;
2163 } else {
2164 xx += 0.25*charheight;
2165 }
2166 }
2167 xx += 0.25*charheight;
2168 }
2169 if ((y0 == y1) && !optionDown && !optionUp) {
2170 if (noExponent) yy += 0.33*charheight;
2171 }
2172 if (n1a == 0)goto L210;
2173 kmod = nbinin/n1a;
2174 if (kmod == 0) kmod=1000000;
2175 if ((nbinin <= n1a) || (j == 1) || (j == nbinin) || ((nbinin > n1a) && (j%kmod == 0))) {
2176 if (labelnumber == 0) {
2177 snprintf(chtemp,256, "1");
2178 } else if (labelnumber == 1) {
2179 snprintf(chtemp,256, "10");
2180 } else {
2181 if (noExponent) {
2182 chtemp = &label[first];
2183 } else {
2184 snprintf(chtemp,256, "10^{%d}", labelnumber);
2185 }
2186 }
2187 if (fNModLabs) {
2188 if (changelablogid == 0) changelablognum = nbinin-j;
2189 changelablogid++;
2190 ChangeLabelAttributes(changelablogid, changelablognum, textaxis, chtemp);
2191 }
2192 typolabel = chtemp;
2193 typolabel.ReplaceAll("-", "#minus");
2194 if (autotoff) {
2195 UInt_t w,h;
2196 textaxis->SetText(0.,0., typolabel.Data());
2197 textaxis->GetBoundingBox(w,h);
2198 double scale=gPad->GetWw()*gPad->GetWNDC();
2199 if (scale>0.0) toffset = TMath::Max(toffset,(double)w/scale);
2200 }
2201 textaxis->PaintLatex(gPad->GetX1() + xx*(gPad->GetX2() - gPad->GetX1()),
2202 gPad->GetY1() + yy*(gPad->GetY2() - gPad->GetY1()),
2203 0, textaxis->GetTextSize(), typolabel.Data());
2204 if (fNModLabs) ResetLabelAttributes(textaxis);
2205 }
2206 labelnumber++;
2207 }
2208L160:
2209 for (k=2;k<10;k++) {
2210
2211// Plot intermediate tick marks
2212 xone = x00+axmul*(TMath::Log10(Double_t(k))+Double_t(decade)-xmnlog);
2213 if (x00 > xone) continue;
2214 if (xone > x11) goto L200;
2215 y = 0;
2216 if (!mside) y -= atick[1];
2217 xtwo = xone;
2218 Rotate(xone,y,cosphi,sinphi,x0,y0,xpl2,ypl2);
2219 Rotate(xtwo,atick[1],cosphi,sinphi,x0,y0,xpl1,ypl1);
2220 if (optionVert) {
2221 if ((x0 != x1) && (y0 != y1)) {
2222 if (mside) {
2223 xpl1 = xpl2;
2224 if (cosphi > 0) ypl1 = ypl2 + atick[1];
2225 else ypl1 = ypl2 - atick[1];
2226 }
2227 else {
2228 xpl1 = 0.5*(xpl1+xpl2);
2229 xpl2 = xpl1;
2230 ypl1 = 0.5*(ypl1+ypl2) + atick[1];
2231 ypl2 = 0.5*(ypl1+ypl2) - atick[1];
2232 }
2233 }
2234 }
2235 idn = n1a*2;
2236 if ((nbinin <= idn) || ((nbinin > idn) && (k == 5))) {
2237 if (!drawGridOnly) PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
2238
2239// Draw the intermediate LOG labels if requested
2240
2241 if (moreLogLabels && !optionUnlab && !drawGridOnly && !overlap) {
2242 if (noExponent) {
2243 rlab = Double_t(k)*TMath::Power(10,labelnumber-1);
2244 snprintf(chtemp,256, "%g", rlab);
2245 } else {
2246 if (labelnumber-1 == 0) {
2247 snprintf(chtemp,256, "%d", k);
2248 } else if (labelnumber-1 == 1) {
2249 snprintf(chtemp,256, "%d", 10*k);
2250 } else {
2251 snprintf(chtemp,256, "%d#times10^{%d}", k, labelnumber-1);
2252 }
2253 }
2254 Rotate (xone,ylabel,cosphi,sinphi,x0,y0,xx,yy);
2255 if ((x0 == x1) && !optionPara) {
2256 if (lside < 0) {
2257 if (mside < 0) {
2258 if (labelnumber == 0) nch=1;
2259 else nch=2;
2260 xx += nch*charheight;
2261 } else {
2262 if (labelnumber >= 0) xx += 0.25*charheight;
2263 else xx += 0.50*charheight;
2264 }
2265 }
2266 xx += 0.25*charheight;
2267 }
2268 if ((y0 == y1) && !optionDown && !optionUp) {
2269 if (noExponent) yy += 0.33*charheight;
2270 }
2271 if (optionVert) {
2272 if ((x0 != x1) && (y0 != y1)) {
2273 Rotate(xone,ylabel,cosphi,sinphi,x0,y0,xx,yy);
2274 if (cosphi > 0) yy += ylabel;
2275 else yy -= ylabel;
2276 }
2277 }
2278 textaxis->SetTitle(chtemp);
2279 Double_t u = gPad->GetX1() + xx*(gPad->GetX2() - gPad->GetX1());
2280 Double_t v = gPad->GetY1() + yy*(gPad->GetY2() - gPad->GetY1());
2281 if (firstintlab) {
2282 textaxis->GetBoundingBox(wi, hi); wi=(UInt_t)(wi*1.3); hi=(UInt_t)(hi*1.3);
2283 xi1 = gPad->XtoAbsPixel(u);
2284 yi1 = gPad->YtoAbsPixel(v);
2285 firstintlab = kFALSE;
2286 typolabel = chtemp;
2287 typolabel.ReplaceAll("-", "#minus");
2288 textaxis->PaintLatex(u,v,0,textaxis->GetTextSize(),typolabel.Data());
2289 } else {
2290 xi2 = gPad->XtoAbsPixel(u);
2291 yi2 = gPad->YtoAbsPixel(v);
2292 xl = TMath::Min(xi1,xi2);
2293 xh = TMath::Max(xi1,xi2);
2294 if ((x0 == x1 && yi1-hi <= yi2) || (y0 == y1 && xl+wi >= xh)){
2295 overlap = kTRUE;
2296 } else {
2297 xi1 = xi2;
2298 yi1 = yi2;
2299 textaxis->GetBoundingBox(wi, hi); wi=(UInt_t)(wi*1.3); hi=(UInt_t)(hi*1.3);
2300 typolabel = chtemp;
2301 typolabel.ReplaceAll("-", "#minus");
2302 textaxis->PaintLatex(u,v,0,textaxis->GetTextSize(),typolabel.Data());
2303 }
2304 }
2305 }
2306
2307// Draw the intermediate LOG grid if only three decades are requested
2308 if (optionGrid && nbinin <= 5 && ndiv > 100) {
2309 Rotate(xone,0,cosphi,sinphi,x0,y0,xpl2, ypl2);
2310 Rotate(xone,grid_side*gridlength,cosphi,sinphi,x0,y0, xpl1,ypl1);
2311 linegrid->PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
2312 }
2313 } //endif ((nbinin <= idn) ||
2314 } //endfor (k=2;k<10;k++)
2315 } //endfor (j=1; j<=nbinin; j++)
2316L200:
2317 Int_t dummy = 0; if (dummy) { }
2318 } //endif (optionLog && ndiv)
2319
2320// Draw axis title if it exists
2321 if (!drawGridOnly && strlen(GetTitle())) {
2322 textaxis->SetTextSize (GetTitleSize());
2323 charheight = GetTitleSize();
2324 if ((GetTextFont() % 10) > 2) {
2325 charheight = charheight/gPad->GetWh();
2326 }
2327 if (x1 == x0) {
2328 if (autotoff) {
2329 if (toffset) ylabel = xlside*charheight+toffset;
2330 else ylabel = xlside*1.6*charheight;
2331 } else {
2332 ylabel = xlside*1.6*charheight*toffset;
2333 }
2334 } else {
2335 ylabel = xlside*1.3*charheight*toffset;
2336 }
2337 if (y1 == y0) {
2338 if (toffset == 0.) toffset = gStyle->GetTitleOffset("X");
2339 ylabel = xlside*1.6*charheight*toffset;
2340 }
2341 Double_t axispos;
2342 if (TestBit(TAxis::kCenterTitle)) axispos = 0.5*axis_length;
2343 else axispos = axis_length;
2345 if (x1 >= x0) {
2346 if (TestBit(TAxis::kCenterTitle)) textaxis->SetTextAlign(22);
2347 else textaxis->SetTextAlign(12);
2348 } else {
2349 if (TestBit(TAxis::kCenterTitle)) textaxis->SetTextAlign(22);
2350 else textaxis->SetTextAlign(32);
2351 }
2352 phil+=kPI;
2353 } else {
2354 if (x1 >= x0) {
2355 if (TestBit(TAxis::kCenterTitle)) textaxis->SetTextAlign(22);
2356 else textaxis->SetTextAlign(32);
2357 } else {
2358 if (TestBit(TAxis::kCenterTitle)) textaxis->SetTextAlign(22);
2359 else textaxis->SetTextAlign(12);
2360 }
2361 }
2362 Rotate(axispos,ylabel,cosphi,sinphi,x0,y0,xpl1,ypl1);
2363 textaxis->SetTextColor(TitleColor);
2364 textaxis->SetTextFont(TitleFont);
2365 textaxis->PaintLatex(gPad->GetX1() + xpl1*(gPad->GetX2() - gPad->GetX1()),
2366 gPad->GetY1() + ypl1*(gPad->GetY2() - gPad->GetY1()),
2367 phil*180/kPI,
2368 GetTitleSize(),
2369 GetTitle());
2370 }
2371
2372L210:
2373 if (optionGrid) delete linegrid;
2374 delete textaxis;
2375}
2376
2377////////////////////////////////////////////////////////////////////////////////
2378/// Internal method for axis labels optimisation. This method adjusts the bining
2379/// of the axis in order to have integer values for the labels.
2380///
2381/// \param[in] A1,A2 Old WMIN,WMAX
2382/// \param[out] binLow,binHigh New WMIN,WMAX
2383/// \param[in] nold Old NDIV (primary divisions)
2384/// \param[out] nbins New NDIV
2385/// \param[out] binWidth Bin width
2386
2388 ,Double_t &binLow, Double_t &binHigh, Int_t &nbins, Double_t &binWidth)
2389{
2390
2391 binWidth = TMath::Abs(A2-A1)/Double_t(nold);
2392 if (binWidth <= 1) { binWidth = 1; binLow = int(A1); }
2393 else {
2394 Int_t width = int(binWidth/5) + 1;
2395 binWidth = 5*width;
2396 binLow = int(A1/binWidth)*binWidth;
2397
2398// We determine binLow to have one tick mark at 0
2399// if there are negative labels.
2400
2401 if (A1 < 0) {
2402 for (Int_t ic=0; ic<1000; ic++) {
2403 Double_t rbl = binLow/binWidth;
2404 Int_t ibl = int(binLow/binWidth);
2405 if ( (rbl-ibl) == 0 || ic > width) { binLow -= 5; break;}
2406 }
2407 }
2408 }
2409 binHigh = int(A2);
2410 nbins = 0;
2411 Double_t xb = binLow;
2412 while (xb <= binHigh) {
2413 xb += binWidth;
2414 nbins++;
2415 }
2416 binHigh = xb - binWidth;
2417}
2418
2419////////////////////////////////////////////////////////////////////////////////
2420/// Internal method to find first and last character of a label.
2421
2422void TGaxis::LabelsLimits(const char *label, Int_t &first, Int_t &last)
2423{
2424 last = strlen(label)-1;
2425 for (Int_t i=0; i<=last; i++) {
2426 if (strchr("1234567890-+.", label[i]) ) { first = i; return; }
2427 }
2428 Error("LabelsLimits", "attempt to draw a blank label");
2429}
2430
2431////////////////////////////////////////////////////////////////////////////////
2432/// Internal method to rotate axis coordinates.
2433
2435 ,Double_t XT, Double_t YT, Double_t &U, Double_t &V)
2436{
2437 U = CFI*X-SFI*Y+XT;
2438 V = SFI*X+CFI*Y+YT;
2439}
2440
2441////////////////////////////////////////////////////////////////////////////////
2442/// Save primitive as a C++ statement(s) on output stream out
2443
2444void TGaxis::SavePrimitive(std::ostream &out, Option_t * /*= ""*/)
2445{
2446 char quote = '"';
2447 if (gROOT->ClassSaved(TGaxis::Class())) {
2448 out<<" ";
2449 } else {
2450 out<<" TGaxis *";
2451 }
2452 out<<"gaxis = new TGaxis("<<fX1<<","<<fY1<<","<<fX2<<","<<fY2
2453 <<","<<fWmin<<","<<fWmax<<","<<fNdiv<<","<<quote<<fChopt.Data()<<quote<<");"<<std::endl;
2454 out<<" gaxis->SetLabelOffset("<<GetLabelOffset()<<");"<<std::endl;
2455 out<<" gaxis->SetLabelSize("<<GetLabelSize()<<");"<<std::endl;
2456 out<<" gaxis->SetTickSize("<<GetTickSize()<<");"<<std::endl;
2457 out<<" gaxis->SetGridLength("<<GetGridLength()<<");"<<std::endl;
2458 out<<" gaxis->SetTitleOffset("<<GetTitleOffset()<<");"<<std::endl;
2459 out<<" gaxis->SetTitleSize("<<GetTitleSize()<<");"<<std::endl;
2460 out<<" gaxis->SetTitleColor("<<GetTextColor()<<");"<<std::endl;
2461 out<<" gaxis->SetTitleFont("<<GetTextFont()<<");"<<std::endl;
2462
2463 if (strlen(GetName())) {
2464 out<<" gaxis->SetName("<<quote<<GetName()<<quote<<");"<<std::endl;
2465 }
2466 if (strlen(GetTitle())) {
2467 out<<" gaxis->SetTitle("<<quote<<GetTitle()<<quote<<");"<<std::endl;
2468 }
2469
2470 if (fLabelColor != 1) {
2471 if (fLabelColor > 228) {
2473 out<<" gaxis->SetLabelColor(ci);" << std::endl;
2474 } else
2475 out<<" gaxis->SetLabelColor("<<GetLabelColor()<<");"<<std::endl;
2476 }
2477 if (fLineColor != 1) {
2478 if (fLineColor > 228) {
2480 out<<" gaxis->SetLineColor(ci);" << std::endl;
2481 } else
2482 out<<" gaxis->SetLineColor("<<GetLineColor()<<");"<<std::endl;
2483 }
2484 if (fLineStyle != 1) {
2485 out<<" gaxis->SetLineStyle("<<GetLineStyle()<<");"<<std::endl;
2486 }
2487 if (fLineWidth != 1) {
2488 out<<" gaxis->SetLineWidth("<<GetLineWidth()<<");"<<std::endl;
2489 }
2490 if (fLabelFont != 62) {
2491 out<<" gaxis->SetLabelFont("<<GetLabelFont()<<");"<<std::endl;
2492 }
2494 out<<" gaxis->SetMoreLogLabels();"<<std::endl;
2495 }
2497 out<<" gaxis->SetNoExponent();"<<std::endl;
2498 }
2499
2500 out<<" gaxis->Draw();"<<std::endl;
2501}
2502
2503////////////////////////////////////////////////////////////////////////////////
2504/// Set the decimals flag. By default, blank characters are stripped, and then the
2505/// label is correctly aligned. The dot, if last character of the string, is also
2506/// stripped, unless this option is specified. One can disable the option by
2507/// calling `axis.SetDecimals(kTRUE)`.
2508/// Note the bit is set in fBits (as opposed to fBits2 in TAxis!)
2509
2511{
2512 if (dot) SetBit(TAxis::kDecimals);
2514}
2515
2516////////////////////////////////////////////////////////////////////////////////
2517/// Specify a function to map the axis values.
2518
2519void TGaxis::SetFunction(const char *funcname)
2520{
2521 fFunctionName = funcname;
2522 if (!funcname[0]) {
2523 fFunction = 0;
2524 return;
2525 }
2526 fFunction = (TF1*)gROOT->GetFunction(funcname);
2527 if (!fFunction) {
2528 Error("SetFunction", "unknown function: %s", funcname);
2529 } else {
2530 fWmin = fFunction->GetXmin();
2531 fWmax = fFunction->GetXmax();
2532 }
2533}
2534
2535////////////////////////////////////////////////////////////////////////////////
2536/// Define new text attributes for the label number "labNum". It allows to do a
2537/// fine tuning of the labels. All the attributes can be changed, even the
2538/// label text itself.
2539///
2540/// \param[in] labNum Number of the label to be changed, negative numbers start from the end
2541/// \param[in] labAngle New angle value
2542/// \param[in] labSize New size (0 erase the label)
2543/// \param[in] labAlign New alignment value
2544/// \param[in] labColor New label color
2545/// \param[in] labFont New label font
2546/// \param[in] labText New label text
2547///
2548/// If an attribute should not be changed just give the value
2549/// "-1".The following macro gives an example:
2550///
2551/// Begin_Macro(source)
2552/// {
2553/// c1 = new TCanvas("c1","Examples of Gaxis",10,10,900,500);
2554/// c1->Range(-6,-0.1,6,0.1);
2555/// TGaxis *axis1 = new TGaxis(-5.5,0.,5.5,0.,0.0,100,510,"");
2556/// axis1->SetName("axis1");
2557/// axis1->SetTitle("Axis Title");
2558/// axis1->SetTitleSize(0.05);
2559/// axis1->SetTitleColor(kBlue);
2560/// axis1->SetTitleFont(42);
2561/// axis1->ChangeLabel(1,-1,-1,-1,2);
2562/// axis1->ChangeLabel(3,-1,0.);
2563/// axis1->ChangeLabel(5,30.,-1,0);
2564/// axis1->ChangeLabel(6,-1,-1,-1,3,-1,"6th label");
2565/// axis1->ChangeLabel(-2,-1,-1,-1,3,-1,"2nd to last label");
2566/// axis1->Draw();
2567/// }
2568/// End_Macro
2569///
2570/// If labnum=0 the list of modified labels is reset.
2571
2572void TGaxis::ChangeLabel(Int_t labNum, Double_t labAngle, Double_t labSize,
2573 Int_t labAlign, Int_t labColor, Int_t labFont,
2574 TString labText)
2575{
2576 fNModLabs++;
2577 if (!fModLabs) fModLabs = new TList();
2578
2579 // Reset the list of modified labels.
2580 if (labNum == 0) {
2581 delete fModLabs;
2582 fModLabs = 0;
2583 fNModLabs = 0;
2584 return;
2585 }
2586
2587 TAxisModLab *ml = new TAxisModLab();
2588 ml->SetLabNum(labNum);
2589 ml->SetAngle(labAngle);
2590 ml->SetSize(labSize);
2591 ml->SetAlign(labAlign);
2592 ml->SetColor(labColor);
2593 ml->SetFont(labFont);
2594 ml->SetText(labText);
2595
2596 fModLabs->Add((TObject*)ml);
2597}
2598
2599////////////////////////////////////////////////////////////////////////////////
2600/// Change the label attributes of label number i. If needed.
2601///
2602/// \param[in] i Current label number to be changed if needed
2603/// \param[in] nlabels Totals number of labels on for this axis (useful when i is counted from the end)
2604/// \param[in] t Original TLatex string holding the label to be changed
2605/// \param[in] c Text string to be drawn
2606
2612
2614{
2615 if (!fModLabs) return;
2616
2617 TIter next(fModLabs);
2618 TAxisModLab *ml;
2619 Int_t labNum;
2620 while ( (ml = (TAxisModLab*)next()) ) {
2626 labNum = ml->GetLabNum();
2627 if (labNum < 0) labNum = nlabels + labNum + 2;
2628 if (i == labNum) {
2629 if (ml->GetAngle()>=0.) t->SetTextAngle(ml->GetAngle());
2630 if (ml->GetSize()>=0.) t->SetTextSize(ml->GetSize());
2631 if (ml->GetAlign()>0) t->SetTextAlign(ml->GetAlign());
2632 if (ml->GetColor()>=0) t->SetTextColor(ml->GetColor());
2633 if (ml->GetFont()>0) t->SetTextFont(ml->GetFont());
2634 if (!(ml->GetText().IsNull())) strncpy(c, (ml->GetText()).Data(), 256);
2635 return;
2636 }
2637 }
2638}
2639
2640////////////////////////////////////////////////////////////////////////////////
2641/// Reset the label attributes to the value they have before the last call to
2642/// ChangeLabelAttributes.
2643
2645{
2651}
2652
2653////////////////////////////////////////////////////////////////////////////////
2654/// Static function to set `fgMaxDigits` for axis.`fgMaxDigits` is
2655/// the maximum number of digits permitted for the axis labels above which the
2656/// notation with 10^N is used.For example, to accept 6 digits number like 900000
2657/// on an axis call `TGaxis::SetMaxDigits(6)`. The default value is 5.
2658/// `fgMaxDigits` must be greater than 0.
2659
2661{
2662 fgMaxDigits = maxd;
2663 if (maxd < 1) fgMaxDigits = 1;
2664}
2665
2666////////////////////////////////////////////////////////////////////////////////
2667/// Change the name of the axis.
2668
2669void TGaxis::SetName(const char *name)
2670{
2671 fName = name;
2672}
2673
2674////////////////////////////////////////////////////////////////////////////////
2675/// Set the kMoreLogLabels bit flag. When this option is selected more labels are
2676/// drawn when in logarithmic scale and there is a small number of decades (less than 3).
2677/// Note that this option is automatically inherited from TAxis
2678
2680{
2681 if (more) SetBit(TAxis::kMoreLogLabels);
2683}
2684
2685////////////////////////////////////////////////////////////////////////////////
2686/// Set the NoExponent flag. By default, an exponent of the form 10^N is used
2687/// when the label values are either all very small or very large. One can disable
2688/// the exponent by calling axis.SetNoExponent(kTRUE).
2689
2691{
2692 if (noExponent) SetBit(TAxis::kNoExponent);
2694}
2695
2696////////////////////////////////////////////////////////////////////////////////
2697/// To set axis options.
2698
2700{
2701 fChopt = option;
2702}
2703
2704////////////////////////////////////////////////////////////////////////////////
2705/// Change the title of the axis.
2706
2707void TGaxis::SetTitle(const char *title)
2708{
2709 fTitle = title;
2710}
2711
2712////////////////////////////////////////////////////////////////////////////////
2713/// Change the format used for time plotting.
2714/// The format string for date and time use the same options as the one used
2715/// in the standard strftime C function, i.e. :
2716///
2717/// for date :
2718///
2719/// - `%a` abbreviated weekday name
2720/// - `%b` abbreviated month name
2721/// - `%d` day of the month (01-31)
2722/// - `%m` month (01-12)
2723/// - `%y` year without century
2724///
2725/// for time :
2726///
2727/// - `%H` hour (24-hour clock)
2728/// - `%I` hour (12-hour clock)
2729/// - `%p` local equivalent of AM or PM
2730/// - `%M` minute (00-59)
2731/// - `%S` seconds (00-61)
2732/// - `%%` %
2733
2734void TGaxis::SetTimeFormat(const char *tformat)
2735{
2736 TString timeformat = tformat;
2737
2738 if (timeformat.Index("%F")>=0 || timeformat.IsNull()) {
2739 fTimeFormat = timeformat;
2740 return;
2741 }
2742
2743 Int_t idF = fTimeFormat.Index("%F");
2744 if (idF>=0) {
2745 Int_t lnF = fTimeFormat.Length();
2746 TString stringtimeoffset = fTimeFormat(idF,lnF);
2747 fTimeFormat = tformat;
2748 fTimeFormat.Append(stringtimeoffset);
2749 } else {
2750 fTimeFormat = tformat;
2752 }
2753}
2754
2755////////////////////////////////////////////////////////////////////////////////
2756/// Change the time offset. If option = "gmt", set display mode to GMT.
2757
2759{
2760 TString opt = option;
2761 opt.ToLower();
2762
2763 char tmp[20];
2764 time_t timeoff;
2765 struct tm* utctis;
2766 Int_t idF = fTimeFormat.Index("%F");
2767 if (idF>=0) fTimeFormat.Remove(idF);
2768 fTimeFormat.Append("%F");
2769
2770 timeoff = (time_t)((Long_t)(toffset));
2771
2772 // offset is always saved in GMT to allow file transport
2773 // to different time zones
2774 utctis = gmtime(&timeoff);
2775
2776 if (utctis != nullptr) {
2777 strftime(tmp, 20,"%Y-%m-%d %H:%M:%S",utctis);
2778 fTimeFormat.Append(tmp);
2779 } else {
2780 fTimeFormat.Append("1970-01-01 00:00:00");
2781 }
2782
2783 // append the decimal part of the time offset
2784 Double_t ds = toffset-(Int_t)toffset;
2785 snprintf(tmp,20,"s%g",ds);
2786 fTimeFormat.Append(tmp);
2787
2788 // add GMT/local option
2789 if (opt.Contains("gmt")) fTimeFormat.Append(" GMT");
2790}
2791
2792////////////////////////////////////////////////////////////////////////////////
2793/// Static function to set X and Y offset of the axis 10^n notation.
2794/// It is in % of the pad size. It can be negative.
2795/// axis specifies which axis ("x","y"), default = "x"
2796/// if axis="xz" set the two axes
2797
2799{
2800 TString opt = axis;
2801 opt.ToLower();
2802
2803 if (opt.Contains("x")) {
2804 fXAxisExpXOffset = xoff;
2805 fXAxisExpYOffset = yoff;
2806 }
2807 if (opt.Contains("y")) {
2808 fYAxisExpXOffset = xoff;
2809 fYAxisExpYOffset = yoff;
2810 }
2811}
2812
2813////////////////////////////////////////////////////////////////////////////////
2814/// Stream an object of class TGaxis.
2815
2816void TGaxis::Streamer(TBuffer &R__b)
2817{
2818 if (R__b.IsReading()) {
2819 UInt_t R__s, R__c;
2820 Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
2821 if (R__v > 3) {
2822 R__b.ReadClassBuffer(TGaxis::Class(), this, R__v, R__s, R__c);
2823 return;
2824 }
2825 //====process old versions before automatic schema evolution
2826 TLine::Streamer(R__b);
2827 TAttText::Streamer(R__b);
2828 R__b >> fNdiv;
2829 R__b >> fWmin;
2830 R__b >> fWmax;
2831 R__b >> fGridLength;
2832 R__b >> fTickSize;
2833 R__b >> fLabelOffset;
2834 R__b >> fLabelSize;
2835 R__b >> fTitleOffset;
2836 R__b >> fTitleSize;
2837 R__b >> fLabelFont;
2838 if (R__v > 2) {
2839 R__b >> fLabelColor;
2840 }
2841 fChopt.Streamer(R__b);
2842 fName.Streamer(R__b);
2843 fTitle.Streamer(R__b);
2844 fTimeFormat.Streamer(R__b);
2845 if (R__v > 1) {
2846 fFunctionName.Streamer(R__b);
2847 fFunction = (TF1*)gROOT->GetFunction(fFunctionName.Data());
2848 }
2849 R__b.CheckByteCount(R__s, R__c, TGaxis::IsA());
2850 //====end of old versions
2851
2852 } else {
2853 R__b.WriteClassBuffer(TGaxis::Class(),this);
2854 }
2855}
#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]
int Int_t
Definition RtypesCore.h:45
short Version_t
Definition RtypesCore.h:65
unsigned int UInt_t
Definition RtypesCore.h:46
const Bool_t kFALSE
Definition RtypesCore.h:92
long Long_t
Definition RtypesCore.h:54
double Double_t
Definition RtypesCore.h:59
float Float_t
Definition RtypesCore.h:57
const Bool_t kTRUE
Definition RtypesCore.h:91
const char Option_t
Definition RtypesCore.h:66
#define BIT(n)
Definition Rtypes.h:85
#define ClassImp(name)
Definition Rtypes.h:364
include TDocParser_001 C image html pict1_TDocParser_001 png width
const Double_t kPI
Definition TEllipse.cxx:24
char name[80]
Definition TGX11.cxx:110
static Double_t SavedTextAngle
Change the label attributes of label number i.
Definition TGaxis.cxx:2607
const Int_t kHori
Definition TGaxis.cxx:43
static Int_t SavedTextFont
Definition TGaxis.cxx:2611
static Int_t SavedTextAlign
Definition TGaxis.cxx:2609
static Double_t SavedTextSize
Definition TGaxis.cxx:2608
static Int_t SavedTextColor
Definition TGaxis.cxx:2610
float xmin
#define hi
float ymin
float xmax
float ymax
double exp(double)
#define gROOT
Definition TROOT.h:406
R__EXTERN TStyle * gStyle
Definition TStyle.h:412
#define gPad
#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.
void SetSize(Double_t s=-1.)
Set modified label size.
void SetFont(Int_t f=-1)
Set modified label font.
Double_t GetSize()
Definition TAxisModLab.h:42
void SetText(TString t="")
Set modified label text.
Int_t GetFont()
Definition TAxisModLab.h:45
void SetAlign(Int_t a=-1)
Set modified label alignment.
Int_t GetColor()
Definition TAxisModLab.h:44
void SetLabNum(Int_t n=0)
Set modified label number.
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.
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:478
@ kTickMinus
Definition TAxis.h:60
@ kLabelsUp
Definition TAxis.h:70
@ kCenterTitle
Definition TAxis.h:62
@ kRotateTitle
Definition TAxis.h:64
@ kNoExponent
Definition TAxis.h:66
@ kMoreLogLabels
Definition TAxis.h:72
@ kTickPlus
Definition TAxis.h:59
@ kLabelsDown
Definition TAxis.h:69
@ kLabelsHori
Definition TAxis.h:67
@ kDecimals
Definition TAxis.h:58
@ kCenterLabels
Definition TAxis.h:63
@ kLabelsVert
Definition TAxis.h:68
const char * GetBinLabel(Int_t bin) const
Return label for bin.
Definition TAxis.cxx:440
Bool_t GetDecimals() const
Definition TAxis.h:116
Int_t GetLast() const
Return last bin on the axis i.e.
Definition TAxis.cxx:469
TList * GetModifiedLabels() const
Definition TAxis.h:118
virtual const char * GetTimeFormat() const
Definition TAxis.h: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:458
THashList * GetLabels() const
Definition TAxis.h:117
Buffer base class used for serializing objects.
Definition TBuffer.h:43
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
virtual Version_t ReadVersion(UInt_t *start=0, UInt_t *bcnt=0, const TClass *cl=0)=0
virtual Int_t CheckByteCount(UInt_t startpos, UInt_t bcnt, const TClass *clss)=0
Bool_t IsReading() const
Definition TBuffer.h:86
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
virtual Int_t GetSize() const
Return the capacity of the collection, i.e.
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:213
virtual Double_t GetXmax() const
Definition TF1.h:556
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:1434
virtual Double_t GetXmin() const
Definition TF1.h:552
The axis painter class.
Definition TGaxis.h:23
static Int_t fgMaxDigits
! Number of digits above which the 10>N notation is used
Definition TGaxis.h:48
virtual void 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:2690
void SetTimeFormat(const char *tformat)
Change the format used for time plotting.
Definition TGaxis.cxx:2734
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:52
TGaxis()
TGaxis default constructor.
Definition TGaxis.cxx:678
void SetTitleOffset(Float_t titleoffset=1)
Definition TGaxis.h:124
Float_t fTitleSize
Size of title in NDC.
Definition TGaxis.h:34
Float_t GetGridLength() const
Definition TGaxis.h:76
void SetLabelFont(Int_t labelfont)
Definition TGaxis.h:105
void SetTitleSize(Float_t titlesize)
Definition TGaxis.h:125
Float_t fTitleOffset
Offset of title wrt axis.
Definition TGaxis.h:33
TString fTitle
Axis title.
Definition TGaxis.h:41
Int_t fLabelFont
Font for labels.
Definition TGaxis.h:37
virtual void SetTitle(const char *title="")
Change the title of the axis.
Definition TGaxis.cxx:2707
Int_t fLabelColor
Color for labels.
Definition TGaxis.h:36
TAxis * fAxis
! Pointer to original TAxis axis (if any)
Definition TGaxis.h:45
TString fTimeFormat
Time format, ex: 09/12/99 12:34:00.
Definition TGaxis.h:42
TString fFunctionName
Name of mapping function pointed by fFunction.
Definition TGaxis.h:43
void SetLabelOffset(Float_t labeloffset)
Definition TGaxis.h:106
virtual void Rotate(Double_t X, Double_t Y, Double_t CFI, Double_t SFI, Double_t XT, Double_t YT, Double_t &U, Double_t &V)
Internal method to rotate axis coordinates.
Definition TGaxis.cxx:2434
virtual const char * GetTitle() const
Returns title of object.
Definition TGaxis.h:86
Float_t GetLabelOffset() const
Definition TGaxis.h:80
void SetTimeOffset(Double_t toffset, Option_t *option="local")
Change the time offset. If option = "gmt", set display mode to GMT.
Definition TGaxis.cxx:2758
Float_t fTickSize
Size of primary tick mark in NDC.
Definition TGaxis.h:30
static Float_t fYAxisExpXOffset
! Exponent X offset for the Y axis
Definition TGaxis.h:51
virtual void ImportAxisAttributes(TAxis *axis)
Internal method to import TAxis attributes to this TGaxis.
Definition TGaxis.cxx:904
Int_t GetLabelFont() const
Definition TGaxis.h:79
Float_t fLabelOffset
Offset of label wrt axis.
Definition TGaxis.h:31
static Float_t fXAxisExpYOffset
! Exponent Y offset for the X axis
Definition TGaxis.h:50
static void SetMaxDigits(Int_t maxd=5)
Static function to set fgMaxDigits for axis.
Definition TGaxis.cxx:2660
Float_t fLabelSize
Size of labels in NDC.
Definition TGaxis.h:32
Int_t fNdiv
Number of divisions.
Definition TGaxis.h:35
void ChangeLabelAttributes(Int_t i, Int_t nlabels, TLatex *t, char *c)
Definition TGaxis.cxx:2613
TString fChopt
Axis options.
Definition TGaxis.h:39
TGaxis & operator=(const TGaxis &)
Assignment operator.
Definition TGaxis.cxx:803
TString fName
Axis name.
Definition TGaxis.h:40
Float_t GetTitleOffset() const
Definition TGaxis.h:82
TF1 * fFunction
! Pointer to function computing axis values
Definition TGaxis.h:44
void ChangeLabel(Int_t labNum=0, Double_t labAngle=-1., Double_t labSize=-1., Int_t labAlign=-1, Int_t labColor=-1, Int_t labFont=-1, TString labText="")
Define new text attributes for the label number "labNum".
Definition TGaxis.cxx:2572
Float_t GetTitleSize() const
Definition TGaxis.h:83
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:104
virtual const char * GetName() const
Returns name of object.
Definition TGaxis.h:84
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:2387
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:29
Int_t GetLabelColor() const
Definition TGaxis.h:78
virtual void SetDecimals(Bool_t dot=kTRUE)
Set the decimals flag.
Definition TGaxis.cxx:2510
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:2444
virtual ~TGaxis()
TGaxis default destructor.
Definition TGaxis.cxx:836
TList * fModLabs
List of modified labels.
Definition TGaxis.h:46
Double_t fWmax
Highest value on the axis.
Definition TGaxis.h:28
Float_t GetTickSize() const
Definition TGaxis.h:91
virtual void SetMoreLogLabels(Bool_t more=kTRUE)
Set the kMoreLogLabels bit flag.
Definition TGaxis.cxx:2679
void SetTickSize(Float_t ticksize)
Definition TGaxis.h:118
Int_t fNModLabs
Number of modified labels.
Definition TGaxis.h:38
Double_t fWmin
Lowest value on the axis.
Definition TGaxis.h:27
void SetLabelSize(Float_t labelsize)
Definition TGaxis.h:107
void SetFunction(const char *funcname="")
Specify a function to map the axis values.
Definition TGaxis.cxx:2519
static Float_t fXAxisExpXOffset
! Exponent X offset for the X axis
Definition TGaxis.h:49
virtual void SetName(const char *name)
Change the name of the axis.
Definition TGaxis.cxx:2669
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:2798
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:2422
void ResetLabelAttributes(TLatex *t)
Reset the label attributes to the value they have before the last call to ChangeLabelAttributes.
Definition TGaxis.cxx:2644
Float_t GetLabelSize() const
Definition TGaxis.h:81
void SetOption(Option_t *option="")
To set axis options.
Definition TGaxis.cxx:2699
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:2522
void GetBoundingBox(UInt_t &w, UInt_t &h, Bool_t angle=kFALSE)
Return text size in pixels.
Definition TLatex.cxx:2553
virtual void PaintLatex(Double_t x, Double_t y, Double_t angle, Double_t size, const char *text)
Main drawing function.
Definition TLatex.cxx:2051
A simple line.
Definition TLine.h:22
Double_t fY1
Y of 1st point.
Definition TLine.h:26
Double_t fX1
X of 1st point.
Definition TLine.h:25
TLine()
Definition TLine.h:38
Double_t fX2
X of 2nd point.
Definition TLine.h:27
TLine & operator=(const TLine &src)
Assignment operator.
Definition TLine.cxx:54
Double_t fY2
Y of 2nd point.
Definition TLine.h:28
virtual void PaintLineNDC(Double_t u1, Double_t v1, Double_t u2, Double_t v2)
Draw this line with new coordinates in NDC.
Definition TLine.cxx: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:107
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition TObject.cxx:696
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:893
void ResetBit(UInt_t f)
Definition TObject.h:186
@ kCanDelete
if object in a list can be deleted
Definition TObject.h:58
Basic string class.
Definition TString.h:136
Ssiz_t Length() const
Definition TString.h:410
void ToLower()
Change string to lower-case.
Definition TString.cxx:1145
const char * Data() const
Definition TString.h:369
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition TString.h:692
Bool_t IsNull() const
Definition TString.h:407
TString & Remove(Ssiz_t pos)
Definition TString.h:673
TString & Append(const char *cs)
Definition TString.h:564
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:624
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition TString.h:639
Double_t GetTimeOffset() const
Definition TStyle.h:260
Color_t GetGridColor() const
Definition TStyle.h:213
Style_t GetGridStyle() const
Definition TStyle.h:214
Float_t GetTitleOffset(Option_t *axis="X") const
Return title offset.
Definition TStyle.cxx:1176
Width_t GetGridWidth() const
Definition TStyle.h:215
Int_t GetStripDecimals() const
Definition TStyle.h:259
virtual void SetText(Double_t x, Double_t y, const char *text)
Definition TText.h:74
static time_t MktimeFromUTC(tm_t *tmstruct)
Equivalent of standard routine "mktime" but using the assumption that tm struct is filled with UTC,...
Double_t y[n]
Definition legend1.C:17
Short_t Max(Short_t a, Short_t b)
Definition TMathBase.h:212
Double_t ATan2(Double_t y, Double_t x)
Definition TMath.h:679
Double_t Sqrt(Double_t x)
Definition TMath.h:691
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition TMath.h:735
Short_t Min(Short_t a, Short_t b)
Definition TMathBase.h:180
Double_t Cos(Double_t)
Definition TMath.h:643
constexpr Double_t Pi()
Definition TMath.h:37
Double_t Sin(Double_t)
Definition TMath.h:639
Double_t Log10(Double_t x)
Definition TMath.h:764
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