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