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