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