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