ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TLatex.cxx
Go to the documentation of this file.
1 // @(#)root/graf:$Id$
2 // Author: Nicolas Brun, Olivier Couet, Oleksandr Grebenyuk
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 "Riostream.h"
13 #include "TROOT.h"
14 #include "TClass.h"
15 #include "TLatex.h"
16 #include "TMathText.h"
17 #include "TMath.h"
18 #include "TVirtualPad.h"
19 #include "TVirtualPS.h"
20 
22 
24 
25 /** \class TLatex
26 \ingroup BasicGraphics
27 
28 To draw Mathematical Formula.
29 
30 TLatex's purpose is to write mathematical equations. The syntax is very similar
31 to the Latex's one. It provides several functionalities:
32 
33 - [Subscripts and Superscripts](#L1)
34 - [Fractions](#L2)
35 - [Splitting Lines](#L3)
36 - [Roots](#L4)
37 - [Mathematical Symbols](#L5)
38 - [Delimiters](#L6)
39 - [Greek Letters](#L7)
40 - [Accents](#L8)
41 - [Changing Style](#L9)
42 - [Alignment Rules](#L10)
43 - [Character Adjustment](#L11)
44 - [Italic and Boldface](#L12)
45 - [Examples](#L13)
46 - [Interface to TMathText](#L14)
47 
48 When the font precision (see `TAttText`) is low (0 or 1), TLatex is
49 painted as a normal TText, the control characters are not interpreted.
50 
51 ## <a name="L1"></a> Subscripts and Superscripts
52 Subscripts and superscripts are made with the `_` and `^`
53 commands. These commands can be combined to make complicated subscript and
54 superscript expressions. You may adjust the display of subscripts and
55 superscripts by using the two functions `SetIndiceSize(Double_t)`,
56 which set relative size of subscripts and superscripts, and
57 `SetLimitIndiceSize(Int_t)`, which set limits for text resizing of
58 subscripts and superscripts.
59 
60 Examples:
61 
62 Begin_Macro(source)
63 {
64  TCanvas *cl = new TCanvas("cl","cl",10,10,700,500);
65  TLatex Tl; Tl.SetTextFont(43); Tl.SetTextSize(20);
66  Double_t dy = 1./7.;
67  Tl.DrawText(.1, dy, "x^{2y} :"); Tl.DrawLatex(.5, dy, "x^{2y}");
68  Tl.DrawText(.1, 2*dy, "x_{2y} :"); Tl.DrawLatex(.5, 2*dy, "x_{2y}");
69  Tl.DrawText(.1, 3*dy, "x^{y^{2}} :"); Tl.DrawLatex(.5, 3*dy, "x^{y^{2}}");
70  Tl.DrawText(.1, 4*dy, "x^{y_{1}} :"); Tl.DrawLatex(.5, 4*dy, "x^{y_{1}}");
71  Tl.DrawText(.1, 5*dy, "x^{y}_{1} :"); Tl.DrawLatex(.5, 5*dy, "x^{y}_{1}");
72  Tl.DrawText(.1, 6*dy, "x_{1}^{y} :"); Tl.DrawLatex(.5, 6*dy, "x_{1}^{y}");
73 }
74 End_Macro
75 
76 The best way to put the subscripts and superscripts before the character and not
77 after, is to use an empty character:
78 
79 Begin_Macro(source)
80 {
81  TCanvas *cl = new TCanvas("cl","cl",10,10,700,100);
82  TLatex Tl; Tl.SetTextFont(43); Tl.SetTextSize(20);
83  Tl.DrawText(.1, .5, "{}^{40}_{20}Ca :"); Tl.DrawLatex(.5, .5, "{}^{40}_{20}Ca");
84 }
85 End_Macro
86 
87 The subscripts and superscripts operators apply not only on one character but
88 on all the "normal text" preceding them. In the following example the second
89 `E` is lower than the first one because the operator `_` is
90 applied on `/f` which has a descending part, and not only on `f`
91 which as no descending part.
92 
93 Begin_Macro(source)
94 {
95  TCanvas *cl = new TCanvas("cl","cl",10,10,700,100);
96  TLatex Tl; Tl.SetTextFont(43); Tl.SetTextSize(20);
97  Tl.DrawText(.1, .5, "f_{E}/f_{E} :"); Tl.DrawLatex(.5, .5, "f_{E}/f_{E}");
98 }
99 End_Macro
100 
101 To make sure the second operator `_` applies only on `f` a dummy operator `^{}`
102 should be introduced to separate the `f` from the `/`.
103 
104 Begin_Macro(source)
105 {
106  TCanvas *cl = new TCanvas("cl","cl",10,10,700,100);
107  TLatex Tl; Tl.SetTextFont(43); Tl.SetTextSize(20);
108  Tl.DrawText(.1, .5, "f_{E}/^{}f_{E} :"); Tl.DrawLatex(.5, .5, "f_{E}/^{}f_{E}");
109 }
110 End_Macro
111 
112 ## <a name="L2"></a> Fractions
113 Fractions denoted by the `/` symbol are made in the obvious way.
114 The `#frac` command is used for large fractions in displayed formula;
115 it has two arguments: the numerator and the denominator.
116 
117 Examples:
118 
119 Begin_Macro(source)
120 {
121  TCanvas *cl = new TCanvas("cl","cl",10,10,700,100);
122  TLatex Tl; Tl.SetTextFont(43); Tl.SetTextSize(20);
123  Tl.DrawText(.1, .5, "x = #frac{y+z/2}{y^{2}+1} :"); Tl.DrawLatex(.5, .5, "x = #frac{y+z/2}{y^{2}+1}");
124 }
125 End_Macro
126 
127 ## <a name="L3"></a> Splitting Lines
128 Text can be split in two lines via the command `#splitline`.
129 
130 Examples:
131 
132 Begin_Macro(source)
133 {
134  TCanvas *cl = new TCanvas("cl","cl",10,10,700,100);
135  TLatex Tl; Tl.SetTextFont(43); Tl.SetTextSize(20);
136  Tl.DrawText(.1, .5, "#splitline{21 April 2003}{14:02:30} :"); Tl.DrawLatex(.6, .5, "#splitline{21 April 2003}{14:02:30}");
137 }
138 End_Macro
139 
140 ## <a name="L4"></a> Roots
141 The `#sqrt` command produces the square root of its argument; it has
142 an optional first argument for other roots.
143 
144 Examples:
145 
146 Begin_Macro(source)
147 {
148  TCanvas *cl = new TCanvas("cl","cl",10,10,700,100);
149  TLatex Tl; Tl.SetTextFont(43); Tl.SetTextSize(20);
150  Tl.DrawText(.1, .5, "#sqrt{10} #sqrt[3]{10} :"); Tl.DrawLatex(.5, .5, "#sqrt{10} #sqrt[3]{10}");
151 }
152 End_Macro
153 
154 ## <a name="L5"></a> Mathematical Symbols
155 TLatex can display dozens of special mathematical symbols. A few of them, such
156 as `+` and `>` , are produced by typing the corresponding
157 keyboard character. Others are obtained with the commands in the following
158 table:
159 
160 Begin_Macro
161 mathsymbols.C
162 End_Macro
163 
164 
165 ## <a name="L6"></a> Delimiters
166 TLatex provides 4 kinds of proportional delimiters:
167 
168  #[]{....} or "a la" Latex #left[.....#right] : big square brackets
169  #{}{....} or #left{.....#right} : big curly brackets
170  #||{....} or #left|.....#right| : big absolute value symbols
171  #(){....} or #left(.....#right) : big parentheses
172 
173 ## <a name="L7"></a> Greek Letters
174 The command to produce a lowercase Greek letter is obtained by adding a
175 `#` to the name of the letter. For an uppercase Greek letter, just
176 capitalize the first letter of the command name. Some letters have two
177 representations. The name of the second one (the "variation") starts with "var".
178 The following table gives the complete list:
179 
180 Begin_Macro
181 greekletters.C
182 End_Macro
183 
184 
185 ## <a name="L8"></a> Accents
186 Several kind of accents are available:
187 
188  #hat = Begin_Latex #hat{a} End_Latex
189  #check = Begin_Latex #check{a} End_Latex
190  #acute = Begin_Latex #acute{a} End_Latex
191  #grave = Begin_Latex #grave{a} End_Latex
192  #dot = Begin_Latex #dot{a} End_Latex
193  #ddot = Begin_Latex #ddot{a} End_Latex
194  #tilde = Begin_Latex #tilde{a} End_Latex
195 
196 
197 The special sign: `#slash` draws a slash on top of the text between brackets:
198 
199 Begin_Macro(source)
200 {
201  TCanvas *cl = new TCanvas("cl","cl",10,10,700,100);
202  TLatex Tl; Tl.SetTextFont(43); Tl.SetTextSize(20);
203  Tl.DrawText(.1, .5, "#slash{E}_{T} :"); Tl.DrawLatex(.5, .5, "#slash{E}_{T}");
204 }
205 End_Macro
206 
207 Bar and vectors sign are done the following way:
208 
209 Begin_Macro(source)
210 {
211  TCanvas *cl = new TCanvas("cl","cl",10,10,700,100);
212  TLatex Tl; Tl.SetTextFont(43); Tl.SetTextSize(20);
213  Tl.DrawText(.1, .5, "#bar{a} and #vec{a} :"); Tl.DrawLatex(.5, .5, "#bar{a} and #vec{a}");
214 }
215 End_Macro
216 
217 ## <a name="L9"></a> Changing Style
218 One can change the font, the text color, or the text size at any time using :
219 `#font[font-number]{...}`, `#color[color-number]{...}`
220 and `#scale[scale-factor]{...}`
221 
222 Examples:
223 
224 Begin_Macro(source)
225 {
226  TCanvas *cl = new TCanvas("cl","cl",10,10,900,300);
227  TLatex Tl; Tl.SetTextFont(43); Tl.SetTextSize(20);
228  Double_t dy = 1./4.;
229  Tl.DrawText(.01, dy, "#font[12]{Times Italic} and #font[22]{Times bold} :"); Tl.DrawLatex(.7, dy, "#font[12]{Times Italic} and #font[22]{Times bold}");
230  Tl.DrawText(.01, 2*dy, "#color[2]{Red} and #color[4]{Blue} :"); Tl.DrawLatex(.7, 2*dy, "#color[2]{Red} and #color[4]{Blue}");
231  Tl.DrawText(.01, 3*dy, "#scale[1.2]{Bigger} and #scale[0.8]{Smaller} :"); Tl.DrawLatex(.7, 3*dy, "#scale[1.2]{Bigger} and #scale[0.8]{Smaller}");
232 }
233 End_Macro
234 
235 ## <a name="L10"></a> Alignment Rules
236 The `TText` alignment rules apply to the `TLatex` objects with one exception
237 concerning the vertical alignment:
238 
239 - if the vertical alignment = 1 , subscripts are not taken into account
240 - if the vertical alignment = 0 , the text is aligned to the box surrounding
241  the full text with sub and superscripts
242 
243 This is illustrated by the following example:
244 
245 Begin_Macro(source)
246 {
247  TCanvas Tlva("Tlva","Tlva",500,500);
248  Tlva.SetGrid();
249  Tlva.DrawFrame(0,0,1,1);
250  const char *longstring = "K_{S}... K^{*0}... #frac{2s}{#pi#alpha^{2}} #frac{d#sigma}{dcos#theta} (e^{+}e^{-} #rightarrow f#bar{f} ) = #left| #frac{1}{1 - #Delta#alpha} #right|^{2} (1+cos^{2}#theta)";
251 
252  TLatex latex;
253  latex.SetTextSize(0.025);
254  latex.SetTextAlign(13); //align at top
255  latex.DrawLatex(.2,.9,"K_{S}");
256  latex.DrawLatex(.3,.9,"K^{*0}");
257  latex.DrawLatex(.2,.8,longstring);
258 
259  latex.SetTextAlign(12); //centered
260  latex.DrawLatex(.2,.6,"K_{S}");
261  latex.DrawLatex(.3,.6,"K^{*0}");
262  latex.DrawLatex(.2,.5,longstring);
263 
264  latex.SetTextAlign(11); //default bottom alignment
265  latex.DrawLatex(.2,.4,"K_{S}");
266  latex.DrawLatex(.3,.4,"K^{*0}");
267  latex.DrawLatex(.2,.3,longstring);
268 
269  latex.SetTextAlign(10); //special bottom alignment
270  latex.DrawLatex(.2,.2,"K_{S}");
271  latex.DrawLatex(.3,.2,"K^{*0}");
272  latex.DrawLatex(.2,.1,longstring);
273 
274  latex.SetTextAlign(12);
275  latex.SetTextFont(72);
276  latex.DrawLatex(.1,.80,"13");
277  latex.DrawLatex(.1,.55,"12");
278  latex.DrawLatex(.1,.35,"11");
279  latex.DrawLatex(.1,.18,"10");
280  return Tlva;
281 }
282 End_Macro
283 
284 
285 ## <a name="L11"></a> Character Adjustment
286 
287 The two commands `#kern` and `#lower` enable a better control
288 over character placement. The command `#kern[(Float_t)dx]{text}` moves
289 the output string horizontally by the fraction `dx` of its length.
290 Similarly, `#lower[(Float_t)dy]{text}` shifts the text up or down by
291 the fraction `dy` of its height.
292 
293 Examples:
294 
295 Begin_Macro(source)
296 {
297  TCanvas *cl = new TCanvas("cl","cl",10,10,900,300);
298  TLatex Tl; Tl.SetTextFont(43); Tl.SetTextSize(20);
299  TLatex Tt; Tt.SetTextFont(43); Tt.SetTextSize(12);
300  Double_t dy = 1./7.;
301  Tl.DrawLatex(.5, dy, "Positive k#kern[0.3]{e}#kern[0.3]{r}#kern[0.3]{n}#kern[0.3]{i}#kern[0.3]{n}#kern[0.3]{g}");
302  Tt.DrawText(.01, 2*dy, "Positive k#kern[0.3]{e}#kern[0.3]{r}#kern[0.3]{n}#kern[0.3]{i}#kern[0.3]{n}#kern[0.3]{g} :");
303  Tl.DrawLatex(.5, 3*dy, "Negative k#kern[-0.3]{e}#kern[-0.3]{r}#kern[-0.3]{n}#kern[-0.3]{i}#kern[-0.3]{n}#kern[-0.3]{g}");
304  Tt.DrawText(.01, 4*dy, "Negative k#kern[-0.3]{e}#kern[-0.3]{r}#kern[-0.3]{n}#kern[-0.3]{i}#kern[-0.3]{n}#kern[-0.3]{g} :");
305  Tl.DrawLatex(.5, 5*dy, "Vertical a#lower[0.2]{d}#lower[0.4]{j}#lower[0.1]{u}#lower[-0.1]{s}#lower[-0.3]{t}#lower[-0.4]{m}#lower[-0.2]{e}#lower[0.1]{n}t");
306  Tt.DrawText(.01, 6*dy, "Vertical a#lower[0.2]{d}#lower[0.4]{j}#lower[0.1]{u}#lower[-0.1]{s}#lower[-0.3]{t}#lower[-0.4]{m}#lower[-0.2]{e}#lower[0.1]{n}t :");
307 
308 }
309 End_Macro
310 
311 ## <a name="L12"></a> Italic and Boldface
312 Text can be turned italic or boldface using the commands
313 `#it` and `#bf`.
314 
315 Examples:
316 
317 Begin_Macro(source)
318 {
319  TCanvas *cl = new TCanvas("cl","cl",10,10,900,300);
320  TLatex Tl; Tl.SetTextFont(43); Tl.SetTextSize(20);
321  Double_t dy = 1./3.;
322  Tl.DrawText(.01, dy, "abc#alpha#beta#gamma, #it{abc#alpha#beta#gamma} :"); Tl.DrawLatex(.7, dy, "abc#alpha#beta#gamma, #it{abc#alpha#beta#gamma}");
323  Tl.DrawText(.01, 2*dy, "#bf{bold}, #it{italic}, #bf{#it{bold italic}}, #bf{#bf{unbold}}} :"); Tl.DrawLatex(.7, 2*dy, "#bf{bold}, #it{italic}, #bf{#it{bold italic}}, #bf{#bf{unbold}}}");
324 }
325 End_Macro
326 
327 ## <a name="L13"></a> Examples
328 
329 Begin_Macro(source)
330 {
331  TCanvas ex1("ex1","Latex",500,600);
332  TLatex Tl;
333  Tl.SetTextAlign(12);
334  Tl.SetTextSize(0.04);
335  Tl.DrawLatex(0.1,0.8,"1) C(x) = d #sqrt{#frac{2}{#lambdaD}} #int^{x}_{0}cos(#frac{#pi}{2}t^{2})dt");
336  Tl.DrawLatex(0.1,0.6,"2) C(x) = d #sqrt{#frac{2}{#lambdaD}} #int^{x}cos(#frac{#pi}{2}t^{2})dt");
337  Tl.DrawLatex(0.1,0.4,"3) R = |A|^{2} = #frac{1}{2}(#[]{#frac{1}{2}+C(V)}^{2}+#[]{#frac{1}{2}+S(V)}^{2})");
338  Tl.DrawLatex(0.1,0.2,"4) F(t) = #sum_{i=-#infty}^{#infty}A(i)cos#[]{#frac{i}{t+i}}");
339  return ex1;
340 }
341 End_Macro
342 Begin_Macro(source)
343 {
344  TCanvas ex2("ex2","Latex",500,300);
345  TLatex Tl;
346  Tl.SetTextAlign(23);
347  Tl.SetTextSize(0.08);
348  Tl.DrawLatex(0.5,0.95,"e^{+}e^{-}#rightarrowZ^{0}#rightarrowI#bar{I}, q#bar{q}");
349  Tl.DrawLatex(0.5,0.75,"|#vec{a}#bullet#vec{b}|=#Sigmaa^{i}_{jk}+b^{bj}_{i}");
350  Tl.DrawLatex(0.5,0.5,"i(#partial_{#mu}#bar{#psi}#gamma^{#mu}+m#bar{#psi}=0#Leftrightarrow(#Box+m^{2})#psi=0");
351  Tl.DrawLatex(0.5,0.3,"L_{em}=eJ^{#mu}_{em}A_{#mu} , J^{#mu}_{em}=#bar{I}#gamma_{#mu}I , M^{j}_{i}=#SigmaA_{#alpha}#tau^{#alphaj}_{i}");
352  return ex2;
353 }
354 End_Macro
355 Begin_Macro(source)
356 {
357  TCanvas ex3("ex3","Latex",500,300);
358  TPaveText pt(.1,.1,.9,.9);
359  pt.AddText("#frac{2s}{#pi#alpha^{2}} #frac{d#sigma}{dcos#theta} (e^{+}e^{-} #rightarrow f#bar{f} ) = ");
360  pt.AddText("#left| #frac{1}{1 - #Delta#alpha} #right|^{2} (1+cos^{2}#theta");
361  pt.AddText("+ 4 Re #left{ #frac{2}{1 - #Delta#alpha} #chi(s) #[]{#hat{g}_{#nu}^{e}#hat{g}_{#nu}^{f} (1 + cos^{2}#theta) + 2 #hat{g}_{a}^{e}#hat{g}_{a}^{f} cos#theta) } #right}");
362  pt.SetLabel("Born equation");
363  pt.Draw();
364  return ex3;
365 }
366 End_Macro
367 
368 
369 ## <a name="L14"></a> Interface to TMathText
370 
371 The class `TMathText` is a TeX math formulae interpreter. It uses plain
372 TeX syntax and uses "\" as control instead of "#". If a piece of text containing
373 "\" is given to `TLatex` then `TMathText` is automatically invoked.
374 Therefore, as histograms' titles, axis titles, labels etc ... are drawn using
375 `TLatex`, the `TMathText` syntax can be used for them also.
376 */
377 
378 ////////////////////////////////////////////////////////////////////////////////
379 /// Default constructor.
380 
381 TLatex::TLatex()
382 {
383  fFactorSize = 1.5;
384  fFactorPos = 0.6;
385  fError = 0;
386  fShow = kFALSE;
387  fPos = 0;
388  fTabMax = 0;
389  fOriginSize = 0.04;
390  fTabSize = 0;
391  fItalic = kFALSE;
392  fLimitFactorSize = 3;
393  SetLineWidth(2);
394 }
395 
396 ////////////////////////////////////////////////////////////////////////////////
397 /// Normal constructor.
398 
399 TLatex::TLatex(Double_t x, Double_t y, const char *text)
400  :TText(x,y,text)
401 {
402  fFactorSize = 1.5;
403  fFactorPos = 0.6;
404  fError = 0;
405  fShow = kFALSE;
406  fPos = 0;
407  fTabMax = 0;
408  fOriginSize = 0.04;
409  fTabSize = 0;
410  fItalic = kFALSE;
411  fLimitFactorSize = 3;
412  SetLineWidth(2);
413 }
414 
415 ////////////////////////////////////////////////////////////////////////////////
416 /// Destructor.
417 
419 {
420 }
421 
422 ////////////////////////////////////////////////////////////////////////////////
423 /// Copy constructor.
424 
425 TLatex::TLatex(const TLatex &text) : TText(text), TAttLine(text)
426 {
427  fFactorSize = 1.5;
428  fFactorPos = 0.6;
429  fError = 0;
430  fShow = kFALSE;
431  fPos = 0;
432  fTabMax = 0;
433  fOriginSize = 0.04;
434  fTabSize = 0;
435  fItalic = kFALSE;
436  fLimitFactorSize = 3;
437  ((TLatex&)text).Copy(*this);
438 }
439 
440 ////////////////////////////////////////////////////////////////////////////////
441 ///assignment operator
442 
443 TLatex& TLatex::operator=(const TLatex& lt)
444 {
445  if(this!=&lt) {
446  TText::operator=(lt);
451  fError=lt.fError;
452  fShow=lt.fShow;
453  fTabSize=lt.fTabSize;
455  fTabSize=lt.fTabSize;
456  fTabSize=lt.fTabSize;
457  fItalic=lt.fItalic;
458  }
459  return *this;
460 }
461 
462 ////////////////////////////////////////////////////////////////////////////////
463 /// Copy this TLatex object to another TLatex.
464 
465 void TLatex::Copy(TObject &obj) const
466 {
467  ((TLatex&)obj).fFactorSize = fFactorSize;
468  ((TLatex&)obj).fFactorPos = fFactorPos;
469  ((TLatex&)obj).fLimitFactorSize = fLimitFactorSize;
470  ((TLatex&)obj).fError = fError;
471  ((TLatex&)obj).fShow = fShow;
472  ((TLatex&)obj).fTabSize = 0;
473  ((TLatex&)obj).fOriginSize = fOriginSize;
474  ((TLatex&)obj).fTabMax = fTabMax;
475  ((TLatex&)obj).fPos = fPos;
476  ((TLatex&)obj).fItalic = fItalic;
477  TText::Copy(obj);
478  TAttLine::Copy(((TAttLine&)obj));
479 }
480 
481 ////////////////////////////////////////////////////////////////////////////////
482 /// Analyse function.
483 
484 TLatex::TLatexFormSize TLatex::Anal1(TextSpec_t spec, const Char_t* t, Int_t length)
485 {
486  return Analyse(0,0,spec,t,length);
487 }
488 
489 ////////////////////////////////////////////////////////////////////////////////
490 /// Analyse and paint the TLatex formula
491 ///
492 /// It is called twice : first for calculating the size of
493 /// each portion of the formula, then to paint the formula.
494 /// When analyse finds an operator or separator, it calls
495 /// itself recursively to analyse the arguments of the operator.
496 /// when the argument is an atom (normal text), it calculates
497 /// the size of it and return it as the result.
498 /// for example : if the operator #frac{arg1}{arg2} is found :
499 /// Analyse(arg1) return the size of arg1 (width, up, down)
500 /// Analyse(arg2) return the size of arg2
501 /// now, we know the size of #frac{arg1}{arg2}:
502 /// ~~~ {.cpp}
503 /// width = max(width_arg1, width_arg2)
504 /// up = up_arg1 + down_arg1
505 /// down = up_arg2 + down_arg2
506 /// ~~~
507 /// so, when the user wants to paint a fraction at position (x,y),
508 /// the rect used for the formula is : (x,y-up,x+width,y+down)
509 ///
510 /// return size of zone occupied by the text/formula
511 /// - t : chain to be analyzed
512 /// - length : number of chars in t.
513 
515 {
516  const char *tab[] = { "alpha","beta","chi","delta","varepsilon","phi","gamma","eta","iota","varphi","kappa","lambda",
517  "mu","nu","omicron","pi","theta","rho","sigma","tau","upsilon","varomega","omega","xi","psi","zeta",
518  "Alpha","Beta","Chi","Delta","Epsilon","Phi","Gamma","Eta","Iota","vartheta",
519  "Kappa","Lambda","Mu","Nu","Omicron","Pi","Theta","Rho","Sigma","Tau",
520  "Upsilon","varsigma","Omega","Xi","Psi","Zeta","varUpsilon","epsilon"};
521 
522  const char *tab2[] = { "leq","/","infty","voidb","club","diamond","heart",
523  "spade","leftrightarrow","leftarrow","uparrow","rightarrow",
524  "downarrow","circ","pm","doublequote","geq","times","propto",
525  "partial","bullet","divide","neq","equiv","approx","3dots",
526  "cbar","topbar","downleftarrow","aleph","Jgothic","Rgothic","voidn",
527  "otimes","oplus","oslash","cap","cup","supset","supseteq",
528  "notsubset","subset","subseteq","in","notin","angle","nabla",
529  "oright","ocopyright","trademark","prod","surd","upoint","corner","wedge",
530  "vee","Leftrightarrow","Leftarrow","Uparrow","Rightarrow",
531  "Downarrow","diamond","LT","void1","copyright","void3","sum",
532  "arctop","lbar","arcbottom","topbar","void8", "bottombar","arcbar",
533  "ltbar","AA","aa","void06","GT","int","forall","exists" };
534 
535  const char *tab3[] = { "bar","vec","dot","hat","ddot","acute","grave","check","tilde","slash"};
536 
537  if (fError != 0) return TLatexFormSize(0,0,0);
538 
539  Int_t nBlancDeb=0,nBlancFin=0,l_nBlancDeb=0,l_nBlancFin=0;
540  Int_t i,k;
541  Int_t min=0, max=0;
542  Bool_t cont = kTRUE;
543  while(cont) {
544  // count leading blanks
545  //while(nBlancDeb+nBlancFin<length && t[nBlancDeb]==' ') nBlancDeb++;
546 
547  if (nBlancDeb==length) return TLatexFormSize(0,0,0); // empty string
548 
549  // count trailing blanks
550  //while(nBlancDeb+nBlancFin<length && t[length-nBlancFin-1]==' ') nBlancFin++;
551 
552  if (nBlancDeb==l_nBlancDeb && nBlancFin==l_nBlancFin) cont = kFALSE;
553 
554  // remove characters { }
555  if (t[nBlancDeb]=='{' && t[length-nBlancFin-1]=='}') {
556  Int_t nBrackets = 0;
557  Bool_t sameBrackets = kTRUE;
558  for(i=nBlancDeb;i<length-nBlancFin;i++) {
559  if (t[i] == '{' && !(i>0 && t[i-1] == '@')) nBrackets++;
560  if (t[i] == '}' && t[i-1]!= '@') nBrackets--;
561  if (nBrackets==0 && i<length-nBlancFin-2) {
562  sameBrackets=kFALSE;
563  break;
564  }
565  }
566 
567  if (sameBrackets) {
568  // begin and end brackets match
569  nBlancDeb++;
570  nBlancFin++;
571  if (nBlancDeb+nBlancFin==length) return TLatexFormSize(0,0,0); // empty string
572  cont = kTRUE;
573  }
574 
575  }
576 
577  l_nBlancDeb = nBlancDeb;
578  l_nBlancFin = nBlancFin;
579  }
580 
581  // make a copy of the current processed chain of characters
582  // removing leading and trailing blanks
583  length -= nBlancFin+nBlancDeb; // length of string without blanks
584  Char_t* text = new Char_t[length+1];
585  strncpy(text,t+nBlancDeb,length);
586  text[length] = 0;
587 
588  // compute size of subscripts and superscripts
589  Double_t indiceSize = spec.fSize/fFactorSize;
591  indiceSize = spec.fSize;
592  // subtract 0.001 because of rounding errors
593  TextSpec_t specNewSize = spec;
594  specNewSize.fSize = indiceSize;
595 
596  // recherche des operateurs
597  Int_t opPower = -1; // Position of first ^ (power)
598  Int_t opUnder = -1; // Position of first _ (indice)
599  Int_t opFrac = -1; // Position of first \frac
600  Int_t opSqrt = -1; // Position of first \sqrt
601  Int_t nBrackets = 0; // Nesting level in { }
602  Int_t nCroch = 0; // Nesting level in [ ]
603  Int_t opCurlyCurly = -1; // Position of first }{
604  Int_t opSquareCurly = -1; // Position of first ]{
605  Int_t opCloseCurly = -2; // Position of first }
606  Int_t opColor = -1; // Position of first #color
607  Int_t opFont = -1; // Position of first #font
608  Int_t opScale = -1; // Position of first #scale
609  Int_t opGreek = -1; // Position of a Greek letter
610  Int_t opSpec = -1; // position of a special character
611  Int_t opAbove = -1; // position of a vector/overline
612  Int_t opSquareBracket = 0 ; // position of a "[]{" operator (#[]{arg})
613  Int_t opBigCurly = 0 ; // position of a "{}{" operator (big curly bracket #{}{arg})
614  Int_t opAbs = 0 ; // position of a "||{" operator (absolute value) (#||{arg})
615  Int_t opParen = 0 ; // position of a "(){" operator (big parenthesis #(){arg})
616  Int_t abovePlace = 0 ; // true if subscripts must be written above and not after
617  Int_t opBox = 0 ; // position of #Box
618  Int_t opPerp = 0; // position of #perp
619  Int_t opOdot = 0; // position of #odot
620  Int_t opHbar = 0; // position of #hbar
621  Int_t opMinus = 0; // position of #minus
622  Int_t opPlus = 0; // position of #plus
623  Int_t opMp = 0; // position of #mp
624  Int_t opBackslash = 0; // position of #backslash
625  Int_t opParallel = 0; // position of #parallel
626  Int_t opSplitLine = -1; // Position of first #splitline
627  Int_t opKern = -1; // Position of first #kern
628  Int_t opLower = -1; // Position of first #lower
629  Int_t opBf = -1; // Position of first #bf
630  Int_t opIt = -1; // Position of first #it
631  Int_t opMbox = -1; // Position of first #mbox
632 
633  Bool_t opFound = kFALSE;
634  Bool_t quote1 = kFALSE, quote2 = kFALSE ;
635 
636  for(i=0;i<length;i++) {
637  switch (text[i]) {
638  case '\'' : quote1 = !quote1 ; break ;
639  case '"' : quote2 = !quote2 ; break ;
640  }
641  //if (quote1 || quote2) continue ;
642  switch (text[i]) {
643  case '{':
644  if (nCroch==0) {
645  if (!(i>0 && text[i-1] == '@')) nBrackets++;
646  }
647  break;
648  case '}':
649  if (nCroch==0) {
650  if (!(i>0 && text[i-1] == '@')) nBrackets--;
651  if (nBrackets==0) {
652  if (i<length-1) if (text[i+1]=='{' && opCurlyCurly==-1) opCurlyCurly=i;
653  if (i<length-2) {
654  if (text[i+1]!='{' && !(text[i+2]=='{' && (text[i+1]=='^' || text[i+1]=='_'))
655  && opCloseCurly==-2) opCloseCurly=i;
656  }
657  else if (i<length-1) {
658  if (text[i+1]!='{' && opCloseCurly==-2) opCloseCurly=i;
659  }
660  else if (opCloseCurly==-2) opCloseCurly=i;
661  }
662  }
663  break;
664  case '[':
665  if (nBrackets==0) {
666  if (!(i>0 && text[i-1] == '@')) nCroch++;
667  }
668  break;
669  case ']':
670  if (nBrackets==0) {
671  if (!(i>0 && text[i-1] == '@')) nCroch--;
672  if (nCroch<0) {
673  // more "]" than "["
674  fError = "Missing \"[\"";
675  delete [] text;
676  return TLatexFormSize(0,0,0);
677  }
678  }
679  break;
680  }
681  if (length>i+1) {
682  Char_t buf[3];
683  strncpy(buf,&text[i],2);
684  if (strncmp(buf,"^{",2)==0) {
685  if (opPower==-1 && nBrackets==0 && nCroch==0) opPower=i;
686  if (i>3) {
687  Char_t buf1[5];
688  strncpy(buf1,&text[i-4],4);
689  if (strncmp(buf1,"#int",4)==0) {
690  abovePlace = 1;
691  if (i>4 && opCloseCurly==-2) opCloseCurly=i-5;
692  }
693  if (strncmp(buf1,"#sum",4)==0) {
694  abovePlace = 2;
695  if (i>4 && opCloseCurly==-2) opCloseCurly=i-5;
696  }
697  }
698  }
699  if (strncmp(buf,"_{",2)==0) {
700  if (opUnder==-1 && nBrackets==0 && nCroch==0) opUnder=i;
701  if (i>3) {
702  Char_t buf2[5];
703  strncpy(buf2,&text[i-4],4);
704  if (strncmp(buf2,"#int",4)==0) {
705  abovePlace = 1;
706  if (i>4 && opCloseCurly==-2) opCloseCurly=i-5;
707  }
708  if (strncmp(buf2,"#sum",4)==0) {
709  abovePlace = 2;
710  if (i>4 && opCloseCurly==-2) opCloseCurly=i-5;
711  }
712  }
713  }
714  if (strncmp(buf,"]{",2)==0)
715  if (opSquareCurly==-1 && nBrackets==0 && nCroch==0) opSquareCurly=i;
716  }
717  // detect other operators
718  if (text[i]=='\\' || (text[i]=='#' && !opFound && nBrackets==0 && nCroch==0)) {
719 
720  if (length>i+10) {
721  Char_t buf[11];
722  strncpy(buf,&text[i+1],10);
723  if (strncmp(buf,"splitline{",10)==0) {
724  opSplitLine=i; opFound = kTRUE;
725  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
726  continue;
727  }
728  }
729  if (length>i+9) {
730  Char_t buf[10];
731  strncpy(buf,&text[i+1],9);
732  if (!opBackslash && strncmp(buf,"backslash",9)==0) {
733  opBackslash=1; opFound = kTRUE;
734  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
735  continue;
736  }
737  }
738  if (length>i+8) {
739  Char_t buf[9];
740  strncpy(buf,&text[i+1],8);
741  if (!opParallel && strncmp(buf,"parallel",8)==0) {
742  opParallel=1; opFound = kTRUE;
743  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
744  continue;
745  }
746  }
747  if (length>i+6) {
748  Char_t buf[7];
749  strncpy(buf,&text[i+1],6);
750  if (strncmp(buf,"lower[",6)==0 || strncmp(buf,"lower{",6)==0) {
751  opLower=i; opFound = kTRUE;
752  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
753  continue ;
754  }
755  if (strncmp(buf,"scale[",6)==0 || strncmp(buf,"scale{",6)==0) {
756  opScale=i; opFound = kTRUE;
757  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
758  continue ;
759  }
760  if (strncmp(buf,"color[",6)==0 || strncmp(buf,"color{",6)==0) {
761  opColor=i; opFound = kTRUE;
762  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
763  continue ;
764  }
765  }
766  if (length>i+5) {
767  Char_t buf[6];
768  strncpy(buf,&text[i+1],5);
769  if (strncmp(buf,"frac{",5)==0) {
770  opFrac=i; opFound = kTRUE;
771  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
772  continue;
773  }
774  if (strncmp(buf,"sqrt{",5)==0 || strncmp(buf,"sqrt[",5)==0) {
775  opSqrt=i; opFound = kTRUE;
776  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
777  continue;
778  }
779  if (strncmp(buf,"font{",5)==0 || strncmp(buf,"font[",5)==0) {
780  opFont=i; opFound = kTRUE;
781  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
782  continue;
783  }
784  if (strncmp(buf,"kern[",5)==0 || strncmp(buf,"kern{",5)==0) {
785  opKern=i; opFound = kTRUE;
786  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
787  continue ;
788  }
789  if (!opMinus && strncmp(buf,"minus",5)==0) {
790  opMinus=1; opFound = kTRUE;
791  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
792  continue;
793  }
794  if (strncmp(buf,"mbox[",5)==0 || strncmp(buf,"mbox{",5)==0) {
795  opMbox=i; opFound = kTRUE;
796  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
797  continue ;
798  }
799  }
800  if (length>i+4) {
801  Char_t buf[5];
802  strncpy(buf,&text[i+1],4);
803  if (!opOdot && strncmp(buf,"odot",4)==0) {
804  opOdot=1; opFound = kTRUE;
805  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
806  continue;
807  }
808  if (!opHbar && strncmp(buf,"hbar",4)==0) {
809  opHbar=1; opFound = kTRUE;
810  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
811  continue;
812  }
813  if (!opPerp && strncmp(buf,"perp",4)==0) {
814  opPerp=1; opFound = kTRUE;
815  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
816  continue;
817  }
818  if (!opPlus && strncmp(buf,"plus",4)==0) {
819  opPlus=1; opFound = kTRUE;
820  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
821  continue;
822  }
823  }
824  if (length>i+3) {
825  Char_t buf[4];
826  strncpy(buf,&text[i+1],3);
827  buf[3] = 0;
828  if (strncmp(buf,"[]{",3)==0) {
829  opSquareBracket=1; opFound = kTRUE;
830  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
831  continue;
832  }
833  if (strncmp(buf,"{}{",3)==0 ) {
834  opBigCurly=1; opFound = kTRUE;
835  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
836  continue;
837  }
838  if (strncmp(buf,"||{",3)==0) {
839  opAbs=1; opFound = kTRUE;
840  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
841  continue;
842  }
843  if (strncmp(buf,"(){",3)==0) {
844  opParen=1; opFound = kTRUE;
845  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
846  continue;
847  }
848  if (!opBox && strncmp(buf,"Box",3)==0) {
849  opBox=1; opFound = kTRUE;
850  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
851  continue;
852  }
853  if (strncmp(buf,"bf[",3)==0 || strncmp(buf,"bf{",3)==0) {
854  opBf=i; opFound = kTRUE;
855  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
856  continue ;
857  }
858  if (strncmp(buf,"it[",3)==0 || strncmp(buf,"it{",3)==0) {
859  opIt=i; opFound = kTRUE;
860  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
861  continue ;
862  }
863  }
864  if (length>i+2) {
865  Char_t buf[3];
866  strncpy(buf,&text[i+1],2);
867  if (!opMp && strncmp(buf,"mp",2)==0) {
868  opMp=1; opFound = kTRUE;
869  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
870  continue;
871  }
872  }
873  for(k=0;k<54;k++) {
874  if (!opFound && UInt_t(length)>i+strlen(tab[k])) {
875  if (strncmp(&text[i+1],tab[k],strlen(tab[k]))==0) {
876  opGreek=k;
877  opFound = kTRUE;
878  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
879  }
880  }
881  }
882  for(k=0;k<10;k++) {
883  if (!opFound && UInt_t(length)>i+strlen(tab3[k])) {
884  if (strncmp(&text[i+1],tab3[k],strlen(tab3[k]))==0) {
885  opAbove=k;
886  opFound = kTRUE;
887  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
888  }
889  }
890  }
891  UInt_t lastsize = 0;
892  if (!opFound)
893  for(k=0;k<82;k++) {
894  if ((opSpec==-1 || strlen(tab2[k])>lastsize) && UInt_t(length)>i+strlen(tab2[k])) {
895  if (strncmp(&text[i+1],tab2[k],strlen(tab2[k]))==0) {
896  lastsize = strlen(tab2[k]);
897  opSpec=k;
898  opFound = kTRUE;
899  if (i>0 && opCloseCurly==-2) opCloseCurly=i-1;
900  }
901  }
902  }
903  }
904  }
905 
906  TLatexFormSize fs1;
907  TLatexFormSize fs2;
908  TLatexFormSize fs3;
909  TLatexFormSize result;
910 
911  // analysis of operators found
912  if (opCloseCurly>-1 && opCloseCurly<length-1) { // separator } found
913  if(!fShow) {
914  fs1 = Anal1(spec,text,opCloseCurly+1);
915  fs2 = Anal1(spec,text+opCloseCurly+1,length-opCloseCurly-1);
916  Savefs(&fs1);
917  } else {
918  fs1 = Readfs();
919  Analyse(x+fs1.Width(),y,spec,text+opCloseCurly+1,length-opCloseCurly-1);
920  Analyse(x,y,spec,text,opCloseCurly+1);
921  }
922  result = fs1+fs2;
923  }
924 
925  else if (opPower>-1 && opUnder>-1) { // ^ and _ found
926  min = TMath::Min(opPower,opUnder);
927  max = TMath::Max(opPower,opUnder);
928  Double_t xfpos = 0. ; //GetHeight()*spec.fSize/5.;
929  Double_t prop=1, propU=1; // scale factor for #sum & #int
930  switch (abovePlace) {
931  case 1 :
932  prop = .8 ; propU = 1.75 ; // Int
933  break;
934  case 2:
935  prop = .9 ; propU = 1.75 ; // Sum
936  break;
937  }
938  // propU acts on upper number
939  // when increasing propU value, the upper indice position is higher
940  // when increasing prop values, the lower indice position is lower
941 
942  if (!fShow) {
943  Int_t ltext = min ;
944  if (min >= 2 && strncmp(&text[min-2],"{}",2)==0) {
945  // upper and lower indice before the character
946  // like with chemical element
947  snprintf(&text[ltext-2],length-(ltext-2),"I ") ;
948  ltext-- ;
949  }
950  fs1 = Anal1(spec,text,ltext);
951  fs2 = Anal1(specNewSize,text+min+1,max-min-1);
952  fs3 = Anal1(specNewSize,text+max+1,length-max-1);
953  Savefs(&fs1);
954  Savefs(&fs2);
955  Savefs(&fs3);
956  } else {
957  fs3 = Readfs();
958  fs2 = Readfs();
959  fs1 = Readfs();
960  Double_t pos = 0;
961  if (!abovePlace) {
962  Double_t addW = fs1.Width()+xfpos, addH1, addH2;
963  if (opPower<opUnder) {
964  addH1 = -fs1.Over()*(fFactorPos)-fs2.Under();
965  addH2 = fs1.Under()+fs3.Over()*(fFactorPos);
966  } else {
967  addH1 = fs1.Under()+fs2.Over()*(fFactorPos);
968  addH2 = -fs1.Over()*(fFactorPos)-fs3.Under();
969  }
970  Analyse(x+addW,y+addH2,specNewSize,text+max+1,length-max-1);
971  Analyse(x+addW,y+addH1,specNewSize,text+min+1,max-min-1);
972  } else {
973  Double_t addW1, addW2, addH1, addH2;
974  Double_t m = TMath::Max(fs1.Width(),TMath::Max(fs2.Width(),fs3.Width()));
975  pos = (m-fs1.Width())/2;
976  if (opPower<opUnder) {
977  addH1 = -fs1.Over()*propU-fs2.Under();
978  addW1 = (m-fs2.Width())/2;
979  addH2 = fs1.Under()*prop+fs3.Over();
980  addW2 = (m-fs3.Width())/2;
981  } else {
982  addH1 = fs1.Under()*prop+fs2.Over();
983  addW1 = (m-fs2.Width())/2;
984  addH2 = -fs1.Over()*propU-fs3.Under();
985  addW2 = (m-fs3.Width())/2;
986  }
987 
988  Analyse(x+addW2,y+addH2,specNewSize,text+max+1,length-max-1);
989  Analyse(x+addW1,y+addH1,specNewSize,text+min+1,max-min-1);
990  }
991 
992  if (min >= 2 && strncmp(&text[min-2],"{}",2)==0) {
993  snprintf(&text[min-2],length-(min-2)," ") ;
994  Analyse(x+pos,y,spec,text,min-1);
995  } else {
996  Analyse(x+pos,y,spec,text,min);
997  }
998  }
999 
1000  if (!abovePlace) {
1001  if (opPower<opUnder) {
1002  result.Set(fs1.Width()+xfpos+TMath::Max(fs2.Width(),fs3.Width()),
1003  fs1.Over()*fFactorPos+fs2.Height(),
1004  fs1.Under()+fs3.Height()-fs3.Over()*(1-fFactorPos));
1005  } else {
1006  result.Set(fs1.Width()+xfpos+TMath::Max(fs2.Width(),fs3.Width()),
1007  fs1.Over()*fFactorPos+fs3.Height(),
1008  fs1.Under()+fs2.Height()-fs2.Over()*(1-fFactorPos));
1009  }
1010  } else {
1011  if (opPower<opUnder) {
1012  result.Set(TMath::Max(fs1.Width(),TMath::Max(fs2.Width(),fs3.Width())),
1013  fs1.Over()*propU+fs2.Height(),fs1.Under()*prop+fs3.Height());
1014  } else {
1015  result.Set(TMath::Max(fs1.Width(),TMath::Max(fs2.Width(),fs3.Width())),
1016  fs1.Over()*propU+fs3.Height(),fs1.Under()*prop+fs2.Height());
1017  }
1018  }
1019  }
1020  else if (opPower>-1) { // ^ found
1021  Double_t prop=1;
1022  Double_t xfpos = 0. ; //GetHeight()*spec.fSize/5. ;
1023  switch (abovePlace) {
1024  case 1 : //int
1025  prop = 1.75 ; break ;
1026  case 2 : // sum
1027  prop = 1.75; break ;
1028  }
1029  // When increasing prop, the upper indice position is higher
1030  if (!fShow) {
1031  Int_t ltext = opPower ;
1032  if (ltext >= 2 && strncmp(&text[ltext-2],"{}",2)==0) {
1033  // upper and lower indice before the character
1034  // like with chemical element
1035  snprintf(&text[ltext-2],length-(ltext-2),"I ") ;
1036  ltext-- ;
1037  }
1038  fs1 = Anal1(spec,text,ltext);
1039  fs2 = Anal1(specNewSize,text+opPower+1,length-opPower-1);
1040  Savefs(&fs1);
1041  Savefs(&fs2);
1042  } else {
1043  fs2 = Readfs();
1044  fs1 = Readfs();
1045  Int_t pos = 0;
1046  if (!abovePlace){
1047  Double_t over = fs1.Over();
1048  if (over <= 0) over = 1.5*fs2.Over();
1049  Analyse(x+fs1.Width()+xfpos,y-over*fFactorPos-fs2.Under(),specNewSize,text+opPower+1,length-opPower-1);
1050  } else {
1051  Int_t pos2=0;
1052  if (fs2.Width()>fs1.Width())
1053  pos=Int_t((fs2.Width()-fs1.Width())/2);
1054  else
1055  pos2=Int_t((fs1.Width()-fs2.Width())/2);
1056 
1057  Analyse(x+pos2,y-fs1.Over()*prop-fs2.Under(),specNewSize,text+opPower+1,length-opPower-1);
1058  }
1059  if (opPower >= 2 && strncmp(&text[opPower-2],"{}",2)==0) {
1060  snprintf(&text[opPower-2],length-(opPower-2)," ") ;
1061  Analyse(x+pos,y,spec,text,opPower-1);
1062  } else {
1063  Analyse(x+pos,y,spec,text,opPower);
1064  }
1065  }
1066 
1067  if (!abovePlace)
1068  result.Set(fs1.Width()+xfpos+fs2.Width(),
1069  fs1.Over()*fFactorPos+fs2.Over(),fs1.Under());
1070  else
1071  result.Set(TMath::Max(fs1.Width(),fs2.Width()),fs1.Over()*prop+fs2.Height(),fs1.Under());
1072 
1073  }
1074  else if (opUnder>-1) { // _ found
1075  Double_t prop = .9; // scale factor for #sum & #frac
1076  Double_t xfpos = 0.;//GetHeight()*spec.fSize/5. ;
1077  Double_t fpos = fFactorPos ;
1078  // When increasing prop, the lower indice position is lower
1079  if(!fShow) {
1080  Int_t ltext = opUnder ;
1081  if (ltext >= 2 && strncmp(&text[ltext-2],"{}",2)==0) {
1082  // upper and lower indice before the character
1083  // like with chemical element
1084  snprintf(&text[ltext-2],length-(ltext-2),"I ") ;
1085  ltext-- ;
1086  }
1087  fs1 = Anal1(spec,text,ltext);
1088  fs2 = Anal1(specNewSize,text+opUnder+1,length-opUnder-1);
1089  Savefs(&fs1);
1090  Savefs(&fs2);
1091  } else {
1092  fs2 = Readfs();
1093  fs1 = Readfs();
1094  Int_t pos = 0;
1095  if (!abovePlace)
1096  Analyse(x+fs1.Width()+xfpos,y+fs1.Under()+fs2.Over()*fpos,specNewSize,text+opUnder+1,length-opUnder-1);
1097  else {
1098  Int_t pos2=0;
1099  if (fs2.Width()>fs1.Width())
1100  pos=Int_t((fs2.Width()-fs1.Width())/2);
1101  else
1102  pos2=Int_t((fs1.Width()-fs2.Width())/2);
1103 
1104  Analyse(x+pos2,y+fs1.Under()*prop+fs2.Over(),specNewSize,text+opUnder+1,length-opUnder-1);
1105  }
1106  if (opUnder >= 2 && strncmp(&text[opUnder-2],"{}",2)==0) {
1107  snprintf(&text[opUnder-2],length-(opUnder-2)," ") ;
1108  Analyse(x+pos,y,spec,text,opUnder-1);
1109  } else {
1110  Analyse(x+pos,y,spec,text,opUnder);
1111  }
1112  }
1113  if (!abovePlace)
1114  result.Set(fs1.Width()+xfpos+fs2.Width(),fs1.Over(),
1115  fs1.Under()+fs2.Under()+fs2.Over()*fpos);
1116  else
1117  result.Set(TMath::Max(fs1.Width(),fs2.Width()),fs1.Over(),fs1.Under()*prop+fs2.Height());
1118  }
1119  else if (opBox) {
1120  Double_t square = GetHeight()*spec.fSize/2;
1121  if (!fShow) {
1122  fs1 = Anal1(spec,text+4,length-4);
1123  } else {
1124  fs1 = Analyse(x+square,y,spec,text+4,length-4);
1125  Double_t adjust = GetHeight()*spec.fSize/20;
1126  Double_t x1 = x+adjust ;
1127  Double_t x2 = x-adjust+square ;
1128  Double_t y1 = y;
1129  Double_t y2 = y-square+adjust;
1130  DrawLine(x1,y1,x2,y1,spec);
1131  DrawLine(x2,y1,x2,y2,spec);
1132  DrawLine(x2,y2,x1,y2,spec);
1133  DrawLine(x1,y2,x1,y1,spec);
1134  }
1135  result = fs1 + TLatexFormSize(square,square,0);
1136  }
1137  else if (opOdot) {
1138  Double_t square = GetHeight()*spec.fSize/2;
1139  if (!fShow) {
1140  fs1 = Anal1(spec,text+5,length-5);
1141  } else {
1142  fs1 = Analyse(x+1.3*square,y,spec,text+5,length-5);
1143  Double_t adjust = GetHeight()*spec.fSize/20;
1144  Double_t r1 = 0.62*square;
1145  Double_t y1 = y-0.3*square-adjust;
1146  DrawCircle(x+0.6*square,y1,r1,spec) ;
1147  DrawCircle(x+0.6*square,y1,r1/100,spec) ;
1148  }
1149  result = fs1 + TLatexFormSize(square,square,0);
1150  }
1151  else if (opHbar) {
1152  Double_t square = GetHeight()*spec.fSize/2;
1153  if (!fShow) {
1154  fs1 = Anal1(spec,text+5,length-5);
1155  } else {
1156  fs1 = Analyse(x+square,y,spec,text+5,length-5);
1157  TText hbar;
1158  hbar.SetTextFont(12);
1159  hbar.SetTextColor(spec.fColor);
1160  hbar.SetTextSize(spec.fSize);
1161  hbar.SetTextAngle(fTextAngle);
1162  hbar.SetTextAlign(11);
1163  Double_t xOrigin = (Double_t)gPad->XtoAbsPixel(fX);
1164  Double_t yOrigin = (Double_t)gPad->YtoAbsPixel(fY);
1165  Double_t angle = kPI*spec.fAngle/180.;
1166  Double_t xx = gPad->AbsPixeltoX(Int_t((x-xOrigin)*TMath::Cos(angle)+(y-yOrigin)*TMath::Sin(angle)+xOrigin));
1167  Double_t yy = gPad->AbsPixeltoY(Int_t((x-xOrigin)*TMath::Sin(-angle)+(y-yOrigin)*TMath::Cos(angle)+yOrigin));
1168  hbar.PaintText(xx,yy,"h");
1169  DrawLine(x,y-0.8*square,x+0.75*square,y-square,spec);
1170  }
1171  result = fs1 + TLatexFormSize(square,square,0);
1172  }
1173  else if (opMinus) {
1174  Double_t square = GetHeight()*spec.fSize/2;
1175  if (!fShow) {
1176  fs1 = Anal1(spec,text+6,length-6);
1177  } else {
1178  fs1 = Analyse(x+square,y,spec,text+6,length-6);
1179  TText minus;
1180  minus.SetTextFont(122);
1181  minus.SetTextColor(spec.fColor);
1182  minus.SetTextSize(spec.fSize);
1183  minus.SetTextAngle(fTextAngle);
1184  minus.SetTextAlign(11);
1185  Double_t xOrigin = (Double_t)gPad->XtoAbsPixel(fX);
1186  Double_t yOrigin = (Double_t)gPad->YtoAbsPixel(fY);
1187  Double_t angle = kPI*spec.fAngle/180.;
1188  Double_t xx = gPad->AbsPixeltoX(Int_t((x-xOrigin)*TMath::Cos(angle)+(y-yOrigin)*TMath::Sin(angle)+xOrigin));
1189  Double_t yy = gPad->AbsPixeltoY(Int_t((x-xOrigin)*TMath::Sin(-angle)+(y-yOrigin)*TMath::Cos(angle)+yOrigin));
1190  minus.PaintText(xx,yy,"-");
1191  }
1192  result = fs1 + TLatexFormSize(square,square,0);
1193  }
1194  else if (opPlus) {
1195  Double_t square = GetHeight()*spec.fSize/2;
1196  if (!fShow) {
1197  fs1 = Anal1(spec,text+5,length-5);
1198  } else {
1199  fs1 = Analyse(x+square,y,spec,text+5,length-5);
1200  TText plus;
1201  plus.SetTextFont(122);
1202  plus.SetTextColor(spec.fColor);
1203  plus.SetTextSize(spec.fSize);
1204  plus.SetTextAngle(fTextAngle);
1205  plus.SetTextAlign(11);
1206  Double_t xOrigin = (Double_t)gPad->XtoAbsPixel(fX);
1207  Double_t yOrigin = (Double_t)gPad->YtoAbsPixel(fY);
1208  Double_t angle = kPI*spec.fAngle/180.;
1209  Double_t xx = gPad->AbsPixeltoX(Int_t((x-xOrigin)*TMath::Cos(angle)+(y-yOrigin)*TMath::Sin(angle)+xOrigin));
1210  Double_t yy = gPad->AbsPixeltoY(Int_t((x-xOrigin)*TMath::Sin(-angle)+(y-yOrigin)*TMath::Cos(angle)+yOrigin));
1211  plus.PaintText(xx,yy,"+");
1212  }
1213  result = fs1 + TLatexFormSize(square,square,0);
1214  }
1215  else if (opMp) {
1216  Double_t square = GetHeight()*spec.fSize/2;
1217  if (!fShow) {
1218  fs1 = Anal1(spec,text+3,length-3);
1219  } else {
1220  fs1 = Analyse(x+square,y,spec,text+3,length-3);
1221  TText mp;
1222  mp.SetTextFont(122);
1223  mp.SetTextColor(spec.fColor);
1224  mp.SetTextSize(spec.fSize);
1225  mp.SetTextAngle(fTextAngle+180);
1226  mp.SetTextAlign(11);
1227  Double_t xOrigin = (Double_t)gPad->XtoAbsPixel(fX);
1228  Double_t yOrigin = (Double_t)gPad->YtoAbsPixel(fY);
1229  Double_t angle = kPI*spec.fAngle/180.;
1230  Double_t xx = gPad->AbsPixeltoX(Int_t((x+square-xOrigin)*TMath::Cos(angle)+(y-1.25*square-yOrigin)*TMath::Sin(angle)+xOrigin));
1231  Double_t yy = gPad->AbsPixeltoY(Int_t((x+square-xOrigin)*TMath::Sin(-angle)+(y-1.25*square-yOrigin)*TMath::Cos(angle)+yOrigin));
1232  mp.PaintText(xx,yy,"\261");
1233  }
1234  result = fs1 + TLatexFormSize(square,square,0);
1235  }
1236  else if (opPerp) {
1237  Double_t square = GetHeight()*spec.fSize/1.4;
1238  if (!fShow) {
1239  fs1 = Anal1(spec,text+5,length-5);
1240  } else {
1241  fs1 = Analyse(x+0.5*square,y,spec,text+5,length-5);
1242  Double_t x0 = x + 0.50*square;
1243  Double_t x1 = x0 - 0.48*square;
1244  Double_t x2 = x0 + 0.48*square;
1245  Double_t y1 = y + 0.6*square;
1246  Double_t y2 = y1 - 1.3*square;
1247  DrawLine(x1,y1,x2,y1,spec);
1248  DrawLine(x0,y1,x0,y2,spec);
1249  }
1250  result = fs1;
1251  }
1252  else if (opBackslash) {
1253  Double_t square = GetHeight()*spec.fSize/2;
1254  if (!fShow) {
1255  fs1 = Anal1(spec,text+10,length-10);
1256  } else {
1257  fs1 = Analyse(x+square,y,spec,text+10,length-10);
1258  TText bs;
1259  bs.SetTextFont(GetTextFont());
1260  bs.SetTextColor(spec.fColor);
1261  bs.SetTextSize(spec.fSize);
1263  bs.SetTextAlign(11);
1264  Double_t xOrigin = (Double_t)gPad->XtoAbsPixel(fX);
1265  Double_t yOrigin = (Double_t)gPad->YtoAbsPixel(fY);
1266  Double_t angle = kPI*spec.fAngle/180.;
1267  Double_t xx = gPad->AbsPixeltoX(Int_t((x-xOrigin)*TMath::Cos(angle)+(y-yOrigin)*TMath::Sin(angle)+xOrigin));
1268  Double_t yy = gPad->AbsPixeltoY(Int_t((x-xOrigin)*TMath::Sin(-angle)+(y-yOrigin)*TMath::Cos(angle)+yOrigin));
1269  bs.PaintText(xx,yy,"\\");
1270  }
1271  result = fs1 + TLatexFormSize(square,square,0);
1272  }
1273  else if (opParallel) {
1274  Double_t square = GetHeight()*spec.fSize/1.4;
1275  if (!fShow) {
1276  fs1 = Anal1(spec,text+9,length-9);
1277  } else {
1278  fs1 = Analyse(x+0.5*square,y,spec,text+9,length-9);
1279  Double_t x1 = x + 0.15*square;
1280  Double_t x2 = x + 0.45*square;
1281  Double_t y1 = y + 0.3*square;
1282  Double_t y2 = y1- 1.3*square;
1283  DrawLine(x1,y1,x1,y2,spec);
1284  DrawLine(x2,y1,x2,y2,spec);
1285  }
1286  result = fs1 + TLatexFormSize(square,square,0);
1287  }
1288  else if (opGreek>-1) {
1289  TextSpec_t newSpec = spec;
1290  newSpec.fFont = fItalic ? 152 : 122;
1291  char letter = 97 + opGreek;
1292  Double_t yoffset = 0.; // Greek letter too low
1293  if (opGreek>25) letter -= 58;
1294  if (opGreek == 52) letter = '\241'; //varUpsilon
1295  if (opGreek == 53) letter = '\316'; //epsilon
1296  if (!fShow) {
1297  fs1 = Anal1(newSpec,&letter,1);
1298  fs2 = Anal1(spec,text+strlen(tab[opGreek])+1,length-strlen(tab[opGreek])-1);
1299  Savefs(&fs1);
1300  } else {
1301  fs1 = Readfs();
1302  Analyse(x+fs1.Width(),y,spec,text+strlen(tab[opGreek])+1,length-strlen(tab[opGreek])-1);
1303  Analyse(x,y-yoffset,newSpec,&letter,1);
1304  }
1305  fs1.AddOver(TLatexFormSize(0,yoffset,0)) ;
1306  result = fs1+fs2;
1307  }
1308 
1309  else if (opSpec>-1) {
1310  TextSpec_t newSpec = spec;
1311  newSpec.fFont = fItalic ? 152 : 122;
1312  char letter = '\243' + opSpec;
1313  if(opSpec == 75 || opSpec == 76) {
1314  newSpec.fFont = GetTextFont();
1315  if (opSpec == 75) letter = '\305'; // AA Angstroem
1316  if (opSpec == 76) letter = '\345'; // aa Angstroem
1317  }
1318  if(opSpec == 80 || opSpec == 81) {
1319  if (opSpec == 80) letter = '\042'; // #forall
1320  if (opSpec == 81) letter = '\044'; // #exists
1321  }
1322  Double_t props, propi;
1323  props = 1.8 ; // scale factor for #sum(66)
1324  propi = 2.3 ; // scale factor for #int(79)
1325 
1326  if (opSpec==66 ) {
1327  newSpec.fSize = spec.fSize*props;
1328  } else if (opSpec==79) {
1329  newSpec.fSize = spec.fSize*propi;
1330  }
1331  if (!fShow) {
1332  fs1 = Anal1(newSpec,&letter,1);
1333  if (opSpec == 79 || opSpec == 66)
1334  fs1.Set(fs1.Width(),fs1.Over()*0.45,fs1.Over()*0.45);
1335 
1336  fs2 = Anal1(spec,text+strlen(tab2[opSpec])+1,length-strlen(tab2[opSpec])-1);
1337  Savefs(&fs1);
1338  } else {
1339  fs1 = Readfs();
1340  Analyse(x+fs1.Width(),y,spec,text+strlen(tab2[opSpec])+1,length-strlen(tab2[opSpec])-1);
1341  if (opSpec!=66 && opSpec!=79)
1342  Analyse(x,y,newSpec,&letter,1);
1343  else {
1344  Analyse(x,y+fs1.Under()/2.,newSpec,&letter,1);
1345  }
1346  }
1347  result = fs1+fs2;
1348  }
1349  else if (opAbove>-1) {
1350  if (!fShow) {
1351  fs1 = Anal1(spec,text+strlen(tab3[opAbove])+1,length-strlen(tab3[opAbove])-1);
1352  Savefs(&fs1);
1353  } else {
1354  fs1 = Readfs();
1355  Analyse(x,y,spec,text+strlen(tab3[opAbove])+1,length-strlen(tab3[opAbove])-1);
1356  Double_t sub = GetHeight()*spec.fSize/14;
1357  Double_t x1 , y1 , x2, y2, x3, x4;
1358  switch(opAbove) {
1359  case 0: // bar
1360  Double_t ypos ;
1361  ypos = y-fs1.Over()-sub ;//-GetHeight()*spec.fSize/4. ;
1362  DrawLine(x,ypos,x+fs1.Width(),ypos,spec);
1363  break;
1364  case 1: // vec
1365  Double_t y0 ;
1366  y0 = y-sub-fs1.Over() ;
1367  y1 = y0-GetHeight()*spec.fSize/8 ;
1368  x1 = x+fs1.Width() ;
1369  DrawLine(x,y1,x1,y1,spec);
1370  DrawLine(x1,y1,x1-GetHeight()*spec.fSize/4,y0-GetHeight()*spec.fSize/4,spec);
1371  DrawLine(x1,y1,x1-GetHeight()*spec.fSize/4,y0,spec);
1372  break;
1373  case 2: // dot
1374  x1 = x+fs1.Width()/2-3*sub/4 ;
1375  x2 = x+fs1.Width()/2+3*sub/4 ;
1376  y1 = y-sub-fs1.Over() ;
1377  DrawLine(x1,y1,x2,y1,spec);
1378  break;
1379  case 3: // hat
1380  x2 = x+fs1.Width()/2 ;
1381  y1 = y -9*sub;
1382  y2 = y1-2*sub;
1383  x1 = x2-fs1.Width()/3 ;
1384  x3 = x2+fs1.Width()/3 ;
1385  DrawLine(x1,y1,x2,y2,spec);
1386  DrawLine(x2,y2,x3,y1,spec);
1387  break;
1388  case 4: // ddot
1389  x1 = x+fs1.Width()/2-9*sub/4 ;
1390  x2 = x+fs1.Width()/2-3*sub/4 ;
1391  x3 = x+fs1.Width()/2+9*sub/4 ;
1392  x4 = x+fs1.Width()/2+3*sub/4 ;
1393  y1 = y-sub-fs1.Over() ;
1394  DrawLine(x1,y1,x2,y1,spec);
1395  DrawLine(x3,y1,x4,y1,spec);
1396  break;
1397  case 5: // acute
1398  x1 = x+fs1.Width()/2;
1399  y1 = y +sub -fs1.Over() ;
1400  x2 = x1 +3*sub;
1401  y2 = y1 -2.5*sub;
1402  DrawLine(x1,y1,x2,y2,spec);
1403  break;
1404  case 6: // grave
1405  x1 = x+fs1.Width()/2-sub;
1406  y1 = y-sub-fs1.Over() ;
1407  x2 = x1 +2*sub;
1408  y2 = y1 +2*sub;
1409  DrawLine(x1,y1,x2,y2,spec);
1410  break;
1411  case 7: // check
1412  x1 = x+fs1.Width()/2 ;
1413  x2 = x1 -2*sub ;
1414  x3 = x1 +2*sub ;
1415  y1 = y-sub-fs1.Over() ;
1416  DrawLine(x2,y-3*sub-fs1.Over(),x1,y1,spec);
1417  DrawLine(x3,y-3*sub-fs1.Over(),x1,y1,spec);
1418  break;
1419  case 8: // tilde
1420  x2 = x+fs1.Width()/2 ;
1421  y2 = y -fs1.Over() ;
1422  {
1423  // tilde must be drawn separately on screen and on PostScript
1424  // because an adjustment is required along Y for PostScript.
1425  TVirtualPS *saveps = gVirtualPS;
1426  if (gVirtualPS) gVirtualPS = 0;
1427  Double_t y22 = y2;
1428  if (gVirtualX->InheritsFrom("TGCocoa")) y2 -= 3*sub;
1429  Double_t sinang = TMath::Sin(spec.fAngle/180*kPI);
1430  Double_t cosang = TMath::Cos(spec.fAngle/180*kPI);
1431  Double_t xOrigin = (Double_t)gPad->XtoAbsPixel(fX);
1432  Double_t yOrigin = (Double_t)gPad->YtoAbsPixel(fY);
1433  Double_t xx = gPad->AbsPixeltoX(Int_t((x2-xOrigin)*cosang+(y2-yOrigin)*sinang+xOrigin));
1434  Double_t yy = gPad->AbsPixeltoY(Int_t((x2-xOrigin)*-sinang+(y2-yOrigin)*cosang+yOrigin));
1435  TText tilde;
1436  tilde.SetTextFont(fTextFont);
1437  tilde.SetTextColor(spec.fColor);
1438  tilde.SetTextSize(0.9*spec.fSize);
1439  tilde.SetTextAlign(22);
1440  tilde.SetTextAngle(fTextAngle);
1441  tilde.PaintText(xx,yy,"~");
1442  if (saveps) {
1443  gVirtualPS = saveps;
1444  if (!strstr(gVirtualPS->GetTitle(),"IMG")) y22 -= 4*sub;
1445  xx = gPad->AbsPixeltoX(Int_t((x2-xOrigin)*cosang+(y22-yOrigin)*sinang+xOrigin));
1446  yy = gPad->AbsPixeltoY(Int_t((x2-xOrigin)*-sinang+(y22-yOrigin)*cosang+yOrigin));
1447  gVirtualPS->SetTextAlign(22);
1448  gVirtualPS->Text(xx, yy, "~");
1449  }
1450  }
1451  break;
1452  case 9: // slash
1453  x1 = x + 0.8*fs1.Width();
1454  y1 = y -fs1.Over() ;
1455  x2 = x + 0.3*fs1.Width();
1456  y2 = y1 + 1.2*fs1.Height();
1457  DrawLine(x1,y1,x2,y2,spec);
1458  break;
1459  }
1460  }
1461  Double_t div = 3;
1462  if (opAbove==1) div=4;
1463  result.Set(fs1.Width(),fs1.Over()+GetHeight()*spec.fSize/div,fs1.Under());
1464  }
1465  else if (opSquareBracket) { // operator #[]{arg}
1466  Double_t l = GetHeight()*spec.fSize/4;
1467  Double_t l2 = l/2 ;
1468  if (!fShow) {
1469  fs1 = Anal1(spec,text+3,length-3);
1470  Savefs(&fs1);
1471  } else {
1472  fs1 = Readfs();
1473  Analyse(x+l2+l,y,spec,text+3,length-3);
1474  DrawLine(x+l2,y-fs1.Over(),x+l2,y+fs1.Under(),spec);
1475  DrawLine(x+l2,y-fs1.Over(),x+l2+l,y-fs1.Over(),spec);
1476  DrawLine(x+l2,y+fs1.Under(),x+l2+l,y+fs1.Under(),spec);
1477  DrawLine(x+l2+fs1.Width()+2*l,y-fs1.Over(),x+l2+fs1.Width()+2*l,y+fs1.Under(),spec);
1478  DrawLine(x+l2+fs1.Width()+2*l,y-fs1.Over(),x+l2+fs1.Width()+l,y-fs1.Over(),spec);
1479  DrawLine(x+l2+fs1.Width()+2*l,y+fs1.Under(),x+l2+fs1.Width()+l,y+fs1.Under(),spec);
1480  }
1481  result.Set(fs1.Width()+3*l,fs1.Over(),fs1.Under());
1482  }
1483  else if (opParen) { // operator #(){arg}
1484  Double_t l = GetHeight()*spec.fSize/4;
1485  Double_t radius2,radius1 , dw, l2 = l/2 ;
1486  Double_t angle = 35 ;
1487  if (!fShow) {
1488  fs1 = Anal1(spec,text+3,length-3);
1489  Savefs(&fs1);
1490  radius2 = fs1.Height() ;
1491  radius1 = radius2 * 2 / 3;
1492  dw = radius1*(1 - TMath::Cos(kPI*angle/180)) ;
1493  } else {
1494  fs1 = Readfs();
1495  radius2 = fs1.Height();
1496  radius1 = radius2 * 2 / 3;
1497  dw = radius1*(1 - TMath::Cos(kPI*angle/180)) ;
1498  Double_t x1 = x+l2+radius1 ;
1499  Double_t x2 = x+5*l2+2*dw+fs1.Width()-radius1 ;
1500  Double_t y1 = y - (fs1.Over() - fs1.Under())/2. ;
1501  DrawParenthesis(x1,y1,radius1,radius2,180-angle,180+angle,spec) ;
1502  DrawParenthesis(x2,y1,radius1,radius2,360-angle,360+angle,spec) ;
1503  Analyse(x+3*l2+dw,y,spec,text+3,length-3);
1504  }
1505  // result = TLatexFormSize(fs1.Width()+3*l,fs1.Over(),fs1.Under());
1506  result.Set(fs1.Width()+3*l+2*dw,fs1.Over(),fs1.Under());
1507  }
1508  else if (opAbs) { // operator #||{arg}
1509  Double_t l = GetHeight()*spec.fSize/4;
1510  Double_t l2 = l/2 ;
1511  if (!fShow) {
1512  fs1 = Anal1(spec,text+3,length-3);
1513  Savefs(&fs1);
1514  } else {
1515  fs1 = Readfs();
1516  Analyse(x+l2+l,y,spec,text+3,length-3);
1517  DrawLine(x+l2,y-fs1.Over(),x+l2,y+fs1.Under(),spec);
1518  DrawLine(x+l2+fs1.Width()+2*l,y-fs1.Over(),x+l2+fs1.Width()+2*l,y+fs1.Under(),spec);
1519  }
1520  result.Set(fs1.Width()+3*l,fs1.Over(),fs1.Under());
1521  }
1522  else if (opBigCurly) { // big curly bracket #{}{arg}
1523  Double_t l = GetHeight()*spec.fSize/4;
1524  Double_t l2 = l/2 ;
1525  Double_t l8 , ltip;
1526 
1527  if (!fShow) {
1528  fs1 = Anal1(spec,text+3,length-3);
1529  l8 = fs1.Height()/8 ;
1530  ltip = TMath::Min(l8,l) ;
1531  l = ltip ;
1532  Savefs(&fs1);
1533  } else {
1534  fs1 = Readfs();
1535  Double_t y2 = y + (fs1.Under()-fs1.Over())/2 ;
1536  l8 = fs1.Height()/8 ;
1537  ltip = TMath::Min(l8,l) ;
1538  l = ltip ;
1539  Analyse(x+l+ltip+l2,y,spec,text+3,length-3);
1540  // Draw open curly bracket
1541  // Vertical lines
1542  DrawLine(x+l2+ltip,y-fs1.Over(),x+l2+ltip,y2-ltip,spec);
1543  DrawLine(x+l2+ltip,y2+ltip,x+l2+ltip,y+fs1.Under(),spec);
1544  // top and bottom lines
1545  DrawLine(x+l2+ltip,y-fs1.Over(),x+l2+ltip+l,y-fs1.Over(),spec);
1546  DrawLine(x+l2+ltip,y+fs1.Under(),x+l2+ltip+l,y+fs1.Under(),spec);
1547  // < sign
1548  DrawLine(x+l2,y2,x+l2+ltip,y2-ltip,spec);
1549  DrawLine(x+l2,y2,x+l2+ltip,y2+ltip,spec);
1550 
1551  // Draw close curly bracket
1552  // vertical lines
1553  DrawLine(x+l2+ltip+fs1.Width()+2*l,y-fs1.Over(),x+l2+ltip+fs1.Width()+2*l,y2-ltip,spec);
1554  DrawLine(x+l2+ltip+fs1.Width()+2*l,y2+ltip,x+l2+ltip+fs1.Width()+2*l,y+fs1.Under(),spec);
1555  // Top and bottom lines
1556  DrawLine(x+l2+fs1.Width()+l+ltip,y-fs1.Over(),x+l2+ltip+fs1.Width()+2*l,y-fs1.Over(),spec);
1557  DrawLine(x+l2+fs1.Width()+l+ltip,y+fs1.Under(),x+l2+ltip+fs1.Width()+2*l,y+fs1.Under(),spec);
1558  // > sign
1559  DrawLine(x+l2+ltip+2*l+fs1.Width(),y2-ltip,x+l2+2*l+2*ltip+fs1.Width(),y2,spec);
1560  DrawLine(x+l2+ltip+2*l+fs1.Width(),y2+ltip,x+l2+2*l+2*ltip+fs1.Width(),y2,spec);
1561  }
1562  result.Set(fs1.Width()+3*l+2*ltip,fs1.Over(),fs1.Under()) ;;
1563  }
1564  else if (opFrac>-1) { // \frac found
1565  if (opCurlyCurly==-1) { // }{ not found
1566  // arguments missing for \frac
1567  fError = "Missing denominator for #frac";
1568  return TLatexFormSize(0,0,0);
1569  }
1570  Double_t height = GetHeight()*spec.fSize/8;
1571  if (!fShow) {
1572  fs1 = Anal1(spec,text+opFrac+6,opCurlyCurly-opFrac-6);
1573  fs2 = Anal1(spec,text+opCurlyCurly+2,length-opCurlyCurly-3);
1574  Savefs(&fs1);
1575  Savefs(&fs2);
1576  } else {
1577  fs2 = Readfs();
1578  fs1 = Readfs();
1579  Double_t addW1,addW2;
1580  if (fs1.Width()<fs2.Width()) {
1581  addW1 = (fs2.Width()-fs1.Width())/2;
1582  addW2 = 0;
1583  } else {
1584  addW1 = 0;
1585  addW2 = (fs1.Width()-fs2.Width())/2;
1586  }
1587  Analyse(x+addW2,y+fs2.Over()-height,spec,text+opCurlyCurly+2,length-opCurlyCurly-3); // denominator
1588  Analyse(x+addW1,y-fs1.Under()-3*height,spec,text+opFrac+6,opCurlyCurly-opFrac-6); //numerator
1589 
1590  DrawLine(x,y-2*height,x+TMath::Max(fs1.Width(),fs2.Width()),y-2*height,spec);
1591  }
1592 
1593  result.Set(TMath::Max(fs1.Width(),fs2.Width()),fs1.Height()+3*height,fs2.Height()-height);
1594 
1595  }
1596  else if (opSplitLine>-1) { // \splitline found
1597  if (opCurlyCurly==-1) { // }{ not found
1598  // arguments missing for \splitline
1599  fError = "Missing second line for #splitline";
1600  return TLatexFormSize(0,0,0);
1601  }
1602  Double_t height = GetHeight()*spec.fSize/8;
1603  if (!fShow) {
1604  fs1 = Anal1(spec,text+opSplitLine+11,opCurlyCurly-opSplitLine-11);
1605  fs2 = Anal1(spec,text+opCurlyCurly+2,length-opCurlyCurly-3);
1606  Savefs(&fs1);
1607  Savefs(&fs2);
1608  } else {
1609  fs2 = Readfs();
1610  fs1 = Readfs();
1611  Analyse(x,y+fs2.Over()-height,spec,text+opCurlyCurly+2,length-opCurlyCurly-3); // second line
1612  Analyse(x,y-fs1.Under()-3*height,spec,text+opSplitLine+11,opCurlyCurly-opSplitLine-11); //first line
1613  }
1614 
1615  result.Set(TMath::Max(fs1.Width(),fs2.Width()),fs1.Height()+3*height,fs2.Height()-height);
1616 
1617  }
1618  else if (opSqrt>-1) { // \sqrt found
1619  if (!fShow) {
1620  if (opSquareCurly>-1) {
1621  // power nth #sqrt[n]{arg}
1622  fs1 = Anal1(specNewSize,text+opSqrt+6,opSquareCurly-opSqrt-6);
1623  fs2 = Anal1(spec,text+opSquareCurly+1,length-opSquareCurly-1);
1624  Savefs(&fs1);
1625  Savefs(&fs2);
1626  result.Set(fs2.Width()+ GetHeight()*spec.fSize/10+TMath::Max(GetHeight()*spec.fSize/2,(Double_t)fs1.Width()),
1627  fs2.Over()+fs1.Height()+GetHeight()*spec.fSize/4,fs2.Under());
1628  } else {
1629  fs1 = Anal1(spec,text+opSqrt+5,length-opSqrt-5);
1630  Savefs(&fs1);
1631  result.Set(fs1.Width()+GetHeight()*spec.fSize/2,fs1.Over()+GetHeight()*spec.fSize/4,fs1.Under());
1632  }
1633  } else {
1634  if (opSquareCurly>-1) { // ]{
1635  fs2 = Readfs();
1636  fs1 = Readfs();
1637  Double_t pas = TMath::Max(GetHeight()*spec.fSize/2,(Double_t)fs1.Width());
1638  Double_t pas2 = pas + GetHeight()*spec.fSize/10;
1639  Double_t y1 = y-fs2.Over() ;
1640  Double_t y2 = y+fs2.Under() ;
1641  Double_t y3 = y1-GetHeight()*spec.fSize/4;
1642  Analyse(x+pas2,y,spec,text+opSquareCurly+1,length-opSquareCurly-1);
1643  Analyse(x,y-fs2.Over()-fs1.Under(),specNewSize,text+opSqrt+6,opSquareCurly-opSqrt-6); // indice
1644  DrawLine(x,y1,x+pas,y2,spec);
1645  DrawLine(x+pas,y2,x+pas,y3,spec);
1646  DrawLine(x+pas,y3,x+pas2+fs2.Width(),y3,spec);
1647  } else {
1648  fs1 = Readfs();
1649  Double_t x1 = x+GetHeight()*spec.fSize*2/5 ;
1650  Double_t x2 = x+GetHeight()*spec.fSize/2+fs1.Width() ;
1651  Double_t y1 = y-fs1.Over() ;
1652  Double_t y2 = y+fs1.Under() ;
1653  Double_t y3 = y1-GetHeight()*spec.fSize/4;
1654 
1655  Analyse(x+GetHeight()*spec.fSize/2,y,spec,text+opSqrt+6,length-opSqrt-7);
1656 
1657  Short_t lineW = GetLineWidth();
1658  SetLineWidth(1);
1659  Double_t dx = (y2-y3)/8;
1660  UInt_t a,d;
1661  GetTextAscentDescent(a, d, text);
1662  if (a>12) SetLineWidth(TMath::Max(2,(Int_t)(dx/2)));
1663  DrawLine(x1-2*dx,y1,x1-dx,y2,spec);
1664  if (a>12) SetLineWidth(TMath::Max(1,(Int_t)(dx/4)));
1665  DrawLine(x1-dx,y2,x1,y3,spec);
1666  DrawLine(x1,y3,x2,y3,spec);
1667  SetLineWidth(lineW);
1668  }
1669  }
1670  }
1671  else if (opColor>-1) { // \color found
1672  if (opSquareCurly==-1) {
1673  // color number is not specified
1674  fError = "Missing color number. Syntax is #color[(Int_t)nb]{ ... }";
1675  return TLatexFormSize(0,0,0);
1676  }
1677  TextSpec_t newSpec = spec;
1678  Char_t *nb = new Char_t[opSquareCurly-opColor-6];
1679  strncpy(nb,text+opColor+7,opSquareCurly-opColor-7);
1680  nb[opSquareCurly-opColor-7] = 0;
1681  if (sscanf(nb,"%d",&newSpec.fColor) < 1) {
1682  delete[] nb;
1683  // color number is invalid
1684  fError = "Invalid color number. Syntax is #color[(Int_t)nb]{ ... }";
1685  return TLatexFormSize(0,0,0);
1686  }
1687  delete[] nb;
1688  if (!fShow) {
1689  result = Anal1(newSpec,text+opSquareCurly+1,length-opSquareCurly-1);
1690  } else {
1691  Analyse(x,y,newSpec,text+opSquareCurly+1,length-opSquareCurly-1);
1692  }
1693  }
1694  else if (opFont>-1) { // \font found
1695  if (opSquareCurly==-1) {
1696  // font number is not specified
1697  fError = "Missing font number. Syntax is #font[nb]{ ... }";
1698  return TLatexFormSize(0,0,0);
1699  }
1700  TextSpec_t newSpec = spec;
1701  Char_t *nb = new Char_t[opSquareCurly-opFont-5];
1702  strncpy(nb,text+opFont+6,opSquareCurly-opFont-6);
1703  nb[opSquareCurly-opFont-6] = 0;
1704  if (sscanf(nb,"%d",&newSpec.fFont) < 1) {
1705  delete[] nb;
1706  // font number is invalid
1707  fError = "Invalid font number. Syntax is #font[(Int_t)nb]{ ... }";
1708  return TLatexFormSize(0,0,0);
1709  }
1710  delete[] nb;
1711  if (!fShow) {
1712  result = Anal1(newSpec,text+opSquareCurly+1,length-opSquareCurly-1);
1713  } else {
1714  Analyse(x,y,newSpec,text+opSquareCurly+1,length-opSquareCurly-1);
1715  }
1716  }
1717  else if (opKern>-1) { // #kern found
1718  if (opSquareCurly==-1) {
1719  // horizontal shift is not specified
1720  fError = "Missing horizontal shift number. Syntax is #kern[dx]{ ... }";
1721  return TLatexFormSize(0,0,0);
1722  }
1723  Char_t *dxc = new Char_t[opSquareCurly-opKern-5];
1724  strncpy(dxc,text+opKern+6,opSquareCurly-opKern-6);
1725  dxc[opSquareCurly-opKern-6] = 0;
1726  Float_t dx = 0;
1727  if (sscanf(dxc,"%f",&dx) < 1) {
1728  delete[] dxc;
1729  // horizontal shift number is invalid
1730  fError = "Invalid horizontal shift number. Syntax is #kern[(Float_t)dx]{ ... }";
1731  return TLatexFormSize(0,0,0);
1732  }
1733  delete[] dxc;
1734  if (!fShow) {
1735  fs1 = Anal1(spec,text+opSquareCurly+1,length-opSquareCurly-1);
1736  Savefs(&fs1);
1737  Double_t ddx = dx * fs1.Width();
1738  result = TLatexFormSize(fs1.Width() + ddx, fs1.Over(), fs1.Under());
1739  } else {
1740  fs1 = Readfs();
1741  Double_t ddx = dx * fs1.Width();
1742  Analyse(x + ddx,y,spec,text+opSquareCurly+1,length-opSquareCurly-1);
1743  }
1744  }
1745  else if (opLower>-1) { // #lower found
1746  if (opSquareCurly==-1) {
1747  // vertical shift is not specified
1748  fError = "Missing vertical shift number. Syntax is #lower[dy]{ ... }";
1749  return TLatexFormSize(0,0,0);
1750  }
1751  Char_t *dyc = new Char_t[opSquareCurly-opLower-6];
1752  strncpy(dyc,text+opLower+7,opSquareCurly-opLower-7);
1753  dyc[opSquareCurly-opLower-7] = 0;
1754  Float_t dy = 0;
1755  if (sscanf(dyc,"%f",&dy) < 1) {
1756  delete[] dyc;
1757  // vertical shift number is invalid
1758  fError = "Invalid vertical shift number. Syntax is #lower[(Float_t)dy]{ ... }";
1759  return TLatexFormSize(0,0,0);
1760  }
1761  delete[] dyc;
1762  if (!fShow) {
1763  fs1 = Anal1(spec,text+opSquareCurly+1,length-opSquareCurly-1);
1764  Savefs(&fs1);
1765  Double_t ddy = dy * (fs1.Over() + fs1.Under());
1766  result = TLatexFormSize(fs1.Width(), fs1.Over() + ddy, fs1.Under() + ddy);
1767  } else {
1768  fs1 = Readfs();
1769  Double_t ddy = dy * (fs1.Over() + fs1.Under());
1770  Analyse(x,y + ddy,spec,text+opSquareCurly+1,length-opSquareCurly-1);
1771  }
1772  }
1773  else if (opScale>-1) { // \scale found
1774  if (opSquareCurly==-1) {
1775  // scale factor is not specified
1776  fError = "Missing scale factor. Syntax is #scale[(Double_t)nb]{ ... }";
1777  return TLatexFormSize(0,0,0);
1778  }
1779  TextSpec_t newSpec = spec;
1780  Char_t *nb = new Char_t[opSquareCurly-opScale-6];
1781  strncpy(nb,text+opScale+7,opSquareCurly-opScale-7);
1782  nb[opSquareCurly-opScale-7] = 0;
1783  if (sscanf(nb,"%lf",&newSpec.fSize) < 1) {
1784  delete[] nb;
1785  // scale factor is invalid
1786  fError = "Invalid scale factor. Syntax is #factor[(Double_t)nb]{ ... }";
1787  return TLatexFormSize(0,0,0);
1788  }
1789  newSpec.fSize *= spec.fSize;
1790  delete[] nb;
1791  if (!fShow) {
1792  result = Anal1(newSpec,text+opSquareCurly+1,length-opSquareCurly-1);
1793  } else {
1794  Analyse(x,y,newSpec,text+opSquareCurly+1,length-opSquareCurly-1);
1795  }
1796  }
1797  else if (opBf>-1) { // operator #bf{arg}
1798  TextSpec_t newSpec = spec;
1799  Int_t lut[] = {3, 13, 1, 6, 7, 4, 5, 10, 11, 8, 9, 12, 2, 14, 15};
1800  Int_t fontId = (newSpec.fFont/10);
1801  if ((fontId >= 1) && (fontId <= (Int_t)(sizeof(lut)/sizeof(lut[0])))) fontId = lut[fontId-1];
1802  newSpec.fFont = fontId*10 + newSpec.fFont%10;
1803  if (!fShow) {
1804  fs1 = Anal1(newSpec,text+3,length-3);
1805  Savefs(&fs1);
1806  } else {
1807  fs1 = Readfs();
1808  Analyse(x,y,newSpec,text+3,length-3);
1809  }
1810  result = fs1;
1811  }
1812  else if (opMbox>-1) { // dummy operator #mbox{arg}
1813  TextSpec_t newSpec = spec;
1814  if (!fShow) {
1815  fs1 = Anal1(newSpec,text+5,length-5);
1816  Savefs(&fs1);
1817  } else {
1818  fs1 = Readfs();
1819  Analyse(x,y,newSpec,text+5,length-5);
1820  }
1821  result = fs1;
1822  }
1823  else if (opIt>-1) { // operator #it{arg}
1824  TextSpec_t newSpec = spec;
1825  Int_t lut[] = {13, 3, 2, 5, 4, 7, 6, 9, 8, 11, 10, 15, 1, 14, 12};
1826  Int_t fontId = (newSpec.fFont/10);
1827  if ((fontId >= 1) && (fontId <= (Int_t)(sizeof(lut)/sizeof(lut[0])))) fontId = lut[fontId-1];
1828  newSpec.fFont = fontId*10 + newSpec.fFont%10;
1829  fItalic = !fItalic;
1830  if (!fShow) {
1831  fs1 = Anal1(newSpec,text+3,length-3);
1832  Savefs(&fs1);
1833  } else {
1834  fs1 = Readfs();
1835  Analyse(x,y,newSpec,text+3,length-3);
1836  }
1837  fItalic = !fItalic;
1838  result = fs1;
1839  }
1840  else { // no operators found, it is a character string
1841  SetTextSize(spec.fSize);
1842  SetTextAngle(spec.fAngle);
1843  SetTextColor(spec.fColor);
1844  SetTextFont(spec.fFont);
1845  SetTextAlign(11);
1846  TAttText::Modify();
1847  UInt_t w=0,h=0;
1848 
1849  Int_t leng = strlen(text) ;
1850 
1851  quote1 = quote2 = kFALSE ;
1852  Char_t *p ;
1853  for (i=0 ; i<leng ; i++) {
1854  switch (text[i]) {
1855  case '\'' : quote1 = !quote1 ; break ; // single quote symbol not correctly interpreted when PostScript
1856  case '"' : quote2 = !quote2 ; break ;
1857  }
1858  //if (quote1 || quote2) continue ;
1859  if (text[i] == '@') { // @ symbol not correctly interpreted when PostScript
1860  p = &text[i] ;
1861  if ( *(p+1) == '{' || *(p+1) == '}' || *(p+1) == '[' || *(p+1) == ']') {
1862  while (*p != 0) {
1863  *p = *(p+1) ; p++ ;
1864  }
1865  leng-- ;
1866  }
1867  }
1868  }
1869  text[leng] = 0 ;
1870 
1871  if (fShow) {
1872  // paint the Latex sub-expression per sub-expression
1873  Double_t xOrigin = (Double_t)gPad->XtoAbsPixel(fX);
1874  Double_t yOrigin = (Double_t)gPad->YtoAbsPixel(fY);
1875  Double_t angle = kPI*spec.fAngle/180.;
1876  Double_t xx = gPad->AbsPixeltoX(Int_t((x-xOrigin)*TMath::Cos(angle)+(y-yOrigin)*TMath::Sin(angle)+xOrigin));
1877  Double_t yy = gPad->AbsPixeltoY(Int_t((x-xOrigin)*TMath::Sin(-angle)+(y-yOrigin)*TMath::Cos(angle)+yOrigin));
1878  gPad->PaintText(xx,yy,text);
1879  } else {
1880  GetTextExtent(w,h,text);
1881  Double_t width = w;
1882  UInt_t a,d;
1883  GetTextAscentDescent(a, d, text);
1884  fs1.Set(width,a,d);
1885  }
1886 
1887  result = fs1;
1888  }
1889 
1890  delete[] text;
1891 
1892  return result;
1893 }
1894 
1895 ////////////////////////////////////////////////////////////////////////////////
1896 /// Make a copy of this object with the new parameters
1897 /// And copy object attributes
1898 
1899 TLatex *TLatex::DrawLatex(Double_t x, Double_t y, const char *text)
1900 {
1901  TLatex *newtext = new TLatex(x, y, text);
1902  TAttText::Copy(*newtext);
1903  newtext->SetBit(kCanDelete);
1904  if (TestBit(kTextNDC)) newtext->SetNDC();
1905  newtext->AppendPad();
1906  return newtext;
1907 }
1908 
1909 ////////////////////////////////////////////////////////////////////////////////
1910 /// Draw this TLatex with new coordinates in NDC.
1911 
1912 TLatex *TLatex::DrawLatexNDC(Double_t x, Double_t y, const char *text)
1913 {
1914  TLatex *newtext = DrawLatex(x, y, text);
1915  newtext->SetNDC();
1916  return newtext;
1917 }
1918 
1919 ////////////////////////////////////////////////////////////////////////////////
1920 /// Draw a line in a Latex formula
1921 
1922 void TLatex::DrawLine(Double_t x1, Double_t y1, Double_t x2, Double_t y2, TextSpec_t spec)
1923 {
1924  Double_t sinang = TMath::Sin(spec.fAngle/180*kPI);
1925  Double_t cosang = TMath::Cos(spec.fAngle/180*kPI);
1926  Double_t xOrigin = (Double_t)gPad->XtoAbsPixel(fX);
1927  Double_t yOrigin = (Double_t)gPad->YtoAbsPixel(fY);
1928  Double_t xx = gPad->AbsPixeltoX(Int_t((x1-xOrigin)*cosang+(y1-yOrigin)*sinang+xOrigin));
1929  Double_t yy = gPad->AbsPixeltoY(Int_t((x1-xOrigin)*-sinang+(y1-yOrigin)*cosang+yOrigin));
1930 
1931  Double_t xx2 = gPad->AbsPixeltoX(Int_t((x2-xOrigin)*cosang+(y2-yOrigin)*sinang+xOrigin));
1932  Double_t yy2 = gPad->AbsPixeltoY(Int_t((x2-xOrigin)*-sinang+(y2-yOrigin)*cosang+yOrigin));
1933 
1934  SetLineColor(spec.fColor);
1935  TAttLine::Modify();
1936  gPad->PaintLine(xx,yy,xx2,yy2);
1937 }
1938 
1939 ////////////////////////////////////////////////////////////////////////////////
1940 /// Draw an arc of ellipse in a Latex formula (right or left parenthesis)
1941 
1942 void TLatex::DrawCircle(Double_t x1, Double_t y1, Double_t r, TextSpec_t spec )
1943 {
1944  if (r < 1) r = 1;
1945  Double_t sinang = TMath::Sin(spec.fAngle/180*kPI);
1946  Double_t cosang = TMath::Cos(spec.fAngle/180*kPI);
1947  Double_t xOrigin = (Double_t)gPad->XtoAbsPixel(fX);
1948  Double_t yOrigin = (Double_t)gPad->YtoAbsPixel(fY);
1949 
1950  const Int_t np = 40;
1951  Double_t dphi = 2*kPI/np;
1952  Double_t x[np+3], y[np+3];
1953  Double_t angle,dx,dy;
1954 
1955  SetLineColor(spec.fColor);
1956  TAttLine::Modify(); //Change line attributes only if necessary
1957 
1958  for (Int_t i=0;i<=np;i++) {
1959  angle = Double_t(i)*dphi;
1960  dx = r*TMath::Cos(angle) +x1 -xOrigin;
1961  dy = r*TMath::Sin(angle) +y1 -yOrigin;
1962  x[i] = gPad->AbsPixeltoX(Int_t( dx*cosang+ dy*sinang +xOrigin));
1963  y[i] = gPad->AbsPixeltoY(Int_t(-dx*sinang+ dy*cosang +yOrigin));
1964  }
1965  gPad->PaintPolyLine(np+1,x,y);
1966 }
1967 
1968 ////////////////////////////////////////////////////////////////////////////////
1969 /// Draw an arc of ellipse in a Latex formula (right or left parenthesis)
1970 
1972  Double_t phimin, Double_t phimax, TextSpec_t spec )
1974  if (r1 < 1) r1 = 1;
1975  if (r2 < 1) r2 = 1;
1976  Double_t sinang = TMath::Sin(spec.fAngle/180*kPI);
1977  Double_t cosang = TMath::Cos(spec.fAngle/180*kPI);
1978  Double_t xOrigin = (Double_t)gPad->XtoAbsPixel(fX);
1979  Double_t yOrigin = (Double_t)gPad->YtoAbsPixel(fY);
1980 
1981  const Int_t np = 40;
1982  Double_t dphi = (phimax-phimin)*kPI/(180*np);
1983  Double_t x[np+3], y[np+3];
1984  Double_t angle,dx,dy ;
1985 
1986  SetLineColor(spec.fColor);
1987  TAttLine::Modify(); //Change line attributes only if necessary
1988 
1989  for (Int_t i=0;i<=np;i++) {
1990  angle = phimin*kPI/180 + Double_t(i)*dphi;
1991  dx = r1*TMath::Cos(angle) +x1 -xOrigin;
1992  dy = r2*TMath::Sin(angle) +y1 -yOrigin;
1993  x[i] = gPad->AbsPixeltoX(Int_t( dx*cosang+dy*sinang +xOrigin));
1994  y[i] = gPad->AbsPixeltoY(Int_t(-dx*sinang+dy*cosang +yOrigin));
1995  }
1996  gPad->PaintPolyLine(np+1,x,y);
1997 }
1998 
1999 ////////////////////////////////////////////////////////////////////////////////
2000 /// Paint.
2001 
2002 void TLatex::Paint(Option_t *)
2003 {
2004  Double_t xsave = fX;
2005  Double_t ysave = fY;
2006  if (TestBit(kTextNDC)) {
2007  fX = gPad->GetX1() + xsave*(gPad->GetX2() - gPad->GetX1());
2008  fY = gPad->GetY1() + ysave*(gPad->GetY2() - gPad->GetY1());
2010  } else {
2011  PaintLatex(gPad->XtoPad(fX),gPad->YtoPad(fY),GetTextAngle(),GetTextSize(),GetTitle());
2012  }
2013  fX = xsave;
2014  fY = ysave;
2015 }
2016 
2017 ////////////////////////////////////////////////////////////////////////////////
2018 /// Main drawing function
2019 ///
2020 /// Warning: Unlike most others "XYZ::PaintXYZ" methods, PaintLatex modifies
2021 /// the TLatex data members.
2022 
2023 void TLatex::PaintLatex(Double_t x, Double_t y, Double_t angle, Double_t size, const Char_t *text1)
2024 {
2025  TAttText::Modify(); // Change text attributes only if necessary.
2026 
2027  TVirtualPS *saveps = gVirtualPS;
2028 
2029  if (gVirtualPS) {
2030  if (gVirtualPS->InheritsFrom("TTeXDump")) {
2031  gVirtualPS->SetTextAngle(angle);
2032  TString t(text1);
2033  if (t.Index("#")>=0 || t.Index("^")>=0 || t.Index("\\")>=0) {
2034  t.ReplaceAll("#LT","\\langle");
2035  t.ReplaceAll("#GT","\\rangle");
2036  t.ReplaceAll("#club","\\clubsuit");
2037  t.ReplaceAll("#spade","\\spadesuit");
2038  t.ReplaceAll("#heart","\\heartsuit");
2039  t.ReplaceAll("#diamond","\\diamondsuit");
2040  t.ReplaceAll("#voidn","\\wp");
2041  t.ReplaceAll("#voidb","f");
2042  t.ReplaceAll("#ocopyright","\\copyright");
2043  t.ReplaceAll("#trademark","TM");
2044  t.ReplaceAll("#void3","TM");
2045  t.ReplaceAll("#oright","R");
2046  t.ReplaceAll("#void1","R");
2047  t.ReplaceAll("#3dots","\\ldots");
2048  t.ReplaceAll("#lbar","\\mid");
2049  t.ReplaceAll("#void8","\\mid");
2050  t.ReplaceAll("#divide","\\div");
2051  t.ReplaceAll("#Jgothic","\\Im");
2052  t.ReplaceAll("#Rgothic","\\Re");
2053  t.ReplaceAll("#doublequote","\"");
2054  t.ReplaceAll("#plus","+");
2055  t.ReplaceAll("#minus","-");
2056  t.ReplaceAll("#/","/");
2057  t.ReplaceAll("#upoint",".");
2058  t.ReplaceAll("#aa","\\mbox{\\aa}");
2059  t.ReplaceAll("#AA","\\mbox{\\AA}");
2060 
2061  t.ReplaceAll("#omicron","o");
2062  t.ReplaceAll("#Alpha","A");
2063  t.ReplaceAll("#Beta","B");
2064  t.ReplaceAll("#Epsilon","E");
2065  t.ReplaceAll("#Zeta","Z");
2066  t.ReplaceAll("#Eta","H");
2067  t.ReplaceAll("#Iota","I");
2068  t.ReplaceAll("#Kappa","K");
2069  t.ReplaceAll("#Mu","M");
2070  t.ReplaceAll("#Nu","N");
2071  t.ReplaceAll("#Omicron","O");
2072  t.ReplaceAll("#Rho","P");
2073  t.ReplaceAll("#Tau","T");
2074  t.ReplaceAll("#Chi","X");
2075  t.ReplaceAll("#varomega","\\varpi");
2076 
2077  t.ReplaceAll("#varUpsilon","?");
2078  t.ReplaceAll("#corner","?");
2079  t.ReplaceAll("#ltbar","?");
2080  t.ReplaceAll("#bottombar","?");
2081  t.ReplaceAll("#notsubset","?");
2082  t.ReplaceAll("#arcbottom","?");
2083  t.ReplaceAll("#cbar","?");
2084  t.ReplaceAll("#arctop","?");
2085  t.ReplaceAll("#topbar","?");
2086  t.ReplaceAll("#arcbar","?");
2087  t.ReplaceAll("#downleftarrow","?");
2088  t.ReplaceAll("#splitline","\\genfrac{}{}{0pt}{}");
2089 
2090  t.ReplaceAll("#","\\");
2091  t.ReplaceAll("%","\\%");
2092  }
2093  gVirtualPS->Text(x,y,t.Data());
2094  } else {
2095  Bool_t saveb = gPad->IsBatch();
2096  gPad->SetBatch(kTRUE);
2097  if (!PaintLatex1( x, y, angle, size, text1)) {
2098  if (saveps) gVirtualPS = saveps;
2099  return;
2100  }
2101  gPad->SetBatch(saveb);
2102  }
2103  gVirtualPS = 0;
2104  }
2105 
2106  if (!gPad->IsBatch()) PaintLatex1( x, y, angle, size, text1);
2107  if (saveps) gVirtualPS = saveps;
2108 }
2109 
2110 ////////////////////////////////////////////////////////////////////////////////
2111 /// Drawing function
2112 
2113 Int_t TLatex::PaintLatex1(Double_t x, Double_t y, Double_t angle, Double_t size, const Char_t *text1)
2114 {
2115  TString newText = text1;
2116  if( newText.Length() == 0) return 0;
2117  newText.ReplaceAll("#hbox","#mbox");
2118 
2119  fError = 0 ;
2120  if (CheckLatexSyntax(newText)) {
2121  std::cout<<"\n*ERROR<TLatex>: "<<fError<<std::endl;
2122  std::cout<<"==> "<<text1<<std::endl;
2123  return 0;
2124  }
2125  fError = 0 ;
2126 
2127  // Do not use Latex if font is low precision.
2128  if (fTextFont%10 < 2) {
2129  if (gVirtualX) gVirtualX->SetTextAngle(angle);
2130  if (gVirtualPS) gVirtualPS->SetTextAngle(angle);
2131  gPad->PaintText(x,y,text1);
2132  return 1;
2133  }
2134 
2135  // Paint the text using TMathText if contains a "\"
2136  if (strstr(text1,"\\")) {
2137  TMathText tm;
2138  tm.SetTextAlign(GetTextAlign());
2139  tm.SetTextFont(GetTextFont());
2140  tm.PaintMathText(x, y, angle, size, text1);
2141  return 1;
2142  }
2143 
2144  Double_t saveSize = size;
2145  Int_t saveFont = fTextFont;
2146  if (fTextFont%10 > 2) {
2147  UInt_t w = TMath::Abs(gPad->XtoAbsPixel(gPad->GetX2()) -
2148  gPad->XtoAbsPixel(gPad->GetX1()));
2149  UInt_t h = TMath::Abs(gPad->YtoAbsPixel(gPad->GetY2()) -
2150  gPad->YtoAbsPixel(gPad->GetY1()));
2151  if (w < h)
2152  size = size/w;
2153  else
2154  size = size/h;
2155  SetTextFont(10*(saveFont/10) + 2);
2156  }
2157 
2158  Int_t length = newText.Length() ;
2159  const Char_t *text = newText.Data() ;
2160 
2161  fX=x;
2162  fY=y;
2163  x = gPad->XtoAbsPixel(x);
2164  y = gPad->YtoAbsPixel(y);
2165  fShow = kFALSE ;
2166  TLatexFormSize fs = FirstParse(angle,size,text);
2167 
2168  fOriginSize = size;
2169 
2170  // Get current line attributes.
2171  Short_t lineW = GetLineWidth();
2172  Int_t lineC = GetLineColor();
2173  if (size<0.04) SetLineWidth(1);
2174 
2175  TextSpec_t spec;
2176  spec.fAngle = angle;
2177  spec.fSize = size;
2178  spec.fColor = GetTextColor();
2179  spec.fFont = GetTextFont();
2180  Short_t halign = fTextAlign/10;
2181  Short_t valign = fTextAlign - 10*halign;
2182  TextSpec_t newSpec = spec;
2183  if (fError != 0) {
2184  std::cout<<"*ERROR<TLatex>: "<<fError<<std::endl;
2185  std::cout<<"==> "<<text<<std::endl;
2186  } else {
2187  fShow = kTRUE;
2188  newSpec.fSize = size;
2189 
2190  switch (valign) {
2191  case 0: y -= fs.Under() ; break;
2192  case 1: break;
2193  case 2: y += fs.Height()*0.5-fs.Under(); y++; break;
2194  case 3: y += fs.Over() ; break;
2195  }
2196  switch (halign) {
2197  case 2: x -= fs.Width()/2 ; break;
2198  case 3: x -= fs.Width() ; break;
2199  }
2200  Analyse(x,y,newSpec,text,length);
2201  }
2202 
2203  SetTextSize(saveSize);
2204  SetTextAngle(angle);
2205  SetTextFont(saveFont);
2206  SetTextColor(spec.fColor);
2207  SetTextAlign(valign+10*halign);
2208  SetLineWidth(lineW);
2209  SetLineColor(lineC);
2210  delete[] fTabSize;
2211  if (fError != 0) return 0;
2212  return 1;
2213 }
2214 
2215 ////////////////////////////////////////////////////////////////////////////////
2216 /// Check if the Latex syntax is correct
2217 
2219 {
2220  const Char_t *kWord1[] = {"{}^{","{}_{","^{","_{","#scale{","#color{","#font{","#sqrt{","#[]{","#{}{","#||{",
2221  "#bar{","#vec{","#dot{","#hat{","#ddot{","#acute{","#grave{","#check{","#tilde{","#slash{","#bf{","#it{","#mbox{",
2222  "\\scale{","\\color{","\\font{","\\sqrt{","\\[]{","\\{}{","\\||{","#(){","\\(){",
2223  "\\bar{","\\vec{","\\dot{","\\hat{","\\ddot{","\\acute{","\\grave{","\\check{","\\bf{","\\it{","\\mbox{"}; // check for }
2224  const Char_t *kWord2[] = {"#scale[","#color[","#font[","#sqrt[","#kern[","#lower[","\\scale[","\\color[","\\font[","\\sqrt[","\\kern[","\\lower["}; // check for ]{ + }
2225  const Char_t *kWord3[] = {"#frac{","\\frac{","#splitline{","\\splitline{"}; // check for }{ then }
2226  const Char_t *kLeft1[] = {"#left[","\\left[","#left{","\\left{","#left|","\\left|","#left(","\\left("};
2227  const Char_t *kLeft2[] = {"#[]{","#[]{","#{}{","#{}{","#||{","#||{","#(){","#(){"};
2228  const Char_t *kRight[] = {"#right]","\\right]","#right}","\\right}","#right|","\\right|","#right)","\\right)"};
2229  const Int_t lkWord1[] = {4,4,2,2,7,7,6,6,4,4,4,5,5,5,5,6,7,7,7,7,7,4,4,6,7,7,6,6,4,4,4,4,4,5,5,5,5,6,7,7,7,4,4,6};
2230  const Int_t lkWord2[] = {7,7,6,6,6,7,7,7,6,6,6,7} ;
2231  const Int_t lkWord3[] = {6,6,11,11} ;
2232  Int_t nkWord1 = 44, nkWord2 = 12, nkWord3 = 4;
2233  Int_t i,k ;
2234  Int_t nLeft1 , nRight , nOfLeft, nOfRight;
2235  Int_t lLeft1 = 6 ;
2236  Int_t lLeft2 = 4 ;
2237  Int_t lRight = 7 ;
2238  nLeft1 = nRight = 8 ;
2239  nOfLeft = nOfRight = 0 ;
2240 
2241  Char_t buf[11] ; for (i=0;i<11;i++) buf[i]=0;
2242  Bool_t opFound ;
2243  Int_t opFrac = 0;
2244  Int_t length = text.Length() ;
2245 
2246  Int_t nOfCurlyBracket, nOfKW1, nOfKW2, nOfKW3, nOfSquareCurly, nOfCurlyCurly ;
2247  Int_t nOfExtraCurly = 0 , nOfExtraSquare = 0;
2248  Int_t nOfSquareBracket = 0 ;
2249  Int_t error = 0 ;
2250  Bool_t quote1 = kFALSE , quote2 = kFALSE;
2251 
2252  // first find and replace all occurrences of "kLeft1" keyword by "kLeft2" keyword,
2253  // and all occurrences of "kRight" keyword by "}".
2254  i = 0 ;
2255  while (i < length) {
2256  // The string in 'buf' does not need to be null terminated,
2257  // we will only check with strncmp.
2258  strncpy(buf,&text[i],TMath::Min(7,length-i));
2259  opFound = kFALSE ;
2260  for (k = 0 ; k < nLeft1 ; k++) {
2261  if (strncmp(buf,kLeft1[k],lLeft1)==0) {
2262  nOfLeft++ ;
2263  i+=lLeft1 ;
2264  opFound = kTRUE ;
2265  break ;
2266  }
2267  }
2268  if (opFound) continue ;
2269 
2270  for(k=0;k<nRight;k++) {
2271  if (strncmp(buf,kRight[k],lRight)==0) {
2272  nOfRight++ ;
2273  i+=lRight ;
2274  opFound = kTRUE ;
2275  break ;
2276  }
2277  }
2278  if (!opFound) i++ ;
2279  }
2280  if (nOfLeft != nOfRight) {
2281  printf(" nOfLeft = %d, nOfRight = %d\n",nOfLeft,nOfRight) ;
2282  error = 1 ;
2283  fError = "Operators \"#left\" and \"#right\" don't match !" ;
2284  goto ERROR_END ;
2285  }
2286 
2287  for (k = 0 ; k < nLeft1 ; k++) {
2288  text.ReplaceAll(kLeft1[k],lLeft1,kLeft2[k],lLeft2) ;
2289  }
2290  for (k = 0 ; k < nRight ; k++) {
2291  text.ReplaceAll(kRight[k],lRight,"}",1) ;
2292  }
2293  length = text.Length() ;
2294 
2295  i = nOfCurlyBracket = nOfKW1 = nOfKW2 = nOfKW3 = nOfSquareCurly = nOfCurlyCurly =0 ;
2296  while (i< length){
2297  switch (text[i]) {
2298  case '"' : quote1 = !quote1 ; break ;
2299  case '\'': quote2 = !quote2 ; break ;
2300  }
2301  // The string in 'buf' does not need to be null terminated,
2302  // we will only check with strncmp
2303  strncpy(buf,&text[i],TMath::Min(11,length-i));
2304  opFound = kFALSE ;
2305 
2306  for(k=0;k<nkWord1;k++) {
2307  if (strncmp(buf,kWord1[k],lkWord1[k])==0) {
2308  nOfKW1++ ;
2309  i+=lkWord1[k] ;
2310  opFound = kTRUE ;
2311  nOfCurlyBracket++ ;
2312  break ;
2313  }
2314  }
2315  if (opFound) continue ;
2316 
2317  for(k=0;k<nkWord2;k++) {
2318  if (strncmp(buf,kWord2[k],lkWord2[k])==0) {
2319  nOfKW2++ ;
2320  i+=lkWord2[k] ;
2321  opFound = kTRUE ;
2322  nOfSquareBracket++;
2323  break ;
2324  }
2325  }
2326  if (opFound) continue ;
2327 
2328  for(k=0;k<nkWord3;k++) {
2329  if (strncmp(buf,kWord3[k],lkWord3[k])==0) {
2330  nOfKW3++ ;
2331  i+=lkWord3[k] ;
2332  opFound = kTRUE ;
2333  opFrac++ ;
2334  nOfCurlyBracket++ ;
2335  break ;
2336  }
2337  }
2338  if (opFound) continue ;
2339  if (strncmp(buf,"}{",2) == 0 && opFrac) {
2340  opFrac-- ;
2341  nOfCurlyCurly++ ;
2342  i+= 2;
2343  }
2344  else if (strncmp(buf,"]{",2) == 0 && nOfSquareBracket) {
2345  nOfSquareCurly++ ;
2346  i+= 2 ;
2347  nOfCurlyBracket++ ;
2348  nOfSquareBracket-- ;
2349  }
2350  else if (strncmp(buf,"@{",2) == 0 || strncmp(buf,"@}",2) == 0) {
2351  i+= 2 ;
2352  }
2353  else if (strncmp(buf,"@[",2) == 0 || strncmp(buf,"@]",2) == 0) {
2354  i+= 2 ;
2355  }
2356  else if (text[i] == ']' ) { // not belonging to a key word, add @ in front
2357  text.Insert(i,"@") ;
2358  length++ ;
2359  i+=2 ;
2360  nOfExtraSquare-- ;
2361  }
2362  else if (text[i] == '[' ) { // not belonging to a key word, add @ in front
2363  text.Insert(i,"@") ;
2364  length++ ;
2365  i+=2 ;
2366  nOfExtraSquare++ ;
2367  }
2368  else if (text[i] == '{' ) { // not belonging to a key word, add @ in front
2369  text.Insert(i,"@") ;
2370  length++ ;
2371  i+=2 ;
2372  nOfExtraCurly++ ;
2373  }
2374  else if (text[i] == '}' ) {
2375  if ( nOfCurlyBracket) {
2376  nOfCurlyBracket-- ;
2377  i++ ;
2378  } else { // extra }, add @ in front
2379  text.Insert(i,"@") ;
2380  length++ ;
2381  i+=2 ;
2382  nOfExtraCurly-- ;
2383  }
2384  } else {
2385  i++ ;
2386  buf[1] = 0 ;
2387  }
2388  }
2389 
2390  if (nOfKW2 != nOfSquareCurly) {
2391  error = 1 ;
2392  fError = "Invalid number of \"]{\"" ;
2393  }
2394  else if (nOfKW3 != nOfCurlyCurly) {
2395  error = 1 ;
2396  fError = "Error in syntax of \"#frac\"" ;
2397  }
2398  else if (nOfCurlyBracket < 0) {
2399  error = 1 ;
2400  fError = "Missing \"{\"" ;
2401  }
2402  else if (nOfCurlyBracket > 0) {
2403  error = 1 ;
2404  fError = "Missing \"}\"" ;
2405  }
2406  else if (nOfSquareBracket < 0) {
2407  error = 1 ;
2408  fError = "Missing \"[\"" ;
2409  }
2410  else if (nOfSquareBracket > 0) {
2411  error = 1 ;
2412  fError = "Missing \"]\"" ;
2413  }
2414 
2415  ERROR_END:
2416  return error ;
2417 }
2418 
2419 ////////////////////////////////////////////////////////////////////////////////
2420 /// First parsing of the analyse sequence
2421 
2423 {
2424  fError = 0;
2425  fTabMax = 100;
2426  fTabSize = new FormSize_t[fTabMax];
2427  // we assume less than 100 parts in one formula
2428  // we will reallocate if necessary.
2429  fPos = 0;
2430  fShow = kFALSE;
2431  fOriginSize = size;
2432 
2433  //get current line attributes
2434  Short_t lineW = GetLineWidth();
2435  Int_t lineC = GetLineColor();
2436 
2437  TextSpec_t spec;
2438  spec.fAngle = angle;
2439  if (fTextFont%10 == 3) {
2440  Double_t hw = TMath::Max((Double_t)gPad->XtoPixel(gPad->GetX2()),
2441  (Double_t)gPad->YtoPixel(gPad->GetY1()));
2442  spec.fSize = size/hw;
2443  } else {
2444  spec.fSize = size;
2445  }
2446  spec.fColor = GetTextColor();
2447  spec.fFont = GetTextFont();
2448  Short_t halign = fTextAlign/10;
2449  Short_t valign = fTextAlign - 10*halign;
2450 
2451  TLatexFormSize fs = Anal1(spec,text,strlen(text));
2452 
2453  SetTextSize(size);
2454  SetTextAngle(angle);
2455  SetTextFont(spec.fFont);
2456  SetTextColor(spec.fColor);
2457  SetTextAlign(valign+10*halign);
2458  SetLineWidth(lineW);
2459  SetLineColor(lineC);
2460  return fs;
2461 }
2462 
2463 ////////////////////////////////////////////////////////////////////////////////
2464 /// Return height of current pad in pixels
2465 
2467 {
2468  Double_t w = gPad->GetAbsWNDC()*Double_t(gPad->GetWw());
2469  Double_t h = gPad->GetAbsHNDC()*Double_t(gPad->GetWh());
2470  if (w < h)
2471  return w;
2472  else
2473  return h;
2474 }
2475 
2476 ////////////////////////////////////////////////////////////////////////////////
2477 /// Return size of the formula along X in pad coordinates
2478 
2480 {
2481  if (!gPad) return 0;
2482  TString newText = GetTitle();
2483  if( newText.Length() == 0) return 0;
2484 
2485  // The text is a TMathText.
2486  if ( newText.Contains("\\") ) {
2487  TMathText tm(0., 0., newText.Data());
2488  return tm.GetXsize();
2489  }
2490 
2491  fError = 0 ;
2492  if (CheckLatexSyntax(newText)) {
2493  std::cout<<"\n*ERROR<TLatex>: "<<fError<<std::endl;
2494  std::cout<<"==> "<<GetTitle()<<std::endl;
2495  return 0;
2496  }
2497  fError = 0 ;
2498 
2499  const Char_t *text = newText.Data() ;
2500  Double_t angle_old = GetTextAngle();
2501  TLatexFormSize fs = FirstParse(0,GetTextSize(),text);
2502  SetTextAngle(angle_old);
2503  delete[] fTabSize;
2504  return TMath::Abs(gPad->AbsPixeltoX(Int_t(fs.Width())) - gPad->AbsPixeltoX(0));
2505 }
2506 
2507 ////////////////////////////////////////////////////////////////////////////////
2508 /// Return text size in pixels
2509 
2510 void TLatex::GetBoundingBox(UInt_t &w, UInt_t &h, Bool_t angle)
2511 {
2512  if (!gPad) return;
2513  TString newText = GetTitle();
2514  if( newText.Length() == 0) return;
2515 
2516  // The text is a TMathText.
2517  if ( newText.Contains("\\") ) {
2518  TMathText tm(0., 0., newText.Data());
2519  tm.GetBoundingBox(w, h);
2520  return;
2521  }
2522 
2523  fError = 0 ;
2524  if (CheckLatexSyntax(newText)) {
2525  std::cout<<"\n*ERROR<TLatex>: "<<fError<<std::endl;
2526  std::cout<<"==> "<<GetTitle()<<std::endl;
2527  return;
2528  }
2529  fError = 0 ;
2530 
2531  if (angle) {
2532  Int_t cBoxX[4], cBoxY[4];
2533  Int_t ptx, pty;
2534  if (TestBit(kTextNDC)) {
2535  ptx = gPad->UtoPixel(fX);
2536  pty = gPad->VtoPixel(fY);
2537  } else {
2538  ptx = gPad->XtoAbsPixel(gPad->XtoPad(fX));
2539  pty = gPad->YtoAbsPixel(gPad->YtoPad(fY));
2540  }
2541  GetControlBox(ptx, pty, fTextAngle, cBoxX, cBoxY);
2542  Int_t x1 = cBoxX[0];
2543  Int_t x2 = cBoxX[0];
2544  Int_t y1 = cBoxY[0];
2545  Int_t y2 = cBoxY[0];
2546  for (Int_t i=1; i<4; i++) {
2547  if (cBoxX[i] < x1) x1 = cBoxX[i];
2548  if (cBoxX[i] > x2) x2 = cBoxX[i];
2549  if (cBoxY[i] < y1) y1 = cBoxY[i];
2550  if (cBoxY[i] > y2) y2 = cBoxY[i];
2551  }
2552  w = x2-x1;
2553  h = y2-y1;
2554  } else {
2555  const Char_t *text = newText.Data() ;
2556  TLatexFormSize fs = FirstParse(GetTextAngle(),GetTextSize(),text);
2557  delete[] fTabSize;
2558  w = (UInt_t)fs.Width();
2559  h = (UInt_t)fs.Height();
2560  }
2561 }
2562 
2563 ////////////////////////////////////////////////////////////////////////////////
2564 /// Return size of the formula along Y in pad coordinates
2565 
2567 {
2568  if (!gPad) return 0;
2569  TString newText = GetTitle();
2570  if( newText.Length() == 0) return 0;
2571 
2572  // The text is a TMathText.
2573  if ( newText.Contains("\\") ) {
2574  TMathText tm(0., 0., newText.Data());
2575  return tm.GetYsize();
2576  }
2577 
2578  fError = 0 ;
2579  if (CheckLatexSyntax(newText)) {
2580  std::cout<<"\n*ERROR<TLatex>: "<<fError<<std::endl;
2581  std::cout<<"==> "<<GetTitle()<<std::endl;
2582  return 0;
2583  }
2584  fError = 0 ;
2585 
2586  const Char_t *text = newText.Data() ;
2587  Double_t angsav = fTextAngle;
2588  TLatexFormSize fs = FirstParse(0,GetTextSize(),text);
2589  fTextAngle = angsav;
2590  delete[] fTabSize;
2591  return TMath::Abs(gPad->AbsPixeltoY(Int_t(fs.Height())) - gPad->AbsPixeltoY(0));
2592 }
2593 
2594 ////////////////////////////////////////////////////////////////////////////////
2595 /// Read fs in fTabSize
2596 
2598 {
2599  fPos--;
2600  TLatexFormSize result(fTabSize[fPos].fWidth,fTabSize[fPos].fOver,fTabSize[fPos].fUnder);
2601  return result;
2602 }
2603 
2604 ////////////////////////////////////////////////////////////////////////////////
2605 /// Save fs values in array fTabSize
2606 
2608 {
2610  fTabSize[fPos].fOver = fs->Over();
2611  fTabSize[fPos].fUnder = fs->Under();
2612  fPos++;
2613  if (fPos>=fTabMax) {
2614  // allocate more memory
2615  FormSize_t *temp = new FormSize_t[fTabMax+100];
2616  // copy array
2617  memcpy(temp,fTabSize,fTabMax*sizeof(FormSize_t));
2618  fTabMax += 100;
2619  // free previous array
2620  delete [] fTabSize;
2621  // swap pointers
2622  fTabSize = temp;
2623  }
2624 }
2625 
2626 ////////////////////////////////////////////////////////////////////////////////
2627 /// Save primitive as a C++ statement(s) on output stream out
2628 
2629 void TLatex::SavePrimitive(std::ostream &out, Option_t * /*= ""*/)
2630 {
2631  char quote = '"';
2632 
2633  if (gROOT->ClassSaved(TLatex::Class())) {
2634  out<<" ";
2635  } else {
2636  out<<" TLatex *";
2637  }
2638 
2639  TString s = GetTitle();
2640 
2641  s.ReplaceAll("\\","\\\\");
2642  s.ReplaceAll("\"","\\\"");
2643  out<<" tex = new TLatex("<<fX<<","<<fY<<","<<quote<<s.Data()<<quote<<");"<<std::endl;
2644  if (TestBit(kTextNDC)) out<<"tex->SetNDC();"<<std::endl;
2645 
2646  SaveTextAttributes(out,"tex",11,0,1,62,0.05);
2647  SaveLineAttributes(out,"tex",1,1,1);
2648 
2649  out<<" tex->Draw();"<<std::endl;
2650 }
2651 
2652 ////////////////////////////////////////////////////////////////////////////////
2653 /// Set relative size of subscripts and superscripts
2654 
2655 void TLatex::SetIndiceSize(Double_t factorSize)
2656 {
2657  fFactorSize = factorSize;
2658 }
2659 
2660 ////////////////////////////////////////////////////////////////////////////////
2661 /// Set limit for text resizing of subscripts and superscripts
2662 
2663 void TLatex::SetLimitIndiceSize(Int_t limitFactorSize)
2664 {
2665  fLimitFactorSize = limitFactorSize;
2666 }
Double_t Width() const
Definition: TLatex.h:80
Double_t Under() const
Definition: TLatex.h:82
virtual const char * GetTitle() const
Returns title of object.
Definition: TNamed.h:52
virtual void SetLineWidth(Width_t lwidth)
Definition: TAttLine.h:57
Int_t fLimitFactorSize
Relative position of subscripts and superscripts.
Definition: TLatex.h:88
Bool_t fShow
error code
Definition: TLatex.h:90
virtual void PaintMathText(Double_t x, Double_t y, Double_t angle, Double_t size, const char *text)
Paint text (used by Paint()).
Definition: TMathText.cxx:584
const Double_t kPI
Definition: TLatex.cxx:21
virtual Float_t GetTextAngle() const
Definition: TAttText.h:47
static Vc_ALWAYS_INLINE int_v min(const int_v &x, const int_v &y)
Definition: vector.h:433
virtual Font_t GetTextFont() const
Definition: TAttText.h:49
Double_t Log(Double_t x)
Definition: TMath.h:526
Double_t fOriginSize
array of values for the different zones
Definition: TLatex.h:92
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:487
Ssiz_t Length() const
Definition: TString.h:390
float Float_t
Definition: RtypesCore.h:53
void Copy(TAttText &atttext) const
Copy this text attributes to a new TAttText.
Definition: TAttText.cxx:291
const char Option_t
Definition: RtypesCore.h:62
TString & ReplaceAll(const TString &s1, const TString &s2)
Definition: TString.h:635
TH1 * h
Definition: legend2.C:5
void Copy(TAttLine &attline) const
Copy this line attributes to a new TAttLine.
Definition: TAttLine.cxx:159
virtual Float_t GetTextSize() const
Definition: TAttText.h:50
#define gROOT
Definition: TROOT.h:344
Basic string class.
Definition: TString.h:137
virtual void GetControlBox(Int_t x, Int_t y, Double_t theta, Int_t cBoxX[4], Int_t cBoxY[4])
Return the text control box.
Definition: TText.cxx:419
TLatex helper struct holding the dimensions of a piece of text.
Definition: TLatex.h:46
Int_t fTabMax
Definition: TLatex.h:93
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:170
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
const Bool_t kFALSE
Definition: Rtypes.h:92
void GetBoundingBox(UInt_t &w, UInt_t &h, Bool_t angle=kFALSE)
Get the text width and height.
Definition: TMathText.cxx:486
virtual void Modify()
Change current line attributes if necessary.
Definition: TAttLine.cxx:229
Double_t fOver
Definition: TLatex.h:47
Double_t fAngle
Definition: TLatex.h:40
TLatex & operator=(const TLatex &)
Currently inside italic operator.
Definition: TLatex.cxx:445
TString & Insert(Ssiz_t pos, const char *s)
Definition: TString.h:592
static T square(T a)
Definition: main.cpp:132
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
TFile * f
void SetBit(UInt_t f, Bool_t set)
Set or unset the user status bits as specified in f.
Definition: TObject.cxx:732
Int_t PaintLatex1(Double_t x, Double_t y, Double_t angle, Double_t size, const char *text)
Drawing function.
Definition: TLatex.cxx:2115
virtual void AppendPad(Option_t *option="")
Append graphics object to current pad.
Definition: TObject.cxx:164
virtual void Modify()
Change current text attributes if necessary.
Definition: TAttText.cxx:303
TLatex * DrawLatexNDC(Double_t x, Double_t y, const char *text)
Draw this TLatex with new coordinates in NDC.
Definition: TLatex.cxx:1914
const char * Data() const
Definition: TString.h:349
virtual void SetTextFont(Font_t tfont=62)
Definition: TAttText.h:59
Double_t fSize
Definition: TLatex.h:40
static const double x2[5]
void DrawLine(Double_t x1, Double_t y1, Double_t x2, Double_t y2, TextSpec_t spec)
Draw a line in a Latex formula.
Definition: TLatex.cxx:1924
Double_t GetYsize()
Return size of the formula along Y in pad coordinates.
Definition: TLatex.cxx:2568
virtual Short_t GetTextAlign() const
Definition: TAttText.h:46
Double_t x[n]
Definition: legend1.C:17
void DrawCircle(Double_t x1, Double_t y1, Double_t r, TextSpec_t spec)
Draw an arc of ellipse in a Latex formula (right or left parenthesis)
Definition: TLatex.cxx:1944
virtual void SaveLineAttributes(std::ostream &out, const char *name, Int_t coldef=1, Int_t stydef=1, Int_t widdef=1)
Save line attributes as C++ statement(s) on output stream out.
Definition: TAttLine.cxx:257
void Class()
Definition: Class.C:29
TLatexFormSize Anal1(TextSpec_t spec, const Char_t *t, Int_t length)
Analyse function.
Definition: TLatex.cxx:486
const char * tab3
Definition: TXMLPlayer.cxx:125
int d
Definition: tornado.py:11
virtual void GetTextAscentDescent(UInt_t &a, UInt_t &d, const char *text) const
Return text ascent and descent for string text.
Definition: TText.cxx:519
To draw Mathematical Formula.
Definition: TLatex.h:33
virtual Color_t GetTextColor() const
Definition: TAttText.h:48
virtual void Paint(Option_t *option="")
Paint.
Definition: TLatex.cxx:2004
static const double x4[22]
TLatex * DrawLatex(Double_t x, Double_t y, const char *text)
Make a copy of this object with the new parameters And copy object attributes.
Definition: TLatex.cxx:1901
if on multiple lines(like in C++).**The" * configuration fragment. * * The "import myobject continue
Parses the configuration file.
Definition: HLFactory.cxx:368
tuple np
Definition: multifit.py:30
Base class for several text objects.
Definition: TText.h:42
virtual void SetNDC(Bool_t isNDC=kTRUE)
Set NDC mode on if isNDC = kTRUE, off otherwise.
Definition: TText.cxx:809
Double_t Over() const
Definition: TLatex.h:81
Int_t CheckLatexSyntax(TString &text)
Check if the Latex syntax is correct.
Definition: TLatex.cxx:2220
virtual void PaintLatex(Double_t x, Double_t y, Double_t angle, Double_t size, const char *text)
Main drawing function.
Definition: TLatex.cxx:2025
Double_t GetHeight() const
Return height of current pad in pixels.
Definition: TLatex.cxx:2468
virtual void Text(Double_t x, Double_t y, const char *string)=0
void DrawParenthesis(Double_t x1, Double_t y1, Double_t r1, Double_t r2, Double_t phimin, Double_t phimax, TextSpec_t spec)
Draw an arc of ellipse in a Latex formula (right or left parenthesis)
Definition: TLatex.cxx:1973
char * out
Definition: TBase64.cxx:29
virtual void SetTextAlign(Short_t align=11)
Definition: TAttText.h:55
const Char_t * fError
Definition: TLatex.h:89
const char * tab2
Definition: TXMLPlayer.cxx:124
void Set(Double_t x, Double_t y1, Double_t y2)
Definition: TLatex.h:71
virtual void SetLineColor(Color_t lcolor)
Definition: TAttLine.h:54
Float_t fTextAngle
Definition: TAttText.h:35
TThread * t[5]
Definition: threadsh1.C:13
virtual void SavePrimitive(std::ostream &out, Option_t *option="")
Save primitive as a C++ statement(s) on output stream out.
Definition: TLatex.cxx:2631
Double_t fWidth
Definition: TLatex.h:47
TNamed & operator=(const TNamed &rhs)
TNamed assignment operator.
Definition: TNamed.cxx:40
ROOT::R::TRInterface & r
Definition: Object.C:4
Double_t length(const TVector2 &v)
Definition: CsgOps.cxx:347
Double_t GetXsize()
Return size of the formula along X in pad coordinates.
Definition: TLatex.cxx:2481
TLatexFormSize Readfs()
Read fs in fTabSize.
Definition: TLatex.cxx:2599
virtual void SaveTextAttributes(std::ostream &out, const char *name, Int_t alidef=12, Float_t angdef=0, Int_t coldef=1, Int_t fondef=61, Float_t sizdef=1)
Save text attributes as C++ statement(s) on output stream out.
Definition: TAttText.cxx:344
void Copy(TObject &text) const
Copy this text to text.
Definition: TText.cxx:107
unsigned int r1[N_CITIES]
Definition: simanTSP.cxx:321
virtual void SetTextAngle(Float_t tangle=0)
Definition: TAttText.h:56
Double_t fY
Definition: TText.h:46
unsigned int UInt_t
Definition: RtypesCore.h:42
Bool_t TestBit(UInt_t f) const
Definition: TObject.h:173
TMarker * m
Definition: textangle.C:8
Double_t fUnder
Definition: TLatex.h:47
TLatexFormSize Analyse(Double_t x, Double_t y, TextSpec_t spec, const Char_t *t, Int_t length)
Analyse and paint the TLatex formula.
Definition: TLatex.cxx:516
tuple w
Definition: qtexample.py:51
short Short_t
Definition: RtypesCore.h:35
virtual void GetTextExtent(UInt_t &w, UInt_t &h, const char *text) const
Return text extent for string text.
Definition: TText.cxx:583
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:2512
f1 SetLineWidth(4)
Font_t fTextFont
Definition: TAttText.h:39
virtual Color_t GetLineColor() const
Definition: TAttLine.h:47
TLatexFormSize FirstParse(Double_t angle, Double_t size, const Char_t *text)
First parsing of the analyse sequence.
Definition: TLatex.cxx:2424
void Savefs(TLatexFormSize *fs)
Save fs values in array fTabSize.
Definition: TLatex.cxx:2609
#define gVirtualX
Definition: TVirtualX.h:362
Double_t Cos(Double_t)
Definition: TMath.h:424
FormSize_t * fTabSize
is true during the second pass (Painting)
Definition: TLatex.h:91
Double_t Pi()
Definition: TMath.h:44
virtual ~TLatex()
Destructor.
Definition: TLatex.cxx:420
Double_t Exp(Double_t x)
Definition: TMath.h:495
static const double x1[5]
Int_t fPos
Maximum allocation for array fTabSize;.
Definition: TLatex.h:94
#define ClassImp(name)
Definition: Rtypes.h:279
TLatex helper class used to compute the size of a portion of a formula.
Definition: TLatex.h:54
double Double_t
Definition: RtypesCore.h:55
TText * text
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
Double_t y[n]
Definition: legend1.C:17
To draw TeX Mathematical Formula.
Definition: TMathText.h:26
static Vc_ALWAYS_INLINE int_v max(const int_v &x, const int_v &y)
Definition: vector.h:440
Double_t fFactorPos
Relative size of subscripts and superscripts.
Definition: TLatex.h:87
void Copy(TObject &text) const
Copy this TLatex object to another TLatex.
Definition: TLatex.cxx:467
Binding & operator=(OUT(*fun)(void))
Mother of all ROOT objects.
Definition: TObject.h:58
TLatex helper struct holding the attributes of a piece of text.
Definition: TLatex.h:39
char Char_t
Definition: RtypesCore.h:29
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:202
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition: TString.h:567
Double_t Sin(Double_t)
Definition: TMath.h:421
Bool_t fItalic
Current position in array fTabSize;.
Definition: TLatex.h:95
R__EXTERN TVirtualPS * gVirtualPS
Definition: TVirtualPS.h:91
#define gPad
Definition: TVirtualPad.h:288
virtual void SetLimitIndiceSize(Int_t limitFactorSize)
Set limit for text resizing of subscripts and superscripts.
Definition: TLatex.cxx:2665
virtual void SetTextColor(Color_t tcolor=1)
Definition: TAttText.h:57
double result[121]
TVirtualPS is an abstract interface to Postscript, PDF, SVG.
Definition: TVirtualPS.h:40
virtual void SetTextSize(Float_t tsize=1)
Definition: TAttText.h:60
const Bool_t kTRUE
Definition: Rtypes.h:91
virtual Width_t GetLineWidth() const
Definition: TAttLine.h:49
TObject * obj
virtual void PaintText(Double_t x, Double_t y, const char *text)
Draw this text with new coordinates.
Definition: TText.cxx:740
Double_t GetYsize(void)
Get Y size.
Definition: TMathText.cxx:527
Double_t GetXsize(void)
Get X size.
Definition: TMathText.cxx:506
Line Attributes class.
Definition: TAttLine.h:32
unsigned int r2[N_CITIES]
Definition: simanTSP.cxx:322
Double_t fFactorSize
Definition: TLatex.h:86
Double_t fX
Definition: TText.h:45
virtual void SetIndiceSize(Double_t factorSize)
Set relative size of subscripts and superscripts.
Definition: TLatex.cxx:2657
Short_t fTextAlign
Definition: TAttText.h:37
static const double x3[11]