Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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 <cstdlib>
233#include <cstring>
234#include <cctype>
235#include <cwchar>
236#include <fstream>
237
238#include "strlcpy.h"
239#include "snprintf.h"
240#include "Byteswap.h"
241#include "TROOT.h"
242#include "TDatime.h"
243#include "TColor.h"
244#include "TVirtualPad.h"
245#include "TPoints.h"
246#include "TPostScript.h"
247#include "TStyle.h"
248#include "TMath.h"
249#include "TText.h"
250#include "TSystem.h"
251#include "TEnv.h"
252
253#include "../../../graf2d/mathtext/inc/fontembed.h"
254
255// to scale fonts to the same size as the old TT version
256const Float_t kScale = 0.93376068;
257
258// Array defining if a font must be embedded or not.
259static Bool_t MustEmbed[32];
260
263
264
265////////////////////////////////////////////////////////////////////////////////
266/// Default PostScript constructor
267
269{
270 fStream = nullptr;
271 fType = 0;
272 gVirtualPS = this;
273 fBlue = 0.;
275 fClear = kFALSE;
276 fClip = 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 fWidth = 0.;
304 fStyle = 1;
305 fX1v = 0.;
306 fX1w = 0.;
307 fX2v = 0.;
308 fX2w = 0.;
309 fXC = 0.;
310 fXVP1 = 0.;
311 fXVP2 = 0.;
312 fXVS1 = 0.;
313 fXVS2 = 0.;
314 fXsize = 0.;
315 fY1v = 0.;
316 fY1w = 0.;
317 fY2v = 0.;
318 fY2w = 0.;
319 fYC = 0.;
320 fYVP1 = 0.;
321 fYVP2 = 0.;
322 fYVS1 = 0.;
323 fYVS2 = 0.;
324 fYsize = 0.;
325 fZone = kFALSE;
326 fFileName = "";
328 Int_t i;
329 for (i=0; i<32; i++) fPatterns[i] = 0;
330 for (i=0; i<32; i++) MustEmbed[i] = kFALSE;
331 SetTitle("PS");
332}
333
334////////////////////////////////////////////////////////////////////////////////
335/// Initialize the PostScript interface
336///
337/// - fname : PostScript file name
338/// - wtype : PostScript workstation type
339///
340///
341/// The possible workstation types are:
342/// - 111 ps Portrait
343/// - 112 ps Landscape
344/// - 113 eps
345
348{
349 fStream = nullptr;
350 SetTitle("PS");
351 Open(fname, wtype);
352}
353
354////////////////////////////////////////////////////////////////////////////////
355/// Open a PostScript file
356
358{
359 if (fStream) {
360 Warning("Open", "postscript file already open");
361 return;
362 }
363
364 fMarkerSizeCur = 0;
365 fRed = -1;
366 fGreen = -1;
367 fBlue = -1;
368 fLenBuffer = 0;
369 fClip = 0;
370 fType = abs(wtype);
371 fClear = kTRUE;
372 fZone = kFALSE;
373 fSave = 0;
379 fMode = fType%10;
381 if (gPad) {
382 Double_t ww = gPad->GetWw();
383 Double_t wh = gPad->GetWh();
384 if (fType == 113) {
385 ww *= gPad->GetWNDC();
386 wh *= gPad->GetHNDC();
387 }
388 Double_t ratio = wh/ww;
389 if (fType == 112) {
390 xrange = fYsize;
391 yrange = xrange*ratio;
392 if (yrange > fXsize) { yrange = fXsize; xrange = yrange/ratio;}
393 } else {
394 xrange = fXsize;
395 yrange = fXsize*ratio;
396 if (yrange > fYsize) { yrange = fYsize; xrange = yrange/ratio;}
397 }
399 }
400
401 // Open OS file
403 fStream = new std::ofstream(fFileName.Data(),std::ios::out);
405 Error("Open", "Cannot open file: %s", fFileName.Data());
406 return;
407 }
408 gVirtualPS = this;
409
410 for (Int_t i=0;i<fSizBuffer;i++) fBuffer[i] = ' ';
411 if( fType == 113) {
413 PrintStr("%!PS-Adobe-2.0 EPSF-2.0@");
414 } else {
416 PrintStr("%!PS-Adobe-2.0@");
417 Initialize();
418 }
419
421 fRange = kFALSE;
422
423 // Set a default range
425
427 if (fType == 113) NewPage();
428}
429
430////////////////////////////////////////////////////////////////////////////////
431/// Default PostScript destructor
432
437
438////////////////////////////////////////////////////////////////////////////////
439/// Close a PostScript file
440
442{
443 if (!gVirtualPS) return;
444 if (!fStream) return;
445 if (gPad) gPad->Update();
446 if( fMode != 3) {
447 SaveRestore(-1);
448 if( fPrinted ) { PrintStr("showpage@"); SaveRestore(-1);}
449 PrintStr("@");
450 PrintStr("%%Trailer@");
451 PrintStr("%%Pages: ");
453 PrintStr("@");
454 while (fSave > 0) { SaveRestore(-1); }
455 } else {
456 PrintStr("@");
457 while (fSave > 0) { SaveRestore(-1); }
458 PrintStr("showpage@");
459 PrintStr("end@");
460 }
461 PrintStr("@");
462 PrintStr("%%EOF@");
463
464 // Embed the fonts previously used by TMathText
465 if (!fFontEmbed) {
466 // Close the file fFileName
467 if (fStream) {
468 PrintStr("@");
469 fStream->close(); delete fStream; fStream = nullptr;
470 }
471
472 // Rename the file fFileName
474 if (gSystem->Rename( fFileName.Data() , tmpname.Data())) {
475 Error("Close", "Cannot open temporary file: %s", tmpname.Data());
476 return;
477 }
478
479 // Reopen the file fFileName
480 fStream = new std::ofstream(fFileName.Data(),std::ios::out);
482 Error("Close", "Cannot open file: %s", fFileName.Data());
483 return;
484 }
485
486 // Embed the fonts at the right place
487 FILE *sg = fopen(tmpname.Data(),"r");
488 if (!sg) {
489 Error("Close", "Cannot open file: %s", tmpname.Data());
490 return;
491 }
492 char line[255];
493 while (fgets(line,255,sg)) {
494 if (strstr(line,"EndComments")) PrintStr("%%DocumentNeededResources: ProcSet (FontSetInit)@");
495 fStream->write(line,strlen(line));
496 if (!fFontEmbed && strstr(line,"m5")) {
497 FontEmbed();
498 PrintStr("@");
499 }
500 }
501 fclose(sg);
502 if (gSystem->Unlink(tmpname.Data())) return;
503 }
504
506
507 // Close file stream
508
509 if (fStream) { fStream->close(); delete fStream; fStream = nullptr;}
510
511 gVirtualPS = nullptr;
512}
513
514////////////////////////////////////////////////////////////////////////////////
515/// Activate an already open PostScript file
516
518{
519 if (!fType) {
520 Error("On", "no postscript file open");
521 Off();
522 return;
523 }
524 gVirtualPS = this;
525}
526
527////////////////////////////////////////////////////////////////////////////////
528/// Deactivate an already open PostScript file
529
531{
532 gVirtualPS = nullptr;
533}
534
535////////////////////////////////////////////////////////////////////////////////
536/// Draw a Cell Array
537///
538/// Drawing a PostScript Cell Array is in fact done thanks to three
539/// procedures: CellArrayBegin, CellArrayFill, and CellArrayEnd.
540///
541/// - CellArrayBegin: Initiate the Cell Array by writing the necessary
542/// PostScript procedures and the initial values of the
543/// required parameters. The input parameters are:
544/// - W: number of boxes along the width.
545/// - H: number of boxes along the height
546/// - x1,x2,y1,y2: First box coordinates.
547/// - CellArrayFill: Is called for each box of the Cell Array. The first
548/// box is the top left one and the last box is the
549/// bottom right one. The input parameters are the Red,
550/// Green, and Blue components of the box colour. These
551/// Levels are between 0 and 255.
552/// - CellArrayEnd: Finishes the Cell Array.
553///
554/// PostScript cannot handle arrays larger than 65535. So the Cell Array
555/// is drawn in several pieces.
556
559{
560 Int_t ix1 = XtoPS(x1);
561 Int_t iy1 = YtoPS(y1);
562
563 Float_t wt = (288/2.54)*gPad->GetAbsWNDC()*
564 fXsize*((x2 - x1)/(gPad->GetX2()-gPad->GetX1()));
565 Float_t ht = (288/2.54)*gPad->GetAbsHNDC()*
566 fYsize*((y2 - y1)/(gPad->GetY2()-gPad->GetY1()));
567
568 fLastCellRed = 300;
569 fLastCellGreen = 300;
570 fLastCellBlue = 300;
572
573 fNbinCT = 0;
574 fNbCellW = W;
575 fNbCellLine = 0;
576 fMaxLines = 40000/(3*fNbCellW);
577
578 // Define some parameters
579 PrintStr("@/WT"); WriteReal(wt) ; PrintStr(" def"); // Cells width
580 PrintStr(" /HT"); WriteReal(ht) ; PrintStr(" def"); // Cells height
581 PrintStr(" /XS"); WriteInteger(ix1) ; PrintStr(" def"); // X start
582 PrintStr(" /YY"); WriteInteger(iy1) ; PrintStr(" def"); // Y start
583 PrintStr(" /NX"); WriteInteger(W) ; PrintStr(" def"); // Number of columns
584 PrintStr(" /NY"); WriteInteger(fMaxLines); PrintStr(" def"); // Number of lines
585
586 // This PS procedure draws one cell.
587 PrintStr(" /DrawCell ");
588 PrintStr( "{WT HT XX YY bf");
589 PrintStr( " /NBBD NBBD 1 add def");
590 PrintStr( " NBBD NBB eq {exit} if");
591 PrintStr( " /XX WT XX add def");
592 PrintStr( " IX NX eq ");
593 PrintStr( "{/YY YY HT sub def");
594 PrintStr( " /XX XS def");
595 PrintStr( " /IX 0 def} if");
596 PrintStr( " /IX IX 1 add def} def");
597
598 // This PS procedure draws fMaxLines line. It takes care of duplicated
599 // colors. Values "n" greater than 300 mean than the previous color
600 // should be duplicated n-300 times.
601 PrintStr(" /DrawCT ");
602 PrintStr( "{/NBB NX NY mul def");
603 PrintStr( " /XX XS def");
604 PrintStr( " /IX 1 def");
605 PrintStr( " /NBBD 0 def");
606 PrintStr( " /RC 0 def /GC 1 def /BC 2 def");
607 PrintStr( " 1 1 NBB ");
608 PrintStr( "{/NB CT RC get def");
609 PrintStr( " NB 301 ge ");
610 PrintStr( "{/NBL NB 300 sub def");
611 PrintStr( " 1 1 NBL ");
612 PrintStr( "{DrawCell}");
613 PrintStr( " for");
614 PrintStr( " /RC RC 1 add def");
615 PrintStr( " /GC RC 1 add def");
616 PrintStr( " /BC RC 2 add def}");
617 PrintStr( "{CT RC get 255 div CT GC get 255 div CT BC get 255 div setrgbcolor");
618 PrintStr( " DrawCell");
619 PrintStr( " /RC RC 3 add def");
620 PrintStr( " /GC GC 3 add def");
621 PrintStr( " /BC BC 3 add def} ifelse NBBD NBB eq {exit} if} for");
622 PrintStr( " /YY YY HT sub def clear} def");
623
624 PrintStr(" /CT [");
625}
626
627////////////////////////////////////////////////////////////////////////////////
628/// Paint the Cell Array
629
631{
632 if (fLastCellRed == r && fLastCellGreen == g && fLastCellBlue == b) {
634 } else {
635 if (fNBSameColorCell != 0 ) {
638 }
642 fLastCellRed = r;
645 }
646
647 fNbinCT++;
648 if (fNbinCT == fNbCellW) {
649 fNbCellLine++;
650 fNbinCT = 0;
651 }
652
653 if (fNbCellLine == fMaxLines) {
655 PrintStr("] def DrawCT /CT [");
656 fNbCellLine = 0;
657 fLastCellRed = 300;
658 fLastCellGreen = 300;
659 fLastCellBlue = 300;
661 fNbinCT = 0;
662 }
663}
664
665////////////////////////////////////////////////////////////////////////////////
666/// End the Cell Array painting
667
669{
671 PrintStr("] def /NY");
673 PrintStr(" def DrawCT ");
674}
675
676////////////////////////////////////////////////////////////////////////////////
677/// Define the markers
678
680{
681 PrintStr("/mp {newpath /y exch def /x exch def} def@");
682 PrintStr("/side {[w .77 mul w .23 mul] .385 w mul sd w 0 l currentpoint t -144 r} def@");
683 PrintStr("/mr {mp x y w2 0 360 arc} def /m24 {mr s} def /m20 {mr f} def@");
684 PrintStr("/mb {mp x y w2 add m w2 neg 0 d 0 w neg d w 0 d 0 w d cl} def@");
685 PrintStr("/mt {mp x y w2 add m w2 neg w neg d w 0 d cl} def@");
686 PrintStr("/w4 {w 4 div} def@");
687 PrintStr("/w6 {w 6 div} def@");
688 PrintStr("/w8 {w 8 div} def@");
689 PrintStr("/m21 {mb f} def /m25 {mb s} def /m22 {mt f} def /m26{mt s} def@");
690 PrintStr("/m23 {mp x y w2 sub m w2 w d w neg 0 d cl f} def@");
691 PrintStr("/m27 {mp x y w2 add m w3 neg w2 neg d w3 w2 neg d w3 w2 d cl s} def@");
692 PrintStr("/m28 {mp x w2 sub y w2 sub w3 add m w3 0 d ");
693 PrintStr(" 0 w3 neg d w3 0 d 0 w3 d w3 0 d ");
694 PrintStr(" 0 w3 d w3 neg 0 d 0 w3 d w3 neg 0 d");
695 PrintStr(" 0 w3 neg d w3 neg 0 d cl s } def@");
696 PrintStr("/m29 {mp gsave x w2 sub y w2 add w3 sub m currentpoint t");
697 PrintStr(" 4 {side} repeat cl fill gr} def@");
698 PrintStr("/m30 {mp gsave x w2 sub y w2 add w3 sub m currentpoint t");
699 PrintStr(" 4 {side} repeat cl s gr} def@");
700 PrintStr("/m31 {mp x y w2 sub m 0 w d x w2 sub y m w 0 d");
701 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");
702 PrintStr(" sub m w 1.44 div w 1.44 div d s} def@");
703 PrintStr("/m32 {mp x y w2 sub m w2 w d w neg 0 d cl s} def@");
704 PrintStr("/m33 {mp x y w2 add m w3 neg w2 neg d w3 w2 neg d w3 w2 d cl f} def@");
705 PrintStr("/m34 {mp x w2 sub y w2 sub w3 add m w3 0 d ");
706 PrintStr(" 0 w3 neg d w3 0 d 0 w3 d w3 0 d ");
707 PrintStr(" 0 w3 d w3 neg 0 d 0 w3 d w3 neg 0 d");
708 PrintStr(" 0 w3 neg d w3 neg 0 d cl f } def@");
709 PrintStr("/m35 {mp x y w2 add m w2 neg w2 neg d w2 w2 neg d w2 w2 d w2 neg w2 d");
710 PrintStr(" x y w2 sub m 0 w d x w2 sub y m w 0 d s} def@");
711 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@");
712 PrintStr("/m37 {mp x y m w4 neg w2 d w4 neg w2 neg d w2 0 d ");
713 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@");
714 PrintStr("/m38 {mp x w4 sub y w2 add m w4 neg w4 neg d 0 w2 neg d w4 w4 neg d");
715 PrintStr(" w2 0 d w4 w4 d 0 w2 d w4 neg w4 d w2 neg 0 d");
716 PrintStr(" x y w2 sub m 0 w d x w2 sub y m w 0 d cl s} def@");
717 PrintStr("/m39 {mp x y m w4 neg w2 d w4 neg w2 neg d w2 0 d ");
718 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@");
719 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");
720 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@");
721 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");
722 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@");
723 PrintStr("/m42 {mp x y w2 add m w8 neg w2 -3 4 div mul d w2 -3 4 div mul w8 neg d");
724 PrintStr(" w2 3 4 div mul w8 neg d w8 w2 -3 4 div mul d");
725 PrintStr(" w8 w2 3 4 div mul d w2 3 4 div mul w8 d");
726 PrintStr(" w2 -3 4 div mul w8 d w8 neg w2 3 4 div mul d cl s} def@");
727 PrintStr("/m43 {mp x y w2 add m w8 neg w2 -3 4 div mul d w2 -3 4 div mul w8 neg d");
728 PrintStr(" w2 3 4 div mul w8 neg d w8 w2 -3 4 div mul d");
729 PrintStr(" w8 w2 3 4 div mul d w2 3 4 div mul w8 d");
730 PrintStr(" w2 -3 4 div mul w8 d w8 neg w2 3 4 div mul d cl f} def@");
731 PrintStr("/m44 {mp x y m w6 neg w2 d w2 2 3 div mul 0 d w6 neg w2 neg d");
732 PrintStr(" w2 w6 d 0 w2 -2 3 div mul d w2 neg w6 d");
733 PrintStr(" w6 w2 neg d w2 -2 3 div mul 0 d w6 w2 d");
734 PrintStr(" w2 neg w6 neg d 0 w2 2 3 div mul d w2 w6 neg d cl s} def@");
735 PrintStr("/m45 {mp x y m w6 neg w2 d w2 2 3 div mul 0 d w6 neg w2 neg d");
736 PrintStr(" w2 w6 d 0 w2 -2 3 div mul d w2 neg w6 d");
737 PrintStr(" w6 w2 neg d w2 -2 3 div mul 0 d w6 w2 d");
738 PrintStr(" w2 neg w6 neg d 0 w2 2 3 div mul d w2 w6 neg d cl f} def@");
739 PrintStr("/m46 {mp x y w4 add m w4 neg w4 d w4 neg w4 neg d ");
740 PrintStr(" w4 w4 neg d w4 neg w4 neg d w4 w4 neg d w4 w4 d");
741 PrintStr(" w4 w4 neg d w4 w4 d w4 neg w4 d w4 w4 d w4 neg w4 d w4 neg w4 neg d cl s} def@");
742 PrintStr("/m47 {mp x y w4 add m w4 neg w4 d w4 neg w4 neg d");
743 PrintStr(" w4 w4 neg d w4 neg w4 neg d w4 w4 neg d w4 w4 d");
744 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@");
745 PrintStr("/m48 {mp x y w4 add m w4 neg w4 d w4 neg w4 neg d w4 w4 neg d ");
746 PrintStr(" w4 neg w4 neg d w4 w4 neg d w4 w4 d w4 w4 neg d w4 w4 d");
747 PrintStr(" w4 neg w4 d w4 w4 d w4 neg w4 d w4 neg w4 neg d ");
748 PrintStr(" w4 w4 neg d w4 neg w4 neg d w4 neg w4 d w4 w4 d cl f} def@");
749 PrintStr("/m49 {mp x w2 sub w3 add y w2 sub w3 add m ");
750 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");
751 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@");
752 PrintStr("/m2 {mp x y w2 sub m 0 w d x w2 sub y m w 0 d s} def@");
753 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@");
754}
755
756////////////////////////////////////////////////////////////////////////////////
757/// Draw a Box
758
760{
761 static Double_t x[4], y[4];
762 Int_t ix1 = XtoPS(x1);
763 Int_t ix2 = XtoPS(x2);
764 Int_t iy1 = YtoPS(y1);
765 Int_t iy2 = YtoPS(y2);
766 Int_t fillis = fFillStyle/1000;
767 Int_t fillsi = fFillStyle%1000;
768
769 if (fillis == 3 || fillis == 2) {
770 if (fillsi > 99) {
771 x[0] = x1; y[0] = y1;
772 x[1] = x2; y[1] = y1;
773 x[2] = x2; y[2] = y2;
774 x[3] = x1; y[3] = y2;
775 return;
776 }
777 if (fillsi > 0 && fillsi < 26) {
778 x[0] = x1; y[0] = y1;
779 x[1] = x2; y[1] = y1;
780 x[2] = x2; y[2] = y2;
781 x[3] = x1; y[3] = y2;
782 DrawPS(-4, &x[0], &y[0]);
783 }
784 if (fillsi == -3) {
785 SetColor(5);
790 PrintFast(3," bf");
791 }
792 }
793 if (fillis == 1) {
799 PrintFast(3," bf");
800 }
801 if (fillis == 0) {
802 if (fLineWidth<=0) return;
808 PrintFast(3," bl");
809 }
810}
811
812////////////////////////////////////////////////////////////////////////////////
813/// Draw a Frame around a box
814///
815/// - mode = -1 box looks as it is behind the screen
816/// - mode = 1 box looks as it is in front of the screen
817/// - border is the border size in already precomputed PostScript units
818/// - dark is the color for the dark part of the frame
819/// - light is the color for the light part of the frame
820
822 Int_t mode, Int_t border, Int_t dark, Int_t light)
823{
824 static Int_t xps[7], yps[7];
825 Int_t i, ixd0, iyd0, idx, idy, ixdi, iydi, ix, iy;
826
827 // Draw top&left part of the box
828 if (mode == -1) SetColor(dark);
829 else SetColor(light);
830 Int_t bordPS = 4*border;
831 xps[0] = XtoPS(xl); yps[0] = YtoPS(yl);
832 xps[1] = xps[0] + bordPS; yps[1] = yps[0] + bordPS;
833 xps[2] = xps[1]; yps[2] = YtoPS(yt) - bordPS;
834 xps[3] = XtoPS(xt) - bordPS; yps[3] = yps[2];
835 xps[4] = XtoPS(xt); yps[4] = YtoPS(yt);
836 xps[5] = xps[0]; yps[5] = yps[4];
837 xps[6] = xps[0]; yps[6] = yps[0];
838
839 ixd0 = xps[0];
840 iyd0 = yps[0];
843
844 PrintFast(2," m");
845 idx = 0;
846 idy = 0;
847 for (i=1;i<7;i++) {
848 ixdi = xps[i];
849 iydi = yps[i];
850 ix = ixdi - ixd0;
851 iy = iydi - iyd0;
852 ixd0 = ixdi;
853 iyd0 = iydi;
854 if( ix && iy) {
855 if( idx ) { MovePS(idx,0); idx = 0; }
856 if( idy ) { MovePS(0,idy); idy = 0; }
857 MovePS(ix,iy);
858 continue;
859 }
860 if ( ix ) {
861 if( idy ) { MovePS(0,idy); idy = 0; }
862 if( !idx ) { idx = ix; continue;}
863 if( ix*idx > 0 ) idx += ix;
864 else { MovePS(idx,0); idx = ix; }
865 continue;
866 }
867 if( iy ) {
868 if( idx ) { MovePS(idx,0); idx = 0; }
869 if( !idy) { idy = iy; continue;}
870 if( iy*idy > 0 ) idy += iy;
871 else { MovePS(0,idy); idy = iy; }
872 }
873 }
874 if( idx ) MovePS(idx,0);
875 if( idy ) MovePS(0,idy);
876 PrintFast(2," f");
877
878 // Draw bottom&right part of the box
879 if (mode == -1) SetColor(light);
880 else SetColor(dark);
881 xps[0] = XtoPS(xl); yps[0] = YtoPS(yl);
882 xps[1] = xps[0] + bordPS; yps[1] = yps[0] + bordPS;
883 xps[2] = XtoPS(xt) - bordPS; yps[2] = yps[1];
884 xps[3] = xps[2]; yps[3] = YtoPS(yt) - bordPS;
885 xps[4] = XtoPS(xt); yps[4] = YtoPS(yt);
886 xps[5] = xps[4]; yps[5] = yps[0];
887 xps[6] = xps[0]; yps[6] = yps[0];
888
889 ixd0 = xps[0];
890 iyd0 = yps[0];
893
894 PrintFast(2," m");
895 idx = 0;
896 idy = 0;
897 for (i=1;i<7;i++) {
898 ixdi = xps[i];
899 iydi = yps[i];
900 ix = ixdi - ixd0;
901 iy = iydi - iyd0;
902 ixd0 = ixdi;
903 iyd0 = iydi;
904 if( ix && iy) {
905 if( idx ) { MovePS(idx,0); idx = 0; }
906 if( idy ) { MovePS(0,idy); idy = 0; }
907 MovePS(ix,iy);
908 continue;
909 }
910 if ( ix ) {
911 if( idy ) { MovePS(0,idy); idy = 0; }
912 if( !idx ) { idx = ix; continue;}
913 if( ix*idx > 0 ) idx += ix;
914 else { MovePS(idx,0); idx = ix; }
915 continue;
916 }
917 if( iy ) {
918 if( idx ) { MovePS(idx,0); idx = 0; }
919 if( !idy) { idy = iy; continue;}
920 if( iy*idy > 0 ) idy += iy;
921 else { MovePS(0,idy); idy = iy; }
922 }
923 }
924 if( idx ) MovePS(idx,0);
925 if( idy ) MovePS(0,idy);
926 PrintFast(2," f");
927}
928
929////////////////////////////////////////////////////////////////////////////////
930/// Draw a PolyLine
931///
932/// Draw a polyline through the points xy.
933/// - If nn=1 moves only to point x,y.
934/// - If nn=0 the x,y are written in the PostScript file
935/// according to the current transformation.
936/// - If nn>0 the line is clipped as a line.
937/// - If nn<0 the line is clipped as a fill area.
938
940{
941 Int_t i, n, ixd0, iyd0, idx, idy, ixdi, iydi, ix, iy;
944 if (nn > 0) {
945 if (fLineWidth<=0) return;
946 n = nn;
950 } else {
951 n = -nn;
952 SetStyle(1);
953 SetWidth(1);
955 }
956
957 ixd0 = XtoPS(xy[0].GetX());
958 iyd0 = YtoPS(xy[0].GetY());
961 if( n <= 1) {
962 if( n == 0) goto END;
963 PrintFast(2," m");
964 goto END;
965 }
966
967 PrintFast(2," m");
968 idx = 0;
969 idy = 0;
970 for (i=1;i<n;i++) {
971 ixdi = XtoPS(xy[i].GetX());
972 iydi = YtoPS(xy[i].GetY());
973 ix = ixdi - ixd0;
974 iy = iydi - iyd0;
975 ixd0 = ixdi;
976 iyd0 = iydi;
977 if( ix && iy) {
978 if( idx ) { MovePS(idx,0); idx = 0; }
979 if( idy ) { MovePS(0,idy); idy = 0; }
980 MovePS(ix,iy);
981 continue;
982 }
983 if ( ix ) {
984 if( idy ) { MovePS(0,idy); idy = 0; }
985 if( !idx ) { idx = ix; continue;}
986 if( ix*idx > 0 ) idx += ix;
987 else { MovePS(idx,0); idx = ix; }
988 continue;
989 }
990 if( iy ) {
991 if( idx ) { MovePS(idx,0); idx = 0; }
992 if( !idy) { idy = iy; continue;}
993 if( iy*idy > 0 ) idy += iy;
994 else { MovePS(0,idy); idy = iy; }
995 }
996 }
997 if( idx ) MovePS(idx,0);
998 if( idy ) MovePS(0,idy);
999
1000 if (nn > 0 ) {
1001 if (xy[0].GetX() == xy[n-1].GetX() && xy[0].GetY() == xy[n-1].GetY()) PrintFast(3," cl");
1002 PrintFast(2," s");
1003 } else {
1004 PrintFast(2," f");
1005 }
1006END:
1007 if (nn < 0) {
1010 }
1011}
1012
1013////////////////////////////////////////////////////////////////////////////////
1014/// Draw a PolyLine in NDC space
1015///
1016/// Draw a polyline through the points xy.
1017/// - If nn=1 moves only to point x,y.
1018/// - If nn=0 the x,y are written in the PostScript file
1019/// according to the current transformation.
1020/// - If nn>0 the line is clipped as a line.
1021/// - If nn<0 the line is clipped as a fill area.
1022
1024{
1025 Int_t i, n, ixd0, iyd0, idx, idy, ixdi, iydi, ix, iy;
1028 if (nn > 0) {
1029 if (fLineWidth<=0) return;
1030 n = nn;
1034 } else {
1035 n = -nn;
1036 SetStyle(1);
1037 SetWidth(1);
1039 }
1040
1041 ixd0 = UtoPS(xy[0].GetX());
1042 iyd0 = VtoPS(xy[0].GetY());
1045 if( n <= 1) {
1046 if( n == 0) goto END;
1047 PrintFast(2," m");
1048 goto END;
1049 }
1050
1051 PrintFast(2," m");
1052 idx = 0;
1053 idy = 0;
1054 for (i=1;i<n;i++) {
1055 ixdi = UtoPS(xy[i].GetX());
1056 iydi = VtoPS(xy[i].GetY());
1057 ix = ixdi - ixd0;
1058 iy = iydi - iyd0;
1059 ixd0 = ixdi;
1060 iyd0 = iydi;
1061 if( ix && iy) {
1062 if( idx ) { MovePS(idx,0); idx = 0; }
1063 if( idy ) { MovePS(0,idy); idy = 0; }
1064 MovePS(ix,iy);
1065 continue;
1066 }
1067 if ( ix ) {
1068 if( idy ) { MovePS(0,idy); idy = 0; }
1069 if( !idx ) { idx = ix; continue;}
1070 if( ix*idx > 0 ) idx += ix;
1071 else { MovePS(idx,0); idx = ix; }
1072 continue;
1073 }
1074 if( iy ) {
1075 if( idx ) { MovePS(idx,0); idx = 0; }
1076 if( !idy) { idy = iy; continue;}
1077 if( iy*idy > 0 ) idy += iy;
1078 else { MovePS(0,idy); idy = iy; }
1079 }
1080 }
1081 if( idx ) MovePS(idx,0);
1082 if( idy ) MovePS(0,idy);
1083
1084 if (nn > 0 ) {
1085 if (xy[0].GetX() == xy[n-1].GetX() && xy[0].GetY() == xy[n-1].GetY()) PrintFast(3," cl");
1086 PrintFast(2," s");
1087 } else {
1088 PrintFast(2," f");
1089 }
1090END:
1091 if (nn < 0) {
1094 }
1095}
1096
1097////////////////////////////////////////////////////////////////////////////////
1098/// Draw markers at the n WC points x, y
1099
1101{
1102 Int_t i, np, markerstyle;
1104 static char chtemp[10];
1105
1106 if (!fMarkerSize) return;
1110 SetStyle(1);
1114 if (markerstyle <= 0) strlcpy(chtemp, " m20",10);
1115 if (markerstyle == 1) strlcpy(chtemp, " m20",10);
1116 if (markerstyle == 2) strlcpy(chtemp, " m2",10);
1117 if (markerstyle == 3) strlcpy(chtemp, " m31",10);
1118 if (markerstyle == 4) strlcpy(chtemp, " m24",10);
1119 if (markerstyle == 5) strlcpy(chtemp, " m5",10);
1120 if (markerstyle >= 6 && markerstyle <= 19) strlcpy(chtemp, " m20",10);
1121 if (markerstyle >= 20 && markerstyle <= 49 ) snprintf(chtemp,10," m%d", markerstyle);
1122 if (markerstyle >= 50) strlcpy(chtemp, " m20",10);
1123
1124 // Set the PostScript marker size
1125 if (markerstyle == 1 || (markerstyle >= 9 && markerstyle <= 19)) {
1126 markersize = 2.;
1127 } else if (markerstyle == 6) {
1128 markersize = 4.;
1129 } else if (markerstyle == 7) {
1130 markersize = 8.;
1131 } else {
1133 const Int_t kBASEMARKER = 8;
1135 Float_t s2x = sbase / Float_t(gPad->GetWw() * gPad->GetAbsWNDC());
1136 markersize = this->UtoPS(s2x) - this->UtoPS(0);
1137 }
1138
1139 if (fMarkerSizeCur != markersize) {
1141 PrintFast(3," /w");
1143 PrintFast(40," def /w2 {w 2 div} def /w3 {w 3 div} def");
1144 }
1145
1146 WriteInteger(XtoPS(x[0]));
1147 WriteInteger(YtoPS(y[0]));
1148 if (n == 1) {
1152 return;
1153 }
1154 np = 1;
1155 for (i=1;i<n;i++) {
1156 WriteInteger(XtoPS(x[i]));
1157 WriteInteger(YtoPS(y[i]));
1158 np++;
1159 if (np == 100 || i == n-1) {
1161 PrintFast(2," {");
1163 PrintFast(3,"} R");
1164 np = 0;
1165 }
1166 }
1169}
1170
1171////////////////////////////////////////////////////////////////////////////////
1172/// Draw markers at the n WC points x, y
1173
1175{
1176 Int_t i, np, markerstyle;
1178 static char chtemp[10];
1179
1180 if (!fMarkerSize) return;
1184 SetStyle(1);
1188 if (markerstyle <= 0) strlcpy(chtemp, " m20",10);
1189 if (markerstyle == 1) strlcpy(chtemp, " m20",10);
1190 if (markerstyle == 2) strlcpy(chtemp, " m2",10);
1191 if (markerstyle == 3) strlcpy(chtemp, " m31",10);
1192 if (markerstyle == 4) strlcpy(chtemp, " m24",10);
1193 if (markerstyle == 5) strlcpy(chtemp, " m5",10);
1194 if (markerstyle >= 6 && markerstyle <= 19) strlcpy(chtemp, " m20",10);
1195 if (markerstyle >= 20 && markerstyle <= 49 ) snprintf(chtemp,10," m%d", markerstyle);
1196 if (markerstyle >= 50) strlcpy(chtemp, " m20",10);
1197
1198 // Set the PostScript marker size
1199 if (markerstyle == 1 || (markerstyle >= 9 && markerstyle <= 19)) {
1200 markersize = 2.;
1201 } else if (markerstyle == 6) {
1202 markersize = 4.;
1203 } else if (markerstyle == 7) {
1204 markersize = 8.;
1205 } else {
1207 const Int_t kBASEMARKER = 8;
1209 Float_t s2x = sbase / Float_t(gPad->GetWw() * gPad->GetAbsWNDC());
1210 markersize = this->UtoPS(s2x) - this->UtoPS(0);
1211 }
1212
1213 if (fMarkerSizeCur != markersize) {
1215 PrintFast(3," /w");
1217 PrintFast(40," def /w2 {w 2 div} def /w3 {w 3 div} def");
1218 }
1219
1220 WriteInteger(XtoPS(x[0]));
1221 WriteInteger(YtoPS(y[0]));
1222 if (n == 1) {
1226 return;
1227 }
1228 np = 1;
1229 for (i=1;i<n;i++) {
1230 WriteInteger(XtoPS(x[i]));
1231 WriteInteger(YtoPS(y[i]));
1232 np++;
1233 if (np == 100 || i == n-1) {
1235 PrintFast(2," {");
1237 PrintFast(3,"} R");
1238 np = 0;
1239 }
1240 }
1243}
1244
1245////////////////////////////////////////////////////////////////////////////////
1246/// Draw a PolyLine
1247///
1248/// Draw a polyline through the points xw,yw.
1249/// - If nn=1 moves only to point xw,yw.
1250/// - If nn=0 the XW(1) and YW(1) are written in the PostScript file
1251/// according to the current NT.
1252/// - If nn>0 the line is clipped as a line.
1253/// - If nn<0 the line is clipped as a fill area.
1254
1256{
1257 static Float_t dyhatch[24] = {.0075,.0075,.0075,.0075,.0075,.0075,.0075,.0075,
1258 .01 ,.01 ,.01 ,.01 ,.01 ,.01 ,.01 ,.01 ,
1259 .015 ,.015 ,.015 ,.015 ,.015 ,.015 ,.015 ,.015};
1260 static Float_t anglehatch[24] = {180, 90,135, 45,150, 30,120, 60,
1261 180, 90,135, 45,150, 30,120, 60,
1262 180, 90,135, 45,150, 30,120, 60};
1263 Int_t i, n, ixd0, iyd0, idx, idy, ixdi, iydi, ix, iy, fais, fasi;
1264 fais = fasi = n = 0;
1265 Int_t jxd0 = XtoPS(xw[0]);
1266 Int_t jyd0 = YtoPS(yw[0]);
1269
1270 if (nn > 0) {
1271 if (fLineWidth<=0) return;
1272 n = nn;
1276 }
1277 if (nn < 0) {
1278 n = -nn;
1279 SetStyle(1);
1280 SetWidth(1);
1282 fais = fFillStyle/1000;
1283 fasi = fFillStyle%1000;
1284 if (fais == 3 || fais == 2) {
1285 if (fasi > 100 && fasi <125) {
1286 DrawHatch(dyhatch[fasi-101],anglehatch[fasi-101], n, xw, yw);
1287 goto END;
1288 }
1289 if (fasi > 0 && fasi < 26) {
1291 }
1292 }
1293 }
1294
1295 ixd0 = jxd0;
1296 iyd0 = jyd0;
1299 if( n <= 1) {
1300 if( n == 0) goto END;
1301 PrintFast(2," m");
1302 goto END;
1303 }
1304
1305 PrintFast(2," m");
1306 idx = idy = 0;
1307 for (i=1;i<n;i++) {
1308 ixdi = XtoPS(xw[i]);
1309 iydi = YtoPS(yw[i]);
1310 ix = ixdi - ixd0;
1311 iy = iydi - iyd0;
1312 ixd0 = ixdi;
1313 iyd0 = iydi;
1314 if( ix && iy) {
1315 if( idx ) { MovePS(idx,0); idx = 0; }
1316 if( idy ) { MovePS(0,idy); idy = 0; }
1317 MovePS(ix,iy);
1318 } else if ( ix ) {
1319 if( idy ) { MovePS(0,idy); idy = 0;}
1320 if( !idx ) { idx = ix;}
1321 else if( TMath::Sign(ix,idx) == ix ) idx += ix;
1322 else { MovePS(idx,0); idx = ix;}
1323 } else if( iy ) {
1324 if( idx ) { MovePS(idx,0); idx = 0;}
1325 if( !idy) { idy = iy;}
1326 else if( TMath::Sign(iy,idy) == iy) idy += iy;
1327 else { MovePS(0,idy); idy = iy;}
1328 }
1329 }
1330 if (idx) MovePS(idx,0);
1331 if (idy) MovePS(0,idy);
1332
1333 if (nn > 0 ) {
1334 if (xw[0] == xw[n-1] && yw[0] == yw[n-1]) PrintFast(3," cl");
1335 PrintFast(2," s");
1336 } else {
1337 if (fais == 0) {PrintFast(5," cl s"); goto END;}
1338 if (fais == 3 || fais == 2) {
1339 if (fasi > 0 && fasi < 26) {
1340 PrintFast(3," FA");
1341 fRed = -1;
1342 fGreen = -1;
1343 fBlue = -1;
1344 }
1345 goto END;
1346 }
1347 PrintFast(2," f");
1348 }
1349END:
1350 if (nn < 0) {
1353 }
1354}
1355
1356////////////////////////////////////////////////////////////////////////////////
1357/// Draw a PolyLine
1358///
1359/// Draw a polyline through the points xw,yw.
1360/// - If nn=1 moves only to point xw,yw.
1361/// - If nn=0 the xw(1) and YW(1) are written in the PostScript file
1362/// --- according to the current NT.
1363/// - If nn>0 the line is clipped as a line.
1364/// - If nn<0 the line is clipped as a fill area.
1365
1367{
1368 static Float_t dyhatch[24] = {.0075,.0075,.0075,.0075,.0075,.0075,.0075,.0075,
1369 .01 ,.01 ,.01 ,.01 ,.01 ,.01 ,.01 ,.01 ,
1370 .015 ,.015 ,.015 ,.015 ,.015 ,.015 ,.015 ,.015};
1371 static Float_t anglehatch[24] = {180, 90,135, 45,150, 30,120, 60,
1372 180, 90,135, 45,150, 30,120, 60,
1373 180, 90,135, 45,150, 30,120, 60};
1374 Int_t i, n, ixd0, iyd0, idx, idy, ixdi, iydi, ix, iy, fais, fasi;
1375 fais = fasi = n = 0;
1376 Int_t jxd0 = XtoPS(xw[0]);
1377 Int_t jyd0 = YtoPS(yw[0]);
1380
1381 if (nn > 0) {
1382 if (fLineWidth<=0) return;
1383 n = nn;
1387 }
1388 if (nn < 0) {
1389 n = -nn;
1390 SetStyle(1);
1391 SetWidth(1);
1393 fais = fFillStyle/1000;
1394 fasi = fFillStyle%1000;
1395 if (fais == 3 || fais == 2) {
1396 if (fasi > 100 && fasi <125) {
1397 DrawHatch(dyhatch[fasi-101],anglehatch[fasi-101], n, xw, yw);
1398 goto END;
1399 }
1400 if (fasi > 0 && fasi < 26) {
1402 }
1403 }
1404 }
1405
1406 ixd0 = jxd0;
1407 iyd0 = jyd0;
1410 if( n <= 1) {
1411 if( n == 0) goto END;
1412 PrintFast(2," m");
1413 goto END;
1414 }
1415
1416 PrintFast(2," m");
1417 idx = idy = 0;
1418 for (i=1;i<n;i++) {
1419 ixdi = XtoPS(xw[i]);
1420 iydi = YtoPS(yw[i]);
1421 ix = ixdi - ixd0;
1422 iy = iydi - iyd0;
1423 ixd0 = ixdi;
1424 iyd0 = iydi;
1425 if( ix && iy) {
1426 if( idx ) { MovePS(idx,0); idx = 0; }
1427 if( idy ) { MovePS(0,idy); idy = 0; }
1428 MovePS(ix,iy);
1429 } else if ( ix ) {
1430 if( idy ) { MovePS(0,idy); idy = 0;}
1431 if( !idx ) { idx = ix;}
1432 else if( TMath::Sign(ix,idx) == ix ) idx += ix;
1433 else { MovePS(idx,0); idx = ix;}
1434 } else if( iy ) {
1435 if( idx ) { MovePS(idx,0); idx = 0;}
1436 if( !idy) { idy = iy;}
1437 else if( TMath::Sign(iy,idy) == iy) idy += iy;
1438 else { MovePS(0,idy); idy = iy;}
1439 }
1440 }
1441 if (idx) MovePS(idx,0);
1442 if (idy) MovePS(0,idy);
1443
1444 if (nn > 0 ) {
1445 if (xw[0] == xw[n-1] && yw[0] == yw[n-1]) PrintFast(3," cl");
1446 PrintFast(2," s");
1447 } else {
1448 if (fais == 0) {PrintFast(5," cl s"); goto END;}
1449 if (fais == 3 || fais == 2) {
1450 if (fasi > 0 && fasi < 26) {
1451 PrintFast(3," FA");
1452 fRed = -1;
1453 fGreen = -1;
1454 fBlue = -1;
1455 }
1456 goto END;
1457 }
1458 PrintFast(2," f");
1459 }
1460END:
1461 if (nn < 0) {
1464 }
1465}
1466
1467////////////////////////////////////////////////////////////////////////////////
1468/// Draw Fill area with hatch styles
1469
1471{
1472 Warning("DrawHatch", "hatch fill style not yet implemented");
1473}
1474
1475////////////////////////////////////////////////////////////////////////////////
1476/// Draw Fill area with hatch styles
1477
1479{
1480 Warning("DrawHatch", "hatch fill style not yet implemented");
1481}
1482
1483////////////////////////////////////////////////////////////////////////////////
1484
1486{
1487 std::ifstream font_file(filename, std::ios::binary);
1488
1489 // We cannot read directly using iostream iterators due to
1490 // signedness
1491 font_file.seekg(0, std::ios::end);
1492
1493 const size_t font_file_length = font_file.tellg();
1494
1495 font_file.seekg(0, std::ios::beg);
1496
1497 std::vector<unsigned char> font_data(font_file_length, '\0');
1498
1499 font_file.read(reinterpret_cast<char *>(&font_data[0]),
1501
1502 std::string font_name;
1503 std::string postscript_string =
1504 mathtext::font_embed_postscript_t::font_embed_type_1(
1505 font_name, font_data);
1506
1507 if (!postscript_string.empty()) {
1509 PrintStr("@");
1510
1511 return true;
1512 }
1513
1514 return false;
1515}
1516
1517////////////////////////////////////////////////////////////////////////////////
1518
1520{
1521 std::ifstream font_file(filename, std::ios::binary);
1522
1523 // We cannot read directly using iostream iterators due to
1524 // signedness
1525 font_file.seekg(0, std::ios::end);
1526
1527 const size_t font_file_length = font_file.tellg();
1528
1529 font_file.seekg(0, std::ios::beg);
1530
1531 std::vector<unsigned char> font_data(font_file_length, '\0');
1532
1533 font_file.read(reinterpret_cast<char *>(&font_data[0]), font_file_length);
1534
1535 std::string font_name;
1536 std::string postscript_string =
1537 mathtext::font_embed_postscript_t::font_embed_type_2(font_name, font_data);
1538
1539 if (!postscript_string.empty()) {
1541 PrintStr("@");
1542
1543 return true;
1544 }
1545
1546 return false;
1547}
1548
1549////////////////////////////////////////////////////////////////////////////////
1550
1552{
1553 std::ifstream font_file(filename, std::ios::binary);
1554
1555 // We cannot read directly using iostream iterators due to signedness
1556
1557 font_file.seekg(0, std::ios::end);
1558
1559 const size_t font_file_length = font_file.tellg();
1560
1561 font_file.seekg(0, std::ios::beg);
1562
1563 std::vector<unsigned char> font_data(font_file_length, '\0');
1564
1565 font_file.read(reinterpret_cast<char *>(&font_data[0]), font_file_length);
1566
1567 std::string font_name;
1568 std::string postscript_string =
1569 mathtext::font_embed_postscript_t::font_embed_type_42(font_name, font_data);
1570
1571 if (!postscript_string.empty()) {
1573 PrintStr("@");
1574
1575 return true;
1576 }
1577 fprintf(stderr, "%s:%d:\n", __FILE__, __LINE__);
1578
1579 return false;
1580}
1581
1582////////////////////////////////////////////////////////////////////////////////
1583/// Embed font in PS file.
1584
1586{
1587 static const char *fonttable[32][2] = {
1588 { "Root.TTFont.0", "FreeSansBold.otf" },
1589 { "Root.TTFont.1", "FreeSerifItalic.otf" },
1590 { "Root.TTFont.2", "FreeSerifBold.otf" },
1591 { "Root.TTFont.3", "FreeSerifBoldItalic.otf" },
1592 { "Root.TTFont.4", "FreeSans.otf" },
1593 { "Root.TTFont.5", "FreeSansOblique.otf" },
1594 { "Root.TTFont.6", "FreeSansBold.otf" },
1595 { "Root.TTFont.7", "FreeSansBoldOblique.otf" },
1596 { "Root.TTFont.8", "FreeMono.otf" },
1597 { "Root.TTFont.9", "FreeMonoOblique.otf" },
1598 { "Root.TTFont.10", "FreeMonoBold.otf" },
1599 { "Root.TTFont.11", "FreeMonoBoldOblique.otf" },
1600 { "Root.TTFont.12", "symbol.ttf" },
1601 { "Root.TTFont.13", "FreeSerif.otf" },
1602 { "Root.TTFont.14", "wingding.ttf" },
1603 { "Root.TTFont.15", "symbol.ttf" },
1604 { "Root.TTFont.STIXGen", "STIXGeneral.otf" },
1605 { "Root.TTFont.STIXGenIt", "STIXGeneralItalic.otf" },
1606 { "Root.TTFont.STIXGenBd", "STIXGeneralBol.otf" },
1607 { "Root.TTFont.STIXGenBdIt", "STIXGeneralBolIta.otf" },
1608 { "Root.TTFont.STIXSiz1Sym", "STIXSiz1Sym.otf" },
1609 { "Root.TTFont.STIXSiz1SymBd", "STIXSiz1SymBol.otf" },
1610 { "Root.TTFont.STIXSiz2Sym", "STIXSiz2Sym.otf" },
1611 { "Root.TTFont.STIXSiz2SymBd", "STIXSiz2SymBol.otf" },
1612 { "Root.TTFont.STIXSiz3Sym", "STIXSiz3Sym.otf" },
1613 { "Root.TTFont.STIXSiz3SymBd", "STIXSiz3SymBol.otf" },
1614 { "Root.TTFont.STIXSiz4Sym", "STIXSiz4Sym.otf" },
1615 { "Root.TTFont.STIXSiz4SymBd", "STIXSiz4SymBol.otf" },
1616 { "Root.TTFont.STIXSiz5Sym", "STIXSiz5Sym.otf" },
1617 { "Root.TTFont.ME", "DroidSansFallback.ttf" },
1618 { "Root.TTFont.CJKMing", "DroidSansFallback.ttf" },
1619 { "Root.TTFont.CJKCothic", "DroidSansFallback.ttf" }
1620 };
1621
1622 PrintStr("%%IncludeResource: ProcSet (FontSetInit)@");
1623
1624 // try to load font (font must be in Root.TTFontPath resource)
1625 const char *ttpath = gEnv->GetValue("Root.TTFontPath",
1627
1628 for (Int_t fontid = 1; fontid < 30; fontid++) {
1629 if (fontid != 15 && MustEmbed[fontid-1]) {
1630 const char *filename = gEnv->GetValue(
1631 fonttable[fontid][0], fonttable[fontid][1]);
1633 if (!ttfont) {
1634 Error("TPostScript::FontEmbed",
1635 "font %d (filename `%s') not found in path",
1636 fontid, filename);
1637 } else {
1638 if (FontEmbedType2(ttfont)) {
1639 // nothing
1640 } else if(FontEmbedType1(ttfont)) {
1641 // nothing
1642 } else if(FontEmbedType42(ttfont)) {
1643 // nothing
1644 } else {
1645 Error("TPostScript::FontEmbed",
1646 "failed to embed font %d (filename `%s')",
1647 fontid, filename);
1648 }
1649 delete [] ttfont;
1650 }
1651 }
1652 }
1653 PrintStr("%%IncludeResource: font Times-Roman@");
1654 PrintStr("%%IncludeResource: font Times-Italic@");
1655 PrintStr("%%IncludeResource: font Times-Bold@");
1656 PrintStr("%%IncludeResource: font Times-BoldItalic@");
1657 PrintStr("%%IncludeResource: font Helvetica@");
1658 PrintStr("%%IncludeResource: font Helvetica-Oblique@");
1659 PrintStr("%%IncludeResource: font Helvetica-Bold@");
1660 PrintStr("%%IncludeResource: font Helvetica-BoldOblique@");
1661 PrintStr("%%IncludeResource: font Courier@");
1662 PrintStr("%%IncludeResource: font Courier-Oblique@");
1663 PrintStr("%%IncludeResource: font Courier-Bold@");
1664 PrintStr("%%IncludeResource: font Courier-BoldOblique@");
1665 PrintStr("%%IncludeResource: font Symbol@");
1666 PrintStr("%%IncludeResource: font ZapfDingbats@");
1667
1668 fFontEmbed = kTRUE;
1669}
1670
1671////////////////////////////////////////////////////////////////////////////////
1672/// Font Re-encoding
1673
1675{
1676 PrintStr("/reEncode ");
1677 PrintStr("{exch findfont");
1678 PrintStr(" dup length dict begin");
1679 PrintStr(" {1 index /FID eq ");
1680 PrintStr(" {pop pop}");
1681 PrintStr(" {def} ifelse");
1682 PrintStr(" } forall");
1683 PrintStr(" /Encoding exch def");
1684 PrintStr(" currentdict end");
1685 PrintStr(" dup /FontName get exch");
1686 PrintStr(" definefont pop");
1687 PrintStr(" } def");
1688 PrintStr(" [/Times-Bold /Times-Italic /Times-BoldItalic /Helvetica");
1689 PrintStr(" /Helvetica-Oblique /Helvetica-Bold /Helvetica-BoldOblique");
1690 PrintStr(" /Courier /Courier-Oblique /Courier-Bold /Courier-BoldOblique");
1691 PrintStr(" /Times-Roman /AvantGarde-Book /AvantGarde-BookOblique");
1692 PrintStr(" /AvantGarde-Demi /AvantGarde-DemiOblique /Bookman-Demi");
1693 PrintStr(" /Bookman-DemiItalic /Bookman-Light /Bookman-LightItalic");
1694 PrintStr(" /Helvetica-Narrow /Helvetica-Narrow-Bold /Helvetica-Narrow-BoldOblique");
1695 PrintStr(" /Helvetica-Narrow-Oblique /NewCenturySchlbk-Roman /NewCenturySchlbk-Bold");
1696 PrintStr(" /NewCenturySchlbk-BoldItalic /NewCenturySchlbk-Italic");
1697 PrintStr(" /Palatino-Bold /Palatino-BoldItalic /Palatino-Italic /Palatino-Roman");
1698 PrintStr(" ] {ISOLatin1Encoding reEncode } forall");
1699}
1700
1701////////////////////////////////////////////////////////////////////////////////
1702/// PostScript Initialisation
1703///
1704/// This method initialize the following PostScript procedures:
1705///
1706/// | Macro Name | Input parameters | Explanation |
1707/// |------------|------------------|-----------------------------------|
1708/// | l | x y | Draw a line to the x y position |
1709/// | m | x y | Move to the position x y |
1710/// | box | dx dy x y | Define a box |
1711/// | bl | dx dy x y | Draw a line box |
1712/// | bf | dx dy x y | Draw a filled box |
1713/// | t | x y | Translate |
1714/// | r | angle | Rotate |
1715/// | rl | i j | Roll the stack |
1716/// | d | x y | Draw a relative line to x y |
1717/// | X | x | Draw a relative line to x (y=0) |
1718/// | Y | y | Draw a relative line to y (x=0) |
1719/// | rm | x y | Move relatively to x y |
1720/// | gr | | Restore the graphic context |
1721/// | lw | lwidth | Set line width to lwidth |
1722/// | sd | [] 0 | Set dash line define by [] |
1723/// | s | | Stroke mode |
1724/// | c | r g b | Set rgb color to r g b |
1725/// | cl | | Close path |
1726/// | f | | Fill the last describe path |
1727/// | mXX | x y | Draw the marker type XX at (x,y) |
1728/// | Zone | ix iy | Define the current zone |
1729/// | black | | The color is black |
1730/// | C | dx dy x y | Clipping on |
1731/// | NC | | Clipping off |
1732/// | R | | repeat |
1733/// | ita | | Used to make the symbols italic |
1734/// | K | | kshow |
1735
1737{
1739 rpxmin = rpymin = width = heigth = 0;
1740 Int_t format;
1741 fNpages=1;
1742 for (Int_t i=0;i<32;i++) fPatterns[i]=0;
1743
1744 // Mode is last digit of PostScript Workstation type
1745 // mode=1,2 for portrait/landscape black and white
1746 // mode=3 for Encapsulated PostScript File
1747 // mode=4 for portrait colour
1748 // mode=5 for lanscape colour
1749 Int_t atype = abs(fType);
1750 fMode = atype%10;
1751 if( fMode <= 0 || fMode > 5) {
1752 Error("Initialize", "invalid file type %d", fMode);
1753 return;
1754 }
1755
1756 // fNXzone (fNYzone) is the total number of windows in x (y)
1757 fNXzone = (atype%1000)/100;
1758 fNYzone = (atype%100)/10;
1759 if( fNXzone <= 0 ) fNXzone = 1;
1760 if( fNYzone <= 0 ) fNYzone = 1;
1761 fIXzone = 1;
1762 fIYzone = 1;
1763
1764 // format = 0-99 is the European page format (A4,A3 ...)
1765 // format = 100 is the US format 8.5x11.0 inch
1766 // format = 200 is the US format 8.5x14.0 inch
1767 // format = 300 is the US format 11.0x17.0 inch
1768 format = atype/1000;
1769 if( format == 0 ) format = 4;
1770 if( format == 99 ) format = 0;
1771
1772 PrintStr("%%Title: ");
1773 const char *pstitle = gStyle->GetTitlePS();
1774 if (gPad && !pstitle[0]) pstitle = gPad->GetMother()->GetTitle();
1775 if (strlen(GetName())<=80) PrintStr(GetName());
1776 if(!pstitle[0] && fMode != 3) {;
1777 PrintFast(2," (");
1778 if ( format <= 99 ) {;
1779 PrintFast(2," A");
1781 PrintFast(1,")");
1782 }
1783 else {
1784 if ( format == 100 ) PrintFast(8," Letter)");
1785 if ( format == 200 ) PrintFast(7," Legal)");
1786 if ( format == 300 ) PrintFast(8," Ledger)");
1787 }
1788 PrintStr("@");
1789 PrintStr("%%Pages: (atend)@");
1790 }
1791 else {
1792 if (!strchr(pstitle,'\n')) {
1793 PrintFast(2,": ");
1795 }
1796 PrintStr("@");
1797 }
1798
1799 PrintFast(24,"%%Creator: ROOT Version ");
1800 PrintStr(gROOT->GetVersion());
1801 PrintStr("@");
1802 PrintFast(16,"%%CreationDate: ");
1803 TDatime t;
1804 PrintStr(t.AsString());
1805 PrintStr("@");
1806
1807 if ( fMode == 1 || fMode == 4) PrintStr("%%Orientation: Portrait@");
1808 if ( fMode == 2 || fMode == 5) PrintStr("%%Orientation: Landscape@");
1809
1810 PrintStr("%%EndComments@");
1811 PrintStr("%%BeginProlog@");
1812
1813 if( fMode == 3)PrintStr("80 dict begin@");
1814
1815 // Initialisation of PostScript procedures
1816 PrintStr("/s {stroke} def /l {lineto} def /m {moveto} def /t {translate} def@");
1817 PrintStr("/r {rotate} def /rl {roll} def /R {repeat} def@");
1818 PrintStr("/d {rlineto} def /rm {rmoveto} def /gr {grestore} def /f {eofill} def@");
1819 if (gStyle->GetColorModelPS()) {
1820 PrintStr("/c {setcmykcolor} def /black {0 0 0 1 setcmykcolor} def /sd {setdash} def@");
1821 } else {
1822 PrintStr("/c {setrgbcolor} def /black {0 setgray} def /sd {setdash} def@");
1823 }
1824 PrintStr("/cl {closepath} def /sf {scalefont setfont} def /lw {setlinewidth} def@");
1825 PrintStr("/box {m dup 0 exch d exch 0 d 0 exch neg d cl} def@");
1826 PrintStr("/NC{systemdict begin initclip end}def/C{NC box clip newpath}def@");
1827 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 @");
1828 PrintStr("/K {{pop pop 0 moveto} exch kshow} bind def@");
1829 PrintStr("/ita {/ang 15 def gsave [1 0 ang dup sin exch cos div 1 0 0] concat} def @");
1830
1831 DefineMarkers();
1832
1833 FontEncode();
1834
1835 // mode=1 for portrait black/white
1836 if (fMode == 1) {
1837 rpxmin = 0.7;
1839 switch (format) {
1840 case 100 :
1841 width = (8.5*2.54)-2.*rpxmin;
1842 heigth = (11.*2.54)-2.*rpymin;
1843 break;
1844 case 200 :
1845 width = (8.5*2.54)-2.*rpxmin;
1846 heigth = (14.*2.54)-2.*rpymin;
1847 break;
1848 case 300 :
1849 width = (11.*2.54)-2.*rpxmin;
1850 heigth = (17.*2.54)-2.*rpymin;
1851 break;
1852 default :
1853 width = 21.0-2.*rpxmin;
1854 heigth = 29.7-2.*rpymin;
1855 };
1856 }
1857
1858 // mode=2 for landscape black/white
1859 if (fMode == 2) {
1860 rpymin = 0.7;
1862 switch (format) {
1863 case 100 :
1864 width = (11.*2.54)-2.*rpxmin;
1865 heigth = (8.5*2.54)-2.*rpymin;
1866 break;
1867 case 200 :
1868 width = (14.*2.54)-2.*rpxmin;
1869 heigth = (8.5*2.54)-2.*rpymin;
1870 break;
1871 case 300 :
1872 width = (17.*2.54)-2.*rpxmin;
1873 heigth = (11.*2.54)-2.*rpymin;
1874 break;
1875 default :
1876 width = 29.7-2.*rpxmin;
1877 heigth = 21-2.*rpymin;
1878 };
1879 }
1880
1881 // mode=3 encapsulated PostScript
1882 if (fMode == 3) {
1883 width = 20;
1884 heigth = 20;
1885 format = 4;
1886 fNXzone = 1;
1887 fNYzone = 1;
1888 }
1889
1890 // mode=4 for portrait colour
1891 if (fMode == 4) {
1892 rpxmin = 0.7;
1893 rpymin = 3.4;
1894 switch (format) {
1895 case 100 :
1896 width = (8.5*2.54)-2.*rpxmin;
1897 heigth = (11.*2.54)-2.*rpymin;
1898 break;
1899 case 200 :
1900 width = (8.5*2.54)-2.*rpxmin;
1901 heigth = (14.*2.54)-2.*rpymin;
1902 break;
1903 case 300 :
1904 width = (11.*2.54)-2.*rpxmin;
1905 heigth = (17.*2.54)-2.*rpymin;
1906 break;
1907 default :
1908 width = (21.0-2*rpxmin);
1909 heigth = (29.7-2.*rpymin);
1910 };
1911 }
1912
1913 // mode=5 for lanscape colour
1914 if (fMode == 5) {
1915 rpxmin = 3.4;
1916 rpymin = 0.7;
1917 switch (format) {
1918 case 100 :
1919 width = (11.*2.54)-2.*rpxmin;
1920 heigth = (8.5*2.54)-2.*rpymin;
1921 break;
1922 case 200 :
1923 width = (14.*2.54)-2.*rpxmin;
1924 heigth = (8.5*2.54)-2.*rpymin;
1925 break;
1926 case 300 :
1927 width = (17.*2.54)-2.*rpxmin;
1928 heigth = (11.*2.54)-2.*rpymin;
1929 break;
1930 default :
1931 width = (29.7-2*rpxmin);
1932 heigth = (21-2.*rpymin);
1933 };
1934 }
1935
1936 Double_t value = 0;
1937 if (format < 100) value = 21*TMath::Power(TMath::Sqrt(2.), 4-format);
1938 else if (format == 100) value = 8.5*2.54;
1939 else if (format == 200) value = 8.5*2.54;
1940 else if (format == 300) value = 11.*2.54;
1941 if (format >= 100) format = 4;
1942
1943 // Compute size (in points) of the window for each picture = f(fNXzone,fNYzone)
1946 Int_t npx = 4*CMtoPS(sizex);
1947 Int_t npy = 4*CMtoPS(sizey);
1948 if (sizex > sizey) fMaxsize = CMtoPS(sizex);
1949 else fMaxsize = CMtoPS(sizey);
1950
1951 // Procedure Zone
1952 if (fMode != 3) {
1953 PrintFast(33,"/Zone {/iy exch def /ix exch def ");
1954 PrintFast(10," ix 1 sub ");
1956 PrintFast(5," mul ");
1958 PrintFast(8," iy sub ");
1960 PrintStr(" mul t} def@");
1961 } else {
1962 PrintStr("@");
1963 }
1964
1965 PrintStr("%%EndProlog@");
1966 PrintStr("%%BeginSetup@");
1967 PrintStr("%%EndSetup@");
1968 PrintFast(8,"newpath ");
1969 SaveRestore(1);
1970 if (fMode == 1 || fMode == 4) {
1973 PrintFast(2," t");
1974 }
1975 if (fMode == 2 || fMode == 5) {
1976 PrintFast(7," 90 r 0");
1978 PrintFast(3," t ");
1981 PrintFast(2," t");
1982 }
1983
1984 PrintFast(15," .25 .25 scale ");
1985 if (fMode != 3) {
1986 SaveRestore(1);
1987 PrintStr("@");
1988 PrintStr("%%Page: 1 1@");
1989 SaveRestore(1);
1990 }
1991
1992 //Check is user has defined a special header in the current style
1994 if (nh) {
1996 if (fMode != 3) SaveRestore(1);
1997 }
1998}
1999
2000////////////////////////////////////////////////////////////////////////////////
2001/// Move to a new position
2002
2004{
2005 if (ix != 0 && iy != 0) {
2006 WriteInteger(ix);
2007 WriteInteger(iy);
2008 PrintFast(2," d");
2009 } else if (ix != 0) {
2010 WriteInteger(ix);
2011 PrintFast(2," X");
2012 } else if (iy != 0) {
2013 WriteInteger(iy);
2014 PrintFast(2," Y");
2015 }
2016}
2017
2018////////////////////////////////////////////////////////////////////////////////
2019/// Move to a new PostScript page
2020
2022{
2023 // Compute pad conversion coefficients
2024 if (gPad) {
2025 // if (!gPad->GetPadPaint()) gPad->Update();
2026 Double_t ww = gPad->GetWw();
2027 Double_t wh = gPad->GetWh();
2028 fYsize = fXsize*wh/ww;
2029 } else fYsize = 27;
2030
2031 if(fType == 113 && !fBoundingBox) {
2033 PrintStr("@%%BoundingBox: ");
2034 Double_t xlow=0, ylow=0, xup=1, yup=1;
2035 if (gPad) {
2036 xlow = gPad->GetAbsXlowNDC();
2037 xup = xlow + gPad->GetAbsWNDC();
2038 ylow = gPad->GetAbsYlowNDC();
2039 yup = ylow + gPad->GetAbsHNDC();
2040 }
2041 WriteInteger(CMtoPS(fXsize*xlow));
2042 WriteInteger(CMtoPS(fYsize*ylow));
2045 PrintStr("@");
2046 Initialize();
2048 fPrinted = psave;
2049 }
2050 if (fPrinted) {
2051 if (fSave) SaveRestore(-1);
2052 fClear = kTRUE;
2053 fPrinted = kFALSE;
2054 }
2055 Zone();
2056}
2057
2058////////////////////////////////////////////////////////////////////////////////
2059/// Set the range for the paper in centimeters
2060
2062{
2063 Float_t xps=0, yps=0, xncm=0, yncm=0, dxwn=0, dywn=0, xwkwn=0, ywkwn=0, xymax=0;
2064
2065 fXsize = xsize;
2066 fYsize = ysize;
2067 if( fType != 113) { xps = fXsize; yps = fYsize; }
2068 else { xps = xsize; yps = ysize; }
2069
2070 if( xsize <= xps && ysize < yps) {
2071 if ( xps > yps ) xymax = xps;
2072 else xymax = yps;
2073 xncm = xsize/xymax;
2074 yncm = ysize/xymax;
2075 dxwn = ((xps/xymax)-xncm)/2;
2076 dywn = ((yps/xymax)-yncm)/2;
2077 } else {
2078 if (xps/yps < 1) xwkwn = xps/yps;
2079 else xwkwn = 1;
2080 if (yps/xps < 1) ywkwn = yps/xps;
2081 else ywkwn = 1;
2082
2083 if (xsize < ysize) {
2085 yncm = ywkwn;
2086 dxwn = (xwkwn-xncm)/2;
2087 dywn = 0;
2088 if( dxwn < 0) {
2089 xncm = xwkwn;
2090 dxwn = 0;
2092 dywn = (ywkwn-yncm)/2;
2093 }
2094 } else {
2095 xncm = xwkwn;
2097 dxwn = 0;
2098 dywn = (ywkwn-yncm)/2;
2099 if( dywn < 0) {
2100 yncm = ywkwn;
2101 dywn = 0;
2103 dxwn = (xwkwn-xncm)/2;
2104 }
2105 }
2106 }
2107 fXVP1 = dxwn;
2108 fXVP2 = xncm+dxwn;
2109 fYVP1 = dywn;
2110 fYVP2 = yncm+dywn;
2111 fRange = kTRUE;
2112}
2113
2114////////////////////////////////////////////////////////////////////////////////
2115/// Compute number of gsaves for restore
2116/// This allows to write the correct number of grestore at the
2117/// end of the PS file.
2118
2120{
2121 if (flag == 1) { PrintFast(7," gsave "); fSave++; }
2122 else { PrintFast(4," gr "); fSave--; }
2123}
2124
2125////////////////////////////////////////////////////////////////////////////////
2126/// Set color index for fill areas
2127
2129{
2131 if (gStyle->GetFillColor() <= 0) cindex = 0;
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()) {
2433 if (colBlack==1) {
2434 WriteReal(0);
2435 WriteReal(0);
2436 WriteReal(0);
2438 } else {
2446 }
2447 } else {
2451 }
2452 }
2453 PrintFast(4,cpat);
2454 PrintFast(9," setcolor");
2455}
2456
2457////////////////////////////////////////////////////////////////////////////////
2458/// Set color index for lines
2459
2464
2465////////////////////////////////////////////////////////////////////////////////
2466/// Set the value of the global parameter TPostScript::fgLineJoin.
2467/// This parameter determines the appearance of joining lines in a PostScript
2468/// output.
2469/// It takes one argument which may be:
2470/// - 0 (miter join)
2471/// - 1 (round join)
2472/// - 2 (bevel join)
2473/// The default value is 0 (miter join).
2474///
2475/// \image html postscript_1.png
2476///
2477/// To change the line join behaviour just do:
2478/// ~~~ {.cpp}
2479/// gStyle->SetJoinLinePS(2); // Set the PS line join to bevel.
2480/// ~~~
2481
2483{
2485 if (fgLineJoin<0) fgLineJoin=0;
2486 if (fgLineJoin>2) fgLineJoin=2;
2487}
2488
2489////////////////////////////////////////////////////////////////////////////////
2490/// Set the value of the global parameter TPostScript::fgLineCap.
2491/// This parameter determines the appearance of line caps in a PostScript
2492/// output.
2493/// It takes one argument which may be:
2494/// - 0 (butt caps)
2495/// - 1 (round caps)
2496/// - 2 (projecting caps)
2497/// The default value is 0 (butt caps).
2498///
2499/// \image html postscript_2.png
2500///
2501/// To change the line cap behaviour just do:
2502/// ~~~ {.cpp}
2503/// gStyle->SetCapLinePS(2); // Set the PS line cap to projecting.
2504/// ~~~
2505
2507{
2509 if (fgLineCap<0) fgLineCap=0;
2510 if (fgLineCap>2) fgLineCap=2;
2511}
2512
2513////////////////////////////////////////////////////////////////////////////////
2514/// Change the line style
2515///
2516/// - linestyle = 2 dashed
2517/// - linestyle = 3 dotted
2518/// - linestyle = 4 dash-dotted
2519/// - linestyle = else = solid
2520///
2521/// See TStyle::SetLineStyleString for style definition
2522
2528
2529////////////////////////////////////////////////////////////////////////////////
2530/// Change the line style in the output file
2531
2533{
2534 if (linestyle == fStyle) return;
2535
2536 fStyle = linestyle;
2537
2538 const char *st = gStyle->GetLineStyleString(fStyle);
2539 PrintFast(1,"[");
2540 Int_t nch = strlen(st);
2541 PrintFast(nch,st);
2542 PrintFast(6,"] 0 sd");
2543}
2544
2545////////////////////////////////////////////////////////////////////////////////
2546/// Change the line width
2547
2553
2554////////////////////////////////////////////////////////////////////////////////
2555/// Change the line width in the output file
2556
2558{
2559 if (linewidth == fWidth) return;
2560
2561 fWidth = linewidth;
2562
2563 if (fWidth!=0) {
2565 PrintFast(3," lw");
2566 }
2567}
2568
2569////////////////////////////////////////////////////////////////////////////////
2570/// Set color index for markers
2571
2576
2577////////////////////////////////////////////////////////////////////////////////
2578/// Set the current color.
2579
2581{
2582 if (color < 0) color = 0;
2583 TColor *col = gROOT->GetColor(color);
2584 if (col)
2585 SetColor(col->GetRed(), col->GetGreen(), col->GetBlue());
2586 else
2587 SetColor(1., 1., 1.);
2588}
2589
2590////////////////////////////////////////////////////////////////////////////////
2591/// Set directly current color (don't go via TColor).
2592
2594{
2595 if (r == fRed && g == fGreen && b == fBlue) return;
2596
2597 fRed = r;
2598 fGreen = g;
2599 fBlue = b;
2600
2601 if (fRed <= 0 && fGreen <= 0 && fBlue <= 0 ) {
2602 PrintFast(6," black");
2603 } else {
2604 if (gStyle->GetColorModelPS()) {
2613 } else {
2614 WriteReal(fRed);
2617 }
2618 PrintFast(2," c");
2619 }
2620}
2621
2622////////////////////////////////////////////////////////////////////////////////
2623/// Set color index for text
2624
2629
2630////////////////////////////////////////////////////////////////////////////////
2631/// Write a string of characters
2632///
2633/// This method writes the string chars into a PostScript file
2634/// at position xx,yy in world coordinates.
2635
2637{
2638 static const char *psfont[31][2] = {
2639 { "Root.PSFont.1", "/Times-Italic" },
2640 { "Root.PSFont.2", "/Times-Bold" },
2641 { "Root.PSFont.3", "/Times-BoldItalic" },
2642 { "Root.PSFont.4", "/Helvetica" },
2643 { "Root.PSFont.5", "/Helvetica-Oblique" },
2644 { "Root.PSFont.6", "/Helvetica-Bold" },
2645 { "Root.PSFont.7", "/Helvetica-BoldOblique" },
2646 { "Root.PSFont.8", "/Courier" },
2647 { "Root.PSFont.9", "/Courier-Oblique" },
2648 { "Root.PSFont.10", "/Courier-Bold" },
2649 { "Root.PSFont.11", "/Courier-BoldOblique" },
2650 { "Root.PSFont.12", "/Symbol" },
2651 { "Root.PSFont.13", "/Times-Roman" },
2652 { "Root.PSFont.14", "/ZapfDingbats" },
2653 { "Root.PSFont.15", "/Symbol" },
2654 { "Root.PSFont.STIXGen", "/STIXGeneral" },
2655 { "Root.PSFont.STIXGenIt", "/STIXGeneral-Italic" },
2656 { "Root.PSFont.STIXGenBd", "/STIXGeneral-Bold" },
2657 { "Root.PSFont.STIXGenBdIt", "/STIXGeneral-BoldItalic" },
2658 { "Root.PSFont.STIXSiz1Sym", "/STIXSize1Symbols" },
2659 { "Root.PSFont.STIXSiz1SymBd", "/STIXSize1Symbols-Bold" },
2660 { "Root.PSFont.STIXSiz2Sym", "/STIXSize2Symbols" },
2661 { "Root.PSFont.STIXSiz2SymBd", "/STIXSize2Symbols-Bold" },
2662 { "Root.PSFont.STIXSiz3Sym", "/STIXSize3Symbols" },
2663 { "Root.PSFont.STIXSiz3SymBd", "/STIXSize3Symbols-Bold" },
2664 { "Root.PSFont.STIXSiz4Sym", "/STIXSize4Symbols" },
2665 { "Root.PSFont.STIXSiz4SymBd", "/STIXSize4Symbols-Bold" },
2666 { "Root.PSFont.STIXSiz5Sym", "/STIXSize5Symbols" },
2667 { "Root.PSFont.ME", "/DroidSansFallback" },
2668 { "Root.PSFont.CJKMing", "/DroidSansFallback" },
2669 { "Root.PSFont.CJKGothic", "/DroidSansFallback" }
2670 };
2671
2672 const Double_t kDEGRAD = TMath::Pi()/180.;
2673 Double_t x = xx;
2674 Double_t y = yy;
2675 if (!gPad) return;
2676
2677 // Compute the font size. Exit if it is 0
2678 // The font size is computed from the TTF size to get exactly the same
2679 // size on the screen and in the PostScript file.
2680 Double_t wh = (Double_t)gPad->XtoPixel(gPad->GetX2());
2681 Double_t hh = (Double_t)gPad->YtoPixel(gPad->GetY1());
2683
2684 if (wh < hh) {
2685 tsize = fTextSize*wh;
2686 Int_t sizeTTF = (Int_t)(tsize*kScale+0.5); // TTF size
2687 ftsize = (sizeTTF*fXsize*gPad->GetAbsWNDC())/wh;
2688 } else {
2689 tsize = fTextSize*hh;
2690 Int_t sizeTTF = (Int_t)(tsize*kScale+0.5); // TTF size
2691 ftsize = (sizeTTF*fYsize*gPad->GetAbsHNDC())/hh;
2692 }
2693 Double_t fontsize = 4*(72*(ftsize)/2.54);
2694 if( fontsize <= 0) return;
2695
2696 Float_t tsizex = gPad->AbsPixeltoX(Int_t(tsize))-gPad->AbsPixeltoX(0);
2697 Float_t tsizey = gPad->AbsPixeltoY(0)-gPad->AbsPixeltoY(Int_t(tsize));
2698
2699 Int_t font = abs(fTextFont)/10;
2700 if( font > 31 || font < 1) font = 1;
2701
2702 // Text color.
2704
2705 // Text alignment.
2706 Int_t txalh = fTextAlign/10;
2707 if (txalh <1) txalh = 1; else if (txalh > 3) txalh = 3;
2708 Int_t txalv = fTextAlign%10;
2709 if (txalv <1) txalv = 1; else if (txalv > 3) txalv = 3;
2710 if (txalv == 3) {
2713 } else if (txalv == 2) {
2716 }
2717
2718 UInt_t w = 0, w0 = 0;
2720 // In order to measure the precise character positions we need to trick
2721 // FreeType into rendering high-resolution characters otherwise it will
2722 // stick to the screen pixel grid, which is far worse than we can achieve
2723 // on print.
2724 const Float_t scale = 16.0;
2725 // Save current text attributes.
2727 saveAttText.TAttText::operator=(*this);
2728 const Int_t len=strlen(chars);
2729 Int_t *charWidthsCumul = nullptr;
2730 TText t;
2735 t.TAttText::Modify();
2736 if (w0-w != 0) kerning = kTRUE;
2737 else kerning = kFALSE;
2738 if (kerning) {
2739 // Calculate the individual character placements.
2740 charWidthsCumul = new Int_t[len];
2741 for (Int_t i = len - 1;i >= 0;i--) {
2742 UInt_t ww = 0;
2743 t.GetTextAdvance(ww, chars + i);
2744 Double_t wwl = (gPad->AbsPixeltoX(ww)-gPad->AbsPixeltoX(0));
2745 charWidthsCumul[i] = (Int_t)((XtoPS(wwl) - XtoPS(0)) / scale);
2746 }
2747 }
2748 // Restore text attributes.
2749 saveAttText.TAttText::Modify();
2750
2751 Double_t charsLength = gPad->AbsPixeltoX(w)-gPad->AbsPixeltoX(0);
2753
2754 // Text angle.
2755 Int_t psangle = Int_t(0.5 + fTextAngle);
2756
2757 // Save context.
2758 PrintStr("@");
2759 SaveRestore(1);
2760
2761 // Clipping
2762 Int_t xc1 = XtoPS(gPad->GetX1());
2763 Int_t xc2 = XtoPS(gPad->GetX2());
2764 Int_t yc1 = YtoPS(gPad->GetY1());
2765 Int_t yc2 = YtoPS(gPad->GetY2());
2766 WriteInteger(xc2 - xc1);
2767 WriteInteger(yc2 - yc1);
2770 PrintStr(" C");
2771
2772 // Output text position and angle. The text position is computed
2773 // using Double_t to avoid precision problems.
2774 Double_t vx = (x - gPad->GetX1())/(gPad->GetX2()-gPad->GetX1());
2775 Double_t cmx = fXsize*(gPad->GetAbsXlowNDC()+vx*gPad->GetAbsWNDC());
2776 WriteReal((288.*cmx)/2.54);
2777 Double_t vy = (y - gPad->GetY1())/(gPad->GetY2()-gPad->GetY1());
2778 Double_t cmy = fYsize*(gPad->GetAbsYlowNDC()+vy*gPad->GetAbsHNDC());
2779 WriteReal((288.*cmy)/2.54);
2780 PrintStr(TString::Format(" t %d r ", psangle));
2781 if(txalh == 2) PrintStr(TString::Format(" %d 0 t ", -psCharsLength/2));
2782 if(txalh == 3) PrintStr(TString::Format(" %d 0 t ", -psCharsLength));
2783 PrintStr(gEnv->GetValue(psfont[font-1][0], psfont[font-1][1]));
2784 if (font != 15) {
2785 PrintStr(TString::Format(" findfont %g sf 0 0 m ",fontsize));
2786 } else {
2787 PrintStr(TString::Format(" findfont %g sf 0 0 m ita ",fontsize));
2788 }
2789
2790 if (kerning) {
2791 PrintStr("@");
2792 // Output individual character placements
2793 for (Int_t i = len-1; i >= 1; i--) {
2795 }
2796 delete [] charWidthsCumul;
2797 PrintStr("@");
2798 }
2799
2800 // Output text.
2801 PrintStr("(");
2802
2803 // Inside a PostScript string, the new line (if needed to break up long lines) must be escaped by a backslash.
2804 const char *crsave = fImplicitCREsc;
2805 fImplicitCREsc = "\\";
2806
2807 char str[8];
2808 for (Int_t i=0; i<len;i++) {
2809 if (chars[i]!='\n') {
2810 if (chars[i]=='(' || chars[i]==')' || chars[i]=='\\') {
2811 snprintf(str,8,"\\%c",chars[i]);
2812 PrintStr(str);
2813 } else if ((chars[i]=='-') && (font != 12)) {
2814 PrintStr("\\255");
2815 } else {
2816 snprintf(str,8,"%c",chars[i]);
2817 PrintFast(1,str);
2818 }
2819 }
2820 }
2821 PrintStr(")");
2823
2824 if (kerning) {
2825 if (font != 15) PrintStr(" K NC");
2826 else PrintStr(" K gr NC");
2827 } else {
2828 if (font != 15) PrintStr(" show NC");
2829 else PrintStr(" show gr NC");
2830 }
2831
2832 SaveRestore(-1);
2833}
2834
2835////////////////////////////////////////////////////////////////////////////////
2836/// Write a string of characters
2837///
2838/// This method writes the string chars into a PostScript file
2839/// at position xx,yy in world coordinates.
2840
2842{
2843 static const char *psfont[31][2] = {
2844 { "Root.PSFont.1", "/FreeSerifItalic" },
2845 { "Root.PSFont.2", "/FreeSerifBold" },
2846 { "Root.PSFont.3", "/FreeSerifBoldItalic" },
2847 { "Root.PSFont.4", "/FreeSans" },
2848 { "Root.PSFont.5", "/FreeSansOblique" },
2849 { "Root.PSFont.6", "/FreeSansBold" },
2850 { "Root.PSFont.7", "/FreeSansBoldOblique" },
2851 { "Root.PSFont.8", "/FreeMono" },
2852 { "Root.PSFont.9", "/FreeMonoOblique" },
2853 { "Root.PSFont.10", "/FreeMonoBold" },
2854 { "Root.PSFont.11", "/FreeMonoBoldOblique" },
2855 { "Root.PSFont.12", "/SymbolMT" },
2856 { "Root.PSFont.13", "/FreeSerif" },
2857 { "Root.PSFont.14", "/Wingdings-Regular" },
2858 { "Root.PSFont.15", "/SymbolMT" },
2859 { "Root.PSFont.STIXGen", "/STIXGeneral" },
2860 { "Root.PSFont.STIXGenIt", "/STIXGeneral-Italic" },
2861 { "Root.PSFont.STIXGenBd", "/STIXGeneral-Bold" },
2862 { "Root.PSFont.STIXGenBdIt", "/STIXGeneral-BoldItalic" },
2863 { "Root.PSFont.STIXSiz1Sym", "/STIXSize1Symbols" },
2864 { "Root.PSFont.STIXSiz1SymBd", "/STIXSize1Symbols-Bold" },
2865 { "Root.PSFont.STIXSiz2Sym", "/STIXSize2Symbols" },
2866 { "Root.PSFont.STIXSiz2SymBd", "/STIXSize2Symbols-Bold" },
2867 { "Root.PSFont.STIXSiz3Sym", "/STIXSize3Symbols" },
2868 { "Root.PSFont.STIXSiz3SymBd", "/STIXSize3Symbols-Bold" },
2869 { "Root.PSFont.STIXSiz4Sym", "/STIXSize4Symbols" },
2870 { "Root.PSFont.STIXSiz4SymBd", "/STIXSize4Symbols-Bold" },
2871 { "Root.PSFont.STIXSiz5Sym", "/STIXSize5Symbols" },
2872 { "Root.PSFont.ME", "/DroidSansFallback" },
2873 { "Root.PSFont.CJKMing", "/DroidSansFallback" },
2874 { "Root.PSFont.CJKGothic", "/DroidSansFallback" }
2875 };
2876
2877 Int_t len = wcslen(chars);
2878 if (len<=0) return;
2879
2880 const Double_t kDEGRAD = TMath::Pi()/180.;
2881 Double_t x = xx;
2882 Double_t y = yy;
2883 if (!gPad) return;
2884
2885 // Compute the font size. Exit if it is 0
2886 // The font size is computed from the TTF size to get exactly the same
2887 // size on the screen and in the PostScript file.
2888 Double_t wh = (Double_t)gPad->XtoPixel(gPad->GetX2());
2889 Double_t hh = (Double_t)gPad->YtoPixel(gPad->GetY1());
2891
2892 if (wh < hh) {
2893 tsize = fTextSize*wh;
2894 Int_t sizeTTF = (Int_t)(tsize*kScale+0.5); // TTF size
2895 ftsize = (sizeTTF*fXsize*gPad->GetAbsWNDC())/wh;
2896 } else {
2897 tsize = fTextSize*hh;
2898 Int_t sizeTTF = (Int_t)(tsize*kScale+0.5); // TTF size
2899 ftsize = (sizeTTF*fYsize*gPad->GetAbsHNDC())/hh;
2900 }
2901 Double_t fontsize = 4*(72*(ftsize)/2.54);
2902 if( fontsize <= 0) return;
2903
2904 Float_t tsizex = gPad->AbsPixeltoX(Int_t(tsize))-gPad->AbsPixeltoX(0);
2905 Float_t tsizey = gPad->AbsPixeltoY(0)-gPad->AbsPixeltoY(Int_t(tsize));
2906
2907 Int_t font = abs(fTextFont)/10;
2908 if( font > 29 || font < 1) font = 1;
2909
2910 // Text color.
2912
2913 // Text alignment.
2914 Int_t txalh = fTextAlign/10;
2915 if (txalh <1) txalh = 1; else if (txalh > 3) txalh = 3;
2916 Int_t txalv = fTextAlign%10;
2917 if (txalv <1) txalv = 1; else if (txalv > 3) txalv = 3;
2918 if (txalv == 3) {
2921 } else if (txalv == 2) {
2924 }
2925 UInt_t w = 0, h = 0;
2926
2927 TText t;
2930 t.GetTextExtent(w, h, chars);
2931 Double_t charsLength = gPad->AbsPixeltoX(w)-gPad->AbsPixeltoX(0);
2933
2934 // Text angle.
2935 Int_t psangle = Int_t(0.5 + fTextAngle);
2936
2937 // Save context.
2938 PrintStr("@");
2939 SaveRestore(1);
2940
2941 // Clipping
2942 Int_t xc1 = XtoPS(gPad->GetX1());
2943 Int_t xc2 = XtoPS(gPad->GetX2());
2944 Int_t yc1 = YtoPS(gPad->GetY1());
2945 Int_t yc2 = YtoPS(gPad->GetY2());
2946 WriteInteger(xc2 - xc1);
2947 WriteInteger(yc2 - yc1);
2950 PrintStr(" C");
2951
2952 // Output text position and angle.
2955 PrintStr(TString::Format(" t %d r ", psangle));
2956 if(txalh == 2) PrintStr(TString::Format(" %d 0 t ", -psCharsLength/2));
2957 if(txalh == 3) PrintStr(TString::Format(" %d 0 t ", -psCharsLength));
2958 MustEmbed[font-1] = kTRUE; // This font will be embedded in the file at EOF time.
2959 PrintStr(gEnv->GetValue(psfont[font-1][0], psfont[font-1][1]));
2960 PrintStr(TString::Format(" findfont %g sf 0 0 m ",fontsize));
2961
2962 // Output text.
2963 if (len > 1) PrintStr(TString::Format("%d ", len));
2964 for(Int_t i = 0; i < len; i++) {
2965 // Adobe Glyph Naming Convention
2966 // http://www.adobe.com/devnet/opentype/archives/glyph.html
2967#include "AdobeGlyphList.h"
2968 const wchar_t *lower = std::lower_bound(
2970 chars[i]);
2972 *lower == chars[i]) {
2973 // Named glyph from AGL 1.2
2974 const unsigned long index =
2977 }
2978 else if((unsigned int)chars[i] < 0xffff) {
2979 // Unicode BMP
2980 PrintStr(TString::Format("/uni%04X ",
2981 (unsigned int)chars[i]));
2982 }
2983 else {
2984 // Unicode supplemental planes
2985 PrintStr(TString::Format("/u%04X ",
2986 (unsigned int)chars[i]));
2987 }
2988 }
2989 if(len > 1) {
2990 PrintStr("{glyphshow} repeat ");
2991 }
2992 else {
2993 PrintStr("glyphshow ");
2994 }
2995
2996 PrintStr("NC");
2997
2998 SaveRestore(-1);
2999}
3000
3001////////////////////////////////////////////////////////////////////////////////
3002/// Draw text with URL. Same as Text.
3003///
3004
3005void TPostScript::TextUrl(Double_t x, Double_t y, const char *chars, const char *)
3006{
3007 Text(x, y, chars);
3008}
3009
3010////////////////////////////////////////////////////////////////////////////////
3011/// Write a string of characters in NDC
3012
3014{
3015 Double_t x = gPad->GetX1() + u*(gPad->GetX2() - gPad->GetX1());
3016 Double_t y = gPad->GetY1() + v*(gPad->GetY2() - gPad->GetY1());
3017 Text(x, y, chars);
3018}
3019
3020////////////////////////////////////////////////////////////////////////////////
3021/// Write a string of characters in NDC
3022
3024{
3025 Double_t x = gPad->GetX1() + u*(gPad->GetX2() - gPad->GetX1());
3026 Double_t y = gPad->GetY1() + v*(gPad->GetY2() - gPad->GetY1());
3027 Text(x, y, chars);
3028}
3029
3030////////////////////////////////////////////////////////////////////////////////
3031/// Convert U from NDC coordinate to PostScript
3032
3034{
3035 Double_t cm = fXsize*(gPad->GetAbsXlowNDC() + u*gPad->GetAbsWNDC());
3036 return Int_t(0.5 + 288*cm/2.54);
3037}
3038
3039////////////////////////////////////////////////////////////////////////////////
3040/// Convert V from NDC coordinate to PostScript
3041
3043{
3044 Double_t cm = fYsize*(gPad->GetAbsYlowNDC() + v*gPad->GetAbsHNDC());
3045 return Int_t(0.5 + 288*cm/2.54);
3046}
3047
3048////////////////////////////////////////////////////////////////////////////////
3049/// Convert X from world coordinate to PostScript
3050
3052{
3053 Double_t u = (x - gPad->GetX1())/(gPad->GetX2() - gPad->GetX1());
3054 return UtoPS(u);
3055}
3056
3057////////////////////////////////////////////////////////////////////////////////
3058/// Convert Y from world coordinate to PostScript
3059
3061{
3062 Double_t v = (y - gPad->GetY1())/(gPad->GetY2() - gPad->GetY1());
3063 return VtoPS(v);
3064}
3065
3066////////////////////////////////////////////////////////////////////////////////
3067/// Initialize the PostScript page in zones
3068
3070{
3071 if( !fClear )return;
3072 fClear = kFALSE;
3073
3074 // When Zone has been called, fZone is TRUE
3075 fZone = kTRUE;
3076
3077 if( fIYzone > fNYzone) {
3078 fIYzone=1;
3079 if( fMode != 3) {
3080 PrintStr("@showpage");
3081 SaveRestore(-1);
3082 fNpages++;
3083 PrintStr("@%%Page:");
3086 PrintStr("@");
3087 } else {
3088 PrintFast(9," showpage");
3089 SaveRestore(-1);
3090 }
3091 }
3092
3093 // No grestore the first time
3094 if( fMode != 3) {
3095 if( fIXzone != 1 || fIYzone != 1) SaveRestore(-1);
3096 SaveRestore(1);
3097 PrintStr("@");
3100 PrintFast(5," Zone");
3101 PrintStr("@");
3102 fIXzone++;
3103 if( fIXzone > fNXzone) { fIXzone=1; fIYzone++; }
3104 }
3105
3106 // Picture Initialisation
3107 SaveRestore(1);
3108 if (fgLineJoin) {
3110 PrintFast(12," setlinejoin");
3111 }
3112 if (fgLineCap) {
3114 PrintFast(11," setlinecap");
3115 }
3116 PrintFast(6," 0 0 t");
3117 fRed = -1;
3118 fGreen = -1;
3119 fBlue = -1;
3120 fPrinted = kFALSE;
3121 fLineColor = -1;
3122 fLineStyle = -1;
3123 fLineWidth = -1;
3124 fFillColor = -1;
3125 fFillStyle = -1;
3126 fMarkerSizeCur = -1;
3127}
static const wchar_t adobe_glyph_ucs[nadobe_glyph]
static const unsigned long nadobe_glyph
static const char * adobe_glyph_name[nadobe_glyph]
#define b(i)
Definition RSha256.hxx:100
#define g(i)
Definition RSha256.hxx:105
#define h(i)
Definition RSha256.hxx:106
short Style_t
Style number (short)
Definition RtypesCore.h:96
bool Bool_t
Boolean (0=false, 1=true) (bool)
Definition RtypesCore.h:77
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
short Color_t
Color number (short)
Definition RtypesCore.h:99
short Width_t
Line width (short)
Definition RtypesCore.h:98
float Float_t
Float 4 bytes (float)
Definition RtypesCore.h:71
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
double Double_t
Double 8 bytes.
Definition RtypesCore.h:73
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
const Float_t kScale
Definition TASImage.cxx:135
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
R__EXTERN TEnv * gEnv
Definition TEnv.h:170
Option_t Option_t cindex
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char filename
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t np
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
Option_t Option_t markerstyle
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t UChar_t len
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint xy
Option_t Option_t TPoint TPoint const char mode
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t format
Option_t Option_t width
Option_t Option_t TPoint TPoint const char y1
const Float_t kScale
static Bool_t MustEmbed[32]
#define gROOT
Definition TROOT.h:426
R__EXTERN TStyle * gStyle
Definition TStyle.h:442
@ kReadPermission
Definition TSystem.h:55
@ kWritePermission
Definition TSystem.h:54
R__EXTERN TSystem * gSystem
Definition TSystem.h:582
R__EXTERN TVirtualPS * gVirtualPS
Definition TVirtualPS.h:84
#define gPad
#define snprintf
Definition civetweb.c:1579
virtual Color_t GetFillColor() const
Return the fill area color.
Definition TAttFill.h:32
Style_t fFillStyle
Fill area style.
Definition TAttFill.h:25
Color_t fFillColor
Fill area color.
Definition TAttFill.h:24
Width_t fLineWidth
Line width.
Definition TAttLine.h:26
Style_t fLineStyle
Line style.
Definition TAttLine.h:25
Color_t fLineColor
Line color.
Definition TAttLine.h:24
Color_t fMarkerColor
Marker color.
Definition TAttMarker.h:24
static Width_t GetMarkerLineWidth(Style_t style)
Internal helper function that returns the line width of the given marker style (0 = filled marker)
Size_t fMarkerSize
Marker size.
Definition TAttMarker.h:26
Style_t fMarkerStyle
Marker style.
Definition TAttMarker.h:25
static Style_t GetMarkerStyleBase(Style_t style)
Internal helper function that returns the corresponding marker style with line width 1 for the given ...
Color_t fTextColor
Text color.
Definition TAttText.h:27
Float_t fTextAngle
Text angle.
Definition TAttText.h:24
virtual void SetTextFont(Font_t tfont=62)
Set the text font.
Definition TAttText.h:50
Font_t fTextFont
Text font.
Definition TAttText.h:28
virtual void SetTextSize(Float_t tsize=1)
Set the text size.
Definition TAttText.h:51
Short_t fTextAlign
Text alignment.
Definition TAttText.h:26
Float_t fTextSize
Text size.
Definition TAttText.h:25
The color creation and management class.
Definition TColor.h:22
Float_t GetRed() const
Definition TColor.h:61
static Int_t GetColor(const char *hexcolor)
Static method returning color number for color specified by hex color string of form: "#rrggbb",...
Definition TColor.cxx:1926
Float_t GetBlue() const
Definition TColor.h:63
Float_t GetGreen() const
Definition TColor.h:62
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:503
virtual void SetTitle(const char *title="")
Set the title of the TNamed.
Definition TNamed.cxx:173
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
virtual void Warning(const char *method, const char *msgfmt,...) const
Issue warning message.
Definition TObject.cxx:1081
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1095
2-D graphics point (world coordinates).
Definition TPoints.h:19
void DefineMarkers()
Define the markers.
Float_t fXsize
Page size along X.
Definition TPostScript.h:45
Int_t fIYzone
Current zone along Y.
Definition TPostScript.h:58
Int_t fType
PostScript workstation type.
Definition TPostScript.h:61
void DrawPolyMarker(Int_t n, Float_t *x, Float_t *y) override
Draw markers at the n WC points x, y.
TPostScript()
Default PostScript constructor.
Int_t VtoPS(Double_t v)
Convert V from NDC coordinate to PostScript.
void CellArrayEnd() override
End the Cell Array painting.
static Int_t fgLineCap
Appearance of line caps.
Definition TPostScript.h:82
Float_t fX2w
Definition TPostScript.h:29
void SetMarkerColor(Color_t cindex=1) override
Set color index for markers.
Int_t fNBSameColorCell
Number of boxes with the same color.
Definition TPostScript.h:77
void SetColor(Int_t color=1)
Set the current color.
Int_t fMode
PostScript mode.
Definition TPostScript.h:62
Int_t fIXzone
Current zone along X.
Definition TPostScript.h:57
void DrawPolyLineNDC(Int_t n, TPoints *uv)
Draw a PolyLine in NDC space.
void FontEncode()
Font Re-encoding.
Float_t fXVP2
Definition TPostScript.h:38
Int_t fLastCellBlue
Last blue value.
Definition TPostScript.h:76
Float_t fY1w
Definition TPostScript.h:28
Int_t fNbCellW
Number of boxes per line.
Definition TPostScript.h:71
Float_t fX1v
X bottom left corner of paper.
Definition TPostScript.h:23
Float_t fWidth
Current line width.
Definition TPostScript.h:51
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
void DrawPS(Int_t n, Float_t *xw, Float_t *yw) override
Draw a PolyLine.
Float_t fXVS1
Definition TPostScript.h:41
Bool_t fClipStatus
Clipping Indicator.
Definition TPostScript.h:66
Float_t fLineScale
Line width scale factor.
Definition TPostScript.h:53
Float_t fMarkerSizeCur
current transformed value of marker size
Definition TPostScript.h:59
Bool_t fZone
Zone indicator.
Definition TPostScript.h:68
void Open(const char *filename, Int_t type=-111) override
Open a PostScript file.
void SetFillPatterns(Int_t ipat, Int_t color)
Patterns definition.
Float_t fDXC
Definition TPostScript.h:31
Int_t fLastCellRed
Last red value.
Definition TPostScript.h:74
void SetLineWidth(Width_t linewidth=1) override
Change the line width.
void Off()
Deactivate an already open PostScript file.
Int_t YtoPS(Double_t y)
Convert Y from world coordinate to PostScript.
Bool_t fRange
True when a range has been defined.
Definition TPostScript.h:67
Float_t fBlue
Per cent of blue.
Definition TPostScript.h:50
void Text(Double_t x, Double_t y, const char *string) override
Write a string of characters.
void MovePS(Int_t x, Int_t y)
Move to a new position.
void Initialize()
PostScript Initialisation.
bool FontEmbedType2(const char *filename)
Int_t fNbCellLine
Number of boxes in the current line.
Definition TPostScript.h:72
Float_t fMaxsize
Largest dimension of X and Y.
Definition TPostScript.h:47
Float_t fYVS1
Definition TPostScript.h:43
void SetStyle(Style_t linestyle=1)
Change the line style in the output file.
void SetFillColor(Color_t cindex=1) override
Set color index for fill areas.
Int_t fClip
Clipping mode.
Definition TPostScript.h:63
TString fFileName
PS file name.
Definition TPostScript.h:78
Int_t fNYzone
Number of zones along Y.
Definition TPostScript.h:56
bool FontEmbedType1(const char *filename)
Float_t fGreen
Per cent of green.
Definition TPostScript.h:49
void SetLineColor(Color_t cindex=1) override
Set color index for lines.
void SetLineStyle(Style_t linestyle=1) override
Change the line style.
void SetLineCap(Int_t linecap=0)
Set the value of the global parameter TPostScript::fgLineCap.
Float_t fYVP2
Definition TPostScript.h:40
void NewPage() override
Move to a new PostScript page.
Float_t fDYC
Definition TPostScript.h:32
void DrawPolyLine(Int_t n, TPoints *xy)
Draw a PolyLine.
void DrawBox(Double_t x1, Double_t y1, Double_t x2, Double_t y2) override
Draw a Box.
void DrawHatch(Float_t dy, Float_t angle, Int_t n, Float_t *x, Float_t *y)
Draw Fill area with hatch styles.
Int_t UtoPS(Double_t u)
Convert U from NDC coordinate to PostScript.
Bool_t fFontEmbed
True is FontEmbed has been called.
Definition TPostScript.h:79
void SetWidth(Width_t linewidth=1)
Change the line width in the output file.
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) override
Draw a Frame around a box.
Float_t fX1w
Definition TPostScript.h:27
Float_t fYVP1
Definition TPostScript.h:39
Int_t fSave
Number of gsave for restore.
Definition TPostScript.h:54
Int_t fStyle
Current line style.
Definition TPostScript.h:52
void CellArrayFill(Int_t r, Int_t g, Int_t b) override
Paint the Cell Array.
void SetTextColor(Color_t cindex=1) override
Set color index for text.
void On()
Activate an already open PostScript file.
Int_t CMtoPS(Double_t u)
Definition TPostScript.h:94
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:55
void TextUrl(Double_t x, Double_t y, const char *string, const char *url) override
Draw text with URL.
void SetLineScale(Float_t scale=3)
Bool_t fBoundingBox
True for Encapsulated PostScript.
Definition TPostScript.h:64
Float_t fY2w
Definition TPostScript.h:30
Int_t XtoPS(Double_t x)
Convert X from world coordinate to PostScript.
~TPostScript() override
Default PostScript destructor.
char fPatterns[32]
Indicate if pattern n is defined.
Definition TPostScript.h:69
Int_t fNbinCT
Number of entries in the current Cell Array.
Definition TPostScript.h:70
static Int_t fgLineJoin
Appearance of joining lines.
Definition TPostScript.h:81
Float_t fY2v
Y top right corner of paper.
Definition TPostScript.h:26
Int_t fLastCellGreen
Last green value.
Definition TPostScript.h:75
Float_t fFY
Definition TPostScript.h:36
Float_t fY1v
Y bottom left corner of paper.
Definition TPostScript.h:24
Float_t fYVS2
Definition TPostScript.h:44
void SaveRestore(Int_t flag)
Compute number of gsaves for restore This allows to write the correct number of grestore at the end o...
Bool_t fClear
True when page must be cleared.
Definition TPostScript.h:65
Int_t fMaxLines
Maximum number of lines in a PS array.
Definition TPostScript.h:73
void Zone()
Initialize the PostScript page in zones.
void SetLineJoin(Int_t linejoin=0)
Set the value of the global parameter TPostScript::fgLineJoin.
void Range(Float_t xrange, Float_t yrange)
Set the range for the paper in centimeters.
void CellArrayBegin(Int_t W, Int_t H, Double_t x1, Double_t x2, Double_t y1, Double_t y2) override
Draw a Cell Array.
Float_t fX2v
X top right corner of paper.
Definition TPostScript.h:25
void FontEmbed()
Embed font in PS file.
void Close(Option_t *opt="") override
Close a PostScript file.
Float_t fXVP1
Definition TPostScript.h:37
Int_t fNpages
number of pages
Definition TPostScript.h:60
Float_t fYsize
Page size along Y.
Definition TPostScript.h:46
static const TString & GetTTFFontDir()
Get the fonts directory in the installation. Static utility function.
Definition TROOT.cxx:3527
Basic string class.
Definition TString.h:138
const char * Data() const
Definition TString.h:384
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2384
Int_t GetJoinLinePS() const
Returns the line join method used for PostScript, PDF and SVG output. See TPostScript::SetLineJoin fo...
Definition TStyle.h:289
Int_t GetColorModelPS() const
Definition TStyle.h:198
const char * GetLineStyleString(Int_t i=1) const
Return line style string (used by PostScript).
Definition TStyle.cxx:1167
const char * GetTitlePS() const
Definition TStyle.h:287
Int_t GetCapLinePS() const
Returns the line cap method used for PostScript, PDF and SVG output. See TPostScript::SetLineCap for ...
Definition TStyle.h:290
void GetPaperSize(Float_t &xsize, Float_t &ysize) const
Set paper size for PostScript output.
Definition TStyle.cxx:1184
const char * GetHeaderPS() const
Definition TStyle.h:286
Float_t GetLineScalePS() const
Definition TStyle.h:291
virtual int GetPid()
Get process id.
Definition TSystem.cxx:716
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:1307
virtual int Rename(const char *from, const char *to)
Rename a file.
Definition TSystem.cxx:1361
virtual char * Which(const char *search, const char *file, EAccessMode mode=kFileExists)
Find location of file in a search path.
Definition TSystem.cxx:1559
virtual int Unlink(const char *name)
Unlink, i.e.
Definition TSystem.cxx:1392
Base class for several text objects.
Definition TText.h:22
virtual void GetTextExtent(UInt_t &w, UInt_t &h, const char *text) const
Return text extent for string text.
Definition TText.cxx:590
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:619
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.
virtual void PrintRaw(Int_t len, const char *str)
Print a raw.
virtual void PrintStr(const char *string="")
Output the string str in the output buffer.
virtual void PrintFast(Int_t nch, const char *string="")
Fast version of Print.
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.
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
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:249
Double_t Floor(Double_t x)
Rounds x downward, returning the largest integral value that is not greater than x.
Definition TMath.h:691
T1 Sign(T1 a, T2 b)
Returns a value with the magnitude of a and the sign of b.
Definition TMathBase.h:174
Double_t Sqrt(Double_t x)
Returns the square root of x.
Definition TMath.h:673
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Returns x raised to the power y.
Definition TMath.h:732
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:197
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Definition TMath.h:605
constexpr Double_t Pi()
Definition TMath.h:40
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:599
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:122