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
404 Error("Open", "Cannot open file: %s", fFileName.Data());
405 return;
406 }
407 gVirtualPS = this;
408
409 for (Int_t i=0;i<fSizBuffer;i++) fBuffer[i] = ' ';
410 if( fType == 113) {
412 PrintStr("%!PS-Adobe-2.0 EPSF-2.0@");
413 } else {
415 PrintStr("%!PS-Adobe-2.0@");
416 Initialize();
417 }
418
420 fRange = kFALSE;
421
422 // Set a default range
424
426 if (fType == 113) NewPage();
427}
428
429////////////////////////////////////////////////////////////////////////////////
430/// Default PostScript destructor
431
436
437////////////////////////////////////////////////////////////////////////////////
438/// Close a PostScript file
439
441{
442 if (!gVirtualPS ||!fStream)
443 return;
444 if (gPad)
445 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 CloseStream();
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
480 Error("Close", "Cannot open file: %s", fFileName.Data());
481 return;
482 }
483
484 // Embed the fonts at the right place
485 FILE *sg = fopen(tmpname.Data(),"r");
486 if (!sg) {
487 Error("Close", "Cannot open file: %s", tmpname.Data());
488 return;
489 }
490 char line[255];
491 while (fgets(line,255,sg)) {
492 if (strstr(line,"EndComments")) PrintStr("%%DocumentNeededResources: ProcSet (FontSetInit)@");
493 fStream->write(line,strlen(line));
494 if (!fFontEmbed && strstr(line,"m5")) {
495 FontEmbed();
496 PrintStr("@");
497 }
498 }
499 fclose(sg);
500 if (gSystem->Unlink(tmpname.Data())) return;
501 }
502
504
505 // Close file stream
506
507 CloseStream();
508
509 gVirtualPS = nullptr;
510}
511
512////////////////////////////////////////////////////////////////////////////////
513/// Activate an already open PostScript file
514
516{
517 if (!fType) {
518 Error("On", "no postscript file open");
519 Off();
520 return;
521 }
522 gVirtualPS = this;
523}
524
525////////////////////////////////////////////////////////////////////////////////
526/// Deactivate an already open PostScript file
527
529{
530 gVirtualPS = nullptr;
531}
532
533////////////////////////////////////////////////////////////////////////////////
534/// Draw a Cell Array
535///
536/// Drawing a PostScript Cell Array is in fact done thanks to three
537/// procedures: CellArrayBegin, CellArrayFill, and CellArrayEnd.
538///
539/// - CellArrayBegin: Initiate the Cell Array by writing the necessary
540/// PostScript procedures and the initial values of the
541/// required parameters. The input parameters are:
542/// - W: number of boxes along the width.
543/// - H: number of boxes along the height
544/// - x1,x2,y1,y2: First box coordinates.
545/// - CellArrayFill: Is called for each box of the Cell Array. The first
546/// box is the top left one and the last box is the
547/// bottom right one. The input parameters are the Red,
548/// Green, and Blue components of the box colour. These
549/// Levels are between 0 and 255.
550/// - CellArrayEnd: Finishes the Cell Array.
551///
552/// PostScript cannot handle arrays larger than 65535. So the Cell Array
553/// is drawn in several pieces.
554
557{
558 Int_t ix1 = XtoPS(x1);
559 Int_t iy1 = YtoPS(y1);
560
561 Float_t wt = (288/2.54)*gPad->GetAbsWNDC()*
562 fXsize*((x2 - x1)/(gPad->GetX2()-gPad->GetX1()));
563 Float_t ht = (288/2.54)*gPad->GetAbsHNDC()*
564 fYsize*((y2 - y1)/(gPad->GetY2()-gPad->GetY1()));
565
566 fLastCellRed = 300;
567 fLastCellGreen = 300;
568 fLastCellBlue = 300;
570
571 fNbinCT = 0;
572 fNbCellW = W;
573 fNbCellLine = 0;
574 fMaxLines = 40000/(3*fNbCellW);
575
576 // Define some parameters
577 PrintStr("@/WT"); WriteReal(wt) ; PrintStr(" def"); // Cells width
578 PrintStr(" /HT"); WriteReal(ht) ; PrintStr(" def"); // Cells height
579 PrintStr(" /XS"); WriteInteger(ix1) ; PrintStr(" def"); // X start
580 PrintStr(" /YY"); WriteInteger(iy1) ; PrintStr(" def"); // Y start
581 PrintStr(" /NX"); WriteInteger(W) ; PrintStr(" def"); // Number of columns
582 PrintStr(" /NY"); WriteInteger(fMaxLines); PrintStr(" def"); // Number of lines
583
584 // This PS procedure draws one cell.
585 PrintStr(" /DrawCell ");
586 PrintStr( "{WT HT XX YY bf");
587 PrintStr( " /NBBD NBBD 1 add def");
588 PrintStr( " NBBD NBB eq {exit} if");
589 PrintStr( " /XX WT XX add def");
590 PrintStr( " IX NX eq ");
591 PrintStr( "{/YY YY HT sub def");
592 PrintStr( " /XX XS def");
593 PrintStr( " /IX 0 def} if");
594 PrintStr( " /IX IX 1 add def} def");
595
596 // This PS procedure draws fMaxLines line. It takes care of duplicated
597 // colors. Values "n" greater than 300 mean than the previous color
598 // should be duplicated n-300 times.
599 PrintStr(" /DrawCT ");
600 PrintStr( "{/NBB NX NY mul def");
601 PrintStr( " /XX XS def");
602 PrintStr( " /IX 1 def");
603 PrintStr( " /NBBD 0 def");
604 PrintStr( " /RC 0 def /GC 1 def /BC 2 def");
605 PrintStr( " 1 1 NBB ");
606 PrintStr( "{/NB CT RC get def");
607 PrintStr( " NB 301 ge ");
608 PrintStr( "{/NBL NB 300 sub def");
609 PrintStr( " 1 1 NBL ");
610 PrintStr( "{DrawCell}");
611 PrintStr( " for");
612 PrintStr( " /RC RC 1 add def");
613 PrintStr( " /GC RC 1 add def");
614 PrintStr( " /BC RC 2 add def}");
615 PrintStr( "{CT RC get 255 div CT GC get 255 div CT BC get 255 div setrgbcolor");
616 PrintStr( " DrawCell");
617 PrintStr( " /RC RC 3 add def");
618 PrintStr( " /GC GC 3 add def");
619 PrintStr( " /BC BC 3 add def} ifelse NBBD NBB eq {exit} if} for");
620 PrintStr( " /YY YY HT sub def clear} def");
621
622 PrintStr(" /CT [");
623}
624
625////////////////////////////////////////////////////////////////////////////////
626/// Paint the Cell Array
627
629{
630 if (fLastCellRed == r && fLastCellGreen == g && fLastCellBlue == b) {
632 } else {
633 if (fNBSameColorCell != 0 ) {
636 }
640 fLastCellRed = r;
643 }
644
645 fNbinCT++;
646 if (fNbinCT == fNbCellW) {
647 fNbCellLine++;
648 fNbinCT = 0;
649 }
650
651 if (fNbCellLine == fMaxLines) {
653 PrintStr("] def DrawCT /CT [");
654 fNbCellLine = 0;
655 fLastCellRed = 300;
656 fLastCellGreen = 300;
657 fLastCellBlue = 300;
659 fNbinCT = 0;
660 }
661}
662
663////////////////////////////////////////////////////////////////////////////////
664/// End the Cell Array painting
665
667{
669 PrintStr("] def /NY");
671 PrintStr(" def DrawCT ");
672}
673
674////////////////////////////////////////////////////////////////////////////////
675/// Define the markers
676
678{
679 PrintStr("/mp {newpath /y exch def /x exch def} def@");
680 PrintStr("/side {[w .77 mul w .23 mul] .385 w mul sd w 0 l currentpoint t -144 r} def@");
681 PrintStr("/mr {mp x y w2 0 360 arc} def /m24 {mr s} def /m20 {mr f} def@");
682 PrintStr("/mb {mp x y w2 add m w2 neg 0 d 0 w neg d w 0 d 0 w d cl} def@");
683 PrintStr("/mt {mp x y w2 add m w2 neg w neg d w 0 d cl} def@");
684 PrintStr("/w4 {w 4 div} def@");
685 PrintStr("/w6 {w 6 div} def@");
686 PrintStr("/w8 {w 8 div} def@");
687 PrintStr("/m21 {mb f} def /m25 {mb s} def /m22 {mt f} def /m26{mt s} def@");
688 PrintStr("/m23 {mp x y w2 sub m w2 w d w neg 0 d cl f} def@");
689 PrintStr("/m27 {mp x y w2 add m w3 neg w2 neg d w3 w2 neg d w3 w2 d cl s} def@");
690 PrintStr("/m28 {mp x w2 sub y w2 sub w3 add m w3 0 d ");
691 PrintStr(" 0 w3 neg d w3 0 d 0 w3 d w3 0 d ");
692 PrintStr(" 0 w3 d w3 neg 0 d 0 w3 d w3 neg 0 d");
693 PrintStr(" 0 w3 neg d w3 neg 0 d cl s } def@");
694 PrintStr("/m29 {mp gsave x w2 sub y w2 add w3 sub m currentpoint t");
695 PrintStr(" 4 {side} repeat cl fill gr} def@");
696 PrintStr("/m30 {mp gsave x w2 sub y w2 add w3 sub m currentpoint t");
697 PrintStr(" 4 {side} repeat cl s gr} def@");
698 PrintStr("/m31 {mp x y w2 sub m 0 w d x w2 sub y m w 0 d");
699 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");
700 PrintStr(" sub m w 1.44 div w 1.44 div d s} def@");
701 PrintStr("/m32 {mp x y w2 sub m w2 w d w neg 0 d cl s} def@");
702 PrintStr("/m33 {mp x y w2 add m w3 neg w2 neg d w3 w2 neg d w3 w2 d cl f} def@");
703 PrintStr("/m34 {mp x w2 sub y w2 sub w3 add m w3 0 d ");
704 PrintStr(" 0 w3 neg d w3 0 d 0 w3 d w3 0 d ");
705 PrintStr(" 0 w3 d w3 neg 0 d 0 w3 d w3 neg 0 d");
706 PrintStr(" 0 w3 neg d w3 neg 0 d cl f } def@");
707 PrintStr("/m35 {mp x y w2 add m w2 neg w2 neg d w2 w2 neg d w2 w2 d w2 neg w2 d");
708 PrintStr(" x y w2 sub m 0 w d x w2 sub y m w 0 d s} def@");
709 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@");
710 PrintStr("/m37 {mp x y m w4 neg w2 d w4 neg w2 neg d w2 0 d ");
711 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@");
712 PrintStr("/m38 {mp x w4 sub y w2 add m w4 neg w4 neg d 0 w2 neg d w4 w4 neg d");
713 PrintStr(" w2 0 d w4 w4 d 0 w2 d w4 neg w4 d w2 neg 0 d");
714 PrintStr(" x y w2 sub m 0 w d x w2 sub y m w 0 d cl s} def@");
715 PrintStr("/m39 {mp x y m w4 neg w2 d w4 neg w2 neg d w2 0 d ");
716 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@");
717 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");
718 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@");
719 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");
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 f} def@");
721 PrintStr("/m42 {mp x y w2 add m w8 neg w2 -3 4 div mul d w2 -3 4 div mul w8 neg d");
722 PrintStr(" w2 3 4 div mul w8 neg d w8 w2 -3 4 div mul d");
723 PrintStr(" w8 w2 3 4 div mul d w2 3 4 div mul w8 d");
724 PrintStr(" w2 -3 4 div mul w8 d w8 neg w2 3 4 div mul d cl s} def@");
725 PrintStr("/m43 {mp x y w2 add m w8 neg w2 -3 4 div mul d w2 -3 4 div mul w8 neg d");
726 PrintStr(" w2 3 4 div mul w8 neg d w8 w2 -3 4 div mul d");
727 PrintStr(" w8 w2 3 4 div mul d w2 3 4 div mul w8 d");
728 PrintStr(" w2 -3 4 div mul w8 d w8 neg w2 3 4 div mul d cl f} def@");
729 PrintStr("/m44 {mp x y m w6 neg w2 d w2 2 3 div mul 0 d w6 neg w2 neg d");
730 PrintStr(" w2 w6 d 0 w2 -2 3 div mul d w2 neg w6 d");
731 PrintStr(" w6 w2 neg d w2 -2 3 div mul 0 d w6 w2 d");
732 PrintStr(" w2 neg w6 neg d 0 w2 2 3 div mul d w2 w6 neg d cl s} def@");
733 PrintStr("/m45 {mp x y m w6 neg w2 d w2 2 3 div mul 0 d w6 neg w2 neg d");
734 PrintStr(" w2 w6 d 0 w2 -2 3 div mul d w2 neg w6 d");
735 PrintStr(" w6 w2 neg d w2 -2 3 div mul 0 d w6 w2 d");
736 PrintStr(" w2 neg w6 neg d 0 w2 2 3 div mul d w2 w6 neg d cl f} def@");
737 PrintStr("/m46 {mp x y w4 add m w4 neg w4 d w4 neg w4 neg d ");
738 PrintStr(" w4 w4 neg d w4 neg w4 neg d w4 w4 neg d w4 w4 d");
739 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@");
740 PrintStr("/m47 {mp x y w4 add m w4 neg w4 d w4 neg w4 neg d");
741 PrintStr(" w4 w4 neg d w4 neg w4 neg d w4 w4 neg d w4 w4 d");
742 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@");
743 PrintStr("/m48 {mp x y w4 add m w4 neg w4 d w4 neg w4 neg d w4 w4 neg d ");
744 PrintStr(" w4 neg w4 neg d w4 w4 neg d w4 w4 d w4 w4 neg d w4 w4 d");
745 PrintStr(" w4 neg w4 d w4 w4 d w4 neg w4 d w4 neg w4 neg d ");
746 PrintStr(" w4 w4 neg d w4 neg w4 neg d w4 neg w4 d w4 w4 d cl f} def@");
747 PrintStr("/m49 {mp x w2 sub w3 add y w2 sub w3 add m ");
748 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");
749 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@");
750 PrintStr("/m2 {mp x y w2 sub m 0 w d x w2 sub y m w 0 d s} def@");
751 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@");
752}
753
754////////////////////////////////////////////////////////////////////////////////
755/// Draw a Box
756
758{
759 static Double_t x[4], y[4];
760 Int_t ix1 = XtoPS(x1);
761 Int_t ix2 = XtoPS(x2);
762 Int_t iy1 = YtoPS(y1);
763 Int_t iy2 = YtoPS(y2);
764 Int_t fillis = fFillStyle/1000;
765 Int_t fillsi = fFillStyle%1000;
766
767 if (fillis == 3 || fillis == 2) {
768 if (fillsi > 99) {
769 x[0] = x1; y[0] = y1;
770 x[1] = x2; y[1] = y1;
771 x[2] = x2; y[2] = y2;
772 x[3] = x1; y[3] = y2;
773 return;
774 }
775 if (fillsi > 0 && fillsi < 26) {
776 x[0] = x1; y[0] = y1;
777 x[1] = x2; y[1] = y1;
778 x[2] = x2; y[2] = y2;
779 x[3] = x1; y[3] = y2;
780 DrawPS(-4, &x[0], &y[0]);
781 }
782 if (fillsi == -3) {
783 SetColor(5);
788 PrintFast(3," bf");
789 }
790 }
791 if (fillis == 1) {
797 PrintFast(3," bf");
798 }
799 if (fillis == 0) {
800 if (fLineWidth<=0) return;
806 PrintFast(3," bl");
807 }
808}
809
810////////////////////////////////////////////////////////////////////////////////
811/// Draw a Frame around a box
812///
813/// - mode = -1 box looks as it is behind the screen
814/// - mode = 1 box looks as it is in front of the screen
815/// - border is the border size in already precomputed PostScript units
816/// - dark is the color for the dark part of the frame
817/// - light is the color for the light part of the frame
818
820 Int_t mode, Int_t border, Int_t dark, Int_t light)
821{
822 static Int_t xps[7], yps[7];
823 Int_t i, ixd0, iyd0, idx, idy, ixdi, iydi, ix, iy;
824
825 // Draw top&left part of the box
826 if (mode == -1) SetColor(dark);
827 else SetColor(light);
828 Int_t bordPS = 4*border;
829 xps[0] = XtoPS(xl); yps[0] = YtoPS(yl);
830 xps[1] = xps[0] + bordPS; yps[1] = yps[0] + bordPS;
831 xps[2] = xps[1]; yps[2] = YtoPS(yt) - bordPS;
832 xps[3] = XtoPS(xt) - bordPS; yps[3] = yps[2];
833 xps[4] = XtoPS(xt); yps[4] = YtoPS(yt);
834 xps[5] = xps[0]; yps[5] = yps[4];
835 xps[6] = xps[0]; yps[6] = yps[0];
836
837 ixd0 = xps[0];
838 iyd0 = yps[0];
841
842 PrintFast(2," m");
843 idx = 0;
844 idy = 0;
845 for (i=1;i<7;i++) {
846 ixdi = xps[i];
847 iydi = yps[i];
848 ix = ixdi - ixd0;
849 iy = iydi - iyd0;
850 ixd0 = ixdi;
851 iyd0 = iydi;
852 if( ix && iy) {
853 if( idx ) { MovePS(idx,0); idx = 0; }
854 if( idy ) { MovePS(0,idy); idy = 0; }
855 MovePS(ix,iy);
856 continue;
857 }
858 if ( ix ) {
859 if( idy ) { MovePS(0,idy); idy = 0; }
860 if( !idx ) { idx = ix; continue;}
861 if( ix*idx > 0 ) idx += ix;
862 else { MovePS(idx,0); idx = ix; }
863 continue;
864 }
865 if( iy ) {
866 if( idx ) { MovePS(idx,0); idx = 0; }
867 if( !idy) { idy = iy; continue;}
868 if( iy*idy > 0 ) idy += iy;
869 else { MovePS(0,idy); idy = iy; }
870 }
871 }
872 if( idx ) MovePS(idx,0);
873 if( idy ) MovePS(0,idy);
874 PrintFast(2," f");
875
876 // Draw bottom&right part of the box
877 if (mode == -1) SetColor(light);
878 else SetColor(dark);
879 xps[0] = XtoPS(xl); yps[0] = YtoPS(yl);
880 xps[1] = xps[0] + bordPS; yps[1] = yps[0] + bordPS;
881 xps[2] = XtoPS(xt) - bordPS; yps[2] = yps[1];
882 xps[3] = xps[2]; yps[3] = YtoPS(yt) - bordPS;
883 xps[4] = XtoPS(xt); yps[4] = YtoPS(yt);
884 xps[5] = xps[4]; yps[5] = yps[0];
885 xps[6] = xps[0]; yps[6] = yps[0];
886
887 ixd0 = xps[0];
888 iyd0 = yps[0];
891
892 PrintFast(2," m");
893 idx = 0;
894 idy = 0;
895 for (i=1;i<7;i++) {
896 ixdi = xps[i];
897 iydi = yps[i];
898 ix = ixdi - ixd0;
899 iy = iydi - iyd0;
900 ixd0 = ixdi;
901 iyd0 = iydi;
902 if( ix && iy) {
903 if( idx ) { MovePS(idx,0); idx = 0; }
904 if( idy ) { MovePS(0,idy); idy = 0; }
905 MovePS(ix,iy);
906 continue;
907 }
908 if ( ix ) {
909 if( idy ) { MovePS(0,idy); idy = 0; }
910 if( !idx ) { idx = ix; continue;}
911 if( ix*idx > 0 ) idx += ix;
912 else { MovePS(idx,0); idx = ix; }
913 continue;
914 }
915 if( iy ) {
916 if( idx ) { MovePS(idx,0); idx = 0; }
917 if( !idy) { idy = iy; continue;}
918 if( iy*idy > 0 ) idy += iy;
919 else { MovePS(0,idy); idy = iy; }
920 }
921 }
922 if( idx ) MovePS(idx,0);
923 if( idy ) MovePS(0,idy);
924 PrintFast(2," f");
925}
926
927////////////////////////////////////////////////////////////////////////////////
928/// Draw a PolyLine
929///
930/// Draw a polyline through the points xy.
931/// - If nn=1 moves only to point x,y.
932/// - If nn=0 the x,y are written in the PostScript file
933/// according to the current transformation.
934/// - If nn>0 the line is clipped as a line.
935/// - If nn<0 the line is clipped as a fill area.
936
938{
939 Int_t i, n, ixd0, iyd0, idx, idy, ixdi, iydi, ix, iy;
942 if (nn > 0) {
943 if (fLineWidth<=0) return;
944 n = nn;
948 } else {
949 n = -nn;
950 SetStyle(1);
951 SetWidth(1);
953 }
954
955 ixd0 = XtoPS(xy[0].GetX());
956 iyd0 = YtoPS(xy[0].GetY());
959 if( n <= 1) {
960 if( n == 0) goto END;
961 PrintFast(2," m");
962 goto END;
963 }
964
965 PrintFast(2," m");
966 idx = 0;
967 idy = 0;
968 for (i=1;i<n;i++) {
969 ixdi = XtoPS(xy[i].GetX());
970 iydi = YtoPS(xy[i].GetY());
971 ix = ixdi - ixd0;
972 iy = iydi - iyd0;
973 ixd0 = ixdi;
974 iyd0 = iydi;
975 if( ix && iy) {
976 if( idx ) { MovePS(idx,0); idx = 0; }
977 if( idy ) { MovePS(0,idy); idy = 0; }
978 MovePS(ix,iy);
979 continue;
980 }
981 if ( ix ) {
982 if( idy ) { MovePS(0,idy); idy = 0; }
983 if( !idx ) { idx = ix; continue;}
984 if( ix*idx > 0 ) idx += ix;
985 else { MovePS(idx,0); idx = ix; }
986 continue;
987 }
988 if( iy ) {
989 if( idx ) { MovePS(idx,0); idx = 0; }
990 if( !idy) { idy = iy; continue;}
991 if( iy*idy > 0 ) idy += iy;
992 else { MovePS(0,idy); idy = iy; }
993 }
994 }
995 if( idx ) MovePS(idx,0);
996 if( idy ) MovePS(0,idy);
997
998 if (nn > 0 ) {
999 if (xy[0].GetX() == xy[n-1].GetX() && xy[0].GetY() == xy[n-1].GetY()) PrintFast(3," cl");
1000 PrintFast(2," s");
1001 } else {
1002 PrintFast(2," f");
1003 }
1004END:
1005 if (nn < 0) {
1008 }
1009}
1010
1011////////////////////////////////////////////////////////////////////////////////
1012/// Draw a PolyLine in NDC space
1013///
1014/// Draw a polyline through the points xy.
1015/// - If nn=1 moves only to point x,y.
1016/// - If nn=0 the x,y are written in the PostScript file
1017/// according to the current transformation.
1018/// - If nn>0 the line is clipped as a line.
1019/// - If nn<0 the line is clipped as a fill area.
1020
1022{
1023 Int_t i, n, ixd0, iyd0, idx, idy, ixdi, iydi, ix, iy;
1026 if (nn > 0) {
1027 if (fLineWidth<=0) return;
1028 n = nn;
1032 } else {
1033 n = -nn;
1034 SetStyle(1);
1035 SetWidth(1);
1037 }
1038
1039 ixd0 = UtoPS(xy[0].GetX());
1040 iyd0 = VtoPS(xy[0].GetY());
1043 if( n <= 1) {
1044 if( n == 0) goto END;
1045 PrintFast(2," m");
1046 goto END;
1047 }
1048
1049 PrintFast(2," m");
1050 idx = 0;
1051 idy = 0;
1052 for (i=1;i<n;i++) {
1053 ixdi = UtoPS(xy[i].GetX());
1054 iydi = VtoPS(xy[i].GetY());
1055 ix = ixdi - ixd0;
1056 iy = iydi - iyd0;
1057 ixd0 = ixdi;
1058 iyd0 = iydi;
1059 if( ix && iy) {
1060 if( idx ) { MovePS(idx,0); idx = 0; }
1061 if( idy ) { MovePS(0,idy); idy = 0; }
1062 MovePS(ix,iy);
1063 continue;
1064 }
1065 if ( ix ) {
1066 if( idy ) { MovePS(0,idy); idy = 0; }
1067 if( !idx ) { idx = ix; continue;}
1068 if( ix*idx > 0 ) idx += ix;
1069 else { MovePS(idx,0); idx = ix; }
1070 continue;
1071 }
1072 if( iy ) {
1073 if( idx ) { MovePS(idx,0); idx = 0; }
1074 if( !idy) { idy = iy; continue;}
1075 if( iy*idy > 0 ) idy += iy;
1076 else { MovePS(0,idy); idy = iy; }
1077 }
1078 }
1079 if( idx ) MovePS(idx,0);
1080 if( idy ) MovePS(0,idy);
1081
1082 if (nn > 0 ) {
1083 if (xy[0].GetX() == xy[n-1].GetX() && xy[0].GetY() == xy[n-1].GetY()) PrintFast(3," cl");
1084 PrintFast(2," s");
1085 } else {
1086 PrintFast(2," f");
1087 }
1088END:
1089 if (nn < 0) {
1092 }
1093}
1094
1095////////////////////////////////////////////////////////////////////////////////
1096/// Draw markers at the n WC points x, y
1097
1099{
1100 Int_t i, np, markerstyle;
1102 static char chtemp[10];
1103
1104 if (!fMarkerSize) return;
1108 SetStyle(1);
1112 if (markerstyle <= 0) strlcpy(chtemp, " m20",10);
1113 if (markerstyle == 1) strlcpy(chtemp, " m20",10);
1114 if (markerstyle == 2) strlcpy(chtemp, " m2",10);
1115 if (markerstyle == 3) strlcpy(chtemp, " m31",10);
1116 if (markerstyle == 4) strlcpy(chtemp, " m24",10);
1117 if (markerstyle == 5) strlcpy(chtemp, " m5",10);
1118 if (markerstyle >= 6 && markerstyle <= 19) strlcpy(chtemp, " m20",10);
1119 if (markerstyle >= 20 && markerstyle <= 49 ) snprintf(chtemp,10," m%d", markerstyle);
1120 if (markerstyle >= 50) strlcpy(chtemp, " m20",10);
1121
1122 // Set the PostScript marker size
1123 if (markerstyle == 1 || (markerstyle >= 9 && markerstyle <= 19)) {
1124 markersize = 2.;
1125 } else if (markerstyle == 6) {
1126 markersize = 4.;
1127 } else if (markerstyle == 7) {
1128 markersize = 8.;
1129 } else {
1131 const Int_t kBASEMARKER = 8;
1133 Float_t s2x = sbase / Float_t(gPad->GetWw() * gPad->GetAbsWNDC());
1134 markersize = this->UtoPS(s2x) - this->UtoPS(0);
1135 }
1136
1137 if (fMarkerSizeCur != markersize) {
1139 PrintFast(3," /w");
1141 PrintFast(40," def /w2 {w 2 div} def /w3 {w 3 div} def");
1142 }
1143
1144 WriteInteger(XtoPS(x[0]));
1145 WriteInteger(YtoPS(y[0]));
1146 if (n == 1) {
1150 return;
1151 }
1152 np = 1;
1153 for (i=1;i<n;i++) {
1154 WriteInteger(XtoPS(x[i]));
1155 WriteInteger(YtoPS(y[i]));
1156 np++;
1157 if (np == 100 || i == n-1) {
1159 PrintFast(2," {");
1161 PrintFast(3,"} R");
1162 np = 0;
1163 }
1164 }
1167}
1168
1169////////////////////////////////////////////////////////////////////////////////
1170/// Draw markers at the n WC points x, y
1171
1173{
1174 Int_t i, np, markerstyle;
1176 static char chtemp[10];
1177
1178 if (!fMarkerSize) return;
1182 SetStyle(1);
1186 if (markerstyle <= 0) strlcpy(chtemp, " m20",10);
1187 if (markerstyle == 1) strlcpy(chtemp, " m20",10);
1188 if (markerstyle == 2) strlcpy(chtemp, " m2",10);
1189 if (markerstyle == 3) strlcpy(chtemp, " m31",10);
1190 if (markerstyle == 4) strlcpy(chtemp, " m24",10);
1191 if (markerstyle == 5) strlcpy(chtemp, " m5",10);
1192 if (markerstyle >= 6 && markerstyle <= 19) strlcpy(chtemp, " m20",10);
1193 if (markerstyle >= 20 && markerstyle <= 49 ) snprintf(chtemp,10," m%d", markerstyle);
1194 if (markerstyle >= 50) strlcpy(chtemp, " m20",10);
1195
1196 // Set the PostScript marker size
1197 if (markerstyle == 1 || (markerstyle >= 9 && markerstyle <= 19)) {
1198 markersize = 2.;
1199 } else if (markerstyle == 6) {
1200 markersize = 4.;
1201 } else if (markerstyle == 7) {
1202 markersize = 8.;
1203 } else {
1205 const Int_t kBASEMARKER = 8;
1207 Float_t s2x = sbase / Float_t(gPad->GetWw() * gPad->GetAbsWNDC());
1208 markersize = this->UtoPS(s2x) - this->UtoPS(0);
1209 }
1210
1211 if (fMarkerSizeCur != markersize) {
1213 PrintFast(3," /w");
1215 PrintFast(40," def /w2 {w 2 div} def /w3 {w 3 div} def");
1216 }
1217
1218 WriteInteger(XtoPS(x[0]));
1219 WriteInteger(YtoPS(y[0]));
1220 if (n == 1) {
1224 return;
1225 }
1226 np = 1;
1227 for (i=1;i<n;i++) {
1228 WriteInteger(XtoPS(x[i]));
1229 WriteInteger(YtoPS(y[i]));
1230 np++;
1231 if (np == 100 || i == n-1) {
1233 PrintFast(2," {");
1235 PrintFast(3,"} R");
1236 np = 0;
1237 }
1238 }
1241}
1242
1243////////////////////////////////////////////////////////////////////////////////
1244/// Draw a PolyLine
1245///
1246/// Draw a polyline through the points xw,yw.
1247/// - If nn=1 moves only to point xw,yw.
1248/// - If nn=0 the XW(1) and YW(1) are written in the PostScript file
1249/// according to the current NT.
1250/// - If nn>0 the line is clipped as a line.
1251/// - If nn<0 the line is clipped as a fill area.
1252
1254{
1255 static Float_t dyhatch[24] = {.0075,.0075,.0075,.0075,.0075,.0075,.0075,.0075,
1256 .01 ,.01 ,.01 ,.01 ,.01 ,.01 ,.01 ,.01 ,
1257 .015 ,.015 ,.015 ,.015 ,.015 ,.015 ,.015 ,.015};
1258 static Float_t anglehatch[24] = {180, 90,135, 45,150, 30,120, 60,
1259 180, 90,135, 45,150, 30,120, 60,
1260 180, 90,135, 45,150, 30,120, 60};
1261 Int_t i, n, ixd0, iyd0, idx, idy, ixdi, iydi, ix, iy, fais, fasi;
1262 fais = fasi = n = 0;
1263 Int_t jxd0 = XtoPS(xw[0]);
1264 Int_t jyd0 = YtoPS(yw[0]);
1267
1268 if (nn > 0) {
1269 if (fLineWidth<=0) return;
1270 n = nn;
1274 }
1275 if (nn < 0) {
1276 n = -nn;
1277 SetStyle(1);
1278 SetWidth(1);
1280 fais = fFillStyle/1000;
1281 fasi = fFillStyle%1000;
1282 if (fais == 3 || fais == 2) {
1283 if (fasi > 100 && fasi <125) {
1284 DrawHatch(dyhatch[fasi-101],anglehatch[fasi-101], n, xw, yw);
1285 goto END;
1286 }
1287 if (fasi > 0 && fasi < 26) {
1289 }
1290 }
1291 }
1292
1293 ixd0 = jxd0;
1294 iyd0 = jyd0;
1297 if( n <= 1) {
1298 if( n == 0) goto END;
1299 PrintFast(2," m");
1300 goto END;
1301 }
1302
1303 PrintFast(2," m");
1304 idx = idy = 0;
1305 for (i=1;i<n;i++) {
1306 ixdi = XtoPS(xw[i]);
1307 iydi = YtoPS(yw[i]);
1308 ix = ixdi - ixd0;
1309 iy = iydi - iyd0;
1310 ixd0 = ixdi;
1311 iyd0 = iydi;
1312 if( ix && iy) {
1313 if( idx ) { MovePS(idx,0); idx = 0; }
1314 if( idy ) { MovePS(0,idy); idy = 0; }
1315 MovePS(ix,iy);
1316 } else if ( ix ) {
1317 if( idy ) { MovePS(0,idy); idy = 0;}
1318 if( !idx ) { idx = ix;}
1319 else if( TMath::Sign(ix,idx) == ix ) idx += ix;
1320 else { MovePS(idx,0); idx = ix;}
1321 } else if( iy ) {
1322 if( idx ) { MovePS(idx,0); idx = 0;}
1323 if( !idy) { idy = iy;}
1324 else if( TMath::Sign(iy,idy) == iy) idy += iy;
1325 else { MovePS(0,idy); idy = iy;}
1326 }
1327 }
1328 if (idx) MovePS(idx,0);
1329 if (idy) MovePS(0,idy);
1330
1331 if (nn > 0 ) {
1332 if (xw[0] == xw[n-1] && yw[0] == yw[n-1]) PrintFast(3," cl");
1333 PrintFast(2," s");
1334 } else {
1335 if (fais == 0) {PrintFast(5," cl s"); goto END;}
1336 if (fais == 3 || fais == 2) {
1337 if (fasi > 0 && fasi < 26) {
1338 PrintFast(3," FA");
1339 fRed = -1;
1340 fGreen = -1;
1341 fBlue = -1;
1342 }
1343 goto END;
1344 }
1345 PrintFast(2," f");
1346 }
1347END:
1348 if (nn < 0) {
1351 }
1352}
1353
1354////////////////////////////////////////////////////////////////////////////////
1355/// Draw a PolyLine
1356///
1357/// Draw a polyline through the points xw,yw.
1358/// - If nn=1 moves only to point xw,yw.
1359/// - If nn=0 the xw(1) and YW(1) are written in the PostScript file
1360/// --- according to the current NT.
1361/// - If nn>0 the line is clipped as a line.
1362/// - If nn<0 the line is clipped as a fill area.
1363
1365{
1366 static Float_t dyhatch[24] = {.0075,.0075,.0075,.0075,.0075,.0075,.0075,.0075,
1367 .01 ,.01 ,.01 ,.01 ,.01 ,.01 ,.01 ,.01 ,
1368 .015 ,.015 ,.015 ,.015 ,.015 ,.015 ,.015 ,.015};
1369 static Float_t anglehatch[24] = {180, 90,135, 45,150, 30,120, 60,
1370 180, 90,135, 45,150, 30,120, 60,
1371 180, 90,135, 45,150, 30,120, 60};
1372 Int_t i, n, ixd0, iyd0, idx, idy, ixdi, iydi, ix, iy, fais, fasi;
1373 fais = fasi = n = 0;
1374 Int_t jxd0 = XtoPS(xw[0]);
1375 Int_t jyd0 = YtoPS(yw[0]);
1378
1379 if (nn > 0) {
1380 if (fLineWidth<=0) return;
1381 n = nn;
1385 }
1386 if (nn < 0) {
1387 n = -nn;
1388 SetStyle(1);
1389 SetWidth(1);
1391 fais = fFillStyle/1000;
1392 fasi = fFillStyle%1000;
1393 if (fais == 3 || fais == 2) {
1394 if (fasi > 100 && fasi <125) {
1395 DrawHatch(dyhatch[fasi-101],anglehatch[fasi-101], n, xw, yw);
1396 goto END;
1397 }
1398 if (fasi > 0 && fasi < 26) {
1400 }
1401 }
1402 }
1403
1404 ixd0 = jxd0;
1405 iyd0 = jyd0;
1408 if( n <= 1) {
1409 if( n == 0) goto END;
1410 PrintFast(2," m");
1411 goto END;
1412 }
1413
1414 PrintFast(2," m");
1415 idx = idy = 0;
1416 for (i=1;i<n;i++) {
1417 ixdi = XtoPS(xw[i]);
1418 iydi = YtoPS(yw[i]);
1419 ix = ixdi - ixd0;
1420 iy = iydi - iyd0;
1421 ixd0 = ixdi;
1422 iyd0 = iydi;
1423 if( ix && iy) {
1424 if( idx ) { MovePS(idx,0); idx = 0; }
1425 if( idy ) { MovePS(0,idy); idy = 0; }
1426 MovePS(ix,iy);
1427 } else if ( ix ) {
1428 if( idy ) { MovePS(0,idy); idy = 0;}
1429 if( !idx ) { idx = ix;}
1430 else if( TMath::Sign(ix,idx) == ix ) idx += ix;
1431 else { MovePS(idx,0); idx = ix;}
1432 } else if( iy ) {
1433 if( idx ) { MovePS(idx,0); idx = 0;}
1434 if( !idy) { idy = iy;}
1435 else if( TMath::Sign(iy,idy) == iy) idy += iy;
1436 else { MovePS(0,idy); idy = iy;}
1437 }
1438 }
1439 if (idx) MovePS(idx,0);
1440 if (idy) MovePS(0,idy);
1441
1442 if (nn > 0 ) {
1443 if (xw[0] == xw[n-1] && yw[0] == yw[n-1]) PrintFast(3," cl");
1444 PrintFast(2," s");
1445 } else {
1446 if (fais == 0) {PrintFast(5," cl s"); goto END;}
1447 if (fais == 3 || fais == 2) {
1448 if (fasi > 0 && fasi < 26) {
1449 PrintFast(3," FA");
1450 fRed = -1;
1451 fGreen = -1;
1452 fBlue = -1;
1453 }
1454 goto END;
1455 }
1456 PrintFast(2," f");
1457 }
1458END:
1459 if (nn < 0) {
1462 }
1463}
1464
1465////////////////////////////////////////////////////////////////////////////////
1466/// Draw Fill area with hatch styles
1467
1469{
1470 Warning("DrawHatch", "hatch fill style not yet implemented");
1471}
1472
1473////////////////////////////////////////////////////////////////////////////////
1474/// Draw Fill area with hatch styles
1475
1477{
1478 Warning("DrawHatch", "hatch fill style not yet implemented");
1479}
1480
1481////////////////////////////////////////////////////////////////////////////////
1482
1484{
1485 std::ifstream font_file(filename, std::ios::binary);
1486
1487 // We cannot read directly using iostream iterators due to
1488 // signedness
1489 font_file.seekg(0, std::ios::end);
1490
1491 const size_t font_file_length = font_file.tellg();
1492
1493 font_file.seekg(0, std::ios::beg);
1494
1495 std::vector<unsigned char> font_data(font_file_length, '\0');
1496
1497 font_file.read(reinterpret_cast<char *>(&font_data[0]),
1499
1500 std::string font_name;
1501 std::string postscript_string =
1502 mathtext::font_embed_postscript_t::font_embed_type_1(
1503 font_name, font_data);
1504
1505 if (!postscript_string.empty()) {
1507 PrintStr("@");
1508
1509 return true;
1510 }
1511
1512 return false;
1513}
1514
1515////////////////////////////////////////////////////////////////////////////////
1516
1518{
1519 std::ifstream font_file(filename, std::ios::binary);
1520
1521 // We cannot read directly using iostream iterators due to
1522 // signedness
1523 font_file.seekg(0, std::ios::end);
1524
1525 const size_t font_file_length = font_file.tellg();
1526
1527 font_file.seekg(0, std::ios::beg);
1528
1529 std::vector<unsigned char> font_data(font_file_length, '\0');
1530
1531 font_file.read(reinterpret_cast<char *>(&font_data[0]), font_file_length);
1532
1533 std::string font_name;
1534 std::string postscript_string =
1535 mathtext::font_embed_postscript_t::font_embed_type_2(font_name, font_data);
1536
1537 if (!postscript_string.empty()) {
1539 PrintStr("@");
1540
1541 return true;
1542 }
1543
1544 return false;
1545}
1546
1547////////////////////////////////////////////////////////////////////////////////
1548
1550{
1551 std::ifstream font_file(filename, std::ios::binary);
1552
1553 // We cannot read directly using iostream iterators due to signedness
1554
1555 font_file.seekg(0, std::ios::end);
1556
1557 const size_t font_file_length = font_file.tellg();
1558
1559 font_file.seekg(0, std::ios::beg);
1560
1561 std::vector<unsigned char> font_data(font_file_length, '\0');
1562
1563 font_file.read(reinterpret_cast<char *>(&font_data[0]), font_file_length);
1564
1565 std::string font_name;
1566 std::string postscript_string =
1567 mathtext::font_embed_postscript_t::font_embed_type_42(font_name, font_data);
1568
1569 if (!postscript_string.empty()) {
1571 PrintStr("@");
1572
1573 return true;
1574 }
1575 fprintf(stderr, "%s:%d:\n", __FILE__, __LINE__);
1576
1577 return false;
1578}
1579
1580////////////////////////////////////////////////////////////////////////////////
1581/// Embed font in PS file.
1582
1584{
1585 static const char *fonttable[32][2] = {
1586 { "Root.TTFont.0", "FreeSansBold.otf" },
1587 { "Root.TTFont.1", "FreeSerifItalic.otf" },
1588 { "Root.TTFont.2", "FreeSerifBold.otf" },
1589 { "Root.TTFont.3", "FreeSerifBoldItalic.otf" },
1590 { "Root.TTFont.4", "FreeSans.otf" },
1591 { "Root.TTFont.5", "FreeSansOblique.otf" },
1592 { "Root.TTFont.6", "FreeSansBold.otf" },
1593 { "Root.TTFont.7", "FreeSansBoldOblique.otf" },
1594 { "Root.TTFont.8", "FreeMono.otf" },
1595 { "Root.TTFont.9", "FreeMonoOblique.otf" },
1596 { "Root.TTFont.10", "FreeMonoBold.otf" },
1597 { "Root.TTFont.11", "FreeMonoBoldOblique.otf" },
1598 { "Root.TTFont.12", "symbol.ttf" },
1599 { "Root.TTFont.13", "FreeSerif.otf" },
1600 { "Root.TTFont.14", "wingding.ttf" },
1601 { "Root.TTFont.15", "symbol.ttf" },
1602 { "Root.TTFont.STIXGen", "STIXGeneral.otf" },
1603 { "Root.TTFont.STIXGenIt", "STIXGeneralItalic.otf" },
1604 { "Root.TTFont.STIXGenBd", "STIXGeneralBol.otf" },
1605 { "Root.TTFont.STIXGenBdIt", "STIXGeneralBolIta.otf" },
1606 { "Root.TTFont.STIXSiz1Sym", "STIXSiz1Sym.otf" },
1607 { "Root.TTFont.STIXSiz1SymBd", "STIXSiz1SymBol.otf" },
1608 { "Root.TTFont.STIXSiz2Sym", "STIXSiz2Sym.otf" },
1609 { "Root.TTFont.STIXSiz2SymBd", "STIXSiz2SymBol.otf" },
1610 { "Root.TTFont.STIXSiz3Sym", "STIXSiz3Sym.otf" },
1611 { "Root.TTFont.STIXSiz3SymBd", "STIXSiz3SymBol.otf" },
1612 { "Root.TTFont.STIXSiz4Sym", "STIXSiz4Sym.otf" },
1613 { "Root.TTFont.STIXSiz4SymBd", "STIXSiz4SymBol.otf" },
1614 { "Root.TTFont.STIXSiz5Sym", "STIXSiz5Sym.otf" },
1615 { "Root.TTFont.ME", "DroidSansFallback.ttf" },
1616 { "Root.TTFont.CJKMing", "DroidSansFallback.ttf" },
1617 { "Root.TTFont.CJKCothic", "DroidSansFallback.ttf" }
1618 };
1619
1620 PrintStr("%%IncludeResource: ProcSet (FontSetInit)@");
1621
1622 // try to load font (font must be in Root.TTFontPath resource)
1623 const char *ttpath = gEnv->GetValue("Root.TTFontPath",
1625
1626 for (Int_t fontid = 1; fontid < 30; fontid++) {
1627 if (fontid != 15 && MustEmbed[fontid-1]) {
1628 const char *filename = gEnv->GetValue(
1629 fonttable[fontid][0], fonttable[fontid][1]);
1631 if (!ttfont) {
1632 Error("TPostScript::FontEmbed",
1633 "font %d (filename `%s') not found in path",
1634 fontid, filename);
1635 } else {
1636 if (FontEmbedType2(ttfont)) {
1637 // nothing
1638 } else if(FontEmbedType1(ttfont)) {
1639 // nothing
1640 } else if(FontEmbedType42(ttfont)) {
1641 // nothing
1642 } else {
1643 Error("TPostScript::FontEmbed",
1644 "failed to embed font %d (filename `%s')",
1645 fontid, filename);
1646 }
1647 delete [] ttfont;
1648 }
1649 }
1650 }
1651 PrintStr("%%IncludeResource: font Times-Roman@");
1652 PrintStr("%%IncludeResource: font Times-Italic@");
1653 PrintStr("%%IncludeResource: font Times-Bold@");
1654 PrintStr("%%IncludeResource: font Times-BoldItalic@");
1655 PrintStr("%%IncludeResource: font Helvetica@");
1656 PrintStr("%%IncludeResource: font Helvetica-Oblique@");
1657 PrintStr("%%IncludeResource: font Helvetica-Bold@");
1658 PrintStr("%%IncludeResource: font Helvetica-BoldOblique@");
1659 PrintStr("%%IncludeResource: font Courier@");
1660 PrintStr("%%IncludeResource: font Courier-Oblique@");
1661 PrintStr("%%IncludeResource: font Courier-Bold@");
1662 PrintStr("%%IncludeResource: font Courier-BoldOblique@");
1663 PrintStr("%%IncludeResource: font Symbol@");
1664 PrintStr("%%IncludeResource: font ZapfDingbats@");
1665
1666 fFontEmbed = kTRUE;
1667}
1668
1669////////////////////////////////////////////////////////////////////////////////
1670/// Font Re-encoding
1671
1673{
1674 PrintStr("/reEncode ");
1675 PrintStr("{exch findfont");
1676 PrintStr(" dup length dict begin");
1677 PrintStr(" {1 index /FID eq ");
1678 PrintStr(" {pop pop}");
1679 PrintStr(" {def} ifelse");
1680 PrintStr(" } forall");
1681 PrintStr(" /Encoding exch def");
1682 PrintStr(" currentdict end");
1683 PrintStr(" dup /FontName get exch");
1684 PrintStr(" definefont pop");
1685 PrintStr(" } def");
1686 PrintStr(" [/Times-Bold /Times-Italic /Times-BoldItalic /Helvetica");
1687 PrintStr(" /Helvetica-Oblique /Helvetica-Bold /Helvetica-BoldOblique");
1688 PrintStr(" /Courier /Courier-Oblique /Courier-Bold /Courier-BoldOblique");
1689 PrintStr(" /Times-Roman /AvantGarde-Book /AvantGarde-BookOblique");
1690 PrintStr(" /AvantGarde-Demi /AvantGarde-DemiOblique /Bookman-Demi");
1691 PrintStr(" /Bookman-DemiItalic /Bookman-Light /Bookman-LightItalic");
1692 PrintStr(" /Helvetica-Narrow /Helvetica-Narrow-Bold /Helvetica-Narrow-BoldOblique");
1693 PrintStr(" /Helvetica-Narrow-Oblique /NewCenturySchlbk-Roman /NewCenturySchlbk-Bold");
1694 PrintStr(" /NewCenturySchlbk-BoldItalic /NewCenturySchlbk-Italic");
1695 PrintStr(" /Palatino-Bold /Palatino-BoldItalic /Palatino-Italic /Palatino-Roman");
1696 PrintStr(" ] {ISOLatin1Encoding reEncode } forall");
1697}
1698
1699////////////////////////////////////////////////////////////////////////////////
1700/// PostScript Initialisation
1701///
1702/// This method initialize the following PostScript procedures:
1703///
1704/// | Macro Name | Input parameters | Explanation |
1705/// |------------|------------------|-----------------------------------|
1706/// | l | x y | Draw a line to the x y position |
1707/// | m | x y | Move to the position x y |
1708/// | box | dx dy x y | Define a box |
1709/// | bl | dx dy x y | Draw a line box |
1710/// | bf | dx dy x y | Draw a filled box |
1711/// | t | x y | Translate |
1712/// | r | angle | Rotate |
1713/// | rl | i j | Roll the stack |
1714/// | d | x y | Draw a relative line to x y |
1715/// | X | x | Draw a relative line to x (y=0) |
1716/// | Y | y | Draw a relative line to y (x=0) |
1717/// | rm | x y | Move relatively to x y |
1718/// | gr | | Restore the graphic context |
1719/// | lw | lwidth | Set line width to lwidth |
1720/// | sd | [] 0 | Set dash line define by [] |
1721/// | s | | Stroke mode |
1722/// | c | r g b | Set rgb color to r g b |
1723/// | cl | | Close path |
1724/// | f | | Fill the last describe path |
1725/// | mXX | x y | Draw the marker type XX at (x,y) |
1726/// | Zone | ix iy | Define the current zone |
1727/// | black | | The color is black |
1728/// | C | dx dy x y | Clipping on |
1729/// | NC | | Clipping off |
1730/// | R | | repeat |
1731/// | ita | | Used to make the symbols italic |
1732/// | K | | kshow |
1733
1735{
1737 rpxmin = rpymin = width = heigth = 0;
1738 Int_t format;
1739 fNpages=1;
1740 for (Int_t i=0;i<32;i++) fPatterns[i]=0;
1741
1742 // Mode is last digit of PostScript Workstation type
1743 // mode=1,2 for portrait/landscape black and white
1744 // mode=3 for Encapsulated PostScript File
1745 // mode=4 for portrait colour
1746 // mode=5 for lanscape colour
1747 Int_t atype = abs(fType);
1748 fMode = atype%10;
1749 if( fMode <= 0 || fMode > 5) {
1750 Error("Initialize", "invalid file type %d", fMode);
1751 return;
1752 }
1753
1754 // fNXzone (fNYzone) is the total number of windows in x (y)
1755 fNXzone = (atype%1000)/100;
1756 fNYzone = (atype%100)/10;
1757 if( fNXzone <= 0 ) fNXzone = 1;
1758 if( fNYzone <= 0 ) fNYzone = 1;
1759 fIXzone = 1;
1760 fIYzone = 1;
1761
1762 // format = 0-99 is the European page format (A4,A3 ...)
1763 // format = 100 is the US format 8.5x11.0 inch
1764 // format = 200 is the US format 8.5x14.0 inch
1765 // format = 300 is the US format 11.0x17.0 inch
1766 format = atype/1000;
1767 if( format == 0 ) format = 4;
1768 if( format == 99 ) format = 0;
1769
1770 PrintStr("%%Title: ");
1771 const char *pstitle = gStyle->GetTitlePS();
1772 if (gPad && !pstitle[0]) pstitle = gPad->GetMother()->GetTitle();
1773 if (strlen(GetName())<=80) PrintStr(GetName());
1774 if(!pstitle[0] && fMode != 3) {;
1775 PrintFast(2," (");
1776 if ( format <= 99 ) {;
1777 PrintFast(2," A");
1779 PrintFast(1,")");
1780 }
1781 else {
1782 if ( format == 100 ) PrintFast(8," Letter)");
1783 if ( format == 200 ) PrintFast(7," Legal)");
1784 if ( format == 300 ) PrintFast(8," Ledger)");
1785 }
1786 PrintStr("@");
1787 PrintStr("%%Pages: (atend)@");
1788 }
1789 else {
1790 if (!strchr(pstitle,'\n')) {
1791 PrintFast(2,": ");
1793 }
1794 PrintStr("@");
1795 }
1796
1797 PrintFast(24,"%%Creator: ROOT Version ");
1798 PrintStr(gROOT->GetVersion());
1799 PrintStr("@");
1800 PrintFast(16,"%%CreationDate: ");
1801 TDatime t;
1802 PrintStr(t.AsString());
1803 PrintStr("@");
1804
1805 if ( fMode == 1 || fMode == 4) PrintStr("%%Orientation: Portrait@");
1806 if ( fMode == 2 || fMode == 5) PrintStr("%%Orientation: Landscape@");
1807
1808 PrintStr("%%EndComments@");
1809 PrintStr("%%BeginProlog@");
1810
1811 if( fMode == 3)PrintStr("80 dict begin@");
1812
1813 // Initialisation of PostScript procedures
1814 PrintStr("/s {stroke} def /l {lineto} def /m {moveto} def /t {translate} def@");
1815 PrintStr("/r {rotate} def /rl {roll} def /R {repeat} def@");
1816 PrintStr("/d {rlineto} def /rm {rmoveto} def /gr {grestore} def /f {eofill} def@");
1817 if (gStyle->GetColorModelPS()) {
1818 PrintStr("/c {setcmykcolor} def /black {0 0 0 1 setcmykcolor} def /sd {setdash} def@");
1819 } else {
1820 PrintStr("/c {setrgbcolor} def /black {0 setgray} def /sd {setdash} def@");
1821 }
1822 PrintStr("/cl {closepath} def /sf {scalefont setfont} def /lw {setlinewidth} def@");
1823 PrintStr("/box {m dup 0 exch d exch 0 d 0 exch neg d cl} def@");
1824 PrintStr("/NC{systemdict begin initclip end}def/C{NC box clip newpath}def@");
1825 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 @");
1826 PrintStr("/K {{pop pop 0 moveto} exch kshow} bind def@");
1827 PrintStr("/ita {/ang 15 def gsave [1 0 ang dup sin exch cos div 1 0 0] concat} def @");
1828
1829 DefineMarkers();
1830
1831 FontEncode();
1832
1833 // mode=1 for portrait black/white
1834 if (fMode == 1) {
1835 rpxmin = 0.7;
1837 switch (format) {
1838 case 100 :
1839 width = (8.5*2.54)-2.*rpxmin;
1840 heigth = (11.*2.54)-2.*rpymin;
1841 break;
1842 case 200 :
1843 width = (8.5*2.54)-2.*rpxmin;
1844 heigth = (14.*2.54)-2.*rpymin;
1845 break;
1846 case 300 :
1847 width = (11.*2.54)-2.*rpxmin;
1848 heigth = (17.*2.54)-2.*rpymin;
1849 break;
1850 default :
1851 width = 21.0-2.*rpxmin;
1852 heigth = 29.7-2.*rpymin;
1853 };
1854 }
1855
1856 // mode=2 for landscape black/white
1857 if (fMode == 2) {
1858 rpymin = 0.7;
1860 switch (format) {
1861 case 100 :
1862 width = (11.*2.54)-2.*rpxmin;
1863 heigth = (8.5*2.54)-2.*rpymin;
1864 break;
1865 case 200 :
1866 width = (14.*2.54)-2.*rpxmin;
1867 heigth = (8.5*2.54)-2.*rpymin;
1868 break;
1869 case 300 :
1870 width = (17.*2.54)-2.*rpxmin;
1871 heigth = (11.*2.54)-2.*rpymin;
1872 break;
1873 default :
1874 width = 29.7-2.*rpxmin;
1875 heigth = 21-2.*rpymin;
1876 };
1877 }
1878
1879 // mode=3 encapsulated PostScript
1880 if (fMode == 3) {
1881 width = 20;
1882 heigth = 20;
1883 format = 4;
1884 fNXzone = 1;
1885 fNYzone = 1;
1886 }
1887
1888 // mode=4 for portrait colour
1889 if (fMode == 4) {
1890 rpxmin = 0.7;
1891 rpymin = 3.4;
1892 switch (format) {
1893 case 100 :
1894 width = (8.5*2.54)-2.*rpxmin;
1895 heigth = (11.*2.54)-2.*rpymin;
1896 break;
1897 case 200 :
1898 width = (8.5*2.54)-2.*rpxmin;
1899 heigth = (14.*2.54)-2.*rpymin;
1900 break;
1901 case 300 :
1902 width = (11.*2.54)-2.*rpxmin;
1903 heigth = (17.*2.54)-2.*rpymin;
1904 break;
1905 default :
1906 width = (21.0-2*rpxmin);
1907 heigth = (29.7-2.*rpymin);
1908 };
1909 }
1910
1911 // mode=5 for lanscape colour
1912 if (fMode == 5) {
1913 rpxmin = 3.4;
1914 rpymin = 0.7;
1915 switch (format) {
1916 case 100 :
1917 width = (11.*2.54)-2.*rpxmin;
1918 heigth = (8.5*2.54)-2.*rpymin;
1919 break;
1920 case 200 :
1921 width = (14.*2.54)-2.*rpxmin;
1922 heigth = (8.5*2.54)-2.*rpymin;
1923 break;
1924 case 300 :
1925 width = (17.*2.54)-2.*rpxmin;
1926 heigth = (11.*2.54)-2.*rpymin;
1927 break;
1928 default :
1929 width = (29.7-2*rpxmin);
1930 heigth = (21-2.*rpymin);
1931 };
1932 }
1933
1934 Double_t value = 0;
1935 if (format < 100) value = 21*TMath::Power(TMath::Sqrt(2.), 4-format);
1936 else if (format == 100) value = 8.5*2.54;
1937 else if (format == 200) value = 8.5*2.54;
1938 else if (format == 300) value = 11.*2.54;
1939 if (format >= 100) format = 4;
1940
1941 // Compute size (in points) of the window for each picture = f(fNXzone,fNYzone)
1944 Int_t npx = 4*CMtoPS(sizex);
1945 Int_t npy = 4*CMtoPS(sizey);
1946 if (sizex > sizey) fMaxsize = CMtoPS(sizex);
1947 else fMaxsize = CMtoPS(sizey);
1948
1949 // Procedure Zone
1950 if (fMode != 3) {
1951 PrintFast(33,"/Zone {/iy exch def /ix exch def ");
1952 PrintFast(10," ix 1 sub ");
1954 PrintFast(5," mul ");
1956 PrintFast(8," iy sub ");
1958 PrintStr(" mul t} def@");
1959 } else {
1960 PrintStr("@");
1961 }
1962
1963 PrintStr("%%EndProlog@");
1964 PrintStr("%%BeginSetup@");
1965 PrintStr("%%EndSetup@");
1966 PrintFast(8,"newpath ");
1967 SaveRestore(1);
1968 if (fMode == 1 || fMode == 4) {
1971 PrintFast(2," t");
1972 }
1973 if (fMode == 2 || fMode == 5) {
1974 PrintFast(7," 90 r 0");
1976 PrintFast(3," t ");
1979 PrintFast(2," t");
1980 }
1981
1982 PrintFast(15," .25 .25 scale ");
1983 if (fMode != 3) {
1984 SaveRestore(1);
1985 PrintStr("@");
1986 PrintStr("%%Page: 1 1@");
1987 SaveRestore(1);
1988 }
1989
1990 //Check is user has defined a special header in the current style
1992 if (nh) {
1994 if (fMode != 3) SaveRestore(1);
1995 }
1996}
1997
1998////////////////////////////////////////////////////////////////////////////////
1999/// Move to a new position
2000
2002{
2003 if (ix != 0 && iy != 0) {
2004 WriteInteger(ix);
2005 WriteInteger(iy);
2006 PrintFast(2," d");
2007 } else if (ix != 0) {
2008 WriteInteger(ix);
2009 PrintFast(2," X");
2010 } else if (iy != 0) {
2011 WriteInteger(iy);
2012 PrintFast(2," Y");
2013 }
2014}
2015
2016////////////////////////////////////////////////////////////////////////////////
2017/// Move to a new PostScript page
2018
2020{
2021 // Compute pad conversion coefficients
2022 if (gPad) {
2023 // if (!gPad->GetPadPaint()) gPad->Update();
2024 Double_t ww = gPad->GetWw();
2025 Double_t wh = gPad->GetWh();
2026 fYsize = fXsize*wh/ww;
2027 } else fYsize = 27;
2028
2029 if(fType == 113 && !fBoundingBox) {
2031 PrintStr("@%%BoundingBox: ");
2032 Double_t xlow=0, ylow=0, xup=1, yup=1;
2033 if (gPad) {
2034 xlow = gPad->GetAbsXlowNDC();
2035 xup = xlow + gPad->GetAbsWNDC();
2036 ylow = gPad->GetAbsYlowNDC();
2037 yup = ylow + gPad->GetAbsHNDC();
2038 }
2039 WriteInteger(CMtoPS(fXsize*xlow));
2040 WriteInteger(CMtoPS(fYsize*ylow));
2043 PrintStr("@");
2044 Initialize();
2046 fPrinted = psave;
2047 }
2048 if (fPrinted) {
2049 if (fSave) SaveRestore(-1);
2050 fClear = kTRUE;
2051 fPrinted = kFALSE;
2052 }
2053 Zone();
2054}
2055
2056////////////////////////////////////////////////////////////////////////////////
2057/// Set the range for the paper in centimeters
2058
2060{
2061 Float_t xps=0, yps=0, xncm=0, yncm=0, dxwn=0, dywn=0, xwkwn=0, ywkwn=0, xymax=0;
2062
2063 fXsize = xsize;
2064 fYsize = ysize;
2065 if( fType != 113) { xps = fXsize; yps = fYsize; }
2066 else { xps = xsize; yps = ysize; }
2067
2068 if( xsize <= xps && ysize < yps) {
2069 if ( xps > yps ) xymax = xps;
2070 else xymax = yps;
2071 xncm = xsize/xymax;
2072 yncm = ysize/xymax;
2073 dxwn = ((xps/xymax)-xncm)/2;
2074 dywn = ((yps/xymax)-yncm)/2;
2075 } else {
2076 if (xps/yps < 1) xwkwn = xps/yps;
2077 else xwkwn = 1;
2078 if (yps/xps < 1) ywkwn = yps/xps;
2079 else ywkwn = 1;
2080
2081 if (xsize < ysize) {
2083 yncm = ywkwn;
2084 dxwn = (xwkwn-xncm)/2;
2085 dywn = 0;
2086 if( dxwn < 0) {
2087 xncm = xwkwn;
2088 dxwn = 0;
2090 dywn = (ywkwn-yncm)/2;
2091 }
2092 } else {
2093 xncm = xwkwn;
2095 dxwn = 0;
2096 dywn = (ywkwn-yncm)/2;
2097 if( dywn < 0) {
2098 yncm = ywkwn;
2099 dywn = 0;
2101 dxwn = (xwkwn-xncm)/2;
2102 }
2103 }
2104 }
2105 fXVP1 = dxwn;
2106 fXVP2 = xncm+dxwn;
2107 fYVP1 = dywn;
2108 fYVP2 = yncm+dywn;
2109 fRange = kTRUE;
2110}
2111
2112////////////////////////////////////////////////////////////////////////////////
2113/// Compute number of gsaves for restore
2114/// This allows to write the correct number of grestore at the
2115/// end of the PS file.
2116
2118{
2119 if (flag == 1) { PrintFast(7," gsave "); fSave++; }
2120 else { PrintFast(4," gr "); fSave--; }
2121}
2122
2123////////////////////////////////////////////////////////////////////////////////
2124/// Set color index for fill areas
2125
2127{
2129 if (gStyle->GetFillColor() <= 0) cindex = 0;
2130}
2131
2132////////////////////////////////////////////////////////////////////////////////
2133/// Patterns definition
2134///
2135/// Define the pattern ipat in the current PS file. ipat can vary from
2136/// 1 to 25. Together with the pattern, the color (color) in which the
2137/// pattern has to be drawn is also required. A pattern is defined in the
2138/// current PS file only the first time it is used. Some level 2
2139/// Postscript functions are used, so on level 1 printers, patterns will
2140/// not work. This is not a big problem because patterns are
2141/// defined only if they are used, so if they are not used a PS level 1
2142/// file will not be polluted by level 2 features, and in any case the old
2143/// patterns used a lot of memory which made them almost unusable on old
2144/// level 1 printers. Finally we should say that level 1 devices are
2145/// becoming very rare. The official PostScript is now level 3 !
2146
2148{
2149 char cdef[28];
2150 char cpat[5];
2151 snprintf(cpat,5," P%2.2d", ipat);
2152
2153 // fPatterns is used as an array of chars. If fPatterns[ipat] != 0 the
2154 // pattern number ipat as already be defined is this file and it
2155 // is not necessary to redefine it. fPatterns is set to zero in Initialize.
2156 // The pattern number 26 allows to know if the macro "cs" has already
2157 // been defined in the current file (see label 200).
2158 if (fPatterns[ipat] == 0) {
2159
2160 // Define the Patterns. Line width must be 1
2161 // Setting fLineWidth to -1 will force the line width definition next time
2162 // TPostScript::SetLineWidth will be called.
2163 fLineWidth = -1;
2164 PrintFast(5," 1 lw");
2165 PrintStr(" << /PatternType 1 /PaintType 2 /TilingType 1");
2166 switch (ipat) {
2167 case 1 :
2168 PrintStr(" /BBox [ 0 0 98 4 ]");
2169 PrintStr(" /XStep 98 /YStep 4");
2170 PrintStr(" /PaintProc { begin gsave");
2171 PrintStr(" [1] 0 sd 2 4 m 99 4 l s 1 3 m 98 3 l s");
2172 PrintStr(" 2 2 m 99 2 l s 1 1 m 98 1 l s");
2173 PrintStr(" gr end } >> [ 4.0 0 0 4.0 0 0 ]");
2174 break;
2175 case 2 :
2176 PrintStr(" /BBox [ 0 0 96 4 ]");
2177 PrintStr(" /XStep 96 /YStep 4");
2178 PrintStr(" /PaintProc { begin gsave");
2179 PrintStr(" [1 3] 0 sd 2 4 m 98 4 l s 0 3 m 96 3 l s");
2180 PrintStr(" 2 2 m 98 2 l s 0 1 m 96 1 l s");
2181 PrintStr(" gr end } >> [ 3.0 0 0 3.0 0 0 ]");
2182 break;
2183 case 3 :
2184 PrintStr(" /BBox [ 0 0 96 16 ]");
2185 PrintStr(" /XStep 96 /YStep 16");
2186 PrintStr(" /PaintProc { begin gsave");
2187 PrintStr(" [1 3] 0 sd 2 13 m 98 13 l s 0 9 m 96 9 l s");
2188 PrintStr(" 2 5 m 98 5 l s 0 1 m 96 1 l s");
2189 PrintStr(" gr end } >> [ 2.0 0 0 2.0 0 0 ]");
2190 break;
2191 case 4 :
2192 PrintStr(" /BBox [ 0 0 100 100 ]");
2193 PrintStr(" /XStep 100 /YStep 100");
2194 PrintStr(" /PaintProc { begin gsave");
2195 PrintStr(" 0 0 m 100 100 l s");
2196 PrintStr(" gr end } >> [ 0.24 0 0 0.24 0 0 ]");
2197 break;
2198 case 5 :
2199 PrintStr(" /BBox [ 0 0 100 100 ]");
2200 PrintStr(" /XStep 100 /YStep 100");
2201 PrintStr(" /PaintProc { begin gsave");
2202 PrintStr(" 0 100 m 100 0 l s");
2203 PrintStr(" gr end } >> [ 0.24 0 0 0.24 0 0 ]");
2204 break;
2205 case 6 :
2206 PrintStr(" /BBox [ 0 0 100 100 ]");
2207 PrintStr(" /XStep 100 /YStep 100");
2208 PrintStr(" /PaintProc { begin gsave");
2209 PrintStr(" 50 0 m 50 100 l s");
2210 PrintStr(" gr end } >> [ 0.12 0 0 0.12 0 0 ]");
2211 break;
2212 case 7 :
2213 PrintStr(" /BBox [ 0 0 100 100 ]");
2214 PrintStr(" /XStep 100 /YStep 100");
2215 PrintStr(" /PaintProc { begin gsave");
2216 PrintStr(" 0 50 m 100 50 l s");
2217 PrintStr(" gr end } >> [ 0.12 0 0 0.12 0 0 ]");
2218 break;
2219 case 8 :
2220 PrintStr(" /BBox [ 0 0 101 101 ]");
2221 PrintStr(" /XStep 100 /YStep 100");
2222 PrintStr(" /PaintProc { begin gsave");
2223 PrintStr(" 0 0 m 0 30 l 30 0 l f 0 70 m 0 100 l 30 100 l f");
2224 PrintStr(" 70 100 m 100 100 l 100 70 l f 70 0 m 100 0 l");
2225 PrintStr(" 100 30 l f 50 20 m 20 50 l 50 80 l 80 50 l f");
2226 PrintStr(" 50 80 m 30 100 l s 20 50 m 0 30 l s 50 20 m");
2227 PrintStr(" 70 0 l s 80 50 m 100 70 l s");
2228 PrintStr(" gr end } >> [ 0.24 0 0 0.24 0 0 ]");
2229 break;
2230 case 9 :
2231 PrintStr(" /BBox [ 0 0 100 100 ]");
2232 PrintStr(" /XStep 100 /YStep 100");
2233 PrintStr(" /PaintProc { begin gsave");
2234 PrintStr(" 0 50 m 50 50 50 180 360 arc");
2235 PrintStr(" 0 50 m 0 100 50 270 360 arc");
2236 PrintStr(" 50 100 m 100 100 50 180 270 arc s");
2237 PrintStr(" gr end } >> [ 0.24 0 0 0.24 0 0 ]");
2238 break;
2239 case 10 :
2240 PrintStr(" /BBox [ 0 0 100 100 ]");
2241 PrintStr(" /XStep 100 /YStep 100");
2242 PrintStr(" /PaintProc { begin gsave");
2243 PrintStr(" 0 50 m 100 50 l 1 1 m 100 1 l");
2244 PrintStr(" 0 0 m 0 50 l 100 0 m 100 50 l");
2245 PrintStr(" 50 50 m 50 100 l s");
2246 PrintStr(" gr end } >> [ 0.24 0 0 0.24 0 0 ]");
2247 break;
2248 case 11 :
2249 PrintStr(" /BBox [ 0 0 100 100 ]");
2250 PrintStr(" /XStep 100 /YStep 100");
2251 PrintStr(" /PaintProc { begin gsave");
2252 PrintStr(" 0 0 m 0 20 l 50 0 m 50 20 l");
2253 PrintStr(" 100 0 m 100 20 l 0 80 m 0 100 l");
2254 PrintStr(" 50 80 m 50 100 l 100 80 m 100 100 l");
2255 PrintStr(" 25 30 m 25 70 l 75 30 m 75 70 l");
2256 PrintStr(" 0 100 m 20 85 l 50 100 m 30 85 l");
2257 PrintStr(" 50 100 m 70 85 l 100 100 m 80 85 l");
2258 PrintStr(" 0 0 m 20 15 l 50 0 m 30 15 l");
2259 PrintStr(" 50 0 m 70 15 l 100 0 m 80 15 l");
2260 PrintStr(" 5 35 m 45 65 l 5 65 m 45 35 l");
2261 PrintStr(" 55 35 m 95 65 l 55 65 m 95 35 l s");
2262 PrintStr(" gr end } >> [ 0.5 0 0 0.5 0 0 ]");
2263 break;
2264 case 12 :
2265 PrintStr(" /BBox [ 0 0 100 100 ]");
2266 PrintStr(" /XStep 100 /YStep 100");
2267 PrintStr(" /PaintProc { begin gsave");
2268 PrintStr(" 0 80 m 0 100 20 270 360 arc");
2269 PrintStr(" 30 100 m 50 100 20 180 360 arc");
2270 PrintStr(" 80 100 m 100 100 20 180 270 arc");
2271 PrintStr(" 20 0 m 0 0 20 0 90 arc");
2272 PrintStr(" 70 0 m 50 0 20 0 180 arc");
2273 PrintStr(" 100 20 m 100 0 20 90 180 arc");
2274 PrintStr(" 45 50 m 25 50 20 0 360 arc");
2275 PrintStr(" 95 50 m 75 50 20 0 360 arc s");
2276 PrintStr(" gr end } >> [ 0.5 0 0 0.5 0 0 ]");
2277 break;
2278 case 13 :
2279 PrintStr(" /BBox [ 0 0 100 100 ]");
2280 PrintStr(" /XStep 100 /YStep 100");
2281 PrintStr(" /PaintProc { begin gsave");
2282 PrintStr(" 0 0 m 100 100 l 0 100 m 100 0 l s");
2283 PrintStr(" gr end } >> [ 0.24 0 0 0.24 0 0 ]");
2284 break;
2285 case 14 :
2286 PrintStr(" /BBox [ 0 0 100 100 ]");
2287 PrintStr(" /XStep 80 /YStep 80");
2288 PrintStr(" /PaintProc { begin gsave");
2289 PrintStr(" 0 20 m 100 20 l 20 0 m 20 100 l");
2290 PrintStr(" 0 80 m 100 80 l 80 0 m 80 100 l");
2291 PrintStr(" 20 40 m 60 40 l 60 20 m 60 60 l");
2292 PrintStr(" 40 40 m 40 80 l 40 60 m 80 60 l s");
2293 PrintStr(" gr end } >> [ 0.60 0 0 0.60 0 0 ]");
2294 break;
2295 case 15 :
2296 PrintStr(" /BBox [ 0 0 60 60 ]");
2297 PrintStr(" /XStep 60 /YStep 60");
2298 PrintStr(" /PaintProc { begin gsave");
2299 PrintStr(" 0 55 m 0 60 5 270 360 arc");
2300 PrintStr(" 25 60 m 30 60 5 180 360 arc");
2301 PrintStr(" 55 60 m 60 60 5 180 270 arc");
2302 PrintStr(" 20 30 m 15 30 5 0 360 arc");
2303 PrintStr(" 50 30 m 45 30 5 0 360");
2304 PrintStr(" arc 5 0 m 0 0 5 0 90 arc");
2305 PrintStr(" 35 0 m 30 0 5 0 180 arc");
2306 PrintStr(" 60 5 m 60 0 5 90 180 arc s");
2307 PrintStr(" gr end } >> [ 0.41 0 0 0.41 0 0 ]");
2308 break;
2309 case 16 :
2310 PrintStr(" /BBox [ 0 0 100 100 ]");
2311 PrintStr(" /XStep 100 /YStep 100");
2312 PrintStr(" /PaintProc { begin gsave");
2313 PrintStr(" 50 50 m 25 50 25 0 180 arc s");
2314 PrintStr(" 50 50 m 75 50 25 180 360 arc s");
2315 PrintStr(" gr end } >> [ 0.4 0 0 0.2 0 0 ]");
2316 break;
2317 case 17 :
2318 PrintStr(" /BBox [ 0 0 100 100 ]");
2319 PrintStr(" /XStep 100 /YStep 100");
2320 PrintStr(" /PaintProc { begin gsave");
2321 PrintStr(" [24] 0 setdash 0 0 m 100 100 l s");
2322 PrintStr(" gr end } >> [ 0.24 0 0 0.24 0 0 ]");
2323 break;
2324 case 18 :
2325 PrintStr(" /BBox [ 0 0 100 100 ]");
2326 PrintStr(" /XStep 100 /YStep 100");
2327 PrintStr(" /PaintProc { begin gsave");
2328 PrintStr(" [24] 0 setdash 0 100 m 100 0 l s");
2329 PrintStr(" gr end } >> [ 0.24 0 0 0.24 0 0 ]");
2330 break;
2331 case 19 :
2332 PrintStr(" /BBox [ 0 0 100 100 ]");
2333 PrintStr(" /XStep 100 /YStep 100");
2334 PrintStr(" /PaintProc { begin gsave");
2335 PrintStr(" 90 50 m 50 50 40 0 360 arc");
2336 PrintStr(" 0 50 m 0 100 50 270 360 arc");
2337 PrintStr(" 50 0 m 0 0 50 0 90 arc");
2338 PrintStr(" 100 50 m 100 0 50 90 180 arc");
2339 PrintStr(" 50 100 m 100 100 50 180 270 arc s");
2340 PrintStr(" gr end } >> [ 0.47 0 0 0.47 0 0 ]");
2341 break;
2342 case 20 :
2343 PrintStr(" /BBox [ 0 0 100 100 ]");
2344 PrintStr(" /XStep 100 /YStep 100");
2345 PrintStr(" /PaintProc { begin gsave");
2346 PrintStr(" 50 50 m 50 75 25 270 450 arc s");
2347 PrintStr(" 50 50 m 50 25 25 90 270 arc s");
2348 PrintStr(" gr end } >> [ 0.2 0 0 0.4 0 0 ]");
2349 break;
2350 case 21 :
2351 PrintStr(" /BBox [ 0 0 101 101 ]");
2352 PrintStr(" /XStep 100 /YStep 100");
2353 PrintStr(" /PaintProc { begin gsave");
2354 PrintStr(" 1 1 m 25 1 l 25 25 l 50 25 l 50 50 l");
2355 PrintStr(" 75 50 l 75 75 l 100 75 l 100 100 l");
2356 PrintStr(" 50 1 m 75 1 l 75 25 l 100 25 l 100 50 l");
2357 PrintStr(" 0 50 m 25 50 l 25 75 l 50 75 l 50 100 l s");
2358 PrintStr(" gr end } >> [ 0.5 0 0 0.5 0 0 ]");
2359 break;
2360 case 22 :
2361 PrintStr(" /BBox [ 0 0 101 101 ]");
2362 PrintStr(" /XStep 100 /YStep 100");
2363 PrintStr(" /PaintProc { begin gsave");
2364 PrintStr(" 1 100 m 25 100 l 25 75 l 50 75 l 50 50 l");
2365 PrintStr(" 75 50 l 75 25 l 100 25 l 100 1 l");
2366 PrintStr(" 50 100 m 75 100 l 75 75 l 100 75 l 100 50 l");
2367 PrintStr(" 0 50 m 25 50 l 25 25 l 50 25 l 50 1 l s");
2368 PrintStr(" gr end } >> [ 0.5 0 0 0.5 0 0 ]");
2369 break;
2370 case 23 :
2371 PrintStr(" /BBox [ 0 0 100 100 ]");
2372 PrintStr(" /XStep 100 /YStep 100");
2373 PrintStr(" /PaintProc { begin gsave");
2374 PrintStr(" [1 7] 0 sd 0 8 50 { dup dup m 2 mul 0 l s } for");
2375 PrintStr(" 0 8 50 { dup dup 2 mul 100 m 50 add exch 50");
2376 PrintStr(" add l s } for 100 0 m 100 100 l 50 50 l f");
2377 PrintStr(" gr end } >> [ 0.24 0 0 0.24 0 0 ]");
2378 break;
2379 case 24 :
2380 PrintStr(" /BBox [ 0 0 100 100 ]");
2381 PrintStr(" /XStep 100 /YStep 100");
2382 PrintStr(" /PaintProc { begin gsave");
2383 PrintStr(" 100 100 m 100 36 l 88 36 l 88 88 l f");
2384 PrintStr(" 100 0 m 100 12 l 56 12 l 50 0 l f");
2385 PrintStr(" 0 0 m 48 0 l 48 48 l 50 48 l 56 60 l");
2386 PrintStr(" 36 60 l 36 12 l 0 12 l f [1 7] 0 sd");
2387 PrintStr(" 61 8 87 { dup dup dup 12 exch m 88 exch l s");
2388 PrintStr(" 16 exch 4 sub m 88 exch 4 sub l s } for");
2389 PrintStr(" 13 8 35 { dup dup dup 0 exch m 36 exch l s");
2390 PrintStr(" 4 exch 4 sub m 36 exch 4 sub l s } for");
2391 PrintStr(" 37 8 59 { dup dup dup 12 exch m 36 exch l s");
2392 PrintStr(" 16 exch 4 sub m 36 exch 4 sub l s } for");
2393 PrintStr(" 13 8 60 { dup dup dup 56 exch m 100 exch l s");
2394 PrintStr(" 60 exch 4 sub m 100 exch 4 sub l s } for");
2395 PrintStr(" gr end } >> [ 0.5 0 0 0.5 0 0 ]");
2396 break;
2397 case 25 :
2398 PrintStr(" /BBox [ 0 0 101 101 ]");
2399 PrintStr(" /XStep 100 /YStep 100");
2400 PrintStr(" /PaintProc { begin gsave");
2401 PrintStr(" 0 0 m 30 30 l 70 30 l 70 70 l 100 100 l 100 0 l");
2402 PrintStr(" f 30 30 m 30 70 l 70 70 l f");
2403 PrintStr(" gr end } >> [ 0.5 0 0 0.5 0 0 ]");
2404 };
2405 snprintf(cdef,28," makepattern /%s exch def",&cpat[1]);
2406 PrintStr(cdef);
2407 fPatterns[ipat] = 1;
2408 }
2409
2410 // Define the macro cs and FA if they are not yet defined.
2411 if (fPatterns[26] == 0) {
2412 if (gStyle->GetColorModelPS()) {
2413 PrintStr(" /cs {[/Pattern /DeviceCMYK] setcolorspace} def");
2414 PrintStr(" /FA {f [/DeviceCMYK] setcolorspace} def");
2415 } else {
2416 PrintStr(" /cs {[/Pattern /DeviceRGB] setcolorspace} def");
2417 PrintStr(" /FA {f [/DeviceRGB] setcolorspace} def");
2418 }
2419 fPatterns[26] = 1;
2420 }
2421
2422 // Activate the pattern.
2423 PrintFast(3," cs");
2424 TColor *col = gROOT->GetColor(color);
2425 if (col) {
2426 Double_t colRed = col->GetRed();
2427 Double_t colGreen = col->GetGreen();
2428 Double_t colBlue = col->GetBlue();
2429 if (gStyle->GetColorModelPS()) {
2431 if (colBlack==1) {
2432 WriteReal(0);
2433 WriteReal(0);
2434 WriteReal(0);
2436 } else {
2444 }
2445 } else {
2449 }
2450 }
2451 PrintFast(4,cpat);
2452 PrintFast(9," setcolor");
2453}
2454
2455////////////////////////////////////////////////////////////////////////////////
2456/// Set color index for lines
2457
2462
2463////////////////////////////////////////////////////////////////////////////////
2464/// Set the value of the global parameter TPostScript::fgLineJoin.
2465/// This parameter determines the appearance of joining lines in a PostScript
2466/// output.
2467/// It takes one argument which may be:
2468/// - 0 (miter join)
2469/// - 1 (round join)
2470/// - 2 (bevel join)
2471/// The default value is 0 (miter join).
2472///
2473/// \image html postscript_1.png
2474///
2475/// To change the line join behaviour just do:
2476/// ~~~ {.cpp}
2477/// gStyle->SetJoinLinePS(2); // Set the PS line join to bevel.
2478/// ~~~
2479
2481{
2483 if (fgLineJoin<0) fgLineJoin=0;
2484 if (fgLineJoin>2) fgLineJoin=2;
2485}
2486
2487////////////////////////////////////////////////////////////////////////////////
2488/// Set the value of the global parameter TPostScript::fgLineCap.
2489/// This parameter determines the appearance of line caps in a PostScript
2490/// output.
2491/// It takes one argument which may be:
2492/// - 0 (butt caps)
2493/// - 1 (round caps)
2494/// - 2 (projecting caps)
2495/// The default value is 0 (butt caps).
2496///
2497/// \image html postscript_2.png
2498///
2499/// To change the line cap behaviour just do:
2500/// ~~~ {.cpp}
2501/// gStyle->SetCapLinePS(2); // Set the PS line cap to projecting.
2502/// ~~~
2503
2505{
2507 if (fgLineCap<0) fgLineCap=0;
2508 if (fgLineCap>2) fgLineCap=2;
2509}
2510
2511////////////////////////////////////////////////////////////////////////////////
2512/// Change the line style
2513///
2514/// - linestyle = 2 dashed
2515/// - linestyle = 3 dotted
2516/// - linestyle = 4 dash-dotted
2517/// - linestyle = else = solid
2518///
2519/// See TStyle::SetLineStyleString for style definition
2520
2526
2527////////////////////////////////////////////////////////////////////////////////
2528/// Change the line style in the output file
2529
2531{
2532 if (linestyle == fStyle) return;
2533
2534 fStyle = linestyle;
2535
2536 const char *st = gStyle->GetLineStyleString(fStyle);
2537 PrintFast(1,"[");
2538 Int_t nch = strlen(st);
2539 PrintFast(nch,st);
2540 PrintFast(6,"] 0 sd");
2541}
2542
2543////////////////////////////////////////////////////////////////////////////////
2544/// Change the line width
2545
2551
2552////////////////////////////////////////////////////////////////////////////////
2553/// Change the line width in the output file
2554
2556{
2557 if (linewidth == fWidth) return;
2558
2559 fWidth = linewidth;
2560
2561 if (fWidth!=0) {
2563 PrintFast(3," lw");
2564 }
2565}
2566
2567////////////////////////////////////////////////////////////////////////////////
2568/// Set color index for markers
2569
2574
2575////////////////////////////////////////////////////////////////////////////////
2576/// Set the current color.
2577
2579{
2580 if (color < 0) color = 0;
2581 TColor *col = gROOT->GetColor(color);
2582 if (col)
2583 SetColor(col->GetRed(), col->GetGreen(), col->GetBlue());
2584 else
2585 SetColor(1., 1., 1.);
2586}
2587
2588////////////////////////////////////////////////////////////////////////////////
2589/// Set directly current color (don't go via TColor).
2590
2592{
2593 if (r == fRed && g == fGreen && b == fBlue) return;
2594
2595 fRed = r;
2596 fGreen = g;
2597 fBlue = b;
2598
2599 if (fRed <= 0 && fGreen <= 0 && fBlue <= 0 ) {
2600 PrintFast(6," black");
2601 } else {
2602 if (gStyle->GetColorModelPS()) {
2611 } else {
2612 WriteReal(fRed);
2615 }
2616 PrintFast(2," c");
2617 }
2618}
2619
2620////////////////////////////////////////////////////////////////////////////////
2621/// Set color index for text
2622
2627
2628////////////////////////////////////////////////////////////////////////////////
2629/// Write a string of characters
2630///
2631/// This method writes the string chars into a PostScript file
2632/// at position xx,yy in world coordinates.
2633
2635{
2636 static const char *psfont[31][2] = {
2637 { "Root.PSFont.1", "/Times-Italic" },
2638 { "Root.PSFont.2", "/Times-Bold" },
2639 { "Root.PSFont.3", "/Times-BoldItalic" },
2640 { "Root.PSFont.4", "/Helvetica" },
2641 { "Root.PSFont.5", "/Helvetica-Oblique" },
2642 { "Root.PSFont.6", "/Helvetica-Bold" },
2643 { "Root.PSFont.7", "/Helvetica-BoldOblique" },
2644 { "Root.PSFont.8", "/Courier" },
2645 { "Root.PSFont.9", "/Courier-Oblique" },
2646 { "Root.PSFont.10", "/Courier-Bold" },
2647 { "Root.PSFont.11", "/Courier-BoldOblique" },
2648 { "Root.PSFont.12", "/Symbol" },
2649 { "Root.PSFont.13", "/Times-Roman" },
2650 { "Root.PSFont.14", "/ZapfDingbats" },
2651 { "Root.PSFont.15", "/Symbol" },
2652 { "Root.PSFont.STIXGen", "/STIXGeneral" },
2653 { "Root.PSFont.STIXGenIt", "/STIXGeneral-Italic" },
2654 { "Root.PSFont.STIXGenBd", "/STIXGeneral-Bold" },
2655 { "Root.PSFont.STIXGenBdIt", "/STIXGeneral-BoldItalic" },
2656 { "Root.PSFont.STIXSiz1Sym", "/STIXSize1Symbols" },
2657 { "Root.PSFont.STIXSiz1SymBd", "/STIXSize1Symbols-Bold" },
2658 { "Root.PSFont.STIXSiz2Sym", "/STIXSize2Symbols" },
2659 { "Root.PSFont.STIXSiz2SymBd", "/STIXSize2Symbols-Bold" },
2660 { "Root.PSFont.STIXSiz3Sym", "/STIXSize3Symbols" },
2661 { "Root.PSFont.STIXSiz3SymBd", "/STIXSize3Symbols-Bold" },
2662 { "Root.PSFont.STIXSiz4Sym", "/STIXSize4Symbols" },
2663 { "Root.PSFont.STIXSiz4SymBd", "/STIXSize4Symbols-Bold" },
2664 { "Root.PSFont.STIXSiz5Sym", "/STIXSize5Symbols" },
2665 { "Root.PSFont.ME", "/DroidSansFallback" },
2666 { "Root.PSFont.CJKMing", "/DroidSansFallback" },
2667 { "Root.PSFont.CJKGothic", "/DroidSansFallback" }
2668 };
2669
2670 const Double_t kDEGRAD = TMath::Pi()/180.;
2671 Double_t x = xx;
2672 Double_t y = yy;
2673 if (!gPad) return;
2674
2675 // Compute the font size. Exit if it is 0
2676 // The font size is computed from the TTF size to get exactly the same
2677 // size on the screen and in the PostScript file.
2678 Double_t wh = (Double_t) gPad->GetPadWidth();
2679 Double_t hh = (Double_t) gPad->GetPadHeight();
2680 if (wh <= 0 || hh <= 0)
2681 return;
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:87
#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:52
Font_t fTextFont
Text font.
Definition TAttText.h:28
virtual void SetTextSize(Float_t tsize=1)
Set the text size.
Definition TAttText.h:53
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:720
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:1311
virtual int Rename(const char *from, const char *to)
Rename a file.
Definition TSystem.cxx:1365
virtual char * Which(const char *search, const char *file, EAccessMode mode=kFileExists)
Find location of file in a search path.
Definition TSystem.cxx:1563
virtual int Unlink(const char *name)
Unlink, i.e.
Definition TSystem.cxx:1396
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:551
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:568
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.
void CloseStream()
Close existing stream.
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
Bool_t OpenStream(const char *fname, Bool_t binary=kFALSE)
Open output stream.
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