Logo ROOT   6.10/09
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  n1aold = nn1old = 0;
787  Int_t ndyn;
788  Int_t nhilab = 0;
789  Int_t idn;
790  Bool_t flexe = 0;
791  Bool_t flexpo,flexne;
792  char *label;
793  char *chtemp;
794  char *coded;
795  char chlabel[256];
796  char kchtemp[256];
797  char chcoded[8];
798  TLine *linegrid;
799  TString timeformat;
800  TString typolabel;
801  time_t timelabel;
802  Double_t timed, wTimeIni;
803  struct tm* utctis;
804  Double_t rangeOffset = 0;
805 
806  Double_t epsilon = 1e-5;
807  const Double_t kPI = TMath::Pi();
808 
809  Double_t rwmi = wmin;
810  Double_t rwma = wmax;
811  chtemp = &kchtemp[0];
812  label = &chlabel[0];
813  linegrid = 0;
814 
815  fFunction = (TF1*)gROOT->GetFunction(fFunctionName.Data());
816 
817  Bool_t noExponent = TestBit(TAxis::kNoExponent);
818 
819 // If moreLogLabels = kTRUE more Log Intermediate Labels are drawn.
820  Bool_t moreLogLabels = TestBit(TAxis::kMoreLogLabels);
821 
822 // the following parameters correspond to the pad range in NDC
823 // and the user's coordinates in the pad
824 
825  Double_t padh = gPad->GetWh()*gPad->GetAbsHNDC();
826  Double_t rwxmin = gPad->GetX1();
827  Double_t rwxmax = gPad->GetX2();
828  Double_t rwymin = gPad->GetY1();
829  Double_t rwymax = gPad->GetY2();
830 
831  if(strchr(chopt,'G')) optionLog = 1; else optionLog = 0;
832  if(strchr(chopt,'B')) optionBlank= 1; else optionBlank= 0;
833  if(strchr(chopt,'V')) optionVert = 1; else optionVert = 0;
834  if(strchr(chopt,'+')) optionPlus = 1; else optionPlus = 0;
835  if(strchr(chopt,'-')) optionMinus= 1; else optionMinus= 0;
836  if(strchr(chopt,'U')) optionUnlab= 1; else optionUnlab= 0;
837  if(strchr(chopt,'P')) optionPara = 1; else optionPara = 0;
838  if(strchr(chopt,'O')) optionDown = 1; else optionDown = 0;
839  if(strchr(chopt,'R')) optionRight= 1; else optionRight= 0;
840  if(strchr(chopt,'L')) optionLeft = 1; else optionLeft = 0;
841  if(strchr(chopt,'C')) optionCent = 1; else optionCent = 0;
842  if(strchr(chopt,'=')) optionEqual= 1; else optionEqual= 0;
843  if(strchr(chopt,'Y')) optionY = 1; else optionY = 0;
844  if(strchr(chopt,'T')) optionText = 1; else optionText = 0;
845  if(strchr(chopt,'W')) optionGrid = 1; else optionGrid = 0;
846  if(strchr(chopt,'S')) optionSize = 1; else optionSize = 0;
847  if(strchr(chopt,'N')) optionNoopt= 1; else optionNoopt= 0;
848  if(strchr(chopt,'I')) optionInt = 1; else optionInt = 0;
849  if(strchr(chopt,'M')) optionM = 1; else optionM = 0;
850  if(strchr(chopt,'0')) optionUp = 1; else optionUp = 0;
851  if(strchr(chopt,'X')) optionX = 1; else optionX = 0;
852  if(strchr(chopt,'t')) optionTime = 1; else optionTime = 0;
853  if(strchr(chopt,'.')) optionDot = 1; else optionDot = 0;
854  if (TestBit(TAxis::kTickPlus)) optionPlus = 2;
855  if (TestBit(TAxis::kTickMinus)) optionMinus = 2;
856  if (TestBit(TAxis::kCenterLabels)) optionM = 1;
857  if (TestBit(TAxis::kDecimals)) optionDecimals = 1;
858  if (!gStyle->GetStripDecimals()) optionDecimals = 1;
859  if (fAxis) {
860  if (fAxis->GetLabels()) {
861  optionM = 1;
862  optionText = 1;
863  ndiv = fAxis->GetLast()-fAxis->GetFirst()+1;
864  }
865  TList *ml = fAxis->GetModifiedLabels();
866  if (ml) {
867  fModLabs = ml;
869  } else {
870  fModLabs = 0;
871  fNModLabs = 0;
872  }
873  }
874  if (ndiv < 0) {
875  Error(where, "Invalid number of divisions: %d",ndiv);
876  return;
877  }
878 
879 // Set the grid length
880 
881  if (optionGrid) {
882  if (gridlength == 0) gridlength = 0.8;
883  linegrid = new TLine();
884  linegrid->SetLineColor(gStyle->GetGridColor());
885  if (linegrid->GetLineColor() == 0) linegrid->SetLineColor(GetLineColor());
886  linegrid->SetLineStyle(gStyle->GetGridStyle());
887  linegrid->SetLineWidth(gStyle->GetGridWidth());
888  }
889 
890 // No labels if the axis label offset is big.
891 // In that case the labels are not visible anyway.
892 
893  if (GetLabelOffset() > 1.1 ) optionUnlab = 1;
894 
895 // Determine time format
896 
897  Int_t idF = fTimeFormat.Index("%F");
898  if (idF>=0) {
899  timeformat = fTimeFormat(0,idF);
900  } else {
901  timeformat = fTimeFormat;
902  }
903 
904  //GMT option
905  if (fTimeFormat.Index("GMT")>=0) optionTime =2;
906 
907  // Determine the time offset and correct for time offset not being integer.
908  Double_t timeoffset =0;
909  if (optionTime) {
910  if (idF>=0) {
911  Int_t lnF = fTimeFormat.Length();
912  TString stringtimeoffset = fTimeFormat(idF+2,lnF);
913  Int_t year, mm, dd, hh, mi, ss;
914  if (sscanf(stringtimeoffset.Data(), "%d-%d-%d %d:%d:%d", &year, &mm, &dd, &hh, &mi, &ss) == 6) {
915  //Get time offset in seconds since EPOCH:
916  struct tm tp;
917  tp.tm_year = year-1900;
918  tp.tm_mon = mm-1;
919  tp.tm_mday = dd;
920  tp.tm_hour = hh;
921  tp.tm_min = mi;
922  tp.tm_sec = ss;
923  tp.tm_isdst = 0; //no DST for UTC (and forced to 0 in MktimeFromUTC function)
924  timeoffset = TTimeStamp::MktimeFromUTC(&tp);
925 
926  // Add the time offset's decimal part if it is there
927  Int_t ids = stringtimeoffset.Index("s");
928  if (ids >= 0) {
929  Float_t dp;
930  Int_t lns = stringtimeoffset.Length();
931  TString sdp = stringtimeoffset(ids+1,lns);
932  sscanf(sdp.Data(),"%g",&dp);
933  timeoffset += dp;
934  }
935  } else {
936  Error(where, "Time offset has not the right format");
937  }
938  } else {
939  timeoffset = gStyle->GetTimeOffset();
940  }
941  wmin += timeoffset - (int)(timeoffset);
942  wmax += timeoffset - (int)(timeoffset);
943 
944  // correct for time offset at a good limit (min, hour, day, month, year)
945  struct tm* tp0;
946  time_t timetp = (time_t)((Long_t)(timeoffset));
947  Double_t range = wmax - wmin;
948  Long_t rangeBase = 60;
949  if (range>60) rangeBase = 60*20; // minutes
950  if (range>3600) rangeBase = 3600*20; // hours
951  if (range>86400) rangeBase = 86400*20; // days
952  if (range>2419200) rangeBase = 31556736; // months (average # days)
953  rangeOffset = (Double_t) ((Long_t)(timeoffset)%rangeBase);
954  if (range>31536000) {
955  tp0 = gmtime(&timetp);
956  tp0->tm_mon = 0;
957  tp0->tm_mday = 1;
958  tp0->tm_hour = 0;
959  tp0->tm_min = 0;
960  tp0->tm_sec = 0;
961  tp0->tm_isdst = 1; // daylight saving time is on.
962  rangeBase = (timetp-mktime(tp0)); // years
963  rangeOffset = (Double_t) (rangeBase);
964  }
965  wmax += rangeOffset;
966  wmin += rangeOffset;
967  }
968 
969 // Determine number of divisions 1, 2 and 3
970  n1a = ndiv%100;
971  n2a = (ndiv%10000 - n1a)/100;
972  n3a = ndiv/10000;
973  nn3 = TMath::Max(n3a,1);
974  nn2 = TMath::Max(n2a,1)*nn3;
975  nn1 = TMath::Max(n1a,1)*nn2+1;
976  nticks= nn1;
977 
978 // Axis bining optimisation is ignored if:
979 // - the first and the last label are equal
980 // - the number of divisions is 0
981 // - less than 1 primary division is requested
982 // - logarithmic scale is requested
983 
984  if (wmin == wmax || ndiv == 0 || n1a <= 1 || optionLog) {
985  optionNoopt = 1;
986  optionInt = 0;
987  }
988 
989 // Axis bining optimisation
990  if ( (wmax-wmin) < 1 && optionInt) {
991  Error(where, "option I not available");
992  optionInt = 0;
993  }
994  if (!optionNoopt || optionInt ) {
995 
996 // Primary divisions optimisation
997 // When integer labelling is required, Optimize is invoked first
998 // and only if the result is not an integer labelling, AdjustBinSize is invoked.
999 
1000  THLimitsFinder::Optimize(wmin,wmax,n1a,binLow,binHigh,nbins,binWidth,fChopt.Data());
1001  if (optionInt) {
1002  if (binLow != Double_t(int(binLow)) || binWidth != Double_t(int(binWidth))) {
1003  AdjustBinSize(wmin,wmax,n1a,binLow,binHigh,nbins,binWidth);
1004  }
1005  }
1006  if ((wmin-binLow) > epsilon) { binLow += binWidth; nbins--; }
1007  if ((binHigh-wmax) > epsilon) { binHigh -= binWidth; nbins--; }
1008  if (xmax == xmin) {
1009  rtyw = (ymax-ymin)/(wmax-wmin);
1010  xxmin = xmin;
1011  xxmax = xmax;
1012  yymin = rtyw*(binLow-wmin) + ymin;
1013  yymax = rtyw*(binHigh-wmin) + ymin;
1014  }
1015  else {
1016  rtxw = (xmax-xmin)/(wmax-wmin);
1017  xxmin = rtxw*(binLow-wmin) + xmin;
1018  xxmax = rtxw*(binHigh-wmin) + xmin;
1019  if (ymax == ymin) {
1020  yymin = ymin;
1021  yymax = ymax;
1022  }
1023  else {
1024  alfa = (ymax-ymin)/(xmax-xmin);
1025  beta = (ymin*xmax-ymax*xmin)/(xmax-xmin);
1026  yymin = alfa*xxmin + beta;
1027  yymax = alfa*xxmax + beta;
1028  }
1029  }
1030  if (fFunction) {
1031  yymin = ymin;
1032  yymax = ymax;
1033  xxmin = xmin;
1034  xxmax = xmax;
1035  } else {
1036  wmin = binLow;
1037  wmax = binHigh;
1038  }
1039 
1040 // Secondary divisions optimisation
1041  nb2 = n2a;
1042  if (!optionNoopt && n2a > 1 && binWidth > 0) {
1043  THLimitsFinder::Optimize(wmin,wmin+binWidth,n2a,binLow2,binHigh2,nb2,binWidth2,fChopt.Data());
1044  }
1045 
1046 // Tertiary divisions optimisation
1047  nb3 = n3a;
1048  if (!optionNoopt && n3a > 1 && binWidth2 > 0) {
1049  THLimitsFinder::Optimize(binLow2,binLow2+binWidth2,n3a,binLow3,binHigh3,nb3,binWidth3,fChopt.Data());
1050  }
1051  n1aold = n1a;
1052  nn1old = nn1;
1053  n1a = nbins;
1054  nn3 = TMath::Max(nb3,1);
1055  nn2 = TMath::Max(nb2,1)*nn3;
1056  nn1 = TMath::Max(n1a,1)*nn2+1;
1057  nticks = nn1;
1058  }
1059 
1060 // Coordinates are normalized
1061 
1062  ratio1 = 1/(rwxmax-rwxmin);
1063  ratio2 = 1/(rwymax-rwymin);
1064  x0 = ratio1*(xmin-rwxmin);
1065  x1 = ratio1*(xmax-rwxmin);
1066  y0 = ratio2*(ymin-rwymin);
1067  y1 = ratio2*(ymax-rwymin);
1068  if (!optionNoopt || optionInt ) {
1069  xx0 = ratio1*(xxmin-rwxmin);
1070  xx1 = ratio1*(xxmax-rwxmin);
1071  yy0 = ratio2*(yymin-rwymin);
1072  yy1 = ratio2*(yymax-rwymin);
1073  }
1074 
1075  if ((x0 == x1) && (y0 == y1)) {
1076  Error(where, "length of axis is 0");
1077  return;
1078  }
1079 
1080 // Title offset. If 0 it is automatically computed
1081  Double_t toffset = GetTitleOffset();
1082  Bool_t autotoff = kFALSE;
1083  if (toffset==0 && x1 == x0) autotoff = kTRUE;
1084 
1085 // Return wmin, wmax and the number of primary divisions
1086  if (optionX) {
1087  ndiv = n1a;
1088  return;
1089  }
1090 
1091  Int_t maxDigits = fgMaxDigits;
1092 
1093  TLatex *textaxis = new TLatex();
1094  SetLineStyle(1); // axis line style
1095  Int_t TitleColor = GetTextColor();
1096  Int_t TitleFont = GetTextFont();
1097 
1098  if (!gPad->IsBatch()) {
1099  gVirtualX->GetCharacterUp(chupxvsav, chupyvsav);
1100  gVirtualX->SetClipOFF(gPad->GetCanvasID());
1101  }
1102 
1103 // Compute length of axis
1104  axis_length = TMath::Sqrt((x1-x0)*(x1-x0)+(y1-y0)*(y1-y0));
1105  if (axis_length == 0) {
1106  Error(where, "length of axis is 0");
1107  goto L210;
1108  }
1109  if (!optionNoopt || optionInt) {
1110  axis_lengthN = TMath::Sqrt((xx1-xx0)*(xx1-xx0)+(yy1-yy0)*(yy1-yy0));
1111  axis_length0 = TMath::Sqrt((xx0-x0)*(xx0-x0)+(yy0-y0)*(yy0-y0));
1112  axis_length1 = TMath::Sqrt((x1-xx1)*(x1-xx1)+(y1-yy1)*(y1-yy1));
1113  if (axis_lengthN < epsilon) {
1114  optionNoopt = 1;
1115  optionInt = 0;
1116  wmin = rwmi;
1117  wmax = rwma;
1118  n1a = n1aold;
1119  nn1 = nn1old;
1120  nticks = nn1;
1121  if (optionTime) {
1122  wmin += timeoffset - (int)(timeoffset) + rangeOffset;
1123  wmax += timeoffset - (int)(timeoffset) + rangeOffset;
1124  }
1125  }
1126  }
1127 
1128  if (x0 == x1) {
1129  if (y1>=y0) phi = 0.5*kPI;
1130  else phi = 1.5*kPI;
1131  phil = phi;
1132  } else {
1133  phi = TMath::ATan2((y1-y0),(x1-x0));
1134  Int_t px0 = gPad->UtoPixel(x0);
1135  Int_t py0 = gPad->VtoPixel(y0);
1136  Int_t px1 = gPad->UtoPixel(x1);
1137  Int_t py1 = gPad->VtoPixel(y1);
1138  if (x0 < x1) phil = TMath::ATan2(Double_t(py0-py1), Double_t(px1-px0));
1139  else phil = TMath::ATan2(Double_t(py1-py0), Double_t(px0-px1));
1140  }
1141  cosphi = TMath::Cos(phi);
1142  sinphi = TMath::Sin(phi);
1143  acosphi = TMath::Abs(cosphi);
1144  asinphi = TMath::Abs(sinphi);
1145  if (acosphi <= epsilon) { acosphi = 0; cosphi = 0; }
1146  if (asinphi <= epsilon) { asinphi = 0; sinphi = 0; }
1147 
1148 // mside positive, tick marks on positive side
1149 // mside negative, tick marks on negative side
1150 // mside zero, tick marks on both sides
1151 // Default is positive except for vertical axis
1152 
1153  mside=1;
1154  if (x0 == x1 && y1 > y0) mside = -1;
1155  if (optionPlus) mside = 1;
1156  if (optionMinus) mside = -1;
1157  if (optionPlus && optionMinus) mside = 0;
1158  xmside = mside;
1159  lside = -mside;
1160  if (optionEqual) lside = mside;
1161  if (optionPlus && optionMinus) {
1162  lside = -1;
1163  if (optionEqual) lside=1;
1164  }
1165  xlside = lside;
1166 
1167 // Tick marks size
1168  if(xmside >= 0) tick_side = 1;
1169  else tick_side = -1;
1170  if (optionSize) atick[0] = tick_side*axis_length*fTickSize;
1171  else atick[0] = tick_side*axis_length*0.03;
1172 
1173  atick[1] = 0.5*atick[0];
1174  atick[2] = 0.5*atick[1];
1175 
1176 // Set the side of the grid
1177  if ((x0 == x1) && (y1 > y0)) grid_side =-1;
1178  else grid_side = 1;
1179 
1180 // Compute Values if Function is given
1181  if(fFunction) {
1182  rwmi = fFunction->Eval(wmin);
1183  rwma = fFunction->Eval(wmax);
1184  if(rwmi > rwma) {
1185  Double_t t = rwma;
1186  rwma = rwmi;
1187  rwmi = t;
1188  }
1189  }
1190 
1191 // Draw the axis if needed...
1192  if (!optionBlank) {
1193  xpl1 = x0;
1194  xpl2 = x1;
1195  ypl1 = y0;
1196  ypl2 = y1;
1197  PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1198  }
1199 
1200 // No bining
1201 
1202  if (ndiv == 0)goto L210;
1203  if (wmin == wmax) {
1204  Error(where, "wmin (%f) == wmax (%f)", wmin, wmax);
1205  goto L210;
1206  }
1207 
1208 // Labels preparation:
1209 // Get character height
1210 // Compute the labels orientation in case of overlaps
1211 // (with alphanumeric labels for horizontal axis).
1212 
1213  charheight = GetLabelSize();
1214  if (optionText && GetLabelFont()%10 != 3) charheight *= 0.66666;
1215  textaxis->SetTextFont(GetLabelFont());
1216  if ((GetLabelFont()%10 < 2) && optionLog) // force TLatex mode in PaintLatex
1217  textaxis->SetTextFont((Int_t)(GetLabelFont()/10)*10+2);
1218  textaxis->SetTextColor(GetLabelColor());
1219  textaxis->SetTextSize (charheight);
1220  textaxis->SetTextAngle(GetTextAngle());
1221  if (GetLabelFont()%10 > 2) {
1222  charheight /= padh;
1223  }
1224  if (!optionUp && !optionDown && !optionY && !optionUnlab) {
1225  if (!drawGridOnly && optionText && ((ymin == ymax) || (xmin == xmax))) {
1226  textaxis->SetTextAlign(32);
1227  optionText = 2;
1228  Int_t nl = fAxis->GetLast()-fAxis->GetFirst()+1;
1229  Double_t angle = 0;
1230  for (i=fAxis->GetFirst(); i<=fAxis->GetLast(); i++) {
1231  textaxis->SetText(0,0,fAxis->GetBinLabel(i));
1232  if (textaxis->GetXsize() < (xmax-xmin)/nl) continue;
1233  angle = -20;
1234  break;
1235  }
1236  for (i=fAxis->GetFirst(); i<=fAxis->GetLast(); i++) {
1237  if ((!strcmp(fAxis->GetName(),"xaxis") && !gPad->TestBit(kHori))
1238  ||(!strcmp(fAxis->GetName(),"yaxis") && gPad->TestBit(kHori))) {
1239  if (nl > 50) angle = 90;
1240  if (fAxis->TestBit(TAxis::kLabelsHori)) angle = 0;
1241  if (fAxis->TestBit(TAxis::kLabelsVert)) angle = 90;
1242  if (fAxis->TestBit(TAxis::kLabelsUp)) angle = 20;
1243  if (fAxis->TestBit(TAxis::kLabelsDown)) angle =-20;
1244  if (angle == 0) textaxis->SetTextAlign(23);
1245  if (angle == -20) textaxis->SetTextAlign(12);
1246  Double_t s = -3;
1247  if (ymin == gPad->GetUymax()) {
1248  if (angle == 0) textaxis->SetTextAlign(21);
1249  s = 3;
1250  }
1251  textaxis->PaintLatex(fAxis->GetBinCenter(i),
1252  ymin + s*fAxis->GetLabelOffset()*(gPad->GetUymax()-gPad->GetUymin()),
1253  angle,
1254  textaxis->GetTextSize(),
1255  fAxis->GetBinLabel(i));
1256  } else if ((!strcmp(fAxis->GetName(),"yaxis") && !gPad->TestBit(kHori))
1257  || (!strcmp(fAxis->GetName(),"xaxis") && gPad->TestBit(kHori))) {
1258  Double_t s = -3;
1259  if (xmin == gPad->GetUxmax()) {
1260  textaxis->SetTextAlign(12);
1261  s = 3;
1262  }
1263  if (autotoff) {
1264  UInt_t w,h;
1265  textaxis->SetText(0.,0., fAxis->GetBinLabel(i));
1266  textaxis->GetBoundingBox(w,h);
1267  toffset = TMath::Max(toffset,(double)w/((double)gPad->GetWw()*gPad->GetWNDC()));
1268  }
1269  textaxis->PaintLatex(xmin + s*fAxis->GetLabelOffset()*(gPad->GetUxmax()-gPad->GetUxmin()),
1270  fAxis->GetBinCenter(i),
1271  0,
1272  textaxis->GetTextSize(),
1273  fAxis->GetBinLabel(i));
1274  } else {
1275  textaxis->PaintLatex(xmin - 3*fAxis->GetLabelOffset()*(gPad->GetUxmax()-gPad->GetUxmin()),
1276  ymin +(i-0.5)*(ymax-ymin)/nl,
1277  0,
1278  textaxis->GetTextSize(),
1279  fAxis->GetBinLabel(i));
1280  }
1281  }
1282  }
1283  }
1284 
1285 // Now determine orientation of labels on axis
1286  if (!gPad->IsBatch()) {
1287  if (cosphi > 0) gVirtualX->SetCharacterUp(-sinphi,cosphi);
1288  else gVirtualX->SetCharacterUp(sinphi,-cosphi);
1289  if (x0 == x1) gVirtualX->SetCharacterUp(0,1);
1290  if (optionVert) gVirtualX->SetCharacterUp(0,1);
1291  if (optionPara) gVirtualX->SetCharacterUp(-sinphi,cosphi);
1292  if (optionDown) gVirtualX->SetCharacterUp(cosphi,sinphi);
1293  }
1294 
1295 // Now determine text alignment
1296  xalign = 2;
1297  yalign = 1;
1298  if (x0 == x1) xalign = 3;
1299  if (y0 != y1) yalign = 2;
1300  if (optionCent) xalign = 2;
1301  if (optionRight) xalign = 3;
1302  if (optionLeft) xalign = 1;
1303  if (TMath::Abs(cosphi) > 0.9) {
1304  xalign = 2;
1305  } else {
1306  if (cosphi*sinphi > 0) xalign = 1;
1307  if (cosphi*sinphi < 0) xalign = 3;
1308  }
1309  textaxis->SetTextAlign(10*xalign+yalign);
1310 
1311 // Position of labels in Y
1312  if (x0 == x1) {
1313  if (optionPlus && !optionMinus) {
1314  if (optionEqual) ylabel = fLabelOffset/2 + atick[0];
1315  else ylabel = -fLabelOffset;
1316  } else {
1317  ylabel = fLabelOffset;
1318  if (lside < 0) ylabel += atick[0];
1319  }
1320  } else if (y0 == y1) {
1321  if (optionMinus && !optionPlus) {
1322  if ((GetLabelFont() % 10) == 3 ) {
1323  ylabel = fLabelOffset+0.5*
1324  ((gPad->AbsPixeltoY(0)-gPad->AbsPixeltoY((Int_t)fLabelSize))/
1325  (gPad->GetY2() - gPad->GetY1()));
1326  } else {
1327  ylabel = fLabelOffset+0.5*fLabelSize;
1328  }
1329  ylabel += TMath::Abs(atick[0]);
1330  } else {
1331  ylabel = -fLabelOffset;
1332  if (mside <= 0) ylabel -= TMath::Abs(atick[0]);
1333  }
1334  if (optionLog) ylabel -= 0.5*charheight;
1335  } else {
1336  if (mside+lside >= 0) ylabel = fLabelOffset;
1337  else ylabel = -fLabelOffset;
1338  }
1339  if (optionText) ylabel /= 2;
1340 
1341 // Draw the linear tick marks if needed...
1342  if (!optionLog) {
1343  if (ndiv) {
1344  if (fFunction) {
1345  dxtick=(binHigh-binLow)/Double_t(nticks-1);
1346  } else {
1347  if (optionNoopt && !optionInt) dxtick=axis_length/Double_t(nticks-1);
1348  else dxtick=axis_lengthN/Double_t(nticks-1);
1349  }
1350  for (k=0;k<nticks; k++) {
1351  ltick = 2;
1352  if (k%nn3 == 0) ltick = 1;
1353  if (k%nn2 == 0) ltick = 0;
1354  if (fFunction) {
1355  Double_t xf = binLow+Double_t(k)*dxtick;
1356  Double_t zz = fFunction->Eval(xf)-rwmi;
1357  xtick = zz* axis_length / TMath::Abs(rwma-rwmi);
1358  } else {
1359  xtick = Double_t(k)*dxtick;
1360  }
1361  ytick = 0;
1362  if (!mside) ytick -= atick[ltick];
1363  if ( optionNoopt && !optionInt) {
1364  Rotate(xtick,ytick,cosphi,sinphi,x0,y0,xpl2,ypl2);
1365  Rotate(xtick,atick[ltick],cosphi,sinphi,x0,y0,xpl1,ypl1);
1366  }
1367  else {
1368  Rotate(xtick,ytick,cosphi,sinphi,xx0,yy0,xpl2,ypl2);
1369  Rotate(xtick,atick[ltick],cosphi,sinphi,xx0,yy0,xpl1,ypl1);
1370  }
1371  if (optionVert) {
1372  if ((x0 != x1) && (y0 != y1)) {
1373  if (mside) {
1374  xpl1 = xpl2;
1375  if (cosphi > 0) ypl1 = ypl2 + atick[ltick];
1376  else ypl1 = ypl2 - atick[ltick];
1377  }
1378  else {
1379  xpl1 = 0.5*(xpl1 + xpl2);
1380  xpl2 = xpl1;
1381  ypl1 = 0.5*(ypl1 + ypl2) + atick[ltick];
1382  ypl2 = 0.5*(ypl1 + ypl2) - atick[ltick];
1383  }
1384  }
1385  }
1386  if (!drawGridOnly) PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1387 
1388  if (optionGrid) {
1389  if (ltick == 0) {
1390  if (optionNoopt && !optionInt) {
1391  Rotate(xtick,0,cosphi,sinphi,x0,y0 ,xpl2,ypl2);
1392  Rotate(xtick,grid_side*gridlength ,cosphi,sinphi,x0,y0 ,xpl1,ypl1);
1393  }
1394  else {
1395  Rotate(xtick,0,cosphi ,sinphi,xx0,yy0 ,xpl2,ypl2);
1396  Rotate(xtick,grid_side*gridlength ,cosphi,sinphi,xx0,yy0 ,xpl1,ypl1);
1397  }
1398  linegrid->PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1399  }
1400  }
1401  }
1402  xtick0 = 0;
1403  xtick1 = xtick;
1404 
1405  if (fFunction) axis_length0 = binLow-wmin;
1406  if ((!optionNoopt || optionInt) && axis_length0) {
1407  nticks0 = Int_t(axis_length0/dxtick);
1408  if (nticks0 > 1000) nticks0 = 1000;
1409  for (k=0; k<=nticks0; k++) {
1410  ltick = 2;
1411  if (k%nn3 == 0) ltick = 1;
1412  if (k%nn2 == 0) ltick = 0;
1413  ytick0 = 0;
1414  if (!mside) ytick0 -= atick[ltick];
1415  if (fFunction) {
1416  xtick0 = (fFunction->Eval(binLow - Double_t(k)*dxtick)-rwmi)
1417  * axis_length / TMath::Abs(rwma-rwmi);
1418  }
1419  Rotate(xtick0,ytick0,cosphi,sinphi,xx0,yy0 ,xpl2,ypl2);
1420  Rotate(xtick0,atick[ltick],cosphi,sinphi,xx0,yy0 ,xpl1,ypl1);
1421  if (optionVert) {
1422  if ((x0 != x1) && (y0 != y1)) {
1423  if (mside) {
1424  xpl1 = xpl2;
1425  if (cosphi > 0) ypl1 = ypl2 + atick[ltick];
1426  else ypl1 = ypl2 - atick[ltick];
1427  }
1428  else {
1429  xpl1 = 0.5*(xpl1 + xpl2);
1430  xpl2 = xpl1;
1431  ypl1 = 0.5*(ypl1 + ypl2) + atick[ltick];
1432  ypl2 = 0.5*(ypl1 + ypl2) - atick[ltick];
1433  }
1434  }
1435  }
1436  if (!drawGridOnly) PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1437 
1438  if (optionGrid) {
1439  if (ltick == 0) {
1440  Rotate(xtick0,0,cosphi,sinphi,xx0,yy0,xpl2,ypl2);
1441  Rotate(xtick0,grid_side*gridlength ,cosphi,sinphi,xx0,yy0 ,xpl1,ypl1);
1442  linegrid->PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1443  }
1444  }
1445  xtick0 -= dxtick;
1446  }
1447  }
1448 
1449  if (fFunction) axis_length1 = wmax-binHigh;
1450  if ((!optionNoopt || optionInt) && axis_length1) {
1451  nticks1 = int(axis_length1/dxtick);
1452  if (nticks1 > 1000) nticks1 = 1000;
1453  for (k=0; k<=nticks1; k++) {
1454  ltick = 2;
1455  if (k%nn3 == 0) ltick = 1;
1456  if (k%nn2 == 0) ltick = 0;
1457  ytick1 = 0;
1458  if (!mside) ytick1 -= atick[ltick];
1459  if (fFunction) {
1460  xtick1 = (fFunction->Eval(binHigh + Double_t(k)*dxtick)-rwmi)
1461  * axis_length / TMath::Abs(rwma-rwmi);
1462  }
1463  Rotate(xtick1,ytick1,cosphi,sinphi,xx0,yy0 ,xpl2,ypl2);
1464  Rotate(xtick1,atick[ltick],cosphi,sinphi,xx0,yy0 ,xpl1,ypl1);
1465  if (optionVert) {
1466  if ((x0 != x1) && (y0 != y1)) {
1467  if (mside) {
1468  xpl1 = xpl2;
1469  if (cosphi > 0) ypl1 = ypl2 + atick[ltick];
1470  else ypl1 = ypl2 - atick[ltick];
1471  }
1472  else {
1473  xpl1 = 0.5*(xpl1 + xpl2);
1474  xpl2 = xpl1;
1475  ypl1 = 0.5*(ypl1 + ypl2) + atick[ltick];
1476  ypl2 = 0.5*(ypl1 + ypl2) - atick[ltick];
1477  }
1478  }
1479  }
1480  if (!drawGridOnly) PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1481  if (optionGrid) {
1482  if (ltick == 0) {
1483  Rotate(xtick1,0,cosphi,sinphi,xx0,yy0 ,xpl2,ypl2);
1484  Rotate(xtick1,grid_side*gridlength,cosphi,sinphi,xx0,yy0,xpl1,ypl1);
1485  linegrid->PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1486  }
1487  }
1488  xtick1 += dxtick;
1489  }
1490  }
1491  }
1492  }
1493 
1494 // Draw the numeric labels if needed...
1495  if (!drawGridOnly && !optionUnlab) {
1496  if (!optionLog) {
1497  if (n1a) {
1498 // Spacing of labels
1499  if ((wmin == wmax) || (ndiv == 0)) {
1500  Error(where, "wmin (%f) == wmax (%f), or ndiv == 0", wmin, wmax);
1501  goto L210;
1502  }
1503  wlabel = wmin;
1504  dwlabel = (wmax-wmin)/Double_t(n1a);
1505  if (optionNoopt && !optionInt) dxlabel = axis_length/Double_t(n1a);
1506  else dxlabel = axis_lengthN/Double_t(n1a);
1507 
1508  if (!optionText && !optionTime) {
1509 
1510 // We have to decide what format to generate
1511 // (for numeric labels only)
1512 // Test the magnitude, decide format
1513  flexe = kFALSE;
1514  nexe = 0;
1515  flexpo = kFALSE;
1516  flexne = kFALSE;
1517  ww = TMath::Max(TMath::Abs(wmin),TMath::Abs(wmax));
1518 
1519 // First case : (wmax-wmin)/n1a less than 0.001
1520 // (0.001 fgMaxDigits of 5 (fgMaxDigits) characters). Then we use x 10 n
1521 // format. If af >=0 x10 n cannot be used
1522  Double_t xmicros = 0.00099;
1523  if (maxDigits) xmicros = TMath::Power(10,-maxDigits);
1524  if (!noExponent && (TMath::Abs(wmax-wmin)/Double_t(n1a)) < xmicros) {
1525  af = TMath::Log10(ww) + epsilon;
1526  if (af < 0) {
1527  flexe = kTRUE;
1528  nexe = int(af);
1529  iexe = TMath::Abs(nexe);
1530  if (iexe%3 == 1) iexe += 2;
1531  else if(iexe%3 == 2) iexe += 1;
1532  if (nexe < 0) nexe = -iexe;
1533  else nexe = iexe;
1534  wlabel = wlabel*TMath::Power(10,iexe);
1535  dwlabel = dwlabel*TMath::Power(10,iexe);
1536  if1 = maxDigits;
1537  if2 = maxDigits-2;
1538  goto L110;
1539  }
1540  }
1541  if (ww >= 1) af = TMath::Log10(ww);
1542  else af = TMath::Log10(ww*0.0001);
1543  af += epsilon;
1544  nf = Int_t(af)+1;
1545  if (!noExponent && nf > maxDigits) flexpo = kTRUE;
1546  if (!noExponent && nf < -maxDigits) flexne = kTRUE;
1547 
1548 // Use x 10 n format. (only powers of 3 allowed)
1549 
1550  if (flexpo) {
1551  flexe = kTRUE;
1552  while (1) {
1553  nexe++;
1554  ww /= 10;
1555  wlabel /= 10;
1556  dwlabel /= 10;
1557  if (nexe%3 == 0 && ww <= TMath::Power(10,maxDigits-1)) break;
1558  }
1559  }
1560 
1561  if (flexne) {
1562  flexe = kTRUE;
1563  rne = 1/TMath::Power(10,maxDigits-2);
1564  while (1) {
1565  nexe--;
1566  ww *= 10;
1567  wlabel *= 10;
1568  dwlabel *= 10;
1569  if (nexe%3 == 0 && ww >= rne) break;
1570  }
1571  }
1572 
1573  na = 0;
1574  for (i=maxDigits-1; i>0; i--) {
1575  if (TMath::Abs(ww) < TMath::Power(10,i)) na = maxDigits-i;
1576  }
1577  ndyn = n1a;
1578  while (ndyn) {
1579  Double_t wdyn = TMath::Abs((wmax-wmin)/ndyn);
1580  if (wdyn <= 0.999 && na < maxDigits-2) {
1581  na++;
1582  ndyn /= 10;
1583  }
1584  else break;
1585  }
1586 // if1 and if2 are the two digits defining the format used to produce the
1587 // labels. The format used will be %[if1].[if2]f .
1588 // if1 and if2 are positive (small) integers.
1589  if2 = na;
1590  if1 = TMath::Max(nf+na,maxDigits)+1;
1591 L110:
1592  if (TMath::Min(wmin,wmax) < 0)if1 = if1+1;
1593  if1 = TMath::Min(if1,32);
1594 
1595 // In some cases, if1 and if2 are too small....
1596  while (dwlabel < TMath::Power(10,-if2)) {
1597  if1++;
1598  if2++;
1599  }
1600  coded = &chcoded[0];
1601  if (if1 > 14) if1=14;
1602  if (if2 > 14) if2=14;
1603  if (if2>0) snprintf(coded,8,"%%%d.%df",if1,if2);
1604  else {
1605  if (if1 < -100) if1 = -100; // Silence a warning with gcc
1606  snprintf(coded,8,"%%%d.%df",if1+1,1);
1607  }
1608 
1609  }
1610 
1611 // We draw labels
1612 
1613  snprintf(chtemp,256,"%g",dwlabel);
1614  Int_t ndecimals = 0;
1615  if (optionDecimals) {
1616  char *dot = strchr(chtemp,'.');
1617  if (dot) {
1618  ndecimals = chtemp + strlen(chtemp) -dot;
1619  } else {
1620  char *exp;
1621  exp = strstr(chtemp,"e-");
1622  if (exp) {
1623  sscanf(&exp[2],"%d",&ndecimals);
1624  ndecimals++;
1625  }
1626  }
1627  }
1628  if (optionM) nlabels = n1a-1;
1629  else nlabels = n1a;
1630  wTimeIni = wlabel;
1631  for ( k=0; k<=nlabels; k++) {
1632  if (fFunction) {
1633  Double_t xf = binLow+Double_t(k*nn2)*dxtick;
1634  Double_t zz = fFunction->Eval(xf)-rwmi;
1635  wlabel = xf;
1636  xlabel = zz* axis_length / TMath::Abs(rwma-rwmi);
1637  } else {
1638  xlabel = dxlabel*k;
1639  }
1640  if (optionM) xlabel += 0.5*dxlabel;
1641 
1642  if (!optionText && !optionTime) {
1643  snprintf(label,256,&chcoded[0],wlabel);
1644  label[28] = 0;
1645  wlabel += dwlabel;
1646 
1647  LabelsLimits(label,first,last); //Eliminate blanks
1648 
1649  if (label[first] == '.') { //check if '.' is preceded by a digit
1650  strncpy(chtemp, "0",256);
1651  strlcat(chtemp, &label[first],256);
1652  strncpy(label, chtemp,256);
1653  first = 1; last = strlen(label);
1654  }
1655  if (label[first] == '-' && label[first+1] == '.') {
1656  strncpy(chtemp, "-0",256);
1657  strlcat(chtemp, &label[first+1],256);
1658  strncpy(label, chtemp, 256);
1659  first = 1; last = strlen(label);
1660  }
1661 
1662 // We eliminate the non significant 0 after '.'
1663  if (ndecimals) {
1664  char *adot = strchr(label,'.');
1665  if (adot) adot[ndecimals] = 0;
1666  } else {
1667  while (label[last] == '0') { label[last] = 0; last--;}
1668  }
1669 
1670 // We eliminate the dot, unless dot is forced.
1671  if (label[last] == '.') {
1672  if (!optionDot) { label[last] = 0; last--;}
1673  }
1674 
1675 // Make sure the label is not "-0"
1676  if (last-first == 1 && label[first] == '-'
1677  && label[last] == '0') {
1678  strncpy(label, "0", 256);
1679  label[last] = 0;
1680  }
1681  }
1682 
1683 // Generate the time labels
1684 
1685  if (optionTime) {
1686  timed = wlabel + (int)(timeoffset) - rangeOffset;
1687  timelabel = (time_t)((Long_t)(timed));
1688  if (optionTime == 1) {
1689  utctis = localtime(&timelabel);
1690  } else {
1691  utctis = gmtime(&timelabel);
1692  }
1693  TString timeformattmp;
1694  if (timeformat.Length() < 220) timeformattmp = timeformat;
1695  else timeformattmp = "#splitline{Format}{too long}";
1696 
1697 // Appends fractional part if seconds displayed
1698  if (dwlabel<0.9) {
1699  double tmpdb;
1700  int tmplast;
1701  snprintf(label, 256, "%%S%7.5f", modf(timed,&tmpdb));
1702  tmplast = strlen(label)-1;
1703 
1704 // We eliminate the non significant 0 after '.'
1705  while (label[tmplast] == '0') {
1706  label[tmplast] = 0; tmplast--;
1707  }
1708 
1709  timeformattmp.ReplaceAll("%S",label);
1710 // replace the "0." at the beginning by "s"
1711  timeformattmp.ReplaceAll("%S0.","%Ss");
1712 
1713  }
1714 
1715  if (utctis != nullptr) {
1716  strftime(label, 256, timeformattmp.Data(), utctis);
1717  } else {
1718  strncpy(label, "invalid", 256);
1719  }
1720  strncpy(chtemp, &label[0], 256);
1721  first = 0; last=strlen(label)-1;
1722  wlabel = wTimeIni + (k+1)*dwlabel;
1723  }
1724 
1725 // We generate labels (numeric or alphanumeric).
1726 
1727  if (optionNoopt && !optionInt)
1728  Rotate (xlabel,ylabel,cosphi,sinphi,x0,y0,xx,yy);
1729  else Rotate (xlabel,ylabel,cosphi,sinphi,xx0,yy0,xx,yy);
1730  if (y0 == y1 && !optionDown && !optionUp) {
1731  yy -= 0.80*charheight;
1732  }
1733  if (optionVert) {
1734  if (x0 != x1 && y0 != y1) {
1735  if (optionNoopt && !optionInt)
1736  Rotate (xlabel,0,cosphi,sinphi,x0,y0,xx,yy);
1737  else Rotate (xlabel,0,cosphi,sinphi,xx0,yy0,xx,yy);
1738  if (cosphi > 0 ) yy += ylabel;
1739  if (cosphi < 0 ) yy -= ylabel;
1740  }
1741  }
1742  if (!optionY || (x0 == x1)) {
1743  if (!optionText) {
1744  if (first > last) strncpy(chtemp, " ", 256);
1745  else strncpy(chtemp, &label[first], 256);
1746  if (fNModLabs) ChangeLabelAttributes(k+1, nlabels, textaxis, chtemp);
1747  typolabel = chtemp;
1748  if (!optionTime) typolabel.ReplaceAll("-", "#minus");
1749  if (autotoff) {
1750  UInt_t w,h;
1751  textaxis->SetText(0.,0., typolabel.Data());
1752  textaxis->GetBoundingBox(w,h);
1753  toffset = TMath::Max(toffset,(double)w/((double)gPad->GetWw()*gPad->GetWNDC()));
1754  }
1755  textaxis->PaintLatex(gPad->GetX1() + xx*(gPad->GetX2() - gPad->GetX1()),
1756  gPad->GetY1() + yy*(gPad->GetY2() - gPad->GetY1()),
1757  textaxis->GetTextAngle(),
1758  textaxis->GetTextSize(),
1759  typolabel.Data());
1760  if (fNModLabs) ResetLabelAttributes(textaxis);
1761  } else {
1762  if (optionText == 1) textaxis->PaintLatex(gPad->GetX1() + xx*(gPad->GetX2() - gPad->GetX1()),
1763  gPad->GetY1() + yy*(gPad->GetY2() - gPad->GetY1()),
1764  0,
1765  textaxis->GetTextSize(),
1766  fAxis->GetBinLabel(k+fAxis->GetFirst()));
1767  }
1768  } else {
1769 
1770 // Text alignment is down
1771  if (!optionText) lnlen = last-first+1;
1772  else {
1773  if (k+1 > nhilab) lnlen = 0;
1774  }
1775  for ( l=1; l<=lnlen; l++) {
1776  if (!optionText) *chtemp = label[first+l-2];
1777  else {
1778  if (lnlen == 0) strncpy(chtemp, " ", 256);
1779  else strncpy(chtemp, "1", 256);
1780  }
1781  typolabel = chtemp;
1782  typolabel.ReplaceAll("-", "#minus");
1783  textaxis->PaintLatex(gPad->GetX1() + xx*(gPad->GetX2() - gPad->GetX1()),
1784  gPad->GetY1() + yy*(gPad->GetY2() - gPad->GetY1()),
1785  0,
1786  textaxis->GetTextSize(),
1787  typolabel.Data());
1788  yy -= charheight*1.3;
1789  }
1790  }
1791  }
1792 
1793 // We use the format x 10 ** n
1794 
1795  if (flexe && !optionText && nexe) {
1796  snprintf(label,256,"#times10^{%d}", nexe);
1797  if (x0 != x1) { xfactor = axis_length+0.1*charheight; yfactor = 0; }
1798  else { xfactor = y1-y0+0.1*charheight; yfactor = 0; }
1799  Rotate (xfactor,yfactor,cosphi,sinphi,x0,y0,xx,yy);
1800  textaxis->SetTextAlign(11);
1801  if (GetLabelFont()%10 < 2) // force TLatex mode in PaintLatex
1802  textaxis->SetTextFont((Int_t)(GetLabelFont()/10)*10+2);
1803  if (fAxis && !strcmp(fAxis->GetName(),"xaxis")) {
1804  xx = xx + fXAxisExpXOffset;
1805  yy = yy + fXAxisExpYOffset;
1806  }
1807  if (fAxis && !strcmp(fAxis->GetName(),"yaxis")) {
1808  xx = xx + fYAxisExpXOffset;
1809  yy = yy + fYAxisExpYOffset;
1810  }
1811  typolabel = label;
1812  typolabel.ReplaceAll("-", "#minus");
1813  textaxis->PaintLatex(gPad->GetX1() + xx*(gPad->GetX2() - gPad->GetX1()),
1814  gPad->GetY1() + yy*(gPad->GetY2() - gPad->GetY1()),
1815  0,
1816  textaxis->GetTextSize(),
1817  typolabel.Data());
1818  }
1819  }
1820  }
1821  }
1822 
1823 // Log axis
1824 
1825  if (optionLog && ndiv) {
1826  UInt_t xi1=0,xi2=0,wi=0,yi1=0,yi2=0,hi=0,xl=0,xh=0;
1827  Bool_t firstintlab = kTRUE, overlap = kFALSE;
1828  if ((wmin == wmax) || (ndiv == 0)) {
1829  Error(where, "wmin (%f) == wmax (%f), or ndiv == 0", wmin, wmax);
1830  goto L210;
1831  }
1832  if (wmin <= 0) {
1833  Error(where, "negative logarithmic axis");
1834  goto L210;
1835  }
1836  if (wmax <= 0) {
1837  Error(where, "negative logarithmic axis");
1838  goto L210;
1839  }
1840  xmnlog = TMath::Log10(wmin);
1841  if (xmnlog > 0) xmnlog += 1.E-6;
1842  else xmnlog -= 1.E-6;
1843  x00 = 0;
1844  x11 = axis_length;
1845  h2 = TMath::Log10(wmax);
1846  h2sav = h2;
1847  if (h2 > 0) h2 += 1.E-6;
1848  else h2 -= 1.E-6;
1849  ih1 = int(xmnlog);
1850  ih2 = 1+int(h2);
1851  nbinin = ih2-ih1+1;
1852  axmul = (x11-x00)/(h2sav-xmnlog);
1853 
1854 // Plot decade and intermediate tick marks
1855  decade = ih1-2;
1856  labelnumber = ih1;
1857  if ( xmnlog > 0 && (xmnlog-Double_t(ih1) > 0) ) labelnumber++;
1858  Int_t changelablogid = 0;
1859  Int_t changelablognum = 0;
1860  for (j=1; j<=nbinin; j++) {
1861 
1862 // Plot decade
1863  firstintlab = kTRUE, overlap = kFALSE;
1864  decade++;
1865  if (x0 == x1 && j == 1) ylabel += charheight*0.33;
1866  if (y0 == y1 && j == 1) ylabel -= charheight*0.65;
1867  xone = x00+axmul*(Double_t(decade)-xmnlog);
1868  //the following statement is a trick to circumvent a gcc bug
1869  if (j < 0) printf("j=%d\n",j);
1870  if (x00 > xone) goto L160;
1871  if ((xone-x11)>epsilon) break;
1872  xtwo = xone;
1873  y = 0;
1874  if (!mside) y -= atick[0];
1875  Rotate(xone,y,cosphi,sinphi,x0,y0,xpl2,ypl2);
1876  Rotate(xtwo,atick[0],cosphi,sinphi,x0,y0,xpl1,ypl1);
1877  if (optionVert) {
1878  if ((x0 != x1) && (y0 != y1)) {
1879  if (mside) {
1880  xpl1=xpl2;
1881  if (cosphi > 0) ypl1 = ypl2 + atick[0];
1882  else ypl1 = ypl2 - atick[0];
1883  }
1884  else {
1885  xpl1 = 0.5*(xpl1 + xpl2);
1886  xpl2 = xpl1;
1887  ypl1 = 0.5*(ypl1 + ypl2) + atick[0];
1888  ypl2 = 0.5*(ypl1 + ypl2) - atick[0];
1889  }
1890  }
1891  }
1892  if (!drawGridOnly) PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1893 
1894  if (optionGrid) {
1895  Rotate(xone,0,cosphi,sinphi,x0,y0,xpl2,ypl2);
1896  Rotate(xone,grid_side*gridlength,cosphi,sinphi,x0,y0,xpl1,ypl1);
1897  linegrid->PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1898  }
1899 
1900  if (!drawGridOnly && !optionUnlab) {
1901 
1902 // We generate labels (numeric only).
1903  if (noExponent) {
1904  rlab = TMath::Power(10,labelnumber);
1905  snprintf(label,256, "%f", rlab);
1906  LabelsLimits(label,first,last);
1907  while (last > first) {
1908  if (label[last] != '0') break;
1909  label[last] = 0;
1910  last--;
1911  }
1912  if (label[last] == '.') {label[last] = 0; last--;}
1913  } else {
1914  snprintf(label,256, "%d", labelnumber);
1915  LabelsLimits(label,first,last);
1916  }
1917  Rotate (xone,ylabel,cosphi,sinphi,x0,y0,xx,yy);
1918  if ((x0 == x1) && !optionPara) {
1919  if (lside < 0) {
1920  if (mside < 0) {
1921  if (labelnumber == 0) nch=1;
1922  else nch=2;
1923  xx += nch*charheight;
1924  } else {
1925  xx += 0.25*charheight;
1926  }
1927  }
1928  xx += 0.25*charheight;
1929  }
1930  if ((y0 == y1) && !optionDown && !optionUp) {
1931  if (noExponent) yy += 0.33*charheight;
1932  }
1933  if (n1a == 0)goto L210;
1934  kmod = nbinin/n1a;
1935  if (kmod == 0) kmod=1000000;
1936  if ((nbinin <= n1a) || (j == 1) || (j == nbinin) || ((nbinin > n1a) && (j%kmod == 0))) {
1937  if (labelnumber == 0) {
1938  snprintf(chtemp,256, "1");
1939  } else if (labelnumber == 1) {
1940  snprintf(chtemp,256, "10");
1941  } else {
1942  if (noExponent) {
1943  chtemp = &label[first];
1944  } else {
1945  snprintf(chtemp,256, "10^{%d}", labelnumber);
1946  }
1947  }
1948  if (fNModLabs) {
1949  if (changelablogid == 0) changelablognum = nbinin-j;
1950  changelablogid++;
1951  ChangeLabelAttributes(changelablogid, changelablognum, textaxis, chtemp);
1952  }
1953  typolabel = chtemp;
1954  typolabel.ReplaceAll("-", "#minus");
1955  if (autotoff) {
1956  UInt_t w,h;
1957  textaxis->SetText(0.,0., typolabel.Data());
1958  textaxis->GetBoundingBox(w,h);
1959  toffset = TMath::Max(toffset,(double)w/((double)gPad->GetWw()*gPad->GetWNDC()));
1960  }
1961  textaxis->PaintLatex(gPad->GetX1() + xx*(gPad->GetX2() - gPad->GetX1()),
1962  gPad->GetY1() + yy*(gPad->GetY2() - gPad->GetY1()),
1963  0, textaxis->GetTextSize(), typolabel.Data());
1964  if (fNModLabs) ResetLabelAttributes(textaxis);
1965  }
1966  labelnumber++;
1967  }
1968 L160:
1969  for (k=2;k<10;k++) {
1970 
1971 // Plot intermediate tick marks
1972  xone = x00+axmul*(TMath::Log10(Double_t(k))+Double_t(decade)-xmnlog);
1973  if (x00 > xone) continue;
1974  if (xone > x11) goto L200;
1975  y = 0;
1976  if (!mside) y -= atick[1];
1977  xtwo = xone;
1978  Rotate(xone,y,cosphi,sinphi,x0,y0,xpl2,ypl2);
1979  Rotate(xtwo,atick[1],cosphi,sinphi,x0,y0,xpl1,ypl1);
1980  if (optionVert) {
1981  if ((x0 != x1) && (y0 != y1)) {
1982  if (mside) {
1983  xpl1 = xpl2;
1984  if (cosphi > 0) ypl1 = ypl2 + atick[1];
1985  else ypl1 = ypl2 - atick[1];
1986  }
1987  else {
1988  xpl1 = 0.5*(xpl1+xpl2);
1989  xpl2 = xpl1;
1990  ypl1 = 0.5*(ypl1+ypl2) + atick[1];
1991  ypl2 = 0.5*(ypl1+ypl2) - atick[1];
1992  }
1993  }
1994  }
1995  idn = n1a*2;
1996  if ((nbinin <= idn) || ((nbinin > idn) && (k == 5))) {
1997  if (!drawGridOnly) PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
1998 
1999 // Draw the intermediate LOG labels if requested
2000 
2001  if (moreLogLabels && !optionUnlab && !drawGridOnly && !overlap) {
2002  if (noExponent) {
2003  rlab = Double_t(k)*TMath::Power(10,labelnumber-1);
2004  snprintf(chtemp,256, "%g", rlab);
2005  } else {
2006  if (labelnumber-1 == 0) {
2007  snprintf(chtemp,256, "%d", k);
2008  } else if (labelnumber-1 == 1) {
2009  snprintf(chtemp,256, "%d", 10*k);
2010  } else {
2011  snprintf(chtemp,256, "%d#times10^{%d}", k, labelnumber-1);
2012  }
2013  }
2014  Rotate (xone,ylabel,cosphi,sinphi,x0,y0,xx,yy);
2015  if ((x0 == x1) && !optionPara) {
2016  if (lside < 0) {
2017  if (mside < 0) {
2018  if (labelnumber == 0) nch=1;
2019  else nch=2;
2020  xx += nch*charheight;
2021  } else {
2022  if (labelnumber >= 0) xx += 0.25*charheight;
2023  else xx += 0.50*charheight;
2024  }
2025  }
2026  xx += 0.25*charheight;
2027  }
2028  if ((y0 == y1) && !optionDown && !optionUp) {
2029  if (noExponent) yy += 0.33*charheight;
2030  }
2031  if (optionVert) {
2032  if ((x0 != x1) && (y0 != y1)) {
2033  Rotate(xone,ylabel,cosphi,sinphi,x0,y0,xx,yy);
2034  if (cosphi > 0) yy += ylabel;
2035  else yy -= ylabel;
2036  }
2037  }
2038  textaxis->SetTitle(chtemp);
2039  Double_t u = gPad->GetX1() + xx*(gPad->GetX2() - gPad->GetX1());
2040  Double_t v = gPad->GetY1() + yy*(gPad->GetY2() - gPad->GetY1());
2041  if (firstintlab) {
2042  textaxis->GetBoundingBox(wi, hi); wi=(UInt_t)(wi*1.3); hi=(UInt_t)(hi*1.3);
2043  xi1 = gPad->XtoAbsPixel(u);
2044  yi1 = gPad->YtoAbsPixel(v);
2045  firstintlab = kFALSE;
2046  typolabel = chtemp;
2047  typolabel.ReplaceAll("-", "#minus");
2048  textaxis->PaintLatex(u,v,0,textaxis->GetTextSize(),typolabel.Data());
2049  } else {
2050  xi2 = gPad->XtoAbsPixel(u);
2051  yi2 = gPad->YtoAbsPixel(v);
2052  xl = TMath::Min(xi1,xi2);
2053  xh = TMath::Max(xi1,xi2);
2054  if ((x0 == x1 && yi1-hi <= yi2) || (y0 == y1 && xl+wi >= xh)){
2055  overlap = kTRUE;
2056  } else {
2057  xi1 = xi2;
2058  yi1 = yi2;
2059  textaxis->GetBoundingBox(wi, hi); wi=(UInt_t)(wi*1.3); hi=(UInt_t)(hi*1.3);
2060  typolabel = chtemp;
2061  typolabel.ReplaceAll("-", "#minus");
2062  textaxis->PaintLatex(u,v,0,textaxis->GetTextSize(),typolabel.Data());
2063  }
2064  }
2065  }
2066 
2067 // Draw the intermediate LOG grid if only three decades are requested
2068  if (optionGrid && nbinin <= 5 && ndiv > 100) {
2069  Rotate(xone,0,cosphi,sinphi,x0,y0,xpl2, ypl2);
2070  Rotate(xone,grid_side*gridlength,cosphi,sinphi,x0,y0, xpl1,ypl1);
2071  linegrid->PaintLineNDC(xpl1, ypl1, xpl2, ypl2);
2072  }
2073  } //endif ((nbinin <= idn) ||
2074  } //endfor (k=2;k<10;k++)
2075  } //endfor (j=1; j<=nbinin; j++)
2076 L200:
2077  Int_t dummy = 0; if (dummy) { }
2078  } //endif (optionLog && ndiv)
2079 
2080 // Draw axis title if it exists
2081  if (!drawGridOnly && strlen(GetTitle())) {
2082  textaxis->SetTextSize (GetTitleSize());
2083  charheight = GetTitleSize();
2084  if ((GetTextFont() % 10) > 2) {
2085  charheight = charheight/gPad->GetWh();
2086  }
2087  if (x1 == x0) {
2088  if (autotoff) {
2089  if (toffset) ylabel = xlside*charheight+toffset;
2090  else ylabel = xlside*1.6*charheight;
2091  } else {
2092  ylabel = xlside*1.6*charheight*toffset;
2093  }
2094  } else {
2095  ylabel = xlside*1.3*charheight*toffset;
2096  }
2097  if (y1 == y0) ylabel = xlside*1.6*charheight*toffset;
2098  Double_t axispos;
2099  if (TestBit(TAxis::kCenterTitle)) axispos = 0.5*axis_length;
2100  else axispos = axis_length;
2102  if (x1 >= x0) {
2103  if (TestBit(TAxis::kCenterTitle)) textaxis->SetTextAlign(22);
2104  else textaxis->SetTextAlign(12);
2105  } else {
2106  if (TestBit(TAxis::kCenterTitle)) textaxis->SetTextAlign(22);
2107  else textaxis->SetTextAlign(32);
2108  }
2109  phil+=kPI;
2110  } else {
2111  if (x1 >= x0) {
2112  if (TestBit(TAxis::kCenterTitle)) textaxis->SetTextAlign(22);
2113  else textaxis->SetTextAlign(32);
2114  } else {
2115  if (TestBit(TAxis::kCenterTitle)) textaxis->SetTextAlign(22);
2116  else textaxis->SetTextAlign(12);
2117  }
2118  }
2119  Rotate(axispos,ylabel,cosphi,sinphi,x0,y0,xpl1,ypl1);
2120  textaxis->SetTextColor(TitleColor);
2121  textaxis->SetTextFont(TitleFont);
2122  textaxis->PaintLatex(gPad->GetX1() + xpl1*(gPad->GetX2() - gPad->GetX1()),
2123  gPad->GetY1() + ypl1*(gPad->GetY2() - gPad->GetY1()),
2124  phil*180/kPI,
2125  GetTitleSize(),
2126  GetTitle());
2127  }
2128 
2129 L210:
2130  if (optionGrid) delete linegrid;
2131  delete textaxis;
2132 }
2133 
2134 ////////////////////////////////////////////////////////////////////////////////
2135 /// Internal method for axis labels optimisation. This method adjusts the bining
2136 /// of the axis in order to have integer values for the labels.
2137 ///
2138 /// \param[in] A1,A2 Old WMIN,WMAX
2139 /// \param[out] binLow,binHigh New WMIN,WMAX
2140 /// \param[in] nold Old NDIV (primary divisions)
2141 /// \param[out] nbins New NDIV
2142 /// \param[out] binWidth Bin width
2143 
2144 void TGaxis::AdjustBinSize(Double_t A1, Double_t A2, Int_t nold
2145  ,Double_t &binLow, Double_t &binHigh, Int_t &nbins, Double_t &binWidth)
2146 {
2147 
2148  binWidth = TMath::Abs(A2-A1)/Double_t(nold);
2149  if (binWidth <= 1) { binWidth = 1; binLow = int(A1); }
2150  else {
2151  Int_t width = int(binWidth/5) + 1;
2152  binWidth = 5*width;
2153  binLow = int(A1/binWidth)*binWidth;
2154 
2155 // We determine binLow to have one tick mark at 0
2156 // if there are negative labels.
2157 
2158  if (A1 < 0) {
2159  for (Int_t ic=0; ic<1000; ic++) {
2160  Double_t rbl = binLow/binWidth;
2161  Int_t ibl = int(binLow/binWidth);
2162  if ( (rbl-ibl) == 0 || ic > width) { binLow -= 5; break;}
2163  }
2164  }
2165  }
2166  binHigh = int(A2);
2167  nbins = 0;
2168  Double_t xb = binLow;
2169  while (xb <= binHigh) {
2170  xb += binWidth;
2171  nbins++;
2172  }
2173  binHigh = xb - binWidth;
2174 }
2175 
2176 ////////////////////////////////////////////////////////////////////////////////
2177 /// Internal method to find first and last character of a label.
2178 
2179 void TGaxis::LabelsLimits(const char *label, Int_t &first, Int_t &last)
2180 {
2181  last = strlen(label)-1;
2182  for (Int_t i=0; i<=last; i++) {
2183  if (strchr("1234567890-+.", label[i]) ) { first = i; return; }
2184  }
2185  Error("LabelsLimits", "attempt to draw a blank label");
2186 }
2187 
2188 ////////////////////////////////////////////////////////////////////////////////
2189 /// Internal method to rotate axis coordinates.
2190 
2192  ,Double_t XT, Double_t YT, Double_t &U, Double_t &V)
2193 {
2194  U = CFI*X-SFI*Y+XT;
2195  V = SFI*X+CFI*Y+YT;
2196 }
2197 
2198 ////////////////////////////////////////////////////////////////////////////////
2199 /// Save primitive as a C++ statement(s) on output stream out
2200 
2201 void TGaxis::SavePrimitive(std::ostream &out, Option_t * /*= ""*/)
2202 {
2203  char quote = '"';
2204  if (gROOT->ClassSaved(TGaxis::Class())) {
2205  out<<" ";
2206  } else {
2207  out<<" TGaxis *";
2208  }
2209  out<<"gaxis = new TGaxis("<<fX1<<","<<fY1<<","<<fX2<<","<<fY2
2210  <<","<<fWmin<<","<<fWmax<<","<<fNdiv<<","<<quote<<fChopt.Data()<<quote<<");"<<std::endl;
2211  out<<" gaxis->SetLabelOffset("<<GetLabelOffset()<<");"<<std::endl;
2212  out<<" gaxis->SetLabelSize("<<GetLabelSize()<<");"<<std::endl;
2213  out<<" gaxis->SetTickSize("<<GetTickSize()<<");"<<std::endl;
2214  out<<" gaxis->SetGridLength("<<GetGridLength()<<");"<<std::endl;
2215  out<<" gaxis->SetTitleOffset("<<GetTitleOffset()<<");"<<std::endl;
2216  out<<" gaxis->SetTitleSize("<<GetTitleSize()<<");"<<std::endl;
2217  out<<" gaxis->SetTitleColor("<<GetTextColor()<<");"<<std::endl;
2218  out<<" gaxis->SetTitleFont("<<GetTextFont()<<");"<<std::endl;
2219 
2220  if (strlen(GetName())) {
2221  out<<" gaxis->SetName("<<quote<<GetName()<<quote<<");"<<std::endl;
2222  }
2223  if (strlen(GetTitle())) {
2224  out<<" gaxis->SetTitle("<<quote<<GetTitle()<<quote<<");"<<std::endl;
2225  }
2226 
2227  if (fLabelColor != 1) {
2228  if (fLabelColor > 228) {
2230  out<<" gaxis->SetLabelColor(ci);" << std::endl;
2231  } else
2232  out<<" gaxis->SetLabelColor("<<GetLabelColor()<<");"<<std::endl;
2233  }
2234  if (fLineColor != 1) {
2235  if (fLineColor > 228) {
2237  out<<" gaxis->SetLineColor(ci);" << std::endl;
2238  } else
2239  out<<" gaxis->SetLineColor("<<GetLineColor()<<");"<<std::endl;
2240  }
2241  if (fLineStyle != 1) {
2242  out<<" gaxis->SetLineStyle("<<GetLineStyle()<<");"<<std::endl;
2243  }
2244  if (fLineWidth != 1) {
2245  out<<" gaxis->SetLineWidth("<<GetLineWidth()<<");"<<std::endl;
2246  }
2247  if (fLabelFont != 62) {
2248  out<<" gaxis->SetLabelFont("<<GetLabelFont()<<");"<<std::endl;
2249  }
2251  out<<" gaxis->SetMoreLogLabels();"<<std::endl;
2252  }
2253  if (TestBit(TAxis::kNoExponent)) {
2254  out<<" gaxis->SetNoExponent();"<<std::endl;
2255  }
2256 
2257  out<<" gaxis->Draw();"<<std::endl;
2258 }
2259 
2260 ////////////////////////////////////////////////////////////////////////////////
2261 /// Set the decimals flag. By default, blank characters are stripped, and then the
2262 /// label is correctly aligned. The dot, if last character of the string, is also
2263 /// stripped, unless this option is specified. One can disable the option by
2264 /// calling `axis.SetDecimals(kTRUE)`.
2265 /// Note the bit is set in fBits (as opposed to fBits2 in TAxis!)
2266 
2267 void TGaxis::SetDecimals(Bool_t dot)
2268 {
2269  if (dot) SetBit(TAxis::kDecimals);
2270  else ResetBit(TAxis::kDecimals);
2272 
2273 ////////////////////////////////////////////////////////////////////////////////
2274 /// Specify a function to map the axis values.
2275 
2276 void TGaxis::SetFunction(const char *funcname)
2277 {
2278  fFunctionName = funcname;
2279  if (!funcname[0]) {
2281  return;
2282  }
2283  fFunction = (TF1*)gROOT->GetFunction(funcname);
2284  if (!fFunction) {
2285  Error("SetFunction", "unknown function: %s", funcname);
2286  } else {
2287  fWmin = fFunction->GetXmin();
2288  fWmax = fFunction->GetXmax();
2289  }
2290 }
2291 
2292 ////////////////////////////////////////////////////////////////////////////////
2293 /// Define new text attributes for the label number "labNum". It allows to do a
2294 /// fine tuning of the labels. All the attributes can be changed, even the
2295 /// label text itself.
2296 ///
2297 /// \param[in] labNum Number of the label to be changed, negative numbers start from the end
2298 /// \param[in] labAngle New angle value
2299 /// \param[in] labSize New size (0 erase the label)
2300 /// \param[in] labAlign New alignment value
2301 /// \param[in] labColor New label color
2302 /// \param[in] labFont New label font
2303 /// \param[in] labText New label text
2304 ///
2305 /// If an attribute should not be changed just give the value
2306 /// "-1".The following macro gives an example:
2307 ///
2308 /// Begin_Macro(source)
2309 /// {
2310 /// c1 = new TCanvas("c1","Examples of Gaxis",10,10,900,500);
2311 /// c1->Range(-6,-0.1,6,0.1);
2312 /// TGaxis *axis1 = new TGaxis(-5.5,0.,5.5,0.,0.0,100,510,"");
2313 /// axis1->SetName("axis1");
2314 /// axis1->SetTitle("Axis Title");
2315 /// axis1->SetTitleSize(0.05);
2316 /// axis1->SetTitleColor(kBlue);
2317 /// axis1->SetTitleFont(42);
2318 /// axis1->ChangeLabel(1,-1,-1,-1,2);
2319 /// axis1->ChangeLabel(3,-1,0.);
2320 /// axis1->ChangeLabel(5,30.,-1,0);
2321 /// axis1->ChangeLabel(6,-1,-1,-1,3,-1,"6th label");
2322 /// axis1->ChangeLabel(-2,-1,-1,-1,3,-1,"2nd to last label");
2323 /// axis1->Draw();
2324 /// }
2325 /// End_Macro
2326 ///
2327 /// If labnum=0 the list of modified labels is reset.
2328 
2329 void TGaxis::ChangeLabel(Int_t labNum, Double_t labAngle, Double_t labSize,
2330  Int_t labAlign, Int_t labColor, Int_t labFont,
2331  TString labText)
2332 {
2334  if (!fModLabs) fModLabs = new TList();
2335 
2336  // Reset the list of modified labels.
2337  if (labNum == 0) {
2338  delete fModLabs;
2339  fModLabs = 0;
2340  fNModLabs = 0;
2341  return;
2342  }
2343 
2344  TAxisModLab *ml = new TAxisModLab();
2345  ml->SetLabNum(labNum);
2346  ml->SetAngle(labAngle);
2347  ml->SetSize(labSize);
2348  ml->SetAlign(labAlign);
2349  ml->SetColor(labColor);
2350  ml->SetFont(labFont);
2351  ml->SetText(labText);
2352 
2353  fModLabs->Add((TObject*)ml);
2354 }
2355 
2356 ////////////////////////////////////////////////////////////////////////////////
2357 /// Change the label attributes of label number i. If needed.
2358 ///
2359 /// \param[in] i Current label number to be changed if needed
2360 /// \param[in] nlabels Totals number of labels on for this axis (useful when i is counted from the end)
2361 /// \param[in] t Original TLatex string holding the label to be changed
2362 /// \param[in] c Text string to be drawn
2363 
2364 static Double_t SavedTextAngle;
2365 static Double_t SavedTextSize;
2366 static Int_t SavedTextAlign;
2367 static Int_t SavedTextColor;
2370 void TGaxis::ChangeLabelAttributes(Int_t i, Int_t nlabels, TLatex* t, char* c)
2372  if (!fModLabs) return;
2373 
2375  TAxisModLab *ml;
2376  Int_t labNum;
2377  while ( (ml = (TAxisModLab*)next()) ) {
2378  SavedTextAngle = t->GetTextAngle();
2379  SavedTextSize = t->GetTextSize();
2380  SavedTextAlign = t->GetTextAlign();
2381  SavedTextColor = t->GetTextColor();
2382  SavedTextFont = t->GetTextFont();
2383  labNum = ml->GetLabNum();
2384  if (labNum < 0) labNum = nlabels + labNum + 2;
2385  if (i == labNum) {
2386  if (ml->GetAngle()>=0.) t->SetTextAngle(ml->GetAngle());
2387  if (ml->GetSize()>=0.) t->SetTextSize(ml->GetSize());
2388  if (ml->GetAlign()>0) t->SetTextAlign(ml->GetAlign());
2389  if (ml->GetColor()>=0) t->SetTextColor(ml->GetColor());
2390  if (ml->GetFont()>0) t->SetTextFont(ml->GetFont());
2391  if (!(ml->GetText().IsNull())) strncpy(c, (ml->GetText()).Data(), 256);
2392  return;
2393  }
2394  }
2395 }
2396 
2397 ////////////////////////////////////////////////////////////////////////////////
2398 /// Reset the label attributes to the value they have before the last call to
2399 /// ChangeLabelAttributes.
2400 
2402 {
2403  t->SetTextAngle(SavedTextAngle);
2404  t->SetTextSize(SavedTextSize);
2405  t->SetTextAlign(SavedTextAlign);
2406  t->SetTextColor(SavedTextColor);
2407  t->SetTextFont(SavedTextFont);
2408 }
2409 
2410 ////////////////////////////////////////////////////////////////////////////////
2411 /// Static function to set `fgMaxDigits` for axis.`fgMaxDigits` is
2412 /// the maximum number of digits permitted for the axis labels above which the
2413 /// notation with 10^N is used.For example, to accept 6 digits number like 900000
2414 /// on an axis call `TGaxis::SetMaxDigits(6)`. The default value is 5.
2415 /// `fgMaxDigits` must be greater than 0.
2416 
2417 void TGaxis::SetMaxDigits(Int_t maxd)
2418 {
2419  fgMaxDigits = maxd;
2420  if (maxd < 1) fgMaxDigits = 1;
2422 
2423 ////////////////////////////////////////////////////////////////////////////////
2424 /// Change the name of the axis.
2425 
2426 void TGaxis::SetName(const char *name)
2427 {
2428  fName = name;
2429 }
2431 ////////////////////////////////////////////////////////////////////////////////
2432 /// Set the kMoreLogLabels bit flag. When this option is selected more labels are
2433 /// drawn when in logarithmic scale and there is a small number of decades (less than 3).
2434 /// Note that this option is automatically inherited from TAxis
2435 
2437 {
2438  if (more) SetBit(TAxis::kMoreLogLabels);
2441 
2442 ////////////////////////////////////////////////////////////////////////////////
2443 /// Set the NoExponent flag. By default, an exponent of the form 10^N is used
2444 /// when the label values are either all very small or very large. One can disable
2445 /// the exponent by calling axis.SetNoExponent(kTRUE).
2446 
2447 void TGaxis::SetNoExponent(Bool_t noExponent)
2448 {
2449  if (noExponent) SetBit(TAxis::kNoExponent);
2452 
2453 ////////////////////////////////////////////////////////////////////////////////
2454 /// To set axis options.
2455 
2456 void TGaxis::SetOption(Option_t *option)
2457 {
2458  fChopt = option;
2459 }
2461 ////////////////////////////////////////////////////////////////////////////////
2462 /// Change the title of the axis.
2463 
2464 void TGaxis::SetTitle(const char *title)
2465 {
2466  fTitle = title;
2467 }
2469 ////////////////////////////////////////////////////////////////////////////////
2470 /// Change the format used for time plotting.
2471 /// The format string for date and time use the same options as the one used
2472 /// in the standard strftime C function, i.e. :
2473 ///
2474 /// for date :
2475 ///
2476 /// - `%a` abbreviated weekday name
2477 /// - `%b` abbreviated month name
2478 /// - `%d` day of the month (01-31)
2479 /// - `%m` month (01-12)
2480 /// - `%y` year without century
2481 ///
2482 /// for time :
2483 ///
2484 /// - `%H` hour (24-hour clock)
2485 /// - `%I` hour (12-hour clock)
2486 /// - `%p` local equivalent of AM or PM
2487 /// - `%M` minute (00-59)
2488 /// - `%S` seconds (00-61)
2489 /// - `%%` %
2490 
2491 void TGaxis::SetTimeFormat(const char *tformat)
2492 {
2493  TString timeformat = tformat;
2494 
2495  if (timeformat.Index("%F")>=0 || timeformat.IsNull()) {
2496  fTimeFormat = timeformat;
2497  return;
2498  }
2499 
2500  Int_t idF = fTimeFormat.Index("%F");
2501  if (idF>=0) {
2502  Int_t lnF = fTimeFormat.Length();
2503  TString stringtimeoffset = fTimeFormat(idF,lnF);
2504  fTimeFormat = tformat;
2505  fTimeFormat.Append(stringtimeoffset);
2506  } else {
2507  fTimeFormat = tformat;
2509  }
2510 }
2511 
2512 ////////////////////////////////////////////////////////////////////////////////
2513 /// Change the time offset. If option = "gmt", set display mode to GMT.
2514 
2515 void TGaxis::SetTimeOffset(Double_t toffset, Option_t *option)
2516 {
2517  TString opt = option;
2518  opt.ToLower();
2520  char tmp[20];
2521  time_t timeoff;
2522  struct tm* utctis;
2523  Int_t idF = fTimeFormat.Index("%F");
2524  if (idF>=0) fTimeFormat.Remove(idF);
2525  fTimeFormat.Append("%F");
2526 
2527  timeoff = (time_t)((Long_t)(toffset));
2528 
2529  // offset is always saved in GMT to allow file transport
2530  // to different time zones
2531  utctis = gmtime(&timeoff);
2532 
2533  if (utctis != nullptr) {
2534  strftime(tmp, 20,"%Y-%m-%d %H:%M:%S",utctis);
2535  fTimeFormat.Append(tmp);
2536  } else {
2537  fTimeFormat.Append("1970-01-01 00:00:00");
2538  }
2539 
2540  // append the decimal part of the time offset
2541  Double_t ds = toffset-(Int_t)toffset;
2542  snprintf(tmp,20,"s%g",ds);
2543  fTimeFormat.Append(tmp);
2544 
2545  // add GMT/local option
2546  if (opt.Contains("gmt")) fTimeFormat.Append(" GMT");
2547 }
2548 
2549 ////////////////////////////////////////////////////////////////////////////////
2550 /// Static function to set X and Y offset of the axis 10^n notation.
2551 /// It is in % of the pad size. It can be negative.
2552 /// axis specifies which axis ("x","y"), default = "x"
2553 /// if axis="xz" set the two axes
2554 
2555 void TGaxis::SetExponentOffset(Float_t xoff, Float_t yoff, Option_t *axis)
2556 {
2557  TString opt = axis;
2558  opt.ToLower();
2560  if (opt.Contains("x")) {
2561  fXAxisExpXOffset = xoff;
2562  fXAxisExpYOffset = yoff;
2563  }
2564  if (opt.Contains("y")) {
2565  fYAxisExpXOffset = xoff;
2566  fYAxisExpYOffset = yoff;
2567  }
2568 }
2569 
2570 ////////////////////////////////////////////////////////////////////////////////
2571 /// Stream an object of class TGaxis.
2572 
2573 void TGaxis::Streamer(TBuffer &R__b)
2574 {
2575  if (R__b.IsReading()) {
2576  UInt_t R__s, R__c;
2577  Version_t R__v = R__b.ReadVersion(&R__s, &R__c);
2578  if (R__v > 3) {
2579  R__b.ReadClassBuffer(TGaxis::Class(), this, R__v, R__s, R__c);
2580  return;
2581  }
2582  //====process old versions before automatic schema evolution
2583  TLine::Streamer(R__b);
2584  TAttText::Streamer(R__b);
2585  R__b >> fNdiv;
2586  R__b >> fWmin;
2587  R__b >> fWmax;
2588  R__b >> fGridLength;
2589  R__b >> fTickSize;
2590  R__b >> fLabelOffset;
2591  R__b >> fLabelSize;
2592  R__b >> fTitleOffset;
2593  R__b >> fTitleSize;
2594  R__b >> fLabelFont;
2595  if (R__v > 2) {
2596  R__b >> fLabelColor;
2597  }
2598  fChopt.Streamer(R__b);
2599  fName.Streamer(R__b);
2600  fTitle.Streamer(R__b);
2601  fTimeFormat.Streamer(R__b);
2602  if (R__v > 1) {
2603  fFunctionName.Streamer(R__b);
2604  fFunction = (TF1*)gROOT->GetFunction(fFunctionName.Data());
2605  }
2606  R__b.CheckByteCount(R__s, R__c, TGaxis::IsA());
2607  //====end of old versions
2608 
2609  } else {
2610  R__b.WriteClassBuffer(TGaxis::Class(),this);
2611  }
2612 }
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:81
virtual Float_t GetTickLength() const
Definition: TAttAxis.h:44
virtual void SetName(const char *name)
Change the name of the axis.
Definition: TGaxis.cxx:2430
virtual void SetMoreLogLabels(Bool_t more=kTRUE)
Set the kMoreLogLabels bit flag.
Definition: TGaxis.cxx:2440
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:640
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:2559
Bool_t TestBit(UInt_t f) const
Definition: TObject.h:159
#define BIT(n)
Definition: Rtypes.h:75
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:2022
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:375
Ssiz_t Index(const char *pat, Ssiz_t i=0, ECaseCompare cmp=kExact) const
Definition: TString.h:587
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:2468
Basic string class.
Definition: TString.h:129
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
int nbins[3]
virtual Float_t GetLabelSize() const
Definition: TAttAxis.h:41
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:2374
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:628
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:687
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:112
static Double_t SavedTextSize
Definition: TGaxis.cxx:2369
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:2205
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:652
void SetOption(Option_t *option="")
To set axis options.
Definition: TGaxis.cxx:2460
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:497
std::vector< std::vector< double > > Data
void SetTimeFormat(const char *tformat)
Change the format used for time plotting.
Definition: TGaxis.cxx:2495
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:250
TString fName
Axis name.
Definition: TGaxis.h:41
Double_t ATan2(Double_t, Double_t)
Definition: TMath.h:581
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:43
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:2421
SVector< double, 2 > v
Definition: Dict.h:5
void SetFunction(const char *funcname="")
Specify a function to map the axis values.
Definition: TGaxis.cxx:2280
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:873
Ssiz_t Length() const
Definition: TString.h:388
A simple line.
Definition: TLine.h:23
static Int_t SavedTextAlign
Definition: TGaxis.cxx:2370
TLine * l
Definition: textangle.C:4
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:506
#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:551
void SetLabNum(Int_t n=0)
Set modified label number.
Definition: TAxisModLab.cxx:42
const Bool_t kFALSE
Definition: RtypesCore.h:92
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:2405
TString & Remove(Ssiz_t pos)
Definition: TString.h:621
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:1197
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:2148
void SetLabelOffset(Float_t labeloffset)
Definition: TGaxis.h:107
static const double x1[5]
#define ClassImp(name)
Definition: Rtypes.h:336
virtual void SetDecimals(Bool_t dot=kTRUE)
Set the decimals flag.
Definition: TGaxis.cxx:2271
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:510
void SetLabelFont(Int_t labelfont)
Definition: TGaxis.h:106
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:572
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:385
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:2195
virtual void Add(TObject *obj)
Definition: TList.h:77
static Int_t SavedTextColor
Definition: TGaxis.cxx:2371
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:150
static Int_t SavedTextFont
Definition: TGaxis.cxx:2372
void SetTimeOffset(Double_t toffset, Option_t *option="local")
Change the time offset. If option = "gmt", set display mode to GMT.
Definition: TGaxis.cxx:2519
Double_t Sin(Double_t)
Definition: TMath.h:548
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:284
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:158
static Double_t SavedTextAngle
Change the label attributes of label number i.
Definition: TGaxis.cxx:2368
Width_t GetGridWidth() const
Definition: TStyle.h:209
Definition: first.py:1
Double_t Sqrt(Double_t x)
Definition: TMath.h:591
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:89
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:155
const Bool_t kTRUE
Definition: RtypesCore.h:91
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:2451
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
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:2183
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:2333
Float_t GetTitleSize() const
Definition: TGaxis.h:84
const char * Data() const
Definition: TString.h:347
Double_t GetTimeOffset() const
Definition: TStyle.h:253