Logo ROOT   6.10/00
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TPostScript.cxx
Go to the documentation of this file.
1 // @(#)root/postscript:$Id$
2 // Author: Rene Brun, Olivier Couet, Pierre Juillot, Oleksandr Grebenyuk, Yue Shi Lai
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 /** \class TPostScript
13 \ingroup PS
14 
15 Interface to PostScript.
16 
17 To generate a Postscript (or encapsulated ps) file corresponding to
18 a single image in a canvas, you can:
19 
20  - Select the <B>Print PostScript</B> item in the canvas <B>File</B> menu.
21  By default, a Postscript file with the name of the canvas.ps is generated.
22  - Click in the canvas area, near the edges, with the right mouse button
23  and select the <B>Print</B> item. You can select the name of the Postscript
24  file. If the file name is xxx.ps, you will generate a Postscript file named
25  xxx.ps. If the file name is xxx.eps, you generate an encapsulated Postscript
26  file instead.
27  - In your program (or macro), you can type:
28 ~~~ {.cpp}
29  c1->Print("xxx.ps")</B> or <B>c1->Print("xxx.eps").
30 ~~~
31  This will generate a file corresponding to the picture in the canvas
32  pointed by `c1`.
33 ~~~ {.cpp}
34  pad1->Print("xxx.ps");
35 ~~~
36  prints only the picture in the pad pointed by `pad1`.
37 
38 The size of the Postscript picture, by default, is computed to keep the aspect
39 ratio of the picture on the screen, where the size along x is always 20cm. You
40 can set the size of the PostScript picture before generating the picture
41 with a command such as:
42 
43 ~~~ {.cpp}
44  TPostScript myps("myfile.ps",111)
45  myps.Range(xsize,ysize);
46  object->Draw();
47  myps.Close();
48 ~~~
49 You can set the default paper size with:
50 ~~~ {.cpp}
51  gStyle->SetPaperSize(xsize,ysize);
52 ~~~
53 You can resume writing again in this file with `myps.Open();`.
54 Note that you may have several Postscript files opened simultaneously.
55 
56  ## Output type
57 
58 The output type allows to define how the PostScript output will looks like.
59 It allows to define the page format (A4, Legal etc..), the orientation
60 (Portrait, Landscape) and the number of images (zones) per page.
61 The output type has the following form:
62 
63 ~~~ {.cpp}
64  [Format][Nx][Ny][Type]
65 ~~~
66 
67 Where:
68 
69  - Format : Is an integer between 0 and 99 defining the page format:
70 ~~~ {.cpp}
71  Format = 3 the paper is in the standard A3 format.
72  Format = n (1<n<98) is an An format.
73  Format = 4 and Format=0 are the same and define an A4 page.
74  The A0 format is selected by Format=99.
75  The US format Letter is selected by Format = 100.
76  The US format Legal is selected by Format = 200.
77  The US format Ledger is selected by Format = 300.
78 ~~~
79  - Nx, Ny : Specify respectively the number of zones on the x and y axis.
80  Nx and Ny are integers between 1 and 9.
81  - Type : Can be equal to:
82  - 1 : Portrait mode with a small margin at the bottom of the page.
83  - 2 : Landscape mode with a small margin at the bottom of the page.
84  - 4 : Portrait mode with a large margin at the bottom of the page.
85  - 5 : Landscape mode with a large margin at the bottom of the page.
86  The large margin is useful for some PostScript printers (very often
87  for the colour printers) as they need more space to grip the paper
88  for mechanical reasons. Note that some PostScript colour printers
89  can also use the so called special A4 format permitting the full
90  usage of the A4 area; in this case larger margins are not necessary
91  and Type=1 or 2 can be used.
92  - 3 : Encapsulated PostScript. This Type permits the generation of files
93  which can be included in other documents, for example in LaTeX files.
94 
95 ## Making several pictures in the same Postscript file: case 1
96 
97 The following macro is an example illustrating how to open a Postscript
98 file and draw several pictures. The generation of a new Postscript page
99 is automatic when `TCanvas::Clear` is called by `object->Draw()`.
100 
101 ~~~ {.cpp}
102  {
103  TFile f("hsimple.root");
104  TCanvas c1("c1","canvas",800,600);
105 
106  // select postscript output type
107  // type = 111 portrait ps
108  // type = 112 landscape ps
109  // type = 113 eps
110  Int_t type = 111;
111 
112  // create a postscript file and set the paper size
113  TPostScript ps("test.ps",type);
114  ps.Range(16,24); //set x,y of printed page
115 
116  // draw 3 histograms from file hsimple.root on separate pages
117  hpx->Draw();
118  c1.Update(); //force drawing in a macro
119  hprof->Draw();
120  c1.Update();
121  hpx->Draw("lego1");
122  c1.Update();
123  ps.Close();
124  }
125 ~~~
126 
127 ## Making several pictures in the same Postscript file: case 2
128 
129 This example shows 2 pages. The canvas is divided.
130 `TPostScript::NewPage` must be called before starting a new
131 picture.`object->Draw` does not clear the canvas in this case
132 because we clear only the pads and not the main canvas.
133 Note that `c1->Update` must be called at the end of the first picture.
134 
135 ~~~ {.cpp}
136  {
137  TFile *f1 = new TFile("hsimple.root");
138  TCanvas *c1 = new TCanvas("c1");
139  TPostScript *ps = new TPostScript("file.ps",112);
140  c1->Divide(2,1);
141  // picture 1
142  ps->NewPage();
143  c1->cd(1);
144  hpx->Draw();
145  c1->cd(2);
146  hprof->Draw();
147  c1->Update();
148 
149  // picture 2
150  ps->NewPage();
151  c1->cd(1);
152  hpxpy->Draw();
153  c1->cd(2);
154  ntuple->Draw("px");
155  c1->Update();
156  ps->Close();
157 
158  // invoke Postscript viewer
159  gSystem->Exec("gs file.ps");
160  }
161 ~~~
162 
163 ## Making several pictures in the same Postscript file: case 3
164 This is the recommended way. If the Postscript file name finishes with
165 "(", the file remains opened (it is not closed). If the Postscript file name
166 finishes with ")" and the file has been opened with "(", the file is closed.
167 
168 Example:
169 ~~~ {.cpp}
170  {
171  TCanvas c1("c1");
172  h1.Draw();
173  c1.Print("c1.ps("); // write canvas and keep the ps file open
174  h2.Draw();
175  c1.Print("c1.ps"); // canvas is added to "c1.ps"
176  h3.Draw();
177  c1.Print("c1.ps)"); // canvas is added to "c1.ps" and ps file is closed
178  }
179 ~~~
180 The `TCanvas::Print("file.ps(")` mechanism is very useful, but it can
181 be a little inconvenient to have the action of opening/closing a file being
182 atomic with printing a page. Particularly if pages are being generated in some
183 loop one needs to detect the special cases of first and last page and then
184 munge the argument to Print() accordingly.
185 The "[" and "]" can be used instead of "(" and ")" as shown below.
186 
187 Example:
188 ~~~ {.cpp}
189  c1.Print("file.ps["); // No actual print, just open file.ps
190 
191  for (int i=0; i<10; ++i) {
192  // fill canvas for context i
193  // ...
194 
195  c1.Print("file.ps"); // Actually print canvas to the file
196  }
197 
198  c1.Print("file.ps]"); // No actual print, just close the file
199 ~~~
200 
201  ## Color Model
202 
203 TPostScript support two color model RGB and CMYK. CMY and CMYK models are
204 subtractive color models unlike RGB which is an additive. They are mainly
205 used for printing purposes. CMY means Cyan Magenta Yellow to convert RGB
206 to CMY it is enough to do: C=1-R, M=1-G and Y=1-B. CMYK has one more
207 component K (black). The conversion from RGB to CMYK is:
208 
209 ~~~ {.cpp}
210  Double_t Black = TMath::Min(TMath::Min(1-Red,1-Green),1-Blue);
211  Double_t Cyan = (1-Red-Black)/(1-Black);
212  Double_t Magenta = (1-Green-Black)/(1-Black);
213  Double_t Yellow = (1-Blue-Black)/(1-Black);
214 ~~~
215 CMYK add the black component which allows to have a better quality for black
216 printing. PostScript support the CMYK model.
217 
218 To change the color model use gStyle->SetColorModelPS(c).
219 
220  - c = 0 means TPostScript will use RGB color model (default)
221  - c = 1 means TPostScript will use CMYK color model
222 */
223 
224 #ifdef WIN32
225 #pragma optimize("",off)
226 #endif
227 
228 #include <stdlib.h>
229 #include <string.h>
230 #include <ctype.h>
231 #include <wchar.h>
232 #include <fontconfig/fontconfig.h>
233 
234 #include "Riostream.h"
235 #include "Byteswap.h"
236 #include "TROOT.h"
237 #include "TColor.h"
238 #include "TVirtualPad.h"
239 #include "TPoints.h"
240 #include "TPostScript.h"
241 #include "TStyle.h"
242 #include "TMath.h"
243 #include "TText.h"
244 #include "TSystem.h"
245 #include "TEnv.h"
246 
247 #include "../../../graf2d/mathtext/inc/fontembed.h"
248 
249 // to scale fonts to the same size as the old TT version
250 const Float_t kScale = 0.93376068;
251 
252 // Array defining if a font must be embedded or not.
253 static Bool_t MustEmbed[32];
254 
256 
258 
259 ////////////////////////////////////////////////////////////////////////////////
260 /// Default PostScript constructor
261 
263 {
264  fStream = 0;
265  fType = 0;
266  gVirtualPS = this;
267  fBlue = 0.;
268  fBoundingBox = kFALSE;
269  fClear = kFALSE;
270  fClip = 0;
271  fClipStatus = kFALSE;
272  fCurrentColor = 0;
273  fDXC = 0.;
274  fDYC = 0.;
275  fFX = 0.;
276  fFY = 0.;
277  fGreen = 0.;
278  fIXzone = 0;
279  fIYzone = 0;
280  fLastCellBlue = 0;
281  fLastCellGreen = 0;
282  fLastCellRed = 0;
283  fLineJoin = 0;
284  fLineScale = 0.;
285  fMarkerSizeCur = 0.;
286  fMaxLines = 0;
287  fMaxsize = 0;
288  fMode = 0;
289  fNBSameColorCell = 0;
290  fNXzone = 0;
291  fNYzone = 0;
292  fNbCellLine = 0;
293  fNbCellW = 0;
294  fNbinCT = 0;
295  fNpages = 0;
296  fRange = kFALSE;
297  fRed = 0.;
298  fSave = 0;
299  fX1v = 0.;
300  fX1w = 0.;
301  fX2v = 0.;
302  fX2w = 0.;
303  fXC = 0.;
304  fXVP1 = 0.;
305  fXVP2 = 0.;
306  fXVS1 = 0.;
307  fXVS2 = 0.;
308  fXsize = 0.;
309  fY1v = 0.;
310  fY1w = 0.;
311  fY2v = 0.;
312  fY2w = 0.;
313  fYC = 0.;
314  fYVP1 = 0.;
315  fYVP2 = 0.;
316  fYVS1 = 0.;
317  fYVS2 = 0.;
318  fYsize = 0.;
319  fZone = kFALSE;
320  fFileName = "";
321  fFontEmbed = kFALSE;
322  Int_t i;
323  for (i=0; i<32; i++) fPatterns[i] = 0;
324  for (i=0; i<32; i++) MustEmbed[i] = kFALSE;
325  SetTitle("PS");
326 }
327 
328 ////////////////////////////////////////////////////////////////////////////////
329 /// Initialize the PostScript interface
330 ///
331 /// - fname : PostScript file name
332 /// - wtype : PostScript workstation type
333 ///
334 ///
335 /// The possible workstation types are:
336 /// - 111 ps Portrait
337 /// - 112 ps Landscape
338 /// - 113 eps
339 
340 TPostScript::TPostScript(const char *fname, Int_t wtype)
341 :TVirtualPS(fname, wtype)
342 {
343  fStream = 0;
344  SetTitle("PS");
345  Open(fname, wtype);
346 }
347 
348 ////////////////////////////////////////////////////////////////////////////////
349 /// Open a PostScript file
350 
351 void TPostScript::Open(const char *fname, Int_t wtype)
352 {
353  if (fStream) {
354  Warning("Open", "postscript file already open");
355  return;
356  }
357 
358  fMarkerSizeCur = 0;
359  fCurrentColor = 0;
360  fRed = -1;
361  fGreen = -1;
362  fBlue = -1;
363  fLenBuffer = 0;
364  fClip = 0;
365  fType = abs(wtype);
366  fClear = kTRUE;
367  fZone = kFALSE;
368  fSave = 0;
369  fFontEmbed = kFALSE;
373  fMode = fType%10;
374  Float_t xrange, yrange;
375  if (gPad) {
376  Double_t ww = gPad->GetWw();
377  Double_t wh = gPad->GetWh();
378  if (fType == 113) {
379  ww *= gPad->GetWNDC();
380  wh *= gPad->GetHNDC();
381  }
382  Double_t ratio = wh/ww;
383  if (fType == 112) {
384  xrange = fYsize;
385  yrange = xrange*ratio;
386  if (yrange > fXsize) { yrange = fXsize; xrange = yrange/ratio;}
387  } else {
388  xrange = fXsize;
389  yrange = fXsize*ratio;
390  if (yrange > fYsize) { yrange = fYsize; xrange = yrange/ratio;}
391  }
392  fXsize = xrange; fYsize = yrange;
393  }
394 
395  // Open OS file
396  fFileName = fname;
397  fStream = new std::ofstream(fFileName.Data(),std::ios::out);
399  printf("ERROR in TPostScript::Open: Cannot open file:%s\n",fFileName.Data());
400  return;
401  }
402  gVirtualPS = this;
403 
404  for (Int_t i=0;i<fSizBuffer;i++) fBuffer[i] = ' ';
405  if( fType == 113) {
407  PrintStr("%!PS-Adobe-2.0 EPSF-2.0@");
408  } else {
410  PrintStr("%!PS-Adobe-2.0@");
411  Initialize();
412  }
413 
415  fRange = kFALSE;
416 
417  // Set a default range
418  Range(fXsize, fYsize);
419 
420  fPrinted = kFALSE;
421  if (fType == 113) NewPage();
422 }
423 
424 ////////////////////////////////////////////////////////////////////////////////
425 /// Default PostScript destructor
426 
428 {
429  Close();
430 }
431 
432 ////////////////////////////////////////////////////////////////////////////////
433 /// Close a PostScript file
434 
436 {
437  if (!gVirtualPS) return;
438  if (!fStream) return;
439  if (gPad) gPad->Update();
440  if( fMode != 3) {
441  SaveRestore(-1);
442  if( fPrinted ) { PrintStr("showpage@"); SaveRestore(-1);}
443  PrintStr("@");
444  PrintStr("%%Trailer@");
445  PrintStr("%%Pages: ");
447  PrintStr("@");
448  while (fSave > 0) { SaveRestore(-1); }
449  } else {
450  PrintStr("@");
451  while (fSave > 0) { SaveRestore(-1); }
452  PrintStr("showpage@");
453  PrintStr("end@");
454  }
455  PrintStr("@");
456  PrintStr("%%EOF@");
457 
458  // Embed the fonts previously used by TMathText
459  if (!fFontEmbed) {
460  // Close the file fFileName
461  if (fStream) {
462  PrintStr("@");
463  fStream->close(); delete fStream; fStream = 0;
464  }
465 
466  // Rename the file fFileName
467  TString tmpname = Form("%s_tmp_%d",fFileName.Data(),gSystem->GetPid());
468  if (gSystem->Rename( fFileName.Data() , tmpname.Data())) {
469  Error("Text", "Cannot open temporary file: %s\n", tmpname.Data());
470  return;
471  }
472 
473  // Reopen the file fFileName
474  fStream = new std::ofstream(fFileName.Data(),std::ios::out);
476  Error("Text", "Cannot open file: %s\n", fFileName.Data());
477  return;
478  }
479 
480  // Embed the fonts at the right place
481  FILE *sg = fopen(tmpname.Data(),"r");
482  if (sg == 0) {
483  Error("Text", "Cannot open file: %s\n", tmpname.Data());
484  return;
485  }
486  char line[255];
487  while (fgets(line,255,sg)) {
488  if (strstr(line,"EndComments")) PrintStr("%%DocumentNeededResources: ProcSet (FontSetInit)@");
489  fStream->write(line,strlen(line));
490  if (!fFontEmbed && strstr(line,"m5")) {
491  FontEmbed();
492  PrintStr("@");
493  }
494  }
495  fclose(sg);
496  if (gSystem->Unlink(tmpname.Data())) return;
497  }
498 
499  fFontEmbed = kFALSE;
500 
501  // Close file stream
502 
503  if (fStream) { fStream->close(); delete fStream; fStream = 0;}
504 
505  gVirtualPS = 0;
506 }
507 
508 ////////////////////////////////////////////////////////////////////////////////
509 /// Activate an already open PostScript file
510 
512 {
513  if (!fType) {
514  Error("On", "no postscript file open");
515  Off();
516  return;
517  }
518  gVirtualPS = this;
519 }
520 
521 ////////////////////////////////////////////////////////////////////////////////
522 /// Deactivate an already open PostScript file
523 
525 {
526  gVirtualPS = 0;
527 }
528 
529 ////////////////////////////////////////////////////////////////////////////////
530 /// Draw a Cell Array
531 ///
532 /// Drawing a PostScript Cell Array is in fact done thanks to three
533 /// procedures: CellArrayBegin, CellArrayFill, and CellArrayEnd.
534 ///
535 /// - CellArrayBegin: Initiate the Cell Array by writing the necessary
536 /// PostScript procedures and the initial values of the
537 /// required parameters. The input parameters are:
538 /// - W: number of boxes along the width.
539 /// - H: number of boxes along the height
540 /// - x1,x2,y1,y2: First box coordinates.
541 /// - CellArrayFill: Is called for each box of the Cell Array. The first
542 /// box is the top left one and the last box is the
543 /// bottom right one. The input parameters are the Red,
544 /// Green, and Blue components of the box colour. These
545 /// Levels are between 0 and 255.
546 /// - CellArrayEnd: Finishes the Cell Array.
547 ///
548 /// PostScript cannot handle arrays larger than 65535. So the Cell Array
549 /// is drawn in several pieces.
550 
552  Double_t y1, Double_t y2)
553 {
554  Int_t ix1 = XtoPS(x1);
555  Int_t iy1 = YtoPS(y1);
556 
557  Float_t wt = (288/2.54)*gPad->GetAbsWNDC()*
558  fXsize*((x2 - x1)/(gPad->GetX2()-gPad->GetX1()));
559  Float_t ht = (288/2.54)*gPad->GetAbsHNDC()*
560  fYsize*((y2 - y1)/(gPad->GetY2()-gPad->GetY1()));
561 
562  fLastCellRed = 300;
563  fLastCellGreen = 300;
564  fLastCellBlue = 300;
565  fNBSameColorCell = 0;
566 
567  fNbinCT = 0;
568  fNbCellW = W;
569  fNbCellLine = 0;
570  fMaxLines = 40000/(3*fNbCellW);
571 
572  // Define some parameters
573  PrintStr("@/WT"); WriteReal(wt) ; PrintStr(" def"); // Cells width
574  PrintStr(" /HT"); WriteReal(ht) ; PrintStr(" def"); // Cells height
575  PrintStr(" /XS"); WriteInteger(ix1) ; PrintStr(" def"); // X start
576  PrintStr(" /YY"); WriteInteger(iy1) ; PrintStr(" def"); // Y start
577  PrintStr(" /NX"); WriteInteger(W) ; PrintStr(" def"); // Number of columns
578  PrintStr(" /NY"); WriteInteger(fMaxLines); PrintStr(" def"); // Number of lines
579 
580  // This PS procedure draws one cell.
581  PrintStr(" /DrawCell ");
582  PrintStr( "{WT HT XX YY bf");
583  PrintStr( " /NBBD NBBD 1 add def");
584  PrintStr( " NBBD NBB eq {exit} if");
585  PrintStr( " /XX WT XX add def");
586  PrintStr( " IX NX eq ");
587  PrintStr( "{/YY YY HT sub def");
588  PrintStr( " /XX XS def");
589  PrintStr( " /IX 0 def} if");
590  PrintStr( " /IX IX 1 add def} def");
591 
592  // This PS procedure draws fMaxLines line. It takes care of duplicated
593  // colors. Values "n" greater than 300 mean than the previous color
594  // should be duplicated n-300 times.
595  PrintStr(" /DrawCT ");
596  PrintStr( "{/NBB NX NY mul def");
597  PrintStr( " /XX XS def");
598  PrintStr( " /IX 1 def");
599  PrintStr( " /NBBD 0 def");
600  PrintStr( " /RC 0 def /GC 1 def /BC 2 def");
601  PrintStr( " 1 1 NBB ");
602  PrintStr( "{/NB CT RC get def");
603  PrintStr( " NB 301 ge ");
604  PrintStr( "{/NBL NB 300 sub def");
605  PrintStr( " 1 1 NBL ");
606  PrintStr( "{DrawCell}");
607  PrintStr( " for");
608  PrintStr( " /RC RC 1 add def");
609  PrintStr( " /GC RC 1 add def");
610  PrintStr( " /BC RC 2 add def}");
611  PrintStr( "{CT RC get 255 div CT GC get 255 div CT BC get 255 div setrgbcolor");
612  PrintStr( " DrawCell");
613  PrintStr( " /RC RC 3 add def");
614  PrintStr( " /GC GC 3 add def");
615  PrintStr( " /BC BC 3 add def} ifelse NBBD NBB eq {exit} if} for");
616  PrintStr( " /YY YY HT sub def clear} def");
617 
618  PrintStr(" /CT [");
619 }
620 
621 ////////////////////////////////////////////////////////////////////////////////
622 /// Paint the Cell Array
623 
625 {
626  if (fLastCellRed == r && fLastCellGreen == g && fLastCellBlue == b) {
628  } else {
629  if (fNBSameColorCell != 0 ) {
631  fNBSameColorCell = 0;
632  }
633  WriteInteger(r);
634  WriteInteger(g);
635  WriteInteger(b);
636  fLastCellRed = r;
637  fLastCellGreen = g;
638  fLastCellBlue = b;
639  }
640 
641  fNbinCT++;
642  if (fNbinCT == fNbCellW) {
643  fNbCellLine++;
644  fNbinCT = 0;
645  }
646 
647  if (fNbCellLine == fMaxLines) {
649  PrintStr("] def DrawCT /CT [");
650  fNbCellLine = 0;
651  fLastCellRed = 300;
652  fLastCellGreen = 300;
653  fLastCellBlue = 300;
654  fNBSameColorCell = 0;
655  fNbinCT = 0;
656  }
657 }
658 
659 ////////////////////////////////////////////////////////////////////////////////
660 /// End the Cell Array painting
661 
663 {
665  PrintStr("] def /NY");
667  PrintStr(" def DrawCT ");
668 }
669 
670 ////////////////////////////////////////////////////////////////////////////////
671 /// Define the markers
672 
674 {
675  PrintStr("/mp {newpath /y exch def /x exch def} def@");
676  PrintStr("/side {[w .77 mul w .23 mul] .385 w mul sd w 0 l currentpoint t -144 r} def@");
677  PrintStr("/mr {mp x y w2 0 360 arc} def /m24 {mr s} def /m20 {mr f} def@");
678  PrintStr("/mb {mp x y w2 add m w2 neg 0 d 0 w neg d w 0 d 0 w d cl} def@");
679  PrintStr("/mt {mp x y w2 add m w2 neg w neg d w 0 d cl} def@");
680  PrintStr("/w4 {w 4 div} def@");
681  PrintStr("/w6 {w 6 div} def@");
682  PrintStr("/w8 {w 8 div} def@");
683  PrintStr("/m21 {mb f} def /m25 {mb s} def /m22 {mt f} def /m26{mt s} def@");
684  PrintStr("/m23 {mp x y w2 sub m w2 w d w neg 0 d cl f} def@");
685  PrintStr("/m27 {mp x y w2 add m w3 neg w2 neg d w3 w2 neg d w3 w2 d cl s} def@");
686  PrintStr("/m28 {mp x w2 sub y w2 sub w3 add m w3 0 d ");
687  PrintStr(" 0 w3 neg d w3 0 d 0 w3 d w3 0 d ");
688  PrintStr(" 0 w3 d w3 neg 0 d 0 w3 d w3 neg 0 d");
689  PrintStr(" 0 w3 neg d w3 neg 0 d cl s } def@");
690  PrintStr("/m29 {mp gsave x w2 sub y w2 add w3 sub m currentpoint t");
691  PrintStr(" 4 {side} repeat cl fill gr} def@");
692  PrintStr("/m30 {mp gsave x w2 sub y w2 add w3 sub m currentpoint t");
693  PrintStr(" 4 {side} repeat cl s gr} def@");
694  PrintStr("/m31 {mp x y w2 sub m 0 w d x w2 sub y m w 0 d");
695  PrintStr(" x w2 sub y w2 add m w w neg d x w2 sub y w2");
696  PrintStr(" sub m w w d s} def@");
697  PrintStr("/m32 {mp x y w2 sub m w2 w d w neg 0 d cl s} def@");
698  PrintStr("/m33 {mp x y w2 add m w3 neg w2 neg d w3 w2 neg d w3 w2 d cl f} def@");
699  PrintStr("/m34 {mp x w2 sub y w2 sub w3 add m w3 0 d ");
700  PrintStr(" 0 w3 neg d w3 0 d 0 w3 d w3 0 d ");
701  PrintStr(" 0 w3 d w3 neg 0 d 0 w3 d w3 neg 0 d");
702  PrintStr(" 0 w3 neg d w3 neg 0 d cl f } def@");
703  PrintStr("/m35 {mp x y w2 add m w2 neg w2 neg d w2 w2 neg d w2 w2 d w2 neg w2 d");
704  PrintStr(" x y w2 sub m 0 w d x w2 sub y m w 0 d s} def@");
705  PrintStr("/m36 {mb x w2 sub y w2 add m w w neg d x w2 sub y w2 sub m w w d s} def@");
706  PrintStr("/m37 {mp x y m w4 neg w2 d w4 neg w2 neg d w2 0 d ");
707  PrintStr(" w4 neg w2 neg d w2 0 d w4 neg w2 d w2 0 d w4 neg w2 d w4 neg w2 neg d cl s} def@");
708  PrintStr("/m38 {mp x w4 sub y w2 add m w4 neg w4 neg d 0 w2 neg d w4 w4 neg d");
709  PrintStr(" w2 0 d w4 w4 d 0 w2 d w4 neg w4 d w2 neg 0 d");
710  PrintStr(" x y w2 sub m 0 w d x w2 sub y m w 0 d cl s} def@");
711  PrintStr("/m39 {mp x y m w4 neg w2 d w4 neg w2 neg d w2 0 d ");
712  PrintStr(" w4 neg w2 neg d w2 0 d w4 neg w2 d w2 0 d w4 neg w2 d w4 neg w2 neg d cl f} def@");
713  PrintStr("/m40 {mp x y m w4 w2 d w4 w4 neg d w2 neg w4 neg d w2 w4 neg d w4 neg w4 neg d");
714  PrintStr(" w4 neg w2 d w4 neg w2 neg d w4 neg w4 d w2 w4 d w2 neg w4 d w4 w4 d w4 w2 neg d cl s} def@");
715  PrintStr("/m41 {mp x y m w4 w2 d w4 w4 neg d w2 neg w4 neg d w2 w4 neg d w4 neg w4 neg d");
716  PrintStr(" w4 neg w2 d w4 neg w2 neg d w4 neg w4 d w2 w4 d w2 neg w4 d w4 w4 d w4 w2 neg d cl f} def@");
717  PrintStr("/m42 {mp x y w2 add m w8 neg w2 -3 4 div mul d w2 -3 4 div mul w8 neg d");
718  PrintStr(" w2 3 4 div mul w8 neg d w8 w2 -3 4 div mul d");
719  PrintStr(" w8 w2 3 4 div mul d w2 3 4 div mul w8 d");
720  PrintStr(" w2 -3 4 div mul w8 d w8 neg w2 3 4 div mul d cl s} def@");
721  PrintStr("/m43 {mp x y w2 add m w8 neg w2 -3 4 div mul d w2 -3 4 div mul w8 neg d");
722  PrintStr(" w2 3 4 div mul w8 neg d w8 w2 -3 4 div mul d");
723  PrintStr(" w8 w2 3 4 div mul d w2 3 4 div mul w8 d");
724  PrintStr(" w2 -3 4 div mul w8 d w8 neg w2 3 4 div mul d cl f} def@");
725  PrintStr("/m44 {mp x y m w6 neg w2 d w2 2 3 div mul 0 d w6 neg w2 neg d");
726  PrintStr(" w2 w6 d 0 w2 -2 3 div mul d w2 neg w6 d");
727  PrintStr(" w6 w2 neg d w2 -2 3 div mul 0 d w6 w2 d");
728  PrintStr(" w2 neg w6 neg d 0 w2 2 3 div mul d w2 w6 neg d cl s} def@");
729  PrintStr("/m45 {mp x y m w6 neg w2 d w2 2 3 div mul 0 d w6 neg w2 neg d");
730  PrintStr(" w2 w6 d 0 w2 -2 3 div mul d w2 neg w6 d");
731  PrintStr(" w6 w2 neg d w2 -2 3 div mul 0 d w6 w2 d");
732  PrintStr(" w2 neg w6 neg d 0 w2 2 3 div mul d w2 w6 neg d cl f} def@");
733  PrintStr("/m46 {mp x y w4 add m w4 neg w4 d w4 neg w4 neg d ");
734  PrintStr(" w4 w4 neg d w4 neg w4 neg d w4 w4 neg d w4 w4 d");
735  PrintStr(" w4 w4 neg d w4 w4 d w4 neg w4 d w4 w4 d w4 neg w4 d w4 neg w4 neg d cl s} def@");
736  PrintStr("/m47 {mp x y w4 add m w4 neg w4 d w4 neg w4 neg d");
737  PrintStr(" w4 w4 neg d w4 neg w4 neg d w4 w4 neg d w4 w4 d");
738  PrintStr(" w4 w4 neg d w4 w4 d w4 neg w4 d w4 w4 d w4 neg w4 d w4 neg w4 neg d cl f} def@");
739  PrintStr("/m48 {mp x y w4 add m w4 neg w4 d w4 neg w4 neg d w4 w4 neg d ");
740  PrintStr(" w4 neg w4 neg d w4 w4 neg d w4 w4 d w4 w4 neg d w4 w4 d");
741  PrintStr(" w4 neg w4 d w4 w4 d w4 neg w4 d w4 neg w4 neg d ");
742  PrintStr(" w4 w4 neg d w4 neg w4 neg d w4 neg w4 d w4 w4 d cl f} def@");
743  PrintStr("/m49 {mp x w2 sub w3 add y w2 sub w3 add m ");
744  PrintStr(" 0 w3 neg d w3 0 d 0 w3 d w3 0 d 0 w3 d w3 neg 0 d 0 w3 d w3 neg 0 d");
745  PrintStr(" 0 w3 neg d w3 neg 0 d 0 w3 neg d w3 0 d 0 w3 d w3 0 d 0 w3 neg d w3 neg 0 d cl f } def@");
746  PrintStr("/m2 {mp x y w2 sub m 0 w d x w2 sub y m w 0 d s} def@");
747  PrintStr("/m5 {mp x w2 sub y w2 sub m w w d x w2 sub y w2 add m w w neg d s} def@");
748 }
749 
750 ////////////////////////////////////////////////////////////////////////////////
751 /// Draw a Box
752 
754 {
755  static Double_t x[4], y[4];
756  Int_t ix1 = XtoPS(x1);
757  Int_t ix2 = XtoPS(x2);
758  Int_t iy1 = YtoPS(y1);
759  Int_t iy2 = YtoPS(y2);
760  Int_t fillis = fFillStyle/1000;
761  Int_t fillsi = fFillStyle%1000;
762 
763  if (fillis == 3 || fillis == 2) {
764  if (fillsi > 99) {
765  x[0] = x1; y[0] = y1;
766  x[1] = x2; y[1] = y1;
767  x[2] = x2; y[2] = y2;
768  x[3] = x1; y[3] = y2;
769  return;
770  }
771  if (fillsi > 0 && fillsi < 26) {
772  x[0] = x1; y[0] = y1;
773  x[1] = x2; y[1] = y1;
774  x[2] = x2; y[2] = y2;
775  x[3] = x1; y[3] = y2;
776  DrawPS(-4, &x[0], &y[0]);
777  }
778  if (fillsi == -3) {
779  SetColor(5);
780  WriteInteger(ix2 - ix1);
781  WriteInteger(iy2 - iy1);
782  WriteInteger(ix1);
783  WriteInteger(iy1);
784  PrintFast(3," bf");
785  }
786  }
787  if (fillis == 1) {
789  WriteInteger(ix2 - ix1);
790  WriteInteger(iy2 - iy1);
791  WriteInteger(ix1);
792  WriteInteger(iy1);
793  PrintFast(3," bf");
794  }
795  if (fillis == 0) {
796  if (fLineWidth<=0) return;
798  WriteInteger(ix2 - ix1);
799  WriteInteger(iy2 - iy1);
800  WriteInteger(ix1);
801  WriteInteger(iy1);
802  PrintFast(3," bl");
803  }
804 }
805 
806 ////////////////////////////////////////////////////////////////////////////////
807 /// Draw a Frame around a box
808 ///
809 /// - mode = -1 box looks as it is behind the screen
810 /// - mode = 1 box looks as it is in front of the screen
811 /// - border is the border size in already precomputed PostScript units
812 /// - dark is the color for the dark part of the frame
813 /// - light is the color for the light part of the frame
814 
816  Int_t mode, Int_t border, Int_t dark, Int_t light)
817 {
818  static Int_t xps[7], yps[7];
819  Int_t i, ixd0, iyd0, idx, idy, ixdi, iydi, ix, iy;
820 
821  // Draw top&left part of the box
822  if (mode == -1) SetColor(dark);
823  else SetColor(light);
824  Int_t bordPS = 4*border;
825  xps[0] = XtoPS(xl); yps[0] = YtoPS(yl);
826  xps[1] = xps[0] + bordPS; yps[1] = yps[0] + bordPS;
827  xps[2] = xps[1]; yps[2] = YtoPS(yt) - bordPS;
828  xps[3] = XtoPS(xt) - bordPS; yps[3] = yps[2];
829  xps[4] = XtoPS(xt); yps[4] = YtoPS(yt);
830  xps[5] = xps[0]; yps[5] = yps[4];
831  xps[6] = xps[0]; yps[6] = yps[0];
832 
833  ixd0 = xps[0];
834  iyd0 = yps[0];
835  WriteInteger(ixd0);
836  WriteInteger(iyd0);
837 
838  PrintFast(2," m");
839  idx = 0;
840  idy = 0;
841  for (i=1;i<7;i++) {
842  ixdi = xps[i];
843  iydi = yps[i];
844  ix = ixdi - ixd0;
845  iy = iydi - iyd0;
846  ixd0 = ixdi;
847  iyd0 = iydi;
848  if( ix && iy) {
849  if( idx ) { MovePS(idx,0); idx = 0; }
850  if( idy ) { MovePS(0,idy); idy = 0; }
851  MovePS(ix,iy);
852  continue;
853  }
854  if ( ix ) {
855  if( idy ) { MovePS(0,idy); idy = 0; }
856  if( !idx ) { idx = ix; continue;}
857  if( ix*idx > 0 ) idx += ix;
858  else { MovePS(idx,0); idx = ix; }
859  continue;
860  }
861  if( iy ) {
862  if( idx ) { MovePS(idx,0); idx = 0; }
863  if( !idy) { idy = iy; continue;}
864  if( iy*idy > 0 ) idy += iy;
865  else { MovePS(0,idy); idy = iy; }
866  }
867  }
868  if( idx ) MovePS(idx,0);
869  if( idy ) MovePS(0,idy);
870  PrintFast(2," f");
871 
872  // Draw bottom&right part of the box
873  if (mode == -1) SetColor(light);
874  else SetColor(dark);
875  xps[0] = XtoPS(xl); yps[0] = YtoPS(yl);
876  xps[1] = xps[0] + bordPS; yps[1] = yps[0] + bordPS;
877  xps[2] = XtoPS(xt) - bordPS; yps[2] = yps[1];
878  xps[3] = xps[2]; yps[3] = YtoPS(yt) - bordPS;
879  xps[4] = XtoPS(xt); yps[4] = YtoPS(yt);
880  xps[5] = xps[4]; yps[5] = yps[0];
881  xps[6] = xps[0]; yps[6] = yps[0];
882 
883  ixd0 = xps[0];
884  iyd0 = yps[0];
885  WriteInteger(ixd0);
886  WriteInteger(iyd0);
887 
888  PrintFast(2," m");
889  idx = 0;
890  idy = 0;
891  for (i=1;i<7;i++) {
892  ixdi = xps[i];
893  iydi = yps[i];
894  ix = ixdi - ixd0;
895  iy = iydi - iyd0;
896  ixd0 = ixdi;
897  iyd0 = iydi;
898  if( ix && iy) {
899  if( idx ) { MovePS(idx,0); idx = 0; }
900  if( idy ) { MovePS(0,idy); idy = 0; }
901  MovePS(ix,iy);
902  continue;
903  }
904  if ( ix ) {
905  if( idy ) { MovePS(0,idy); idy = 0; }
906  if( !idx ) { idx = ix; continue;}
907  if( ix*idx > 0 ) idx += ix;
908  else { MovePS(idx,0); idx = ix; }
909  continue;
910  }
911  if( iy ) {
912  if( idx ) { MovePS(idx,0); idx = 0; }
913  if( !idy) { idy = iy; continue;}
914  if( iy*idy > 0 ) idy += iy;
915  else { MovePS(0,idy); idy = iy; }
916  }
917  }
918  if( idx ) MovePS(idx,0);
919  if( idy ) MovePS(0,idy);
920  PrintFast(2," f");
921 }
922 
923 ////////////////////////////////////////////////////////////////////////////////
924 /// Draw a PolyLine
925 ///
926 /// Draw a polyline through the points xy.
927 /// - If nn=1 moves only to point x,y.
928 /// - If nn=0 the x,y are written in the PostScript file
929 /// according to the current transformation.
930 /// - If nn>0 the line is clipped as a line.
931 /// - If nn<0 the line is clipped as a fill area.
932 
934 {
935  Int_t i, n, ixd0, iyd0, idx, idy, ixdi, iydi, ix, iy;
936  Style_t linestylesav = fLineStyle;
937  Width_t linewidthsav = fLineWidth;
938  if (nn > 0) {
939  if (fLineWidth<=0) return;
940  n = nn;
944  } else {
945  n = -nn;
946  SetLineStyle(1);
947  SetLineWidth(1);
949  }
950 
951  ixd0 = XtoPS(xy[0].GetX());
952  iyd0 = YtoPS(xy[0].GetY());
953  WriteInteger(ixd0);
954  WriteInteger(iyd0);
955  if( n <= 1) {
956  if( n == 0) goto END;
957  PrintFast(2," m");
958  goto END;
959  }
960 
961  PrintFast(2," m");
962  idx = 0;
963  idy = 0;
964  for (i=1;i<n;i++) {
965  ixdi = XtoPS(xy[i].GetX());
966  iydi = YtoPS(xy[i].GetY());
967  ix = ixdi - ixd0;
968  iy = iydi - iyd0;
969  ixd0 = ixdi;
970  iyd0 = iydi;
971  if( ix && iy) {
972  if( idx ) { MovePS(idx,0); idx = 0; }
973  if( idy ) { MovePS(0,idy); idy = 0; }
974  MovePS(ix,iy);
975  continue;
976  }
977  if ( ix ) {
978  if( idy ) { MovePS(0,idy); idy = 0; }
979  if( !idx ) { idx = ix; continue;}
980  if( ix*idx > 0 ) idx += ix;
981  else { MovePS(idx,0); idx = ix; }
982  continue;
983  }
984  if( iy ) {
985  if( idx ) { MovePS(idx,0); idx = 0; }
986  if( !idy) { idy = iy; continue;}
987  if( iy*idy > 0 ) idy += iy;
988  else { MovePS(0,idy); idy = iy; }
989  }
990  }
991  if( idx ) MovePS(idx,0);
992  if( idy ) MovePS(0,idy);
993 
994  if (nn > 0 ) {
995  if (xy[0].GetX() == xy[n-1].GetX() && xy[0].GetY() == xy[n-1].GetY()) PrintFast(3," cl");
996  PrintFast(2," s");
997  } else {
998  PrintFast(2," f");
999  }
1000 END:
1001  if (nn < 0) {
1002  SetLineStyle(linestylesav);
1003  SetLineWidth(linewidthsav);
1004  }
1005 }
1006 
1007 ////////////////////////////////////////////////////////////////////////////////
1008 /// Draw a PolyLine in NDC space
1009 ///
1010 /// Draw a polyline through the points xy.
1011 /// - If nn=1 moves only to point x,y.
1012 /// - If nn=0 the x,y are written in the PostScript file
1013 /// according to the current transformation.
1014 /// - If nn>0 the line is clipped as a line.
1015 /// - If nn<0 the line is clipped as a fill area.
1016 
1018 {
1019  Int_t i, n, ixd0, iyd0, idx, idy, ixdi, iydi, ix, iy;
1020  Style_t linestylesav = fLineStyle;
1021  Width_t linewidthsav = fLineWidth;
1022  if (nn > 0) {
1023  if (fLineWidth<=0) return;
1024  n = nn;
1028  } else {
1029  n = -nn;
1030  SetLineStyle(1);
1031  SetLineWidth(1);
1033  }
1034 
1035  ixd0 = UtoPS(xy[0].GetX());
1036  iyd0 = VtoPS(xy[0].GetY());
1037  WriteInteger(ixd0);
1038  WriteInteger(iyd0);
1039  if( n <= 1) {
1040  if( n == 0) goto END;
1041  PrintFast(2," m");
1042  goto END;
1043  }
1044 
1045  PrintFast(2," m");
1046  idx = 0;
1047  idy = 0;
1048  for (i=1;i<n;i++) {
1049  ixdi = UtoPS(xy[i].GetX());
1050  iydi = VtoPS(xy[i].GetY());
1051  ix = ixdi - ixd0;
1052  iy = iydi - iyd0;
1053  ixd0 = ixdi;
1054  iyd0 = iydi;
1055  if( ix && iy) {
1056  if( idx ) { MovePS(idx,0); idx = 0; }
1057  if( idy ) { MovePS(0,idy); idy = 0; }
1058  MovePS(ix,iy);
1059  continue;
1060  }
1061  if ( ix ) {
1062  if( idy ) { MovePS(0,idy); idy = 0; }
1063  if( !idx ) { idx = ix; continue;}
1064  if( ix*idx > 0 ) idx += ix;
1065  else { MovePS(idx,0); idx = ix; }
1066  continue;
1067  }
1068  if( iy ) {
1069  if( idx ) { MovePS(idx,0); idx = 0; }
1070  if( !idy) { idy = iy; continue;}
1071  if( iy*idy > 0 ) idy += iy;
1072  else { MovePS(0,idy); idy = iy; }
1073  }
1074  }
1075  if( idx ) MovePS(idx,0);
1076  if( idy ) MovePS(0,idy);
1077 
1078  if (nn > 0 ) {
1079  if (xy[0].GetX() == xy[n-1].GetX() && xy[0].GetY() == xy[n-1].GetY()) PrintFast(3," cl");
1080  PrintFast(2," s");
1081  } else {
1082  PrintFast(2," f");
1083  }
1084 END:
1085  if (nn < 0) {
1086  SetLineStyle(linestylesav);
1087  SetLineWidth(linewidthsav);
1088  }
1089 }
1090 
1091 ////////////////////////////////////////////////////////////////////////////////
1092 /// Draw markers at the n WC points x, y
1093 
1095 {
1096  Int_t i, np, markerstyle;
1097  Float_t markersize;
1098  static char chtemp[10];
1099 
1100  if (!fMarkerSize) return;
1101  Style_t linestylesav = fLineStyle;
1102  Width_t linewidthsav = fLineWidth;
1103  SetLineStyle(1);
1104  SetLineWidth(1);
1106  markerstyle = abs(fMarkerStyle);
1107  if (markerstyle <= 0) strlcpy(chtemp, " m20",10);
1108  if (markerstyle == 1) strlcpy(chtemp, " m20",10);
1109  if (markerstyle == 2) strlcpy(chtemp, " m2",10);
1110  if (markerstyle == 3) strlcpy(chtemp, " m31",10);
1111  if (markerstyle == 4) strlcpy(chtemp, " m24",10);
1112  if (markerstyle == 5) strlcpy(chtemp, " m5",10);
1113  if (markerstyle >= 6 && markerstyle <= 19) strlcpy(chtemp, " m20",10);
1114  if (markerstyle >= 20 && markerstyle <= 49 ) snprintf(chtemp,10," m%d", markerstyle);
1115  if (markerstyle >= 50) strlcpy(chtemp, " m20",10);
1116 
1117  // Set the PostScript marker size
1118  if (markerstyle == 1) {
1119  markersize = 2.;
1120  } else if (markerstyle == 6) {
1121  markersize = 4.;
1122  } else if (markerstyle == 7) {
1123  markersize = 8.;
1124  } else {
1125  Float_t symbolsize = fMarkerSize;
1126  const Int_t kBASEMARKER = 8;
1127  Float_t sbase = symbolsize*kBASEMARKER;
1128  Float_t s2x = sbase / Float_t(gPad->GetWw() * gPad->GetAbsWNDC());
1129  markersize = this->UtoPS(s2x) - this->UtoPS(0);
1130  }
1131 
1132  if (fMarkerSizeCur != markersize) {
1133  fMarkerSizeCur = markersize;
1134  PrintFast(3," /w");
1135  WriteInteger(Int_t(markersize+0.5));
1136  PrintFast(40," def /w2 {w 2 div} def /w3 {w 3 div} def");
1137  }
1138 
1139  WriteInteger(XtoPS(x[0]));
1140  WriteInteger(YtoPS(y[0]));
1141  if (n == 1) {
1142  PrintStr(chtemp);
1143  SetLineStyle(linestylesav);
1144  SetLineWidth(linewidthsav);
1145  return;
1146  }
1147  np = 1;
1148  for (i=1;i<n;i++) {
1149  WriteInteger(XtoPS(x[i]));
1150  WriteInteger(YtoPS(y[i]));
1151  np++;
1152  if (np == 100 || i == n-1) {
1153  WriteInteger(np);
1154  PrintFast(2," {");
1155  PrintStr(chtemp);
1156  PrintFast(3,"} R");
1157  np = 0;
1158  }
1159  }
1160  SetLineStyle(linestylesav);
1161  SetLineWidth(linewidthsav);
1162 }
1163 
1164 ////////////////////////////////////////////////////////////////////////////////
1165 /// Draw markers at the n WC points x, y
1166 
1168 {
1169  Int_t i, np, markerstyle;
1170  Float_t markersize;
1171  static char chtemp[10];
1172 
1173  if (!fMarkerSize) return;
1174  Style_t linestylesav = fLineStyle;
1175  Width_t linewidthsav = fLineWidth;
1176  SetLineStyle(1);
1177  SetLineWidth(1);
1179  markerstyle = abs(fMarkerStyle);
1180  if (markerstyle <= 0) strlcpy(chtemp, " m20",10);
1181  if (markerstyle == 1) strlcpy(chtemp, " m20",10);
1182  if (markerstyle == 2) strlcpy(chtemp, " m2",10);
1183  if (markerstyle == 3) strlcpy(chtemp, " m31",10);
1184  if (markerstyle == 4) strlcpy(chtemp, " m24",10);
1185  if (markerstyle == 5) strlcpy(chtemp, " m5",10);
1186  if (markerstyle >= 6 && markerstyle <= 19) strlcpy(chtemp, " m20",10);
1187  if (markerstyle >= 20 && markerstyle <= 49 ) snprintf(chtemp,10," m%d", markerstyle);
1188  if (markerstyle >= 50) strlcpy(chtemp, " m20",10);
1189 
1190  // Set the PostScript marker size
1191  if (markerstyle == 1) {
1192  markersize = 2.;
1193  } else if (markerstyle == 6) {
1194  markersize = 4.;
1195  } else if (markerstyle == 7) {
1196  markersize = 8.;
1197  } else {
1198  Float_t symbolsize = fMarkerSize;
1199  const Int_t kBASEMARKER = 8;
1200  Float_t sbase = symbolsize*kBASEMARKER;
1201  Float_t s2x = sbase / Float_t(gPad->GetWw() * gPad->GetAbsWNDC());
1202  markersize = this->UtoPS(s2x) - this->UtoPS(0);
1203  }
1204 
1205  if (fMarkerSizeCur != markersize) {
1206  fMarkerSizeCur = markersize;
1207  PrintFast(3," /w");
1208  WriteInteger(Int_t(markersize+0.5));
1209  PrintFast(40," def /w2 {w 2 div} def /w3 {w 3 div} def");
1210  }
1211 
1212  WriteInteger(XtoPS(x[0]));
1213  WriteInteger(YtoPS(y[0]));
1214  if (n == 1) {
1215  PrintStr(chtemp);
1216  SetLineStyle(linestylesav);
1217  SetLineWidth(linewidthsav);
1218  return;
1219  }
1220  np = 1;
1221  for (i=1;i<n;i++) {
1222  WriteInteger(XtoPS(x[i]));
1223  WriteInteger(YtoPS(y[i]));
1224  np++;
1225  if (np == 100 || i == n-1) {
1226  WriteInteger(np);
1227  PrintFast(2," {");
1228  PrintStr(chtemp);
1229  PrintFast(3,"} R");
1230  np = 0;
1231  }
1232  }
1233  SetLineStyle(linestylesav);
1234  SetLineWidth(linewidthsav);
1235 }
1236 
1237 ////////////////////////////////////////////////////////////////////////////////
1238 /// Draw a PolyLine
1239 ///
1240 /// Draw a polyline through the points xw,yw.
1241 /// - If nn=1 moves only to point xw,yw.
1242 /// - If nn=0 the XW(1) and YW(1) are written in the PostScript file
1243 /// according to the current NT.
1244 /// - If nn>0 the line is clipped as a line.
1245 /// - If nn<0 the line is clipped as a fill area.
1246 
1248 {
1249  static Float_t dyhatch[24] = {.0075,.0075,.0075,.0075,.0075,.0075,.0075,.0075,
1250  .01 ,.01 ,.01 ,.01 ,.01 ,.01 ,.01 ,.01 ,
1251  .015 ,.015 ,.015 ,.015 ,.015 ,.015 ,.015 ,.015};
1252  static Float_t anglehatch[24] = {180, 90,135, 45,150, 30,120, 60,
1253  180, 90,135, 45,150, 30,120, 60,
1254  180, 90,135, 45,150, 30,120, 60};
1255  Int_t i, n, ixd0, iyd0, idx, idy, ixdi, iydi, ix, iy, fais, fasi;
1256  fais = fasi = n = 0;
1257  Int_t jxd0 = XtoPS(xw[0]);
1258  Int_t jyd0 = YtoPS(yw[0]);
1259  Style_t linestylesav = fLineStyle;
1260  Width_t linewidthsav = fLineWidth;
1261 
1262  if (nn > 0) {
1263  if (fLineWidth<=0) return;
1264  n = nn;
1268  }
1269  if (nn < 0) {
1270  n = -nn;
1271  SetLineStyle(1);
1272  SetLineWidth(1);
1274  fais = fFillStyle/1000;
1275  fasi = fFillStyle%1000;
1276  if (fais == 3 || fais == 2) {
1277  if (fasi > 100 && fasi <125) {
1278  DrawHatch(dyhatch[fasi-101],anglehatch[fasi-101], n, xw, yw);
1279  goto END;
1280  }
1281  if (fasi > 0 && fasi < 26) {
1283  }
1284  }
1285  }
1286 
1287  ixd0 = jxd0;
1288  iyd0 = jyd0;
1289  WriteInteger(ixd0);
1290  WriteInteger(iyd0);
1291  if( n <= 1) {
1292  if( n == 0) goto END;
1293  PrintFast(2," m");
1294  goto END;
1295  }
1296 
1297  PrintFast(2," m");
1298  idx = idy = 0;
1299  for (i=1;i<n;i++) {
1300  ixdi = XtoPS(xw[i]);
1301  iydi = YtoPS(yw[i]);
1302  ix = ixdi - ixd0;
1303  iy = iydi - iyd0;
1304  ixd0 = ixdi;
1305  iyd0 = iydi;
1306  if( ix && iy) {
1307  if( idx ) { MovePS(idx,0); idx = 0; }
1308  if( idy ) { MovePS(0,idy); idy = 0; }
1309  MovePS(ix,iy);
1310  } else if ( ix ) {
1311  if( idy ) { MovePS(0,idy); idy = 0;}
1312  if( !idx ) { idx = ix;}
1313  else if( TMath::Sign(ix,idx) == ix ) idx += ix;
1314  else { MovePS(idx,0); idx = ix;}
1315  } else if( iy ) {
1316  if( idx ) { MovePS(idx,0); idx = 0;}
1317  if( !idy) { idy = iy;}
1318  else if( TMath::Sign(iy,idy) == iy) idy += iy;
1319  else { MovePS(0,idy); idy = iy;}
1320  }
1321  }
1322  if (idx) MovePS(idx,0);
1323  if (idy) MovePS(0,idy);
1324 
1325  if (nn > 0 ) {
1326  if (xw[0] == xw[n-1] && yw[0] == yw[n-1]) PrintFast(3," cl");
1327  PrintFast(2," s");
1328  } else {
1329  if (fais == 0) {PrintFast(5," cl s"); goto END;}
1330  if (fais == 3 || fais == 2) {
1331  if (fasi > 0 && fasi < 26) {
1332  PrintFast(3," FA");
1333  fRed = -1;
1334  fGreen = -1;
1335  fBlue = -1;
1336  }
1337  goto END;
1338  }
1339  PrintFast(2," f");
1340  }
1341 END:
1342  if (nn < 0) {
1343  SetLineStyle(linestylesav);
1344  SetLineWidth(linewidthsav);
1345  }
1346 }
1347 
1348 ////////////////////////////////////////////////////////////////////////////////
1349 /// Draw a PolyLine
1350 ///
1351 /// Draw a polyline through the points xw,yw.
1352 /// - If nn=1 moves only to point xw,yw.
1353 /// - If nn=0 the xw(1) and YW(1) are written in the PostScript file
1354 /// --- according to the current NT.
1355 /// - If nn>0 the line is clipped as a line.
1356 /// - If nn<0 the line is clipped as a fill area.
1357 
1359 {
1360  static Float_t dyhatch[24] = {.0075,.0075,.0075,.0075,.0075,.0075,.0075,.0075,
1361  .01 ,.01 ,.01 ,.01 ,.01 ,.01 ,.01 ,.01 ,
1362  .015 ,.015 ,.015 ,.015 ,.015 ,.015 ,.015 ,.015};
1363  static Float_t anglehatch[24] = {180, 90,135, 45,150, 30,120, 60,
1364  180, 90,135, 45,150, 30,120, 60,
1365  180, 90,135, 45,150, 30,120, 60};
1366  Int_t i, n, ixd0, iyd0, idx, idy, ixdi, iydi, ix, iy, fais, fasi;
1367  fais = fasi = n = 0;
1368  Int_t jxd0 = XtoPS(xw[0]);
1369  Int_t jyd0 = YtoPS(yw[0]);
1370  Style_t linestylesav = fLineStyle;
1371  Width_t linewidthsav = fLineWidth;
1372 
1373  if (nn > 0) {
1374  if (fLineWidth<=0) return;
1375  n = nn;
1379  }
1380  if (nn < 0) {
1381  n = -nn;
1382  SetLineStyle(1);
1383  SetLineWidth(1);
1385  fais = fFillStyle/1000;
1386  fasi = fFillStyle%1000;
1387  if (fais == 3 || fais == 2) {
1388  if (fasi > 100 && fasi <125) {
1389  DrawHatch(dyhatch[fasi-101],anglehatch[fasi-101], n, xw, yw);
1390  goto END;
1391  }
1392  if (fasi > 0 && fasi < 26) {
1394  }
1395  }
1396  }
1397 
1398  ixd0 = jxd0;
1399  iyd0 = jyd0;
1400  WriteInteger(ixd0);
1401  WriteInteger(iyd0);
1402  if( n <= 1) {
1403  if( n == 0) goto END;
1404  PrintFast(2," m");
1405  goto END;
1406  }
1407 
1408  PrintFast(2," m");
1409  idx = idy = 0;
1410  for (i=1;i<n;i++) {
1411  ixdi = XtoPS(xw[i]);
1412  iydi = YtoPS(yw[i]);
1413  ix = ixdi - ixd0;
1414  iy = iydi - iyd0;
1415  ixd0 = ixdi;
1416  iyd0 = iydi;
1417  if( ix && iy) {
1418  if( idx ) { MovePS(idx,0); idx = 0; }
1419  if( idy ) { MovePS(0,idy); idy = 0; }
1420  MovePS(ix,iy);
1421  } else if ( ix ) {
1422  if( idy ) { MovePS(0,idy); idy = 0;}
1423  if( !idx ) { idx = ix;}
1424  else if( TMath::Sign(ix,idx) == ix ) idx += ix;
1425  else { MovePS(idx,0); idx = ix;}
1426  } else if( iy ) {
1427  if( idx ) { MovePS(idx,0); idx = 0;}
1428  if( !idy) { idy = iy;}
1429  else if( TMath::Sign(iy,idy) == iy) idy += iy;
1430  else { MovePS(0,idy); idy = iy;}
1431  }
1432  }
1433  if (idx) MovePS(idx,0);
1434  if (idy) MovePS(0,idy);
1435 
1436  if (nn > 0 ) {
1437  if (xw[0] == xw[n-1] && yw[0] == yw[n-1]) PrintFast(3," cl");
1438  PrintFast(2," s");
1439  } else {
1440  if (fais == 0) {PrintFast(5," cl s"); goto END;}
1441  if (fais == 3 || fais == 2) {
1442  if (fasi > 0 && fasi < 26) {
1443  PrintFast(3," FA");
1444  fRed = -1;
1445  fGreen = -1;
1446  fBlue = -1;
1447  }
1448  goto END;
1449  }
1450  PrintFast(2," f");
1451  }
1452 END:
1453  if (nn < 0) {
1454  SetLineStyle(linestylesav);
1455  SetLineWidth(linewidthsav);
1456  }
1457 }
1458 
1459 ////////////////////////////////////////////////////////////////////////////////
1460 /// Draw Fill area with hatch styles
1461 
1463 {
1464  Warning("DrawHatch", "hatch fill style not yet implemented");
1465 }
1466 
1467 ////////////////////////////////////////////////////////////////////////////////
1468 /// Draw Fill area with hatch styles
1469 
1471 {
1472  Warning("DrawHatch", "hatch fill style not yet implemented");
1473 }
1474 
1475 ////////////////////////////////////////////////////////////////////////////////
1476 
1478 {
1479  std::ifstream font_file(filename, std::ios::binary);
1480 
1481  // We cannot read directly using iostream iterators due to
1482  // signedness
1483  font_file.seekg(0, std::ios::end);
1484 
1485  const size_t font_file_length = font_file.tellg();
1486 
1487  font_file.seekg(0, std::ios::beg);
1488 
1489  std::vector<unsigned char> font_data(font_file_length, '\0');
1490 
1491  font_file.read(reinterpret_cast<char *>(&font_data[0]),
1492  font_file_length);
1493 
1494  std::string font_name;
1495  std::string postscript_string =
1496  mathtext::font_embed_postscript_t::font_embed_type_1(
1497  font_name, font_data);
1498 
1499  if (!postscript_string.empty()) {
1500  PrintRaw(postscript_string.size(), postscript_string.data());
1501  PrintStr("@");
1502 
1503  return true;
1504  }
1505 
1506  return false;
1507 }
1508 
1509 ////////////////////////////////////////////////////////////////////////////////
1510 
1512 {
1513  std::ifstream font_file(filename, std::ios::binary);
1514 
1515  // We cannot read directly using iostream iterators due to
1516  // signedness
1517  font_file.seekg(0, std::ios::end);
1518 
1519  const size_t font_file_length = font_file.tellg();
1520 
1521  font_file.seekg(0, std::ios::beg);
1522 
1523  std::vector<unsigned char> font_data(font_file_length, '\0');
1524 
1525  font_file.read(reinterpret_cast<char *>(&font_data[0]), font_file_length);
1526 
1527  std::string font_name;
1528  std::string postscript_string =
1529  mathtext::font_embed_postscript_t::font_embed_type_2(font_name, font_data);
1530 
1531  if (!postscript_string.empty()) {
1532  PrintRaw(postscript_string.size(), postscript_string.data());
1533  PrintStr("@");
1534 
1535  return true;
1536  }
1537 
1538  return false;
1539 }
1540 
1541 ////////////////////////////////////////////////////////////////////////////////
1542 
1544 {
1545  std::ifstream font_file(filename, std::ios::binary);
1546 
1547  // We cannot read directly using iostream iterators due to signedness
1548 
1549  font_file.seekg(0, std::ios::end);
1550 
1551  const size_t font_file_length = font_file.tellg();
1552 
1553  font_file.seekg(0, std::ios::beg);
1554 
1555  std::vector<unsigned char> font_data(font_file_length, '\0');
1556 
1557  font_file.read(reinterpret_cast<char *>(&font_data[0]), font_file_length);
1558 
1559  std::string font_name;
1560  std::string postscript_string =
1561  mathtext::font_embed_postscript_t::font_embed_type_42(font_name, font_data);
1562 
1563  if (!postscript_string.empty()) {
1564  PrintRaw(postscript_string.size(), postscript_string.data());
1565  PrintStr("@");
1566 
1567  return true;
1568  }
1569  fprintf(stderr, "%s:%d:\n", __FILE__, __LINE__);
1570 
1571  return false;
1572 }
1573 
1574 ////////////////////////////////////////////////////////////////////////////////
1575 /// Embed font in PS file.
1576 
1578 {
1579  PrintStr("%%IncludeResource: ProcSet (FontSetInit)@");
1580 
1581  for (Int_t fontid = 1; fontid < 30; fontid++) {
1582  if (fontid != 15 && MustEmbed[fontid-1]) {
1583 
1584  char *ttfont = NULL;
1585 
1586  FcPattern *pat, *match;
1587  FcCharSet *set = NULL;
1588  FcResult result;
1589 
1590  pat = FcPatternCreate ();
1591 
1592  switch (fontid) {
1593  case 1:
1594  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"freeserif");
1595  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_REGULAR);
1596  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ITALIC);
1597  break;
1598  case 2:
1599  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"freeserif");
1600  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_BOLD);
1601  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ROMAN);
1602  break;
1603  case 3:
1604  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"freeserif");
1605  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_BOLD);
1606  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ITALIC);
1607  break;
1608  case 4:
1609  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"freesans");
1610  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_REGULAR);
1611  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ROMAN);
1612  break;
1613  case 5:
1614  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"freesans");
1615  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_REGULAR);
1616  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ITALIC);
1617  break;
1618  case 6:
1619  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"freesans");
1620  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_BOLD);
1621  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ROMAN);
1622  break;
1623  case 7:
1624  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"freesans");
1625  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_BOLD);
1626  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ITALIC);
1627  break;
1628  case 8:
1629  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"freemono");
1630  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_REGULAR);
1631  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ROMAN);
1632  break;
1633  case 9:
1634  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"freemono");
1635  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_REGULAR);
1636  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ITALIC);
1637  break;
1638  case 10:
1639  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"freemono");
1640  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_BOLD);
1641  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ROMAN);
1642  break;
1643  case 11:
1644  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"freemono");
1645  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_BOLD);
1646  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ITALIC);
1647  break;
1648  case 12:
1649  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"symbol");
1650  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_REGULAR);
1651  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ROMAN);
1652  break;
1653  case 13:
1654  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"freeserif");
1655  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_REGULAR);
1656  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ROMAN);
1657  break;
1658  case 14:
1659  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"dingbats");
1660  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_REGULAR);
1661  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ROMAN);
1662  break;
1663  case 16:
1664  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"stixgeneral");
1665  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_REGULAR);
1666  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ROMAN);
1667  break;
1668  case 17:
1669  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"stixgeneral");
1670  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_REGULAR);
1671  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ITALIC);
1672  break;
1673  case 18:
1674  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"stixgeneral");
1675  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_BOLD);
1676  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ROMAN);
1677  break;
1678  case 19:
1679  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"stixgeneral");
1680  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_BOLD);
1681  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ITALIC);
1682  break;
1683  case 20:
1684  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"stixsize1");
1685  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_REGULAR);
1686  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ROMAN);
1687  break;
1688  case 21:
1689  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"stixsize1");
1690  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_BOLD);
1691  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ROMAN);
1692  break;
1693  case 22:
1694  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"stixsize2");
1695  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_REGULAR);
1696  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ROMAN);
1697  break;
1698  case 23:
1699  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"stixsize2");
1700  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_BOLD);
1701  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ROMAN);
1702  break;
1703  case 24:
1704  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"stixsize3");
1705  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_REGULAR);
1706  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ROMAN);
1707  break;
1708  case 25:
1709  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"stixsize3");
1710  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_BOLD);
1711  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ROMAN);
1712  break;
1713  case 26:
1714  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"stixsize4");
1715  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_REGULAR);
1716  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ROMAN);
1717  break;
1718  case 27:
1719  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"stixsize4");
1720  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_BOLD);
1721  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ROMAN);
1722  break;
1723  case 28:
1724  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"stixsize5");
1725  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_REGULAR);
1726  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ROMAN);
1727  break;
1728  case 29:
1729  FcPatternAddString (pat, FC_FAMILY, (const FcChar8*)"droidsansfallback");
1730  FcPatternAddInteger (pat, FC_WEIGHT, FC_WEIGHT_REGULAR);
1731  FcPatternAddInteger (pat, FC_SLANT, FC_SLANT_ROMAN);
1732  set = FcCharSetCreate ();
1733  FcCharSetAddChar (set, 0x0410); // Cyrillic
1734  FcCharSetAddChar (set, 0x4e00); // CJK
1735  FcPatternAddCharSet (pat, FC_CHARSET, set);
1736  break;
1737  }
1738 
1739  FcConfigSubstitute (NULL, pat, FcMatchPattern);
1740  FcDefaultSubstitute (pat);
1741  match = FcFontMatch (NULL, pat, &result);
1742  char *ttfnt;
1743  FcPatternGetString (match, FC_FILE, 0, (FcChar8**)&ttfnt);
1744  ttfont = StrDup(ttfnt);
1745  FcPatternDestroy (match);
1746  FcPatternDestroy (pat);
1747  if (set) FcCharSetDestroy (set);
1748 
1749  if (!ttfont) {
1750  Error("TPostScript::FontEmbed", "font %d not found in path",
1751  fontid);
1752  } else {
1753  if (FontEmbedType2(ttfont)) {
1754  // nothing
1755  } else if(FontEmbedType1(ttfont)) {
1756  // nothing
1757  } else if(FontEmbedType42(ttfont)) {
1758  // nothing
1759  } else {
1760  Error("TPostScript::FontEmbed", "failed to embed font %d)",
1761  fontid);
1762  }
1763  delete [] ttfont;
1764  }
1765  }
1766  }
1767  PrintStr("%%IncludeResource: font Times-Roman@");
1768  PrintStr("%%IncludeResource: font Times-Italic@");
1769  PrintStr("%%IncludeResource: font Times-Bold@");
1770  PrintStr("%%IncludeResource: font Times-BoldItalic@");
1771  PrintStr("%%IncludeResource: font Helvetica@");
1772  PrintStr("%%IncludeResource: font Helvetica-Oblique@");
1773  PrintStr("%%IncludeResource: font Helvetica-Bold@");
1774  PrintStr("%%IncludeResource: font Helvetica-BoldOblique@");
1775  PrintStr("%%IncludeResource: font Courier@");
1776  PrintStr("%%IncludeResource: font Courier-Oblique@");
1777  PrintStr("%%IncludeResource: font Courier-Bold@");
1778  PrintStr("%%IncludeResource: font Courier-BoldOblique@");
1779  PrintStr("%%IncludeResource: font Symbol@");
1780  PrintStr("%%IncludeResource: font ZapfDingbats@");
1781 
1782  fFontEmbed = kTRUE;
1783 }
1784 
1785 ////////////////////////////////////////////////////////////////////////////////
1786 /// Font Re-encoding
1787 
1789 {
1790  PrintStr("/reEncode ");
1791  PrintStr("{exch findfont");
1792  PrintStr(" dup length dict begin");
1793  PrintStr(" {1 index /FID eq ");
1794  PrintStr(" {pop pop}");
1795  PrintStr(" {def} ifelse");
1796  PrintStr(" } forall");
1797  PrintStr(" /Encoding exch def");
1798  PrintStr(" currentdict end");
1799  PrintStr(" dup /FontName get exch");
1800  PrintStr(" definefont pop");
1801  PrintStr(" } def");
1802  PrintStr(" [/Times-Bold /Times-Italic /Times-BoldItalic /Helvetica");
1803  PrintStr(" /Helvetica-Oblique /Helvetica-Bold /Helvetica-BoldOblique");
1804  PrintStr(" /Courier /Courier-Oblique /Courier-Bold /Courier-BoldOblique");
1805  PrintStr(" /Times-Roman /AvantGarde-Book /AvantGarde-BookOblique");
1806  PrintStr(" /AvantGarde-Demi /AvantGarde-DemiOblique /Bookman-Demi");
1807  PrintStr(" /Bookman-DemiItalic /Bookman-Light /Bookman-LightItalic");
1808  PrintStr(" /Helvetica-Narrow /Helvetica-Narrow-Bold /Helvetica-Narrow-BoldOblique");
1809  PrintStr(" /Helvetica-Narrow-Oblique /NewCenturySchlbk-Roman /NewCenturySchlbk-Bold");
1810  PrintStr(" /NewCenturySchlbk-BoldItalic /NewCenturySchlbk-Italic");
1811  PrintStr(" /Palatino-Bold /Palatino-BoldItalic /Palatino-Italic /Palatino-Roman");
1812  PrintStr(" ] {ISOLatin1Encoding reEncode } forall");
1813 }
1814 
1815 ////////////////////////////////////////////////////////////////////////////////
1816 /// PostScript Initialisation
1817 ///
1818 /// This method initialize the following PostScript procedures:
1819 ///
1820 /// | Macro Name | Input parameters | Explanation |
1821 /// |------------|------------------|-----------------------------------|
1822 /// | l | x y | Draw a line to the x y position |
1823 /// | m | x y | Move to the position x y |
1824 /// | box | dx dy x y | Define a box |
1825 /// | bl | dx dy x y | Draw a line box |
1826 /// | bf | dx dy x y | Draw a filled box |
1827 /// | t | x y | Translate |
1828 /// | r | angle | Rotate |
1829 /// | rl | i j | Roll the stack |
1830 /// | d | x y | Draw a relative line to x y |
1831 /// | X | x | Draw a relative line to x (y=0) |
1832 /// | Y | y | Draw a relative line to y (x=0) |
1833 /// | rm | x y | Move relatively to x y |
1834 /// | gr | | Restore the graphic context |
1835 /// | lw | lwidth | Set line width to lwidth |
1836 /// | sd | [] 0 | Set dash line define by [] |
1837 /// | s | | Stroke mode |
1838 /// | c | r g b | Set rgb color to r g b |
1839 /// | cl | | Close path |
1840 /// | f | | Fill the last describe path |
1841 /// | mXX | x y | Draw the marker type XX at (x,y) |
1842 /// | Zone | ix iy | Define the current zone |
1843 /// | black | | The color is black |
1844 /// | C | dx dy x y | Clipping on |
1845 /// | NC | | Clipping off |
1846 /// | R | | repeat |
1847 /// | ita | | Used to make the symbols italic |
1848 /// | K | | kshow |
1849 
1851 {
1852  Double_t rpxmin, rpymin, width, heigth;
1853  rpxmin = rpymin = width = heigth = 0;
1854  Int_t format;
1855  fNpages=1;
1856  for (Int_t i=0;i<32;i++) fPatterns[i]=0;
1857 
1858  // Mode is last digit of PostScript Workstation type
1859  // mode=1,2 for portrait/landscape black and white
1860  // mode=3 for Encapsulated PostScript File
1861  // mode=4 for portrait colour
1862  // mode=5 for lanscape colour
1863  Int_t atype = abs(fType);
1864  fMode = atype%10;
1865  if( fMode <= 0 || fMode > 5) {
1866  Error("Initialize", "invalid file type %d", fMode);
1867  return;
1868  }
1869 
1870  // fNXzone (fNYzone) is the total number of windows in x (y)
1871  fNXzone = (atype%1000)/100;
1872  fNYzone = (atype%100)/10;
1873  if( fNXzone <= 0 ) fNXzone = 1;
1874  if( fNYzone <= 0 ) fNYzone = 1;
1875  fIXzone = 1;
1876  fIYzone = 1;
1877 
1878  // format = 0-99 is the European page format (A4,A3 ...)
1879  // format = 100 is the US format 8.5x11.0 inch
1880  // format = 200 is the US format 8.5x14.0 inch
1881  // format = 300 is the US format 11.0x17.0 inch
1882  format = atype/1000;
1883  if( format == 0 ) format = 4;
1884  if( format == 99 ) format = 0;
1885 
1886  PrintStr("%%Title: ");
1887  const char *pstitle = gStyle->GetTitlePS();
1888  if (gPad && !pstitle[0]) pstitle = gPad->GetMother()->GetTitle();
1889  if (strlen(GetName())<=80) PrintStr(GetName());
1890  if(!pstitle[0] && fMode != 3) {;
1891  PrintFast(2," (");
1892  if ( format <= 99 ) {;
1893  PrintFast(2," A");
1894  WriteInteger(format);
1895  PrintFast(1,")");
1896  }
1897  else {
1898  if ( format == 100 ) PrintFast(8," Letter)");
1899  if ( format == 200 ) PrintFast(7," Legal)");
1900  if ( format == 300 ) PrintFast(8," Ledger)");
1901  }
1902  PrintStr("@");
1903  PrintStr("%%Pages: (atend)@");
1904  }
1905  else {
1906  if (!strchr(pstitle,'\n')) {
1907  PrintFast(2,": ");
1908  PrintStr(pstitle);
1909  }
1910  PrintStr("@");
1911  }
1912 
1913  PrintFast(24,"%%Creator: ROOT Version ");
1914  PrintStr(gROOT->GetVersion());
1915  PrintStr("@");
1916  PrintFast(16,"%%CreationDate: ");
1917  TDatime t;
1918  PrintStr(t.AsString());
1919  PrintStr("@");
1920 
1921  if ( fMode == 1 || fMode == 4) PrintStr("%%Orientation: Portrait@");
1922  if ( fMode == 2 || fMode == 5) PrintStr("%%Orientation: Landscape@");
1923 
1924  PrintStr("%%EndComments@");
1925  PrintStr("%%BeginProlog@");
1926 
1927  if( fMode == 3)PrintStr("80 dict begin@");
1928 
1929  // Initialisation of PostScript procedures
1930  PrintStr("/s {stroke} def /l {lineto} def /m {moveto} def /t {translate} def@");
1931  PrintStr("/r {rotate} def /rl {roll} def /R {repeat} def@");
1932  PrintStr("/d {rlineto} def /rm {rmoveto} def /gr {grestore} def /f {eofill} def@");
1933  if (gStyle->GetColorModelPS()) {
1934  PrintStr("/c {setcmykcolor} def /black {0 0 0 1 setcmykcolor} def /sd {setdash} def@");
1935  } else {
1936  PrintStr("/c {setrgbcolor} def /black {0 setgray} def /sd {setdash} def@");
1937  }
1938  PrintStr("/cl {closepath} def /sf {scalefont setfont} def /lw {setlinewidth} def@");
1939  PrintStr("/box {m dup 0 exch d exch 0 d 0 exch neg d cl} def@");
1940  PrintStr("/NC{systemdict begin initclip end}def/C{NC box clip newpath}def@");
1941  PrintStr("/bl {box s} def /bf {gsave box gsave f grestore 1 lw [] 0 sd s grestore} def /Y { 0 exch d} def /X { 0 d} def @");
1942  PrintStr("/K {{pop pop 0 moveto} exch kshow} bind def@");
1943  PrintStr("/ita {/ang 15 def gsave [1 0 ang dup sin exch cos div 1 0 0] concat} def @");
1944 
1945  DefineMarkers();
1946 
1947  FontEncode();
1948 
1949  // mode=1 for portrait black/white
1950  if (fMode == 1) {
1951  rpxmin = 0.7;
1952  rpymin = TMath::Sqrt(2.)*rpxmin;
1953  switch (format) {
1954  case 100 :
1955  width = (8.5*2.54)-2.*rpxmin;
1956  heigth = (11.*2.54)-2.*rpymin;
1957  break;
1958  case 200 :
1959  width = (8.5*2.54)-2.*rpxmin;
1960  heigth = (14.*2.54)-2.*rpymin;
1961  break;
1962  case 300 :
1963  width = (11.*2.54)-2.*rpxmin;
1964  heigth = (17.*2.54)-2.*rpymin;
1965  break;
1966  default :
1967  width = 21.0-2.*rpxmin;
1968  heigth = 29.7-2.*rpymin;
1969  };
1970  }
1971 
1972  // mode=2 for landscape black/white
1973  if (fMode == 2) {
1974  rpymin = 0.7;
1975  rpxmin = TMath::Sqrt(2.)*rpymin;
1976  switch (format) {
1977  case 100 :
1978  width = (11.*2.54)-2.*rpxmin;
1979  heigth = (8.5*2.54)-2.*rpymin;
1980  break;
1981  case 200 :
1982  width = (14.*2.54)-2.*rpxmin;
1983  heigth = (8.5*2.54)-2.*rpymin;
1984  break;
1985  case 300 :
1986  width = (17.*2.54)-2.*rpxmin;
1987  heigth = (11.*2.54)-2.*rpymin;
1988  break;
1989  default :
1990  width = 29.7-2.*rpxmin;
1991  heigth = 21-2.*rpymin;
1992  };
1993  }
1994 
1995  // mode=3 encapsulated PostScript
1996  if (fMode == 3) {
1997  width = 20;
1998  heigth = 20;
1999  format = 4;
2000  fNXzone = 1;
2001  fNYzone = 1;
2002  }
2003 
2004  // mode=4 for portrait colour
2005  if (fMode == 4) {
2006  rpxmin = 0.7;
2007  rpymin = 3.4;
2008  switch (format) {
2009  case 100 :
2010  width = (8.5*2.54)-2.*rpxmin;
2011  heigth = (11.*2.54)-2.*rpymin;
2012  break;
2013  case 200 :
2014  width = (8.5*2.54)-2.*rpxmin;
2015  heigth = (14.*2.54)-2.*rpymin;
2016  break;
2017  case 300 :
2018  width = (11.*2.54)-2.*rpxmin;
2019  heigth = (17.*2.54)-2.*rpymin;
2020  break;
2021  default :
2022  width = (21.0-2*rpxmin);
2023  heigth = (29.7-2.*rpymin);
2024  };
2025  }
2026 
2027  // mode=5 for lanscape colour
2028  if (fMode == 5) {
2029  rpxmin = 3.4;
2030  rpymin = 0.7;
2031  switch (format) {
2032  case 100 :
2033  width = (11.*2.54)-2.*rpxmin;
2034  heigth = (8.5*2.54)-2.*rpymin;
2035  break;
2036  case 200 :
2037  width = (14.*2.54)-2.*rpxmin;
2038  heigth = (8.5*2.54)-2.*rpymin;
2039  break;
2040  case 300 :
2041  width = (17.*2.54)-2.*rpxmin;
2042  heigth = (11.*2.54)-2.*rpymin;
2043  break;
2044  default :
2045  width = (29.7-2*rpxmin);
2046  heigth = (21-2.*rpymin);
2047  };
2048  }
2049 
2050  Double_t value = 0;
2051  if (format < 100) value = 21*TMath::Power(TMath::Sqrt(2.), 4-format);
2052  else if (format == 100) value = 8.5*2.54;
2053  else if (format == 200) value = 8.5*2.54;
2054  else if (format == 300) value = 11.*2.54;
2055  if (format >= 100) format = 4;
2056 
2057  // Compute size (in points) of the window for each picture = f(fNXzone,fNYzone)
2058  Double_t sizex = width/Double_t(fNXzone)*TMath::Power(TMath::Sqrt(2.), 4-format);
2059  Double_t sizey = heigth/Double_t(fNYzone)*TMath::Power(TMath::Sqrt(2.), 4-format);
2060  Int_t npx = 4*CMtoPS(sizex);
2061  Int_t npy = 4*CMtoPS(sizey);
2062  if (sizex > sizey) fMaxsize = CMtoPS(sizex);
2063  else fMaxsize = CMtoPS(sizey);
2064 
2065  // Procedure Zone
2066  if (fMode != 3) {
2067  PrintFast(33,"/Zone {/iy exch def /ix exch def ");
2068  PrintFast(10," ix 1 sub ");
2069  WriteInteger(npx);
2070  PrintFast(5," mul ");
2072  PrintFast(8," iy sub ");
2073  WriteInteger(npy);
2074  PrintStr(" mul t} def@");
2075  } else {
2076  PrintStr("@");
2077  }
2078 
2079  PrintStr("%%EndProlog@");
2080  PrintStr("%%BeginSetup@");
2081  PrintStr("%%EndSetup@");
2082  PrintFast(8,"newpath ");
2083  SaveRestore(1);
2084  if (fMode == 1 || fMode == 4) {
2085  WriteInteger(CMtoPS(rpxmin));
2086  WriteInteger(CMtoPS(rpymin));
2087  PrintFast(2," t");
2088  }
2089  if (fMode == 2 || fMode == 5) {
2090  PrintFast(7," 90 r 0");
2091  WriteInteger(CMtoPS(-value));
2092  PrintFast(3," t ");
2093  WriteInteger(CMtoPS(rpxmin));
2094  WriteInteger(CMtoPS(rpymin));
2095  PrintFast(2," t");
2096  }
2097 
2098  PrintFast(15," .25 .25 scale ");
2099  if (fMode != 3) {
2100  SaveRestore(1);
2101  PrintStr("@");
2102  PrintStr("%%Page: 1 1@");
2103  SaveRestore(1);
2104  }
2105 
2106  //Check is user has defined a special header in the current style
2107  Int_t nh = strlen(gStyle->GetHeaderPS());
2108  if (nh) {
2109  PrintFast(nh,gStyle->GetHeaderPS());
2110  if (fMode != 3) SaveRestore(1);
2111  }
2112 }
2113 
2114 ////////////////////////////////////////////////////////////////////////////////
2115 /// Move to a new position
2116 
2118 {
2119  if (ix != 0 && iy != 0) {
2120  WriteInteger(ix);
2121  WriteInteger(iy);
2122  PrintFast(2," d");
2123  } else if (ix != 0) {
2124  WriteInteger(ix);
2125  PrintFast(2," X");
2126  } else if (iy != 0) {
2127  WriteInteger(iy);
2128  PrintFast(2," Y");
2129  }
2130 }
2131 
2132 ////////////////////////////////////////////////////////////////////////////////
2133 /// Move to a new PostScript page
2134 
2136 {
2137  // Compute pad conversion coefficients
2138  if (gPad) {
2139  // if (!gPad->GetPadPaint()) gPad->Update();
2140  Double_t ww = gPad->GetWw();
2141  Double_t wh = gPad->GetWh();
2142  fYsize = fXsize*wh/ww;
2143  } else fYsize = 27;
2144 
2145  if(fType == 113 && !fBoundingBox) {
2146  Bool_t psave = fPrinted;
2147  PrintStr("@%%BoundingBox: ");
2148  Double_t xlow=0, ylow=0, xup=1, yup=1;
2149  if (gPad) {
2150  xlow = gPad->GetAbsXlowNDC();
2151  xup = xlow + gPad->GetAbsWNDC();
2152  ylow = gPad->GetAbsYlowNDC();
2153  yup = ylow + gPad->GetAbsHNDC();
2154  }
2155  WriteInteger(CMtoPS(fXsize*xlow));
2156  WriteInteger(CMtoPS(fYsize*ylow));
2157  WriteInteger(CMtoPS(fXsize*xup));
2158  WriteInteger(CMtoPS(fYsize*yup));
2159  PrintStr("@");
2160  Initialize();
2161  fBoundingBox = kTRUE;
2162  fPrinted = psave;
2163  }
2164  if (fPrinted) {
2165  if (fSave) SaveRestore(-1);
2166  fClear = kTRUE;
2167  fPrinted = kFALSE;
2168  }
2169  Zone();
2170 }
2171 
2172 ////////////////////////////////////////////////////////////////////////////////
2173 /// Set the range for the paper in centimeters
2174 
2176 {
2177  Float_t xps=0, yps=0, xncm=0, yncm=0, dxwn=0, dywn=0, xwkwn=0, ywkwn=0, xymax=0;
2178 
2179  fXsize = xsize;
2180  fYsize = ysize;
2181  if( fType != 113) { xps = fXsize; yps = fYsize; }
2182  else { xps = xsize; yps = ysize; }
2183 
2184  if( xsize <= xps && ysize < yps) {
2185  if ( xps > yps ) xymax = xps;
2186  else xymax = yps;
2187  xncm = xsize/xymax;
2188  yncm = ysize/xymax;
2189  dxwn = ((xps/xymax)-xncm)/2;
2190  dywn = ((yps/xymax)-yncm)/2;
2191  } else {
2192  if (xps/yps < 1) xwkwn = xps/yps;
2193  else xwkwn = 1;
2194  if (yps/xps < 1) ywkwn = yps/xps;
2195  else ywkwn = 1;
2196 
2197  if (xsize < ysize) {
2198  xncm = ywkwn*xsize/ysize;
2199  yncm = ywkwn;
2200  dxwn = (xwkwn-xncm)/2;
2201  dywn = 0;
2202  if( dxwn < 0) {
2203  xncm = xwkwn;
2204  dxwn = 0;
2205  yncm = xwkwn*ysize/xsize;
2206  dywn = (ywkwn-yncm)/2;
2207  }
2208  } else {
2209  xncm = xwkwn;
2210  yncm = xwkwn*ysize/xsize;
2211  dxwn = 0;
2212  dywn = (ywkwn-yncm)/2;
2213  if( dywn < 0) {
2214  yncm = ywkwn;
2215  dywn = 0;
2216  xncm = ywkwn*xsize/ysize;
2217  dxwn = (xwkwn-xncm)/2;
2218  }
2219  }
2220  }
2221  fXVP1 = dxwn;
2222  fXVP2 = xncm+dxwn;
2223  fYVP1 = dywn;
2224  fYVP2 = yncm+dywn;
2225  fRange = kTRUE;
2226 }
2227 
2228 ////////////////////////////////////////////////////////////////////////////////
2229 /// Compute number of gsaves for restore
2230 /// This allows to write the correct number of grestore at the
2231 /// end of the PS file.
2232 
2234 {
2235  if (flag == 1) { PrintFast(7," gsave "); fSave++; }
2236  else { PrintFast(4," gr "); fSave--; }
2237 }
2238 
2239 ////////////////////////////////////////////////////////////////////////////////
2240 /// Set color index for fill areas
2241 
2243 {
2244  fFillColor = cindex;
2245  if (gStyle->GetFillColor() <= 0) cindex = 0;
2246  SetColor(Int_t(cindex));
2247 }
2248 
2249 ////////////////////////////////////////////////////////////////////////////////
2250 /// Patterns definition
2251 ///
2252 /// Define the pattern ipat in the current PS file. ipat can vary from
2253 /// 1 to 25. Together with the pattern, the color (color) in which the
2254 /// pattern has to be drawn is also required. A pattern is defined in the
2255 /// current PS file only the first time it is used. Some level 2
2256 /// Postscript functions are used, so on level 1 printers, patterns will
2257 /// not work. This is not a big problem because patterns are
2258 /// defined only if they are used, so if they are not used a PS level 1
2259 /// file will not be polluted by level 2 features, and in any case the old
2260 /// patterns used a lot of memory which made them almost unusable on old
2261 /// level 1 printers. Finally we should say that level 1 devices are
2262 /// becoming very rare. The official PostScript is now level 3 !
2263 
2265 {
2266  char cdef[28];
2267  char cpat[5];
2268  snprintf(cpat,5," P%2.2d", ipat);
2269 
2270  // fPatterns is used as an array of chars. If fPatterns[ipat] != 0 the
2271  // pattern number ipat as already be defined is this file and it
2272  // is not necessary to redefine it. fPatterns is set to zero in Initialize.
2273  // The pattern number 26 allows to know if the macro "cs" has already
2274  // been defined in the current file (see label 200).
2275  if (fPatterns[ipat] == 0) {
2276 
2277  // Define the Patterns. Line width must be 1
2278  // Setting fLineWidth to -1 will force the line width definition next time
2279  // TPostScript::SetLineWidth will be called.
2280  fLineWidth = -1;
2281  PrintFast(5," 1 lw");
2282  PrintStr(" << /PatternType 1 /PaintType 2 /TilingType 1");
2283  switch (ipat) {
2284  case 1 :
2285  PrintStr(" /BBox [ 0 0 98 4 ]");
2286  PrintStr(" /XStep 98 /YStep 4");
2287  PrintStr(" /PaintProc { begin gsave");
2288  PrintStr(" [1] 0 sd 2 4 m 99 4 l s 1 3 m 98 3 l s");
2289  PrintStr(" 2 2 m 99 2 l s 1 1 m 98 1 l s");
2290  PrintStr(" gr end } >> [ 4.0 0 0 4.0 0 0 ]");
2291  break;
2292  case 2 :
2293  PrintStr(" /BBox [ 0 0 96 4 ]");
2294  PrintStr(" /XStep 96 /YStep 4");
2295  PrintStr(" /PaintProc { begin gsave");
2296  PrintStr(" [1 3] 0 sd 2 4 m 98 4 l s 0 3 m 96 3 l s");
2297  PrintStr(" 2 2 m 98 2 l s 0 1 m 96 1 l s");
2298  PrintStr(" gr end } >> [ 3.0 0 0 3.0 0 0 ]");
2299  break;
2300  case 3 :
2301  PrintStr(" /BBox [ 0 0 96 16 ]");
2302  PrintStr(" /XStep 96 /YStep 16");
2303  PrintStr(" /PaintProc { begin gsave");
2304  PrintStr(" [1 3] 0 sd 2 13 m 98 13 l s 0 9 m 96 9 l s");
2305  PrintStr(" 2 5 m 98 5 l s 0 1 m 96 1 l s");
2306  PrintStr(" gr end } >> [ 2.0 0 0 2.0 0 0 ]");
2307  break;
2308  case 4 :
2309  PrintStr(" /BBox [ 0 0 100 100 ]");
2310  PrintStr(" /XStep 100 /YStep 100");
2311  PrintStr(" /PaintProc { begin gsave");
2312  PrintStr(" 0 0 m 100 100 l s");
2313  PrintStr(" gr end } >> [ 0.24 0 0 0.24 0 0 ]");
2314  break;
2315  case 5 :
2316  PrintStr(" /BBox [ 0 0 100 100 ]");
2317  PrintStr(" /XStep 100 /YStep 100");
2318  PrintStr(" /PaintProc { begin gsave");
2319  PrintStr(" 0 100 m 100 0 l s");
2320  PrintStr(" gr end } >> [ 0.24 0 0 0.24 0 0 ]");
2321  break;
2322  case 6 :
2323  PrintStr(" /BBox [ 0 0 100 100 ]");
2324  PrintStr(" /XStep 100 /YStep 100");
2325  PrintStr(" /PaintProc { begin gsave");
2326  PrintStr(" 50 0 m 50 100 l s");
2327  PrintStr(" gr end } >> [ 0.12 0 0 0.12 0 0 ]");
2328  break;
2329  case 7 :
2330  PrintStr(" /BBox [ 0 0 100 100 ]");
2331  PrintStr(" /XStep 100 /YStep 100");
2332  PrintStr(" /PaintProc { begin gsave");
2333  PrintStr(" 0 50 m 100 50 l s");
2334  PrintStr(" gr end } >> [ 0.12 0 0 0.12 0 0 ]");
2335  break;
2336  case 8 :
2337  PrintStr(" /BBox [ 0 0 101 101 ]");
2338  PrintStr(" /XStep 100 /YStep 100");
2339  PrintStr(" /PaintProc { begin gsave");
2340  PrintStr(" 0 0 m 0 30 l 30 0 l f 0 70 m 0 100 l 30 100 l f");
2341  PrintStr(" 70 100 m 100 100 l 100 70 l f 70 0 m 100 0 l");
2342  PrintStr(" 100 30 l f 50 20 m 20 50 l 50 80 l 80 50 l f");
2343  PrintStr(" 50 80 m 30 100 l s 20 50 m 0 30 l s 50 20 m");
2344  PrintStr(" 70 0 l s 80 50 m 100 70 l s");
2345  PrintStr(" gr end } >> [ 0.24 0 0 0.24 0 0 ]");
2346  break;
2347  case 9 :
2348  PrintStr(" /BBox [ 0 0 100 100 ]");
2349  PrintStr(" /XStep 100 /YStep 100");
2350  PrintStr(" /PaintProc { begin gsave");
2351  PrintStr(" 0 50 m 50 50 50 180 360 arc");
2352  PrintStr(" 0 50 m 0 100 50 270 360 arc");
2353  PrintStr(" 50 100 m 100 100 50 180 270 arc s");
2354  PrintStr(" gr end } >> [ 0.24 0 0 0.24 0 0 ]");
2355  break;
2356  case 10 :
2357  PrintStr(" /BBox [ 0 0 100 100 ]");
2358  PrintStr(" /XStep 100 /YStep 100");
2359  PrintStr(" /PaintProc { begin gsave");
2360  PrintStr(" 0 50 m 100 50 l 1 1 m 100 1 l");
2361  PrintStr(" 0 0 m 0 50 l 100 0 m 100 50 l");
2362  PrintStr(" 50 50 m 50 100 l s");
2363  PrintStr(" gr end } >> [ 0.24 0 0 0.24 0 0 ]");
2364  break;
2365  case 11 :
2366  PrintStr(" /BBox [ 0 0 100 100 ]");
2367  PrintStr(" /XStep 100 /YStep 100");
2368  PrintStr(" /PaintProc { begin gsave");
2369  PrintStr(" 0 0 m 0 20 l 50 0 m 50 20 l");
2370  PrintStr(" 100 0 m 100 20 l 0 80 m 0 100 l");
2371  PrintStr(" 50 80 m 50 100 l 100 80 m 100 100 l");
2372  PrintStr(" 25 30 m 25 70 l 75 30 m 75 70 l");
2373  PrintStr(" 0 100 m 20 85 l 50 100 m 30 85 l");
2374  PrintStr(" 50 100 m 70 85 l 100 100 m 80 85 l");
2375  PrintStr(" 0 0 m 20 15 l 50 0 m 30 15 l");
2376  PrintStr(" 50 0 m 70 15 l 100 0 m 80 15 l");
2377  PrintStr(" 5 35 m 45 65 l 5 65 m 45 35 l");
2378  PrintStr(" 55 35 m 95 65 l 55 65 m 95 35 l s");
2379  PrintStr(" gr end } >> [ 0.5 0 0 0.5 0 0 ]");
2380  break;
2381  case 12 :
2382  PrintStr(" /BBox [ 0 0 100 100 ]");
2383  PrintStr(" /XStep 100 /YStep 100");
2384  PrintStr(" /PaintProc { begin gsave");
2385  PrintStr(" 0 80 m 0 100 20 270 360 arc");
2386  PrintStr(" 30 100 m 50 100 20 180 360 arc");
2387  PrintStr(" 80 100 m 100 100 20 180 270 arc");
2388  PrintStr(" 20 0 m 0 0 20 0 90 arc");
2389  PrintStr(" 70 0 m 50 0 20 0 180 arc");
2390  PrintStr(" 100 20 m 100 0 20 90 180 arc");
2391  PrintStr(" 45 50 m 25 50 20 0 360 arc");
2392  PrintStr(" 95 50 m 75 50 20 0 360 arc s");
2393  PrintStr(" gr end } >> [ 0.5 0 0 0.5 0 0 ]");
2394  break;
2395  case 13 :
2396  PrintStr(" /BBox [ 0 0 100 100 ]");
2397  PrintStr(" /XStep 100 /YStep 100");
2398  PrintStr(" /PaintProc { begin gsave");
2399  PrintStr(" 0 0 m 100 100 l 0 100 m 100 0 l s");
2400  PrintStr(" gr end } >> [ 0.24 0 0 0.24 0 0 ]");
2401  break;
2402  case 14 :
2403  PrintStr(" /BBox [ 0 0 100 100 ]");
2404  PrintStr(" /XStep 80 /YStep 80");
2405  PrintStr(" /PaintProc { begin gsave");
2406  PrintStr(" 0 20 m 100 20 l 20 0 m 20 100 l");
2407  PrintStr(" 0 80 m 100 80 l 80 0 m 80 100 l");
2408  PrintStr(" 20 40 m 60 40 l 60 20 m 60 60 l");
2409  PrintStr(" 40 40 m 40 80 l 40 60 m 80 60 l s");
2410  PrintStr(" gr end } >> [ 0.60 0 0 0.60 0 0 ]");
2411  break;
2412  case 15 :
2413  PrintStr(" /BBox [ 0 0 60 60 ]");
2414  PrintStr(" /XStep 60 /YStep 60");
2415  PrintStr(" /PaintProc { begin gsave");
2416  PrintStr(" 0 55 m 0 60 5 270 360 arc");
2417  PrintStr(" 25 60 m 30 60 5 180 360 arc");
2418  PrintStr(" 55 60 m 60 60 5 180 270 arc");
2419  PrintStr(" 20 30 m 15 30 5 0 360 arc");
2420  PrintStr(" 50 30 m 45 30 5 0 360");
2421  PrintStr(" arc 5 0 m 0 0 5 0 90 arc");
2422  PrintStr(" 35 0 m 30 0 5 0 180 arc");
2423  PrintStr(" 60 5 m 60 0 5 90 180 arc s");
2424  PrintStr(" gr end } >> [ 0.41 0 0 0.41 0 0 ]");
2425  break;
2426  case 16 :
2427  PrintStr(" /BBox [ 0 0 100 100 ]");
2428  PrintStr(" /XStep 100 /YStep 100");
2429  PrintStr(" /PaintProc { begin gsave");
2430  PrintStr(" 50 50 m 25 50 25 0 180 arc s");
2431  PrintStr(" 50 50 m 75 50 25 180 360 arc s");
2432  PrintStr(" gr end } >> [ 0.4 0 0 0.2 0 0 ]");
2433  break;
2434  case 17 :
2435  PrintStr(" /BBox [ 0 0 100 100 ]");
2436  PrintStr(" /XStep 100 /YStep 100");
2437  PrintStr(" /PaintProc { begin gsave");
2438  PrintStr(" [24] 0 setdash 0 0 m 100 100 l s");
2439  PrintStr(" gr end } >> [ 0.24 0 0 0.24 0 0 ]");
2440  break;
2441  case 18 :
2442  PrintStr(" /BBox [ 0 0 100 100 ]");
2443  PrintStr(" /XStep 100 /YStep 100");
2444  PrintStr(" /PaintProc { begin gsave");
2445  PrintStr(" [24] 0 setdash 0 100 m 100 0 l s");
2446  PrintStr(" gr end } >> [ 0.24 0 0 0.24 0 0 ]");
2447  break;
2448  case 19 :
2449  PrintStr(" /BBox [ 0 0 100 100 ]");
2450  PrintStr(" /XStep 100 /YStep 100");
2451  PrintStr(" /PaintProc { begin gsave");
2452  PrintStr(" 90 50 m 50 50 40 0 360 arc");
2453  PrintStr(" 0 50 m 0 100 50 270 360 arc");
2454  PrintStr(" 50 0 m 0 0 50 0 90 arc");
2455  PrintStr(" 100 50 m 100 0 50 90 180 arc");
2456  PrintStr(" 50 100 m 100 100 50 180 270 arc s");
2457  PrintStr(" gr end } >> [ 0.47 0 0 0.47 0 0 ]");
2458  break;
2459  case 20 :
2460  PrintStr(" /BBox [ 0 0 100 100 ]");
2461  PrintStr(" /XStep 100 /YStep 100");
2462  PrintStr(" /PaintProc { begin gsave");
2463  PrintStr(" 50 50 m 50 75 25 270 450 arc s");
2464  PrintStr(" 50 50 m 50 25 25 90 270 arc s");
2465  PrintStr(" gr end } >> [ 0.2 0 0 0.4 0 0 ]");
2466  break;
2467  case 21 :
2468  PrintStr(" /BBox [ 0 0 101 101 ]");
2469  PrintStr(" /XStep 100 /YStep 100");
2470  PrintStr(" /PaintProc { begin gsave");
2471  PrintStr(" 1 1 m 25 1 l 25 25 l 50 25 l 50 50 l");
2472  PrintStr(" 75 50 l 75 75 l 100 75 l 100 100 l");
2473  PrintStr(" 50 1 m 75 1 l 75 25 l 100 25 l 100 50 l");
2474  PrintStr(" 0 50 m 25 50 l 25 75 l 50 75 l 50 100 l s");
2475  PrintStr(" gr end } >> [ 0.5 0 0 0.5 0 0 ]");
2476  break;
2477  case 22 :
2478  PrintStr(" /BBox [ 0 0 101 101 ]");
2479  PrintStr(" /XStep 100 /YStep 100");
2480  PrintStr(" /PaintProc { begin gsave");
2481  PrintStr(" 1 100 m 25 100 l 25 75 l 50 75 l 50 50 l");
2482  PrintStr(" 75 50 l 75 25 l 100 25 l 100 1 l");
2483  PrintStr(" 50 100 m 75 100 l 75 75 l 100 75 l 100 50 l");
2484  PrintStr(" 0 50 m 25 50 l 25 25 l 50 25 l 50 1 l s");
2485  PrintStr(" gr end } >> [ 0.5 0 0 0.5 0 0 ]");
2486  break;
2487  case 23 :
2488  PrintStr(" /BBox [ 0 0 100 100 ]");
2489  PrintStr(" /XStep 100 /YStep 100");
2490  PrintStr(" /PaintProc { begin gsave");
2491  PrintStr(" [1 7] 0 sd 0 8 50 { dup dup m 2 mul 0 l s } for");
2492  PrintStr(" 0 8 50 { dup dup 2 mul 100 m 50 add exch 50");
2493  PrintStr(" add l s } for 100 0 m 100 100 l 50 50 l f");
2494  PrintStr(" gr end } >> [ 0.24 0 0 0.24 0 0 ]");
2495  break;
2496  case 24 :
2497  PrintStr(" /BBox [ 0 0 100 100 ]");
2498  PrintStr(" /XStep 100 /YStep 100");
2499  PrintStr(" /PaintProc { begin gsave");
2500  PrintStr(" 100 100 m 100 36 l 88 36 l 88 88 l f");
2501  PrintStr(" 100 0 m 100 12 l 56 12 l 50 0 l f");
2502  PrintStr(" 0 0 m 48 0 l 48 48 l 50 48 l 56 60 l");
2503  PrintStr(" 36 60 l 36 12 l 0 12 l f [1 7] 0 sd");
2504  PrintStr(" 61 8 87 { dup dup dup 12 exch m 88 exch l s");
2505  PrintStr(" 16 exch 4 sub m 88 exch 4 sub l s } for");
2506  PrintStr(" 13 8 35 { dup dup dup 0 exch m 36 exch l s");
2507  PrintStr(" 4 exch 4 sub m 36 exch 4 sub l s } for");
2508  PrintStr(" 37 8 59 { dup dup dup 12 exch m 36 exch l s");
2509  PrintStr(" 16 exch 4 sub m 36 exch 4 sub l s } for");
2510  PrintStr(" 13 8 60 { dup dup dup 56 exch m 100 exch l s");
2511  PrintStr(" 60 exch 4 sub m 100 exch 4 sub l s } for");
2512  PrintStr(" gr end } >> [ 0.5 0 0 0.5 0 0 ]");
2513  break;
2514  case 25 :
2515  PrintStr(" /BBox [ 0 0 101 101 ]");
2516  PrintStr(" /XStep 100 /YStep 100");
2517  PrintStr(" /PaintProc { begin gsave");
2518  PrintStr(" 0 0 m 30 30 l 70 30 l 70 70 l 100 100 l 100 0 l");
2519  PrintStr(" f 30 30 m 30 70 l 70 70 l f");
2520  PrintStr(" gr end } >> [ 0.5 0 0 0.5 0 0 ]");
2521  };
2522  snprintf(cdef,28," makepattern /%s exch def",&cpat[1]);
2523  PrintStr(cdef);
2524  fPatterns[ipat] = 1;
2525  }
2526 
2527  // Define the macro cs and FA if they are not yet defined.
2528  if (fPatterns[26] == 0) {
2529  if (gStyle->GetColorModelPS()) {
2530  PrintStr(" /cs {[/Pattern /DeviceCMYK] setcolorspace} def");
2531  PrintStr(" /FA {f [/DeviceCMYK] setcolorspace} def");
2532  } else {
2533  PrintStr(" /cs {[/Pattern /DeviceRGB] setcolorspace} def");
2534  PrintStr(" /FA {f [/DeviceRGB] setcolorspace} def");
2535  }
2536  fPatterns[26] = 1;
2537  }
2538 
2539  // Activate the pattern.
2540  PrintFast(3," cs");
2541  TColor *col = gROOT->GetColor(color);
2542  if (col) {
2543  Double_t colRed = col->GetRed();
2544  Double_t colGreen = col->GetGreen();
2545  Double_t colBlue = col->GetBlue();
2546  if (gStyle->GetColorModelPS()) {
2547  Double_t colBlack = TMath::Min(TMath::Min(1-colRed,1-colGreen),1-colBlue);
2548  if (colBlack==1) {
2549  WriteReal(0);
2550  WriteReal(0);
2551  WriteReal(0);
2552  WriteReal(colBlack);
2553  } else {
2554  Double_t colCyan = (1-colRed-colBlack)/(1-colBlack);
2555  Double_t colMagenta = (1-colGreen-colBlack)/(1-colBlack);
2556  Double_t colYellow = (1-colBlue-colBlack)/(1-colBlack);
2557  WriteReal(colCyan);
2558  WriteReal(colMagenta);
2559  WriteReal(colYellow);
2560  WriteReal(colBlack);
2561  }
2562  } else {
2563  WriteReal(colRed);
2564  WriteReal(colGreen);
2565  WriteReal(colBlue);
2566  }
2567  }
2568  PrintFast(4,cpat);
2569  PrintFast(9," setcolor");
2570 }
2571 
2572 ////////////////////////////////////////////////////////////////////////////////
2573 /// Set color index for lines
2574 
2576 {
2577  fLineColor = cindex;
2578  SetColor(Int_t(cindex));
2579 }
2580 
2581 ////////////////////////////////////////////////////////////////////////////////
2582 /// Set the value of the global parameter TPostScript::fgLineJoin.
2583 /// This parameter determines the appearance of joining lines in a PostScript
2584 /// output.
2585 /// It takes one argument which may be:
2586 /// - 0 (miter join)
2587 /// - 1 (round join)
2588 /// - 2 (bevel join)
2589 /// The default value is 0 (miter join).
2590 ///
2591 /// \image html postscript_1.png
2592 ///
2593 /// To change the line join behaviour just do:
2594 /// ~~~ {.cpp}
2595 /// gStyle->SetJoinLinePS(2); // Set the PS line join to bevel.
2596 /// ~~~
2597 
2599 {
2600  fgLineJoin = linejoin;
2601  if (fgLineJoin<0) fgLineJoin=0;
2602  if (fgLineJoin>2) fgLineJoin=2;
2603 }
2604 
2605 ////////////////////////////////////////////////////////////////////////////////
2606 /// Change the line style
2607 ///
2608 /// - linestyle = 2 dashed
2609 /// - linestyle = 3 dotted
2610 /// - linestyle = 4 dash-dotted
2611 /// - linestyle = else = solid
2612 ///
2613 /// See TStyle::SetLineStyleString for style definition
2614 
2616 {
2617  if ( linestyle == fLineStyle) return;
2618  fLineStyle = linestyle;
2619  const char *st = gStyle->GetLineStyleString(linestyle);
2620  PrintFast(1,"[");
2621  Int_t nch = strlen(st);
2622  PrintFast(nch,st);
2623  PrintFast(6,"] 0 sd");
2624 }
2625 
2626 ////////////////////////////////////////////////////////////////////////////////
2627 /// Change the line width
2628 
2630 {
2631  if ( linewidth == fLineWidth) return;
2632  fLineWidth = linewidth;
2633  if (fLineWidth!=0) {
2635  PrintFast(3," lw");
2636  }
2637 }
2638 
2639 ////////////////////////////////////////////////////////////////////////////////
2640 /// Set color index for markers
2641 
2643 {
2644  fMarkerColor = cindex;
2645  SetColor(Int_t(cindex));
2646 }
2647 
2648 ////////////////////////////////////////////////////////////////////////////////
2649 /// Set the current color.
2650 
2652 {
2653  if (color < 0) color = 0;
2654  fCurrentColor = color;
2655  TColor *col = gROOT->GetColor(color);
2656  if (col)
2657  SetColor(col->GetRed(), col->GetGreen(), col->GetBlue());
2658  else
2659  SetColor(1., 1., 1.);
2660 }
2661 
2662 ////////////////////////////////////////////////////////////////////////////////
2663 /// Set directly current color (don't go via TColor).
2664 
2666 {
2667  if (r == fRed && g == fGreen && b == fBlue) return;
2668 
2669  fRed = r;
2670  fGreen = g;
2671  fBlue = b;
2672 
2673  if (fRed <= 0 && fGreen <= 0 && fBlue <= 0 ) {
2674  PrintFast(6," black");
2675  } else {
2676  if (gStyle->GetColorModelPS()) {
2677  Double_t colBlack = TMath::Min(TMath::Min(1-fRed,1-fGreen),1-fBlue);
2678  Double_t colCyan = (1-fRed-colBlack)/(1-colBlack);
2679  Double_t colMagenta = (1-fGreen-colBlack)/(1-colBlack);
2680  Double_t colYellow = (1-fBlue-colBlack)/(1-colBlack);
2681  WriteReal(colCyan);
2682  WriteReal(colMagenta);
2683  WriteReal(colYellow);
2684  WriteReal(colBlack);
2685  } else {
2686  WriteReal(fRed);
2687  WriteReal(fGreen);
2688  WriteReal(fBlue);
2689  }
2690  PrintFast(2," c");
2691  }
2692 }
2693 
2694 ////////////////////////////////////////////////////////////////////////////////
2695 /// Set color index for text
2696 
2698 {
2699  fTextColor = cindex;
2700 
2701  SetColor( Int_t(cindex) );
2702 }
2703 
2704 ////////////////////////////////////////////////////////////////////////////////
2705 /// Write a string of characters
2706 ///
2707 /// This method writes the string chars into a PostScript file
2708 /// at position xx,yy in world coordinates.
2709 
2710 void TPostScript::Text(Double_t xx, Double_t yy, const char *chars)
2711 {
2712  static const char *psfont[31][2] = {
2713  { "Root.PSFont.1", "/Times-Italic" },
2714  { "Root.PSFont.2", "/Times-Bold" },
2715  { "Root.PSFont.3", "/Times-BoldItalic" },
2716  { "Root.PSFont.4", "/Helvetica" },
2717  { "Root.PSFont.5", "/Helvetica-Oblique" },
2718  { "Root.PSFont.6", "/Helvetica-Bold" },
2719  { "Root.PSFont.7", "/Helvetica-BoldOblique" },
2720  { "Root.PSFont.8", "/Courier" },
2721  { "Root.PSFont.9", "/Courier-Oblique" },
2722  { "Root.PSFont.10", "/Courier-Bold" },
2723  { "Root.PSFont.11", "/Courier-BoldOblique" },
2724  { "Root.PSFont.12", "/Symbol" },
2725  { "Root.PSFont.13", "/Times-Roman" },
2726  { "Root.PSFont.14", "/ZapfDingbats" },
2727  { "Root.PSFont.15", "/Symbol" },
2728  { "Root.PSFont.STIXGen", "/STIXGeneral" },
2729  { "Root.PSFont.STIXGenIt", "/STIXGeneral-Italic" },
2730  { "Root.PSFont.STIXGenBd", "/STIXGeneral-Bold" },
2731  { "Root.PSFont.STIXGenBdIt", "/STIXGeneral-BoldItalic" },
2732  { "Root.PSFont.STIXSiz1Sym", "/STIXSize1Symbols" },
2733  { "Root.PSFont.STIXSiz1SymBd", "/STIXSize1Symbols-Bold" },
2734  { "Root.PSFont.STIXSiz2Sym", "/STIXSize2Symbols" },
2735  { "Root.PSFont.STIXSiz2SymBd", "/STIXSize2Symbols-Bold" },
2736  { "Root.PSFont.STIXSiz3Sym", "/STIXSize3Symbols" },
2737  { "Root.PSFont.STIXSiz3SymBd", "/STIXSize3Symbols-Bold" },
2738  { "Root.PSFont.STIXSiz4Sym", "/STIXSize4Symbols" },
2739  { "Root.PSFont.STIXSiz4SymBd", "/STIXSize4Symbols-Bold" },
2740  { "Root.PSFont.STIXSiz5Sym", "/STIXSize5Symbols" },
2741  { "Root.PSFont.ME", "/DroidSansFallback" },
2742  { "Root.PSFont.CJKMing", "/DroidSansFallback" },
2743  { "Root.PSFont.CJKGothic", "/DroidSansFallback" }
2744  };
2745 
2746  const Double_t kDEGRAD = TMath::Pi()/180.;
2747  Double_t x = xx;
2748  Double_t y = yy;
2749  if (!gPad) return;
2750 
2751  // Compute the font size. Exit if it is 0
2752  // The font size is computed from the TTF size to get exactly the same
2753  // size on the screen and in the PostScript file.
2754  Double_t wh = (Double_t)gPad->XtoPixel(gPad->GetX2());
2755  Double_t hh = (Double_t)gPad->YtoPixel(gPad->GetY1());
2756  Float_t tsize, ftsize;
2757 
2758  if (wh < hh) {
2759  tsize = fTextSize*wh;
2760  Int_t sizeTTF = (Int_t)(tsize*kScale+0.5); // TTF size
2761  ftsize = (sizeTTF*fXsize*gPad->GetAbsWNDC())/wh;
2762  } else {
2763  tsize = fTextSize*hh;
2764  Int_t sizeTTF = (Int_t)(tsize*kScale+0.5); // TTF size
2765  ftsize = (sizeTTF*fYsize*gPad->GetAbsHNDC())/hh;
2766  }
2767  Double_t fontsize = 4*(72*(ftsize)/2.54);
2768  if( fontsize <= 0) return;
2769 
2770  Float_t tsizex = gPad->AbsPixeltoX(Int_t(tsize))-gPad->AbsPixeltoX(0);
2771  Float_t tsizey = gPad->AbsPixeltoY(0)-gPad->AbsPixeltoY(Int_t(tsize));
2772 
2773  Int_t font = abs(fTextFont)/10;
2774  if( font > 31 || font < 1) font = 1;
2775 
2776  // Text color.
2778 
2779  // Text alignment.
2780  Int_t txalh = fTextAlign/10;
2781  if (txalh <1) txalh = 1; else if (txalh > 3) txalh = 3;
2782  Int_t txalv = fTextAlign%10;
2783  if (txalv <1) txalv = 1; else if (txalv > 3) txalv = 3;
2784  if (txalv == 3) {
2785  y -= 0.8*tsizey*TMath::Cos(kDEGRAD*fTextAngle);
2786  x += 0.8*tsizex*TMath::Sin(kDEGRAD*fTextAngle);
2787  } else if (txalv == 2) {
2788  y -= 0.4*tsizey*TMath::Cos(kDEGRAD*fTextAngle);
2789  x += 0.4*tsizex*TMath::Sin(kDEGRAD*fTextAngle);
2790  }
2791 
2792  UInt_t w = 0, w0 = 0;
2793  Bool_t kerning;
2794  // In order to measure the precise character positions we need to trick
2795  // FreeType into rendering high-resolution characters otherwise it will
2796  // stick to the screen pixel grid, which is far worse than we can achieve
2797  // on print.
2798  const Float_t scale = 16.0;
2799  // Save current text attributes.
2800  TText saveAttText;
2801  saveAttText.TAttText::operator=(*this);
2802  const Int_t len=strlen(chars);
2803  Int_t *charWidthsCumul = 0;
2804  TText t;
2805  t.SetTextSize(fTextSize * scale);
2807  t.GetTextAdvance(w, chars);
2808  t.GetTextAdvance(w0, chars, kFALSE);
2809  t.TAttText::Modify();
2810  if (w0-w != 0) kerning = kTRUE;
2811  else kerning = kFALSE;
2812  if (kerning) {
2813  // Calculate the individual character placements.
2814  charWidthsCumul = new Int_t[len];
2815  for (Int_t i = len - 1;i >= 0;i--) {
2816  UInt_t ww = 0;
2817  t.GetTextAdvance(ww, chars + i);
2818  Double_t wwl = (gPad->AbsPixeltoX(ww)-gPad->AbsPixeltoX(0));
2819  charWidthsCumul[i] = (Int_t)((XtoPS(wwl) - XtoPS(0)) / scale);
2820  }
2821  }
2822  // Restore text attributes.
2823  saveAttText.TAttText::Modify();
2824 
2825  Double_t charsLength = gPad->AbsPixeltoX(w)-gPad->AbsPixeltoX(0);
2826  Int_t psCharsLength = (Int_t)((XtoPS(charsLength)-XtoPS(0)) / scale);
2827 
2828  // Text angle.
2829  Int_t psangle = Int_t(0.5 + fTextAngle);
2830 
2831  // Save context.
2832  PrintStr("@");
2833  SaveRestore(1);
2834 
2835  // Clipping
2836  Int_t xc1 = XtoPS(gPad->GetX1());
2837  Int_t xc2 = XtoPS(gPad->GetX2());
2838  Int_t yc1 = YtoPS(gPad->GetY1());
2839  Int_t yc2 = YtoPS(gPad->GetY2());
2840  WriteInteger(xc2 - xc1);
2841  WriteInteger(yc2 - yc1);
2842  WriteInteger(xc1);
2843  WriteInteger(yc1);
2844  PrintStr(" C");
2845 
2846  // Output text position and angle. The text position is computed
2847  // using Double_t to avoid precision problems.
2848  Double_t vx = (x - gPad->GetX1())/(gPad->GetX2()-gPad->GetX1());
2849  Double_t cmx = fXsize*(gPad->GetAbsXlowNDC()+vx*gPad->GetAbsWNDC());
2850  WriteReal((288.*cmx)/2.54);
2851  Double_t vy = (y - gPad->GetY1())/(gPad->GetY2()-gPad->GetY1());
2852  Double_t cmy = fYsize*(gPad->GetAbsYlowNDC()+vy*gPad->GetAbsHNDC());
2853  WriteReal((288.*cmy)/2.54);
2854  PrintStr(Form(" t %d r ", psangle));
2855  if(txalh == 2) PrintStr(Form(" %d 0 t ", -psCharsLength/2));
2856  if(txalh == 3) PrintStr(Form(" %d 0 t ", -psCharsLength));
2857  PrintStr(gEnv->GetValue(psfont[font-1][0], psfont[font-1][1]));
2858  if (font != 15) {
2859  PrintStr(Form(" findfont %g sf 0 0 m ",fontsize));
2860  } else {
2861  PrintStr(Form(" findfont %g sf 0 0 m ita ",fontsize));
2862  }
2863 
2864  if (kerning) {
2865  PrintStr("@");
2866  // Output individual character placements
2867  for (Int_t i = len-1; i >= 1; i--) {
2868  WriteInteger(charWidthsCumul[0] - charWidthsCumul[i]);
2869  }
2870  delete [] charWidthsCumul;
2871  PrintStr("@");
2872  }
2873 
2874  // Output text.
2875  PrintStr("(");
2876 
2877  // Inside a PostScript string, the new line (if needed to break up long lines) must be escaped by a backslash.
2878  const char *crsave = fImplicitCREsc;
2879  fImplicitCREsc = "\\";
2880 
2881  char str[8];
2882  for (Int_t i=0; i<len;i++) {
2883  if (chars[i]!='\n') {
2884  if (chars[i]=='(' || chars[i]==')' || chars[i]=='\\') {
2885  snprintf(str,8,"\\%c",chars[i]);
2886  PrintStr(str);
2887  } else if ((chars[i]=='-') && (font != 12)) {
2888  PrintStr("\\255");
2889  } else {
2890  snprintf(str,8,"%c",chars[i]);
2891  PrintFast(1,str);
2892  }
2893  }
2894  }
2895  PrintStr(")");
2896  fImplicitCREsc = crsave;
2897 
2898  if (kerning) {
2899  if (font != 15) PrintStr(" K NC");
2900  else PrintStr(" K gr NC");
2901  } else {
2902  if (font != 15) PrintStr(" show NC");
2903  else PrintStr(" show gr NC");
2904  }
2905 
2906  SaveRestore(-1);
2907 }
2908 
2909 ////////////////////////////////////////////////////////////////////////////////
2910 /// Write a string of characters
2911 ///
2912 /// This method writes the string chars into a PostScript file
2913 /// at position xx,yy in world coordinates.
2914 
2915 void TPostScript::Text(Double_t xx, Double_t yy, const wchar_t *chars)
2916 {
2917  static const char *psfont[31][2] = {
2918  { "Root.PSFont.1", "/FreeSerifItalic" },
2919  { "Root.PSFont.2", "/FreeSerifBold" },
2920  { "Root.PSFont.3", "/FreeSerifBoldItalic" },
2921  { "Root.PSFont.4", "/FreeSans" },
2922  { "Root.PSFont.5", "/FreeSansOblique" },
2923  { "Root.PSFont.6", "/FreeSansBold" },
2924  { "Root.PSFont.7", "/FreeSansBoldOblique" },
2925  { "Root.PSFont.8", "/FreeMono" },
2926  { "Root.PSFont.9", "/FreeMonoOblique" },
2927  { "Root.PSFont.10", "/FreeMonoBold" },
2928  { "Root.PSFont.11", "/FreeMonoBoldOblique" },
2929  { "Root.PSFont.12", "/StandardSymbolsL" },
2930  { "Root.PSFont.13", "/FreeSerif" },
2931  { "Root.PSFont.14", "/Dingbats" },
2932  { "Root.PSFont.15", "/StandardSymbolsL" },
2933  { "Root.PSFont.STIXGen", "/STIXGeneral" },
2934  { "Root.PSFont.STIXGenIt", "/STIXGeneral-Italic" },
2935  { "Root.PSFont.STIXGenBd", "/STIXGeneral-Bold" },
2936  { "Root.PSFont.STIXGenBdIt", "/STIXGeneral-BoldItalic" },
2937  { "Root.PSFont.STIXSiz1Sym", "/STIXSize1Symbols" },
2938  { "Root.PSFont.STIXSiz1SymBd", "/STIXSize1Symbols-Bold" },
2939  { "Root.PSFont.STIXSiz2Sym", "/STIXSize2Symbols" },
2940  { "Root.PSFont.STIXSiz2SymBd", "/STIXSize2Symbols-Bold" },
2941  { "Root.PSFont.STIXSiz3Sym", "/STIXSize3Symbols" },
2942  { "Root.PSFont.STIXSiz3SymBd", "/STIXSize3Symbols-Bold" },
2943  { "Root.PSFont.STIXSiz4Sym", "/STIXSize4Symbols" },
2944  { "Root.PSFont.STIXSiz4SymBd", "/STIXSize4Symbols-Bold" },
2945  { "Root.PSFont.STIXSiz5Sym", "/STIXSize5Symbols" },
2946  { "Root.PSFont.ME", "/DroidSansFallback" },
2947  { "Root.PSFont.CJKMing", "/DroidSansFallback" },
2948  { "Root.PSFont.CJKGothic", "/DroidSansFallback" }
2949  };
2950 
2951  Int_t len = wcslen(chars);
2952  if (len<=0) return;
2953 
2954  const Double_t kDEGRAD = TMath::Pi()/180.;
2955  Double_t x = xx;
2956  Double_t y = yy;
2957  if (!gPad) return;
2958 
2959  // Compute the font size. Exit if it is 0
2960  // The font size is computed from the TTF size to get exactly the same
2961  // size on the screen and in the PostScript file.
2962  Double_t wh = (Double_t)gPad->XtoPixel(gPad->GetX2());
2963  Double_t hh = (Double_t)gPad->YtoPixel(gPad->GetY1());
2964  Float_t tsize, ftsize;
2965 
2966  if (wh < hh) {
2967  tsize = fTextSize*wh;
2968  Int_t sizeTTF = (Int_t)(tsize*kScale+0.5); // TTF size
2969  ftsize = (sizeTTF*fXsize*gPad->GetAbsWNDC())/wh;
2970  } else {
2971  tsize = fTextSize*hh;
2972  Int_t sizeTTF = (Int_t)(tsize*kScale+0.5); // TTF size
2973  ftsize = (sizeTTF*fYsize*gPad->GetAbsHNDC())/hh;
2974  }
2975  Double_t fontsize = 4*(72*(ftsize)/2.54);
2976  if( fontsize <= 0) return;
2977 
2978  Float_t tsizex = gPad->AbsPixeltoX(Int_t(tsize))-gPad->AbsPixeltoX(0);
2979  Float_t tsizey = gPad->AbsPixeltoY(0)-gPad->AbsPixeltoY(Int_t(tsize));
2980 
2981  Int_t font = abs(fTextFont)/10;
2982  if( font > 29 || font < 1) font = 1;
2983 
2984  // Text color.
2986 
2987  // Text alignment.
2988  Int_t txalh = fTextAlign/10;
2989  if (txalh <1) txalh = 1; else if (txalh > 3) txalh = 3;
2990  Int_t txalv = fTextAlign%10;
2991  if (txalv <1) txalv = 1; else if (txalv > 3) txalv = 3;
2992  if (txalv == 3) {
2993  y -= 0.8*tsizey*TMath::Cos(kDEGRAD*fTextAngle);
2994  x += 0.8*tsizex*TMath::Sin(kDEGRAD*fTextAngle);
2995  } else if (txalv == 2) {
2996  y -= 0.4*tsizey*TMath::Cos(kDEGRAD*fTextAngle);
2997  x += 0.4*tsizex*TMath::Sin(kDEGRAD*fTextAngle);
2998  }
2999  UInt_t w = 0, h = 0;
3000 
3001  TText t;
3004  t.GetTextExtent(w, h, chars);
3005  Double_t charsLength = gPad->AbsPixeltoX(w)-gPad->AbsPixeltoX(0);
3006  Int_t psCharsLength = XtoPS(charsLength)-XtoPS(0);
3007 
3008  // Text angle.
3009  Int_t psangle = Int_t(0.5 + fTextAngle);
3010 
3011  // Save context.
3012  PrintStr("@");
3013  SaveRestore(1);
3014 
3015  // Clipping
3016  Int_t xc1 = XtoPS(gPad->GetX1());
3017  Int_t xc2 = XtoPS(gPad->GetX2());
3018  Int_t yc1 = YtoPS(gPad->GetY1());
3019  Int_t yc2 = YtoPS(gPad->GetY2());
3020  WriteInteger(xc2 - xc1);
3021  WriteInteger(yc2 - yc1);
3022  WriteInteger(xc1);
3023  WriteInteger(yc1);
3024  PrintStr(" C");
3025 
3026  // Output text position and angle.
3027  WriteInteger(XtoPS(x));
3028  WriteInteger(YtoPS(y));
3029  PrintStr(Form(" t %d r ", psangle));
3030  if(txalh == 2) PrintStr(Form(" %d 0 t ", -psCharsLength/2));
3031  if(txalh == 3) PrintStr(Form(" %d 0 t ", -psCharsLength));
3032  MustEmbed[font-1] = kTRUE; // This font will be embedded in the file at EOF time.
3033  PrintStr(gEnv->GetValue(psfont[font-1][0], psfont[font-1][1]));
3034  PrintStr(Form(" findfont %g sf 0 0 m ",fontsize));
3035 
3036  // Output text.
3037  if (len > 1) PrintStr(Form("%d ", len));
3038  for(Int_t i = 0; i < len; i++) {
3039  // Adobe Glyph Naming Convention
3040  // http://www.adobe.com/devnet/opentype/archives/glyph.html
3041 #include "AdobeGlyphList.h"
3042  const wchar_t *lower = std::lower_bound(
3044  chars[i]);
3045  if(lower < adobe_glyph_ucs + nadobe_glyph &&
3046  *lower == chars[i]) {
3047  // Named glyph from AGL 1.2
3048  const unsigned long index =
3049  lower - adobe_glyph_ucs;
3050  PrintStr(Form("/%s ", adobe_glyph_name[index]));
3051  }
3052  else if((unsigned int)chars[i] < 0xffff) {
3053  // Unicode BMP
3054  PrintStr(Form("/uni%04X ",
3055  (unsigned int)chars[i]));
3056  }
3057  else {
3058  // Unicode supplemental planes
3059  PrintStr(Form("/u%04X ",
3060  (unsigned int)chars[i]));
3061  }
3062  }
3063  if(len > 1) {
3064  PrintStr("{glyphshow} repeat ");
3065  }
3066  else {
3067  PrintStr("glyphshow ");
3068  }
3069 
3070  PrintStr("NC");
3071 
3072  SaveRestore(-1);
3073 }
3074 
3075 ////////////////////////////////////////////////////////////////////////////////
3076 /// Write a string of characters in NDC
3077 
3078 void TPostScript::TextNDC(Double_t u, Double_t v, const char *chars)
3079 {
3080  Double_t x = gPad->GetX1() + u*(gPad->GetX2() - gPad->GetX1());
3081  Double_t y = gPad->GetY1() + v*(gPad->GetY2() - gPad->GetY1());
3082  Text(x, y, chars);
3083 }
3084 
3085 ////////////////////////////////////////////////////////////////////////////////
3086 /// Write a string of characters in NDC
3087 
3088 void TPostScript::TextNDC(Double_t u, Double_t v, const wchar_t *chars)
3089 {
3090  Double_t x = gPad->GetX1() + u*(gPad->GetX2() - gPad->GetX1());
3091  Double_t y = gPad->GetY1() + v*(gPad->GetY2() - gPad->GetY1());
3092  Text(x, y, chars);
3093 }
3094 
3095 ////////////////////////////////////////////////////////////////////////////////
3096 /// Convert U from NDC coordinate to PostScript
3097 
3099 {
3100  Double_t cm = fXsize*(gPad->GetAbsXlowNDC() + u*gPad->GetAbsWNDC());
3101  return Int_t(0.5 + 288*cm/2.54);
3102 }
3103 
3104 ////////////////////////////////////////////////////////////////////////////////
3105 /// Convert V from NDC coordinate to PostScript
3106 
3108 {
3109  Double_t cm = fYsize*(gPad->GetAbsYlowNDC() + v*gPad->GetAbsHNDC());
3110  return Int_t(0.5 + 288*cm/2.54);
3111 }
3112 
3113 ////////////////////////////////////////////////////////////////////////////////
3114 /// Convert X from world coordinate to PostScript
3115 
3117 {
3118  Double_t u = (x - gPad->GetX1())/(gPad->GetX2() - gPad->GetX1());
3119  return UtoPS(u);
3120 }
3121 
3122 ////////////////////////////////////////////////////////////////////////////////
3123 /// Convert Y from world coordinate to PostScript
3124 
3126 {
3127  Double_t v = (y - gPad->GetY1())/(gPad->GetY2() - gPad->GetY1());
3128  return VtoPS(v);
3129 }
3130 
3131 ////////////////////////////////////////////////////////////////////////////////
3132 /// Initialize the PostScript page in zones
3133 
3135 {
3136  if( !fClear )return;
3137  fClear = kFALSE;
3138 
3139  // When Zone has been called, fZone is TRUE
3140  fZone = kTRUE;
3141 
3142  if( fIYzone > fNYzone) {
3143  fIYzone=1;
3144  if( fMode != 3) {
3145  PrintStr("@showpage");
3146  SaveRestore(-1);
3147  fNpages++;
3148  PrintStr("@%%Page:");
3151  PrintStr("@");
3152  } else {
3153  PrintFast(9," showpage");
3154  SaveRestore(-1);
3155  }
3156  }
3157 
3158  // No grestore the first time
3159  if( fMode != 3) {
3160  if( fIXzone != 1 || fIYzone != 1) SaveRestore(-1);
3161  SaveRestore(1);
3162  PrintStr("@");
3165  PrintFast(5," Zone");
3166  PrintStr("@");
3167  fIXzone++;
3168  if( fIXzone > fNXzone) { fIXzone=1; fIYzone++; }
3169  }
3170 
3171  // Picture Initialisation
3172  SaveRestore(1);
3173  if (fgLineJoin) {
3175  PrintFast(12," setlinejoin");
3176  }
3177  PrintFast(6," 0 0 t");
3178  fRed = -1;
3179  fGreen = -1;
3180  fBlue = -1;
3181  fPrinted = kFALSE;
3182  fLineColor = -1;
3183  fLineStyle = -1;
3184  fLineWidth = -1;
3185  fFillColor = -1;
3186  fFillStyle = -1;
3187  fMarkerSizeCur = -1;
3188 }
Interface to PostScript.
Definition: TPostScript.h:20
virtual ~TPostScript()
Default PostScript destructor.
const char * GetLineStyleString(Int_t i=1) const
Return line style string (used by PostScript).
Definition: TStyle.cxx:792
void DrawPolyLine(Int_t n, TPoints *xy)
Draw a PolyLine.
virtual Bool_t AccessPathName(const char *path, EAccessMode mode=kFileExists)
Returns FALSE if one can access a file using the specified access mode.
Definition: TSystem.cxx:1272
static Bool_t MustEmbed[32]
const char * GetHeaderPS() const
Definition: TStyle.h:268
Int_t fNBSameColorCell
Number of boxes with the same color.
Definition: TPostScript.h:77
char * fBuffer
Definition: TVirtualPS.h:42
void Text(Double_t x, Double_t y, const char *string)
Write a string of characters.
virtual int GetPid()
Get process id.
Definition: TSystem.cxx:714
Float_t fXsize
Page size along X.
Definition: TPostScript.h:45
Float_t fYVP2
Definition: TPostScript.h:40
T1 Sign(T1 a, T2 b)
Definition: TMathBase.h:153
short Style_t
Definition: RtypesCore.h:76
Bool_t fClear
True when page must be cleared.
Definition: TPostScript.h:65
TLine * line
void On()
Activate an already open PostScript file.
float Float_t
Definition: RtypesCore.h:53
void DrawPolyLineNDC(Int_t n, TPoints *uv)
Draw a PolyLine in NDC space.
const char Option_t
Definition: RtypesCore.h:62
Int_t CMtoPS(Double_t u)
Definition: TPostScript.h:93
void DefineMarkers()
Define the markers.
R__EXTERN TStyle * gStyle
Definition: TStyle.h:402
void CellArrayFill(Int_t r, Int_t g, Int_t b)
Paint the Cell Array.
Float_t fYsize
Page size along Y.
Definition: TPostScript.h:46
TH1 * h
Definition: legend2.C:5
Float_t fXVP2
Definition: TPostScript.h:38
const char * fImplicitCREsc
Definition: TVirtualPS.h:43
Size_t fMarkerSize
Marker size.
Definition: TAttMarker.h:24
Bool_t fBoundingBox
True for Encapsulated PostScript.
Definition: TPostScript.h:64
Int_t fCurrentColor
current Postscript color index
Definition: TPostScript.h:59
#define gROOT
Definition: TROOT.h:375
Basic string class.
Definition: TString.h:129
virtual void WriteInteger(Int_t i, Bool_t space=kTRUE)
Write one Integer to the file.
Definition: TVirtualPS.cxx:170
void GetPaperSize(Float_t &xsize, Float_t &ysize) const
Set paper size for PostScript output.
Definition: TStyle.cxx:810
Float_t GetLineScalePS() const
Definition: TStyle.h:272
Short_t Min(Short_t a, Short_t b)
Definition: TMathBase.h:168
static const unsigned long nadobe_glyph
Definition: AdobeGlyphList.h:1
Int_t fLastCellRed
Last red value.
Definition: TPostScript.h:74
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
Float_t fMarkerSizeCur
current transformed value of marker size
Definition: TPostScript.h:58
void DrawPS(Int_t n, Float_t *xw, Float_t *yw)
Draw a PolyLine.
Float_t GetGreen() const
Definition: TColor.h:57
Int_t fNXzone
Number of zones along X.
Definition: TPostScript.h:54
char fPatterns[32]
Indicate if pattern n is defined.
Definition: TPostScript.h:69
Bool_t fFontEmbed
True is FontEmbed has been called.
Definition: TPostScript.h:79
void DrawBox(Double_t x1, Double_t y1, Double_t x2, Double_t y2)
Draw a Box.
void FontEmbed()
Embed font in PS file.
bool FontEmbedType42(const char *filename)
#define NULL
Definition: RtypesCore.h:88
static std::string format(double x, double y, int digits, int width)
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:628
Int_t GetJoinLinePS() const
Definition: TStyle.h:271
virtual void WriteReal(Float_t r, Bool_t space=kTRUE)
Write a Real number to the file.
Definition: TVirtualPS.cxx:185
virtual int Rename(const char *from, const char *to)
Rename a file.
Definition: TSystem.cxx:1326
Int_t fType
PostScript workstation type.
Definition: TPostScript.h:61
void SetLineStyle(Style_t linestyle=1)
Change the line style.
const char * Data() const
Definition: TString.h:344
virtual void SetTextFont(Font_t tfont=62)
Set the text font.
Definition: TAttText.h:45
Int_t fSave
Number of gsave for restore.
Definition: TPostScript.h:53
virtual int Unlink(const char *name)
Unlink, i.e. remove, a file.
Definition: TSystem.cxx:1353
static const double x2[5]
Float_t fMaxsize
Largest dimension of X and Y.
Definition: TPostScript.h:47
Int_t fNpages
number of pages
Definition: TPostScript.h:60
Double_t x[n]
Definition: legend1.C:17
TPostScript()
Default PostScript constructor.
void SetLineWidth(Width_t linewidth=1)
Change the line width.
void Open(const char *filename, Int_t type=-111)
Open a PostScript file.
static Int_t fgLineJoin
Appearance of joining lines.
Definition: TPostScript.h:81
virtual void PrintStr(const char *string="")
Output the string str in the output buffer.
Definition: TVirtualPS.cxx:72
TString fFileName
PS file name.
Definition: TPostScript.h:78
Int_t fNYzone
Number of zones along Y.
Definition: TPostScript.h:55
bool FontEmbedType2(const char *filename)
Float_t GetBlue() const
Definition: TColor.h:58
Int_t fNbCellW
Number of boxes per line.
Definition: TPostScript.h:71
Bool_t fClipStatus
Clipping Indicator.
Definition: TPostScript.h:66
Base class for several text objects.
Definition: TText.h:23
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:873
constexpr Double_t Pi()
Definition: TMath.h:40
const Float_t kScale
Definition: TASImage.cxx:130
Int_t fIYzone
Current zone along Y.
Definition: TPostScript.h:57
Int_t fNbCellLine
Number of boxes in the current line.
Definition: TPostScript.h:72
short Color_t
Definition: RtypesCore.h:79
Int_t GetColorModelPS() const
Definition: TStyle.h:182
void DrawFrame(Double_t xl, Double_t yl, Double_t xt, Double_t yt, Int_t mode, Int_t border, Int_t dark, Int_t light)
Draw a Frame around a box.
Style_t fMarkerStyle
Marker style.
Definition: TAttMarker.h:23
TCanvas * kerning()
Definition: kerning.C:1
void CellArrayEnd()
End the Cell Array painting.
void MovePS(Int_t x, Int_t y)
Move to a new position.
Int_t UtoPS(Double_t u)
Convert U from NDC coordinate to PostScript.
Float_t fTextAngle
Text angle.
Definition: TAttText.h:21
Style_t fLineStyle
Line style.
Definition: TAttLine.h:22
TRandom2 r(17)
R__EXTERN TSystem * gSystem
Definition: TSystem.h:539
void Range(Float_t xrange, Float_t yrange)
Set the range for the paper in centimeters.
SVector< double, 2 > v
Definition: Dict.h:5
XPoint xy[kMAXMK]
Definition: TGX11.cxx:122
Int_t fSizBuffer
Definition: TVirtualPS.h:39
Color_t fLineColor
Line color.
Definition: TAttLine.h:21
virtual Int_t GetValue(const char *name, Int_t dflt)
Returns the integer value for a resource.
Definition: TEnv.cxx:482
void DrawPolyMarker(Int_t n, Float_t *x, Float_t *y)
Draw markers at the n WC points x, y.
virtual Color_t GetFillColor() const
Return the fill area color.
Definition: TAttFill.h:30
void TextNDC(Double_t u, Double_t v, const char *string)
Write a string of characters in NDC.
unsigned int UInt_t
Definition: RtypesCore.h:42
Width_t fLineWidth
Line width.
Definition: TAttLine.h:23
Float_t fXVP1
Definition: TPostScript.h:37
char * Form(const char *fmt,...)
Int_t fLastCellBlue
Last blue value.
Definition: TPostScript.h:76
void DrawHatch(Float_t dy, Float_t angle, Int_t n, Float_t *x, Float_t *y)
Draw Fill area with hatch styles.
void SetMarkerColor(Color_t cindex=1)
Set color index for markers.
virtual void GetTextExtent(UInt_t &w, UInt_t &h, const char *text) const
Return text extent for string text.
Definition: TText.cxx:588
Bool_t fPrinted
Definition: TVirtualPS.h:40
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
void SetLineJoin(Int_t linejoin=0)
Set the value of the global parameter TPostScript::fgLineJoin.
std::ofstream * fStream
Definition: TVirtualPS.h:41
void FontEncode()
Font Re-encoding.
Float_t fRed
Per cent of red.
Definition: TPostScript.h:48
Font_t fTextFont
Text font.
Definition: TAttText.h:25
Bool_t fRange
True when a range has been defined.
Definition: TPostScript.h:67
Double_t Cos(Double_t)
Definition: TMath.h:551
short Width_t
Definition: RtypesCore.h:78
char * StrDup(const char *str)
Duplicate the string str.
Definition: TString.cxx:2524
void CellArrayBegin(Int_t W, Int_t H, Double_t x1, Double_t x2, Double_t y1, Double_t y2)
Draw a Cell Array.
const Bool_t kFALSE
Definition: RtypesCore.h:92
PyObject * fType
2-D graphics point (world coordinates).
Definition: TPoints.h:19
void SaveRestore(Int_t flag)
Compute number of gsaves for restore This allows to write the correct number of grestore at the end o...
Int_t YtoPS(Double_t y)
Convert Y from world coordinate to PostScript.
virtual void PrintRaw(Int_t len, const char *str)
Print a raw.
Definition: TVirtualPS.cxx:200
Int_t fNbinCT
Number of entries in the current Cell Array.
Definition: TPostScript.h:70
static const double x1[5]
void Off()
Deactivate an already open PostScript file.
#define ClassImp(name)
Definition: Rtypes.h:336
void SetFillPatterns(Int_t ipat, Int_t color)
Patterns definition.
double Double_t
Definition: RtypesCore.h:55
virtual void PrintFast(Int_t nch, const char *string="")
Fast version of Print.
Definition: TVirtualPS.cxx:103
R__EXTERN TEnv * gEnv
Definition: TEnv.h:170
Bool_t fZone
Zone indicator.
Definition: TPostScript.h:68
bool FontEmbedType1(const char *filename)
Double_t y[n]
Definition: legend1.C:17
Int_t fIXzone
Current zone along X.
Definition: TPostScript.h:56
The color creation and management class.
Definition: TColor.h:19
Float_t GetRed() const
Definition: TColor.h:56
Float_t fTextSize
Text size.
Definition: TAttText.h:22
void Zone()
Initialize the PostScript page in zones.
void SetColor(Int_t color=1)
Set the current color.
Int_t fLenBuffer
Definition: TVirtualPS.h:38
static const char * adobe_glyph_name[nadobe_glyph]
void SetLineScale(Float_t scale=3)
Definition: TPostScript.h:125
Int_t VtoPS(Double_t v)
Convert V from NDC coordinate to PostScript.
void Close(Option_t *opt="")
Close a PostScript file.
void SetFillColor(Color_t cindex=1)
Set color index for fill areas.
Color_t fFillColor
Fill area color.
Definition: TAttFill.h:22
Double_t Sin(Double_t)
Definition: TMath.h:548
Int_t fClip
Clipping mode.
Definition: TPostScript.h:63
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
Definition: TRolke.cxx:630
void Initialize()
PostScript Initialisation.
#define snprintf
Definition: civetweb.c:822
R__EXTERN TVirtualPS * gVirtualPS
Definition: TVirtualPS.h:81
#define gPad
Definition: TVirtualPad.h:284
Int_t fMode
PostScript mode.
Definition: TPostScript.h:62
Int_t fLastCellGreen
Last green value.
Definition: TPostScript.h:75
void NewPage()
Move to a new PostScript page.
virtual void GetTextAdvance(UInt_t &a, const char *text, const Bool_t kern=kTRUE) const
Return text advance for string text if kern is true (default) kerning is taken into account...
Definition: TText.cxx:616
Float_t fLineScale
Line width scale factor.
Definition: TPostScript.h:51
double result[121]
TVirtualPS is an abstract interface to Postscript, PDF, SVG.
Definition: TVirtualPS.h:30
const int nn
Double_t Sqrt(Double_t x)
Definition: TMath.h:591
virtual void SetTextSize(Float_t tsize=1)
Set the text size.
Definition: TAttText.h:46
Int_t XtoPS(Double_t x)
Convert X from world coordinate to PostScript.
const char * AsString() const
Return the date &amp; time as a string (ctime() format).
Definition: TDatime.cxx:101
Color_t fMarkerColor
Marker color.
Definition: TAttMarker.h:22
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition: TNamed.cxx:155
const char * GetTitlePS() const
Definition: TStyle.h:269
const Bool_t kTRUE
Definition: RtypesCore.h:91
Color_t fTextColor
Text color.
Definition: TAttText.h:24
void SetTextColor(Color_t cindex=1)
Set color index for text.
const Int_t n
Definition: legend1.C:16
void SetLineColor(Color_t cindex=1)
Set color index for lines.
Float_t fGreen
Per cent of green.
Definition: TPostScript.h:49
Float_t fYVP1
Definition: TPostScript.h:39
Int_t fMaxLines
Maximum number of lines in a PS array.
Definition: TPostScript.h:73
Float_t fBlue
Per cent of blue.
Definition: TPostScript.h:50
This class stores the date and time with a precision of one second in an unsigned 32 bit word (950130...
Definition: TDatime.h:37
Style_t fFillStyle
Fill area style.
Definition: TAttFill.h:23
Short_t fTextAlign
Text alignment.
Definition: TAttText.h:23
static const wchar_t adobe_glyph_ucs[nadobe_glyph]
Definition: AdobeGlyphList.h:2
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition: TObject.cxx:859