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