Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TSpectrum2Painter.cxx
Go to the documentation of this file.
1// @(#)root/spectrumpainter:$Id: TSpectrum2Painter.cxx,v 1.00
2// Author: Miroslav Morhac 29/09/06
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
13/** \class TSpectrum2Painter
14 \ingroup Spectrumpainter
15
16\legacy{TSpectrum2Painter}
17
18Two-dimensional graphics function
19
20TSpectrum2Painter is a set of graphical functions developed by Miroslav
21Morhac to paint 2D-histograms in three dimensions. This package is accessed
22via THistPainter in a transparent way. For the ROOT user it is enough to use
23the "SPEC" option to draw a 2D-Histogram. This option offers many
24functionalities detailed in the header of the PaintSpectrum function.
25
26Reference:
27
28Morhac M., Kliman J., Matousek V., Turzo I.: Sophisticated visualization
29algorithms for analysis of multidimensional experimental nuclear data. Acta
30Pysica Slovaca Vol. 54/ 4 (2004), pp. 385-400.
31*/
32
33#include <climits>
34
35#include "TROOT.h"
36#include "TColor.h"
37#include "TMath.h"
38#include "TLine.h"
39#include "TEllipse.h"
40#include "TVirtualPad.h"
41#include "TBox.h"
42#include "TF1.h"
43#include "TH2.h"
44#include "TGaxis.h"
45#include "THLimitsFinder.h"
46#include "TSpectrum2Painter.h"
47#include "strlcpy.h"
48
50
51
52////////////////////////////////////////////////////////////////////////////////
53/// TSpectrum2Painter normal constructor
54
56 : TNamed ("Spectrum Painter2","Miroslav Morhac Painter")
57{
58 int i, j;
59 double val;
60 gPad->Range(0, 0, 1 ,1);
61 fXmin = 0;
62 fXmax = h2->GetNbinsX() - 1;
63 fYmin = 0;
64 fYmax = h2->GetNbinsY() - 1;
65 fZmin = 0, fZmax = 0;
66 fMaximumXScreenResolution = bs;
67
68 for (i = 0;i <= fXmax; i++) {
69 for (j = 0;j <= fYmax; j++) {
70 val = h2->GetBinContent(i + 1,j + 1);
71 if (val > fZmax) fZmax = val;
72 }
73 }
74
75 fBx1 = gPad->XtoPixel(0.1); //axis positions
76 fBx2 = gPad->XtoPixel(0.99);
77 fBy1 = gPad->YtoPixel(0.99);
78 fBy2 = gPad->YtoPixel(0.05);
79
80 fModeGroup = kModeGroupLightHeight;
81
82 fDisplayMode = kDisplayModeSurface;
83
84 fZscale = kZScaleLinear; // Scale linear, log.
85
86 fNodesx = fXmax-fXmin+1; // Number of nodes in x dimension of grid
87 fNodesy = fYmax-fYmin+1; // Number of nodes in y dimension of grid
88
89 fContWidth = 50; // Width between contours,
90 // applies only for contours display mode.
91 fAlpha = 20; // Angles of display,alfa+beta must be less or equal to 90,
92 // alpha- angle between base line of Canvas and left lower
93 // edge of picture picture base plane
94 fBeta = 60; // Angle between base line of Canvas and right lower edge
95 // of picture base plane
96 fViewAngle = 0; // Rotation angle of the view,
97 // it can be 0, 90, 180, 270 degrees.
98
99 fLevels = 256; // Number of color levels for rainbowed display modes,
100 // It does not apply for simple display modes
101 // algorithm group
102 fRainbow1Step = 1; // Determines the first component step for
103 // neighbouring color levels, applies only for
104 // rainbowed display modes, it does not apply for
105 // simple display modes algorithm group.
106 fRainbow2Step = 1; // Determines the second component step for
107 // neighbouring color levels, applies only for
108 // rainbowed display modes, it does not apply for
109 // simple display modes algorithm group.
110 fRainbow3Step = 1; // Determines the third component step for
111 // neighbouring color levels, applies only for
112 // rainbowed display modes, it does not apply for
113 // simple display modes algorithm group.
114
115 fColorAlg = kColorAlgRgbSmooth; // Applies only for rainbowed display modes
116 // (rgb smooth algorithm, rgb modulo color
117 // component, cmy smooth algorithm, cmy
118 // modulo color component, cie smooth
119 // algorithm, cie modulo color component,
120 // yiq smooth algorithm, yiq modulo color
121 // component, hsv smooth algorithm, hsv
122 // modulo color component, it does not
123 // apply for simple display modes
124 // algorithm group.
125
126 fLHweight = 0.5; // Weight between shading according to fictive light
127 // source and according to channels counts, applies only
128 // for kModeGroupLightHeight modes group.
129
130 fXlight = 1000; // X position of fictive light source, applies only for
131 // rainbowed display modes with shading according to light.
132 fYlight = 1000; // Y position of fictive light source, applies only for
133 // rainbowed display modes with shading according to light.
134 fZlight = 1000; // Z position of fictive light source, applies only for
135 // rainbowed display modes with shading according to light.
136
137 fShadow = kShadowsNotPainted; // Determines whether shadow will be drawn
138 // (no shadow, shadow), for rainbowed
139 // display modes with shading according to
140 // light.
141
142 fShading = kShaded; // Determines whether the picture will shaded,
143 // smoothed (no shading, shading), for rainbowed
144 // display modes only.
145
146 fBezier = kNoBezierInterpol; // Determines Bezier interpolation (applies
147 // only for simple display modes group for
148 // grid, x_lines, y_lines display modes).
149
150 fPenColor = kBlack; // Color of spectrum.
151 fPenWidth = 1; // Width of line.
152 fPenDash = kPenStyleSolid; // Style of pen.
153
154 fChanmarkEnDis = kChannelMarksNotDrawn; // Decides whether the channel
155 // marks are shown.
156 fChanmarkColor = kBlue; // Color of channel marks.
157 fChanmarkWidth = 8; // Width of channel marks.
158 fChanmarkHeight = 8; // Height of channel marks.
159 fChanmarkStyle = kChannelMarksStyleDot; // Style of channel marks.
160
161 fChanlineEnDis = kChannelGridNotDrawn; // Decides whether the channel lines
162 // (grid) are shown.
163 fChanlineColor = kRed; // Color of channel marks.
164 fNewColor = nullptr;
165 fEnvelope = new Short_t [fMaximumXScreenResolution];
166 fEnvelopeContour = new Short_t [fMaximumXScreenResolution];
167 for (i=0;i<fMaximumXScreenResolution;i++) {
168 fEnvelope[i] = fBy2;
169 fEnvelopeContour[i] = fBy2;
170 }
171 fH2 = h2;
172}
173
174
175////////////////////////////////////////////////////////////////////////////////
176/// TSpectrum2Painter destructor
177
179{
180 TColor* col;
181 for (int i=0; i<256; i++) {
182 col = gROOT->GetColor(250+i);
183 if (col) delete col;
184 }
185 if (fEnvelope) delete [] fEnvelope;
186 if (fEnvelopeContour) delete [] fEnvelopeContour;
187}
188
189
190////////////////////////////////////////////////////////////////////////////////
191/// Reads out the value from histogram and calculates screen coordinates
192///
193/// Parameters:
194/// - it - node in x- direction
195/// - jt - node in y- direction
196/// - zmt - control variable
197
199{
200 Int_t lxt,lyt,ix,iy;
201 Double_t zf = 0;
202 Double_t p1,p2;
203 p1 = fXmin+fKx*(Double_t)it;
204 p2 = fYmin+fKy*(Double_t)jt;
205 ix = (Int_t)p1;
206 iy = (Int_t)p2;
207 fDxspline = p1;
208 fDyspline = p2;
209 if ((zmt==0)||(zmt==-3)||(zmt==-4)) {
210 zf = fH2->GetBinContent(ix+1,iy+1);
211 } else if (zmt==-2) zf = fZPresetValue;
212 if (zf<fZmin) zf = fZmin;
213 fZeq = zf;
214 switch (fZscale) {
215 case kZScaleLog:
216 if (zf>=1.0) zf = log(zf);
217 else zf = 0;
218 break;
219 case kZScaleSqrt:
220 if (zf>0) zf = sqrt(zf);
221 else zf = 0;
222 break;
223 }
224 lxt = (Int_t)(fTxx*(Double_t)it+fTxy*(Double_t)jt+fVx);
225 lyt = (Int_t)(fTyx*(Double_t)it+fTyy*(Double_t)jt+fTyz*zf+fVy);
226 if (lxt<fBx1) lxt = fBx1;
227 if (lxt>fBx2) lxt = fBx2;
228 if (lyt<fBy1) lyt = fBy1;
229 if (lyt>fBy2) lyt = fBy2;
230 fXt = lxt;
231 fYt = lyt;
232 fZ = zf;
233 return;
234}
235
236
237////////////////////////////////////////////////////////////////////////////////
238/// Calculates and returns color value for the surface triangle
239/// given by function parameters:
240/// -dx1,dy1,z1 coordinates of the first point in 3d space
241/// -dx2,dy2,z2 coordinates of the second point in 3d space
242/// -dx3,dy3,z3 coordinates of the third point in 3d space
243
245 Double_t dx1, Double_t dy1, Double_t z1,
246 Double_t dx2, Double_t dy2, Double_t z2,
247 Double_t dx3, Double_t dy3, Double_t z3)
248{
249 Double_t da,db,dc=0,dd,dl,dm,dn,xtaz,ytaz,ztaz,v=0,v1;
250 Double_t pi=3.1415927;
251 Int_t i;
252 switch (fZscale) {
253 case kZScaleLog:
254 if (z1>900) z1 = 900;
255 z1 = exp(z1);
256 if (z2>900) z2 = 900;
257 z2 = exp(z2);
258 if (z3>900) z3 = 900;
259 z3 = exp(z3);
260 break;
261 case kZScaleSqrt:
262 z1 = z1*z1;
263 z2 = z2*z2;
264 z3 = z3*z3;
265 break;
266 }
267 i = fViewAngle;
268 i = i/90;
269 if ((i==1)||(i==3)) {
270 da = dx1;
271 dx1 = dx2;
272 dx2 = da;
273 da = dy1;
274 dy1 = dy2;
275 dy2 = da;
276 da = z1;
277 z1 = z2;
278 z2 = da;
279 }
280 xtaz = (dx1+dx2+dx3)/3;
281 ytaz = (dy1+dy2+dy3)/3;
282 ztaz = (z1+z2+z3)/3;
284 dn = (Double_t)fZlight-ztaz;
285 dm = (Double_t)fYlight-ytaz;
286 dl = (Double_t)fXlight-xtaz;
287 da = (dy2-dy1)*(z3-z1)-(dy3-dy1)*(z2-z1);
288 db = (z2-z1)*(dx3-dx1)-(z3-z1)*(dx2-dx1);
289 dc = (dx2-dx1)*(dy3-dy1)-(dx3-dx1)*(dy2-dy1);
290 dd = (da*da+db*db+dc*dc)*(dl*dl+dm*dm+dn*dn);
291 dd = sqrt(dd);
292 if (dd!=0) v = (da*dl+db*dm+dc*dn)/dd;
293 else v = 0;
294 if (v<-1) v=-1;
295 if (v>1) v=1;
296 v = asin(v);
297 v = v+pi/2;
298 v = v/pi;
299 } else if (fModeGroup==kModeGroupHeight) {
300 da = fZmax-fZmin;
301 if (ztaz<fZmin) ztaz=fZmin;
302 if (ztaz>=fZmax) ztaz=fZmax-1;
303 db = ztaz-fZmin;
304 if (da!=0) {
305 switch (fZscale) {
306 case kZScaleLinear:
307 dc = db/da;
308 break;
309 case kZScaleLog:
310 if (da>=1) da=log(da);
311 if (db>=1) db=log(db);
312 if (da!=0) dc=db/da;
313 else dc=0;
314 break;
315 case kZScaleSqrt:
316 da = sqrt(da);
317 db = sqrt(db);
318 dc = db/da;
319 break;
320 }
321 } else {
322 dc=0;
323 }
324 i = (Int_t)dc;
325 v = dc-i;
326 } else if (fModeGroup==kModeGroupLightHeight) {
327 dn = (Double_t)fZlight-ztaz;
328 dm = (Double_t)fYlight-ytaz;
329 dl = (Double_t)fXlight-xtaz;
330 da = (dy2-dy1)*(z3-z1)-(dy3-dy1)*(z2-z1);
331 db = (z2-z1)*(dx3-dx1)-(z3-z1)*(dx2-dx1);
332 dc = (dx2-dx1)*(dy3-dy1)-(dx3-dx1)*(dy2-dy1);
333 dd = (da*da+db*db+dc*dc)*(dl*dl+dm*dm+dn*dn);
334 dd = sqrt(dd);
335 if (dd!=0) v = (da*dl+db*dm+dc*dn)/dd;
336 else v = 0;
337 if (v<-1) v=-1;
338 if (v>1) v=1;
339 v = asin(v);
340 v = v+pi/2;
341 v = v/pi;
342 da = fZmax-fZmin;
343 if (ztaz<fZmin) ztaz = fZmin;
344 if (ztaz>=fZmax) ztaz = fZmax-1;
345 db = ztaz-fZmin;
346 if (da!=0) {
347 switch (fZscale) {
348 case kZScaleLinear:
349 dc = db/da;
350 break;
351 case kZScaleLog:
352 if (da>=1) da = log(da);
353 if (db>=1) db = log(db);
354 if (da!=0) dc = db/da;
355 else dc = 0;
356 break;
357 case kZScaleSqrt:
358 da = sqrt(da);
359 db = sqrt(db);
360 dc = db/da;
361 break;
362 }
363 } else {
364 dc = 0;
365 }
366 i = (Int_t)dc;
367 v1 = dc-i;
368 v = fLHweight*v+(1-fLHweight)*v1;
369 }
371 da = 1.0/(Double_t)fLevels;
372 if (v<da) v = da;
373 } else {
374 da = 2.0/(Double_t)fLevels;
375 if (v<da) v = da;
376 }
377 return(v);
378}
379
380
381////////////////////////////////////////////////////////////////////////////////
382/// Determines whether the center of the triangle in 3-d space
383/// given by function parameters:
384/// - xtaz,ytaz,ztaz
385/// is in shadow or not. If yes it return 1 otherwise it returns 0.
386
388 Double_t ztaz,
389 Double_t shad_noise)
390{
391 Int_t sx2,sy2,sz1,sz2,skrokx,skroky,xmax,ymax;
392 Double_t sx1,sy1;
393 Double_t pom1,pom2,sdx1=0,sdx2=0,sdy1,sdy2,spriz;
394 switch (fZscale) {
395 case kZScaleLog:
396 if (ztaz>900) ztaz = 900;
397 ztaz = exp(ztaz);
398 if (ztaz>32767) ztaz = 32767;
399 break;
400 case kZScaleSqrt:
401 ztaz = ztaz*ztaz;
402 break;
403 }
404 spriz = 0;
405 sx1 = xtaz;
406 sy1 = ytaz;
407 sz1 = (Int_t)ztaz;
408 sx2 = fXlight;
409 sy2 = fYlight;
410 sz2 = fZlight;
411 xmax = fXmax;
412 ymax = fYmax;
413 if (sx1!=sx2) {
414 if (sx1<sx2) skrokx = 1;
415 else skrokx = -1;
416 if (sy1<sy2) skroky = 1;
417 else skroky = -1;
418 pom1 = sx2-sx1;
419 pom2 = sy2-sy1;
420 if (TMath::Abs(pom1)>0.0000001) sdx1 = pom2/pom1;
421 pom1 = sx1;
422 pom2 = sy1;
423 sdy1 = pom2-sdx1*pom1;
424 pom1 = sx2-sx1;
425 pom2 = sz2-sz1;
426 if (TMath::Abs(pom1)>0.0000001) sdx2 = pom2/pom1;
427 pom1 = sx1;
428 pom2 = sz1;
429 sdy2 = pom2-sdx2*pom1;
430 spriz = 0;
431 pom1 = sx1;
432 pom2 = pom1*sdx1+sdy1;
433 sy1 = pom2;
434 for (;(sx1>(fXmin-skrokx)) && (sx1<(xmax-skrokx)) &&
435 (sy1>(fYmin-skroky)) && (sy1<(ymax-skroky)) &&
436 (spriz==0);sx1+=skrokx) {
437 pom1 = sx1;
438 pom2 = pom1*sdx1+sdy1;
439 sy1 = pom2+skroky;
440 if ((sy1>=fYmin)&&(sy1<=fYmax)) {
441 sz1 = (Int_t)(fH2->GetBinContent((Int_t)sx1+1,(Int_t)sy1+1));
442 pom2 = pom1*sdx2+sdy2;
443 sz2 = (Int_t)(pom2+shad_noise);
444 if (sz1>sz2) spriz = 1;
445 }
446 }
447 } else if (sy1!=sy2) {
448 if (sy1<sy2) skroky = 1;
449 else skroky = -1;
450 pom1 = sy2-sy1;
451 pom2 = sz2-sz1;
452 if (TMath::Abs(pom1)>0.0000001) sdx2 = pom2/pom1;
453 pom1 = sy1;
454 pom2 = sz1;
455 sdy2 = pom2-sdx2*pom1;
456 spriz = 0;
457 for (;(sy1>(fYmin-skroky)) && (sy1<(ymax-skroky)) &&
458 (spriz==0);sy1+=skroky) {
459 sz1 = (Int_t)(fH2->GetBinContent((Int_t)sx1+1,(Int_t)sy1+1));
460 pom1 = sy1;
461 pom2 = pom1*sdx2+sdy2;
462 sz2 = (Int_t)(pom2+shad_noise);
463 if (sz1>sz2) spriz=1;
464 }
465 }
466 return(spriz);
467}
468
469
470////////////////////////////////////////////////////////////////////////////////
471/// This function calculates color for one palette entry
472/// given by function parameter ui. Other parameters
473/// -ui1,ui2,ui3
474/// represent r, g, b color components of the basic pen color.
475
476void TSpectrum2Painter::ColorModel(unsigned ui, unsigned ui1, unsigned ui2,
477 unsigned ui3)
478{
479 unsigned long uinc1=0,uinc2=0,uinc3=0,upom,i;
480 Double_t a,b,c,d,h,v,s,f;
481 Int_t j,iv=ui;
482 Double_t red=0,green=0,blue=0;
483 if (iv<0) iv = 0;
484 else if (iv>255) iv = 255;
485 if (gROOT->GetColor(250+iv)) {
486 fNewColorIndex = 250+iv;
487 return;
488 }
489 if (fColorAlg%2==0) {
491 a = ui*a;
492 a = ui1+a;
493 if (a >= UINT_MAX) uinc1 = UINT_MAX;
494 else uinc1 = (unsigned)a;
495 upom = uinc1%256;
496 i = (uinc1-upom)/256;
497 if ((i%2)==0) uinc1 = upom;
498 else uinc1 = 255-upom;
500 b = ui*b;
501 b = ui2+b;
502 uinc2 = (Int_t)b;
503 upom = uinc2%256;
504 i = (uinc2-upom)/256;
505 if ((i%2)==0) uinc2 = upom;
506 else uinc2 = 255-upom;
508 c = ui*c;
509 c = ui3+c;
510 uinc3 = (Int_t)c;
511 upom = uinc3%256;
512 i = (uinc3-upom)/256;
513 if ((i%2)==0) uinc3 = upom;
514 else uinc3 = 255-upom;
516 uinc1 = 255-uinc1;
517 uinc2 = 255-uinc2;
518 uinc3 = 255-uinc3;
519 } else if (fColorAlg==kColorAlgCieSmooth) {
520 a = uinc1;
521 b = uinc2;
522 c = uinc3;
523 d = a+b+c;
524 if (d!=0) {
525 a = a/d;
526 b = b/d;
527 c = c/d;
528 }
529 red = a*255;
530 green = b*255;
531 blue = c*255;
532 uinc1 = (Int_t)red;
533 uinc2 = (Int_t)green;
534 uinc3 = (Int_t)blue;
535 } else if (fColorAlg==kColorAlgYiqSmooth) {
536 a = uinc1;
537 b = uinc2;
538 c = uinc3;
539 a = a/256;
540 b = b/256;
541 c = c/256;
542 red = a+0.956*b+0.62*c;
543 green = a-0.272*b-0.647*c;
544 blue = a-1.108*b+1.705*c;
545 if (red>=2) red = red-2;
546 else if (red>=1) red = 2-red;
547 if (green<0) green = -green;
548 if (blue>=2) blue = blue-2;
549 else if (blue>=1) blue = 2-blue;
550 else if (blue<-1) blue = 2+blue;
551 else if (blue<0) blue = -blue;
552 red = red*255;
553 green = green*255;
554 blue = blue*255;
555 uinc1 = (Int_t)red;
556 uinc2 = (Int_t)green;
557 uinc3 = (Int_t)blue;
558 } else if (fColorAlg==kColorAlgHvsSmooth) {
559 h = uinc1;
560 v = uinc2;
561 s = uinc3;
562 h = h/256;
563 v = v/256;
564 s = s/256;
565 if (s==0) {
566 red = v;
567 green = v;
568 blue = v;
569 } else {
570 if (h==1.0) h=0;
571 h = 6.0*h;
572 j = (Int_t)h;
573 f = h-j;
574 a = v*(1-s);
575 b = v*(1-s*f);
576 c = v*(1-s*(1-f));
577 switch (j) {
578 case 0:
579 red = v;
580 green = c;
581 blue = a;
582 break;
583 case 1:
584 red = b;
585 green = v;
586 blue = a;
587 break;
588 case 2:
589 red = a;
590 green = v;
591 blue = c;
592 break;
593 case 3:
594 red = a;
595 green = b;
596 blue = v;
597 break;
598 case 4:
599 red = c;
600 green = a;
601 blue = v;
602 break;
603 case 5:
604 red = v;
605 green = a;
606 blue = b;
607 break;
608 }
609 }
610 red = red*255;
611 green = green*255;
612 blue = blue*255;
613 uinc1 = (Int_t)red;
614 uinc2 = (Int_t)green;
615 uinc3 = (Int_t)blue;
616 }
617 ui = uinc1+uinc2*256+uinc3*65536;
618 } else if (fColorAlg%2==1) {
620 a = ui*a;
621 a = ui1/2+a;
622 uinc1 = (Int_t)a;
623 uinc1 = uinc1%256;
625 b = ui*b;
626 b = ui2/2+b;
627 uinc2 = (Int_t)b;
628 uinc2 = uinc2%256;
630 c = ui*c;
631 c = ui3/2+c;
632 uinc3 = (Int_t)c;
633 uinc3 = uinc3%256;
635 uinc1 = 255-uinc1;
636 uinc2 = 255-uinc2;
637 uinc3 = 255-uinc3;
638 } else if (fColorAlg==kColorAlgCieModulo) {
639 a = uinc1;
640 b = uinc2;
641 c = uinc3;
642 d = a+b+c;
643 if (d!=0) {
644 a = a/d;
645 b = b/d;
646 c = c/d;
647 }
648 red = a*255;
649 green = b*255;
650 blue = c*255;
651 uinc1 = (Int_t)red;
652 uinc2 = (Int_t)green;
653 uinc3 = (Int_t)blue;
654 } else if (fColorAlg==kColorAlgYiqModulo) {
655 a = uinc1;
656 b = uinc2;
657 c = uinc3;
658 a = a/256;
659 b = b/256;
660 c = c/256;
661 red = a+0.956*b+0.62*c;
662 green = a-0.272*b-0.647*c;
663 blue = a-1.108*b+1.705*c;
664 if (red>=2) red = red-2;
665 else if (red>=1) red = red-1;
666 if (green<0) green = 1+green;
667 if (blue>=2) blue = blue-2;
668 else if (blue>=1) blue = blue-1;
669 else if (blue<-1) blue = 2+blue;
670 else if (blue<0) blue = 1+blue;
671 red = red*255;
672 green = green*255;
673 blue = blue*255;
674 uinc1 = (Int_t)red;
675 uinc2 = (Int_t)green;
676 uinc3 = (Int_t)blue;
677 } else if (fColorAlg==kColorAlgHvsModulo) {
678 h = uinc1;
679 v = uinc2;
680 s = uinc3;
681 h = h/256;
682 v = v/256;
683 s = s/256;
684 if (s==0) {
685 red = v;
686 green = v;
687 blue = v;
688 } else {
689 if (h==1.0) h = 0;
690 h = 6.0*h;
691 j = (Int_t)h;
692 f = h-j;
693 a = v*(1-s);
694 b = v*(1-s*f);
695 c = v*(1-s*(1-f));
696 switch (j) {
697 case 0:
698 red = v;
699 green = c;
700 blue = a;
701 break;
702 case 1:
703 red = b;
704 green = v;
705 blue = a;
706 break;
707 case 2:
708 red = a;
709 green = v;
710 blue = c;
711 break;
712 case 3:
713 red = a;
714 green = b;
715 blue = v;
716 break;
717 case 4:
718 red = c;
719 green = a;
720 blue = v;
721 break;
722 case 5:
723 red = v;
724 green = a;
725 blue = b;
726 break;
727 }
728 }
729 red = red*255;
730 green = green*255;
731 blue = blue*255;
732 uinc1 = (Int_t)red;
733 uinc2 = (Int_t)green;
734 uinc3 = (Int_t)blue;
735 }
736 ui = uinc1+uinc2*256+uinc3*65536;
737 }
738 red = uinc1;
739 green = uinc2;
740 blue = uinc3;
741 red = red/255.0;
742 green = green/255.0;
743 blue = blue/255.0;
744 fNewColor = new TColor(250+iv,red,green,blue);
745 fNewColorIndex = 250+iv;
746 return;
747}
748
749
750////////////////////////////////////////////////////////////////////////////////
751/// This function is called from BezierBlend function.
752
754{
755 Int_t j,a;
756 a = 1;
757 for (j=i+1;j<=3;j++) a = a*j;
758 for (j=1;j<=3-i;j++) a = a/j;
759 return a;
760}
761
762
763////////////////////////////////////////////////////////////////////////////////
764/// This function calculates Bezier approximation.
765
767{
768 Int_t j;
769 Double_t v;
770 v = BezC(i);
771 for (j=1;j<=i;j++) v = v*bezf;
772 for (j=1;j<=3-i;j++) v = v*(1-bezf);
773 return v;
774}
775
776
777////////////////////////////////////////////////////////////////////////////////
778/// Calculates screen coordinates of the smoothed point.
779/// Parameter bezf changes within the interval 0 to 1 in 0.1 steps.
780
782{
783 Int_t i;
784 Double_t b;
785 fGbezx = 0;
786 fGbezy = 0;
787 for (i=0;i<4;i++) {
788 b = BezierBlend(i,bezf);
789 fGbezx += fBzX[i]*b;
790 fGbezy += fBzY[i]*b;
791 }
792 return;
793}
794
795
796////////////////////////////////////////////////////////////////////////////////
797/// Ensures hidden surface removal.
798
800{
801 Int_t x,y,krok,xold=0,yold=0,prvy,yprv=0;
802 Double_t fx,fy,fx1,fy1;
803 if (y1<fBy1) y1 = fBy1;
804 if (y2<fBy1) y2 = fBy1;
805 if (x1==x2) {
806 if ((y1>=fEnvelope[x1]) && (y2>=fEnvelope[x1])) {
807 if (x1>0) {
808 if (y1<=fEnvelope[x1-1]||y2<=fEnvelope[x1-1]) {
809 if (y1>fEnvelope[x1-1]) y1 = fEnvelope[x1-1];
810 if (y2>fEnvelope[x1-1]) y2 = fEnvelope[x1-1];
811 fLine = 2;
812 fXs = x1;
813 fYs = y1;
814 fXe = x2;
815 fYe = y2;
816 return;
817 }
818 }
819 if (x1<fBx2) {
820 if (y1<=fEnvelope[x1+1]||y2<=fEnvelope[x1+1]) {
821 if (y1>fEnvelope[x1+1]) y1 = fEnvelope[x1+1];
822 if (y2>fEnvelope[x1+1]) y2 = fEnvelope[x1+1];
823 fLine = 2;
824 fXs = x1;
825 fYs = y1;
826 fXe = x2;
827 fYe = y2;
828 return;
829 }
830 }
831 fLine=0;
832 return;
833 }
834 if ((y1<fEnvelope[x1]) && (y2<fEnvelope[x1])) {
835 fLine = 2;
836 fXs = x1;
837 fYs = y1;
838 fXe = x2;
839 fYe = y2;
840 if (y1<y2) fEnvelope[x1] = y1;
841 else fEnvelope[x1] = y2;
842 return;
843 }
844 if (y1<y2) {
845 fLine = 2;
846 fXs = x1;
847 fYs = y1;
848 fXe = x1;
849 fYe = fEnvelope[x1];
850 fEnvelope[x1] = y1;
851 return;
852 } else {
853 fLine = 2;
854 fXs = x1;
855 fYs = y2;
856 fXe = x1;
857 fYe = fEnvelope[x1];
858 fEnvelope[x1] = y2;
859 return;
860 }
861 }
862 krok = (x1<x2)? 1:-1;
863 fLine = 0;
864 prvy = 0;
865 x = x1;
866 y = y1;
867l1:
868 if (y<=fEnvelope[x]) {
869 xold = x;
870 yold = y;
871 if (fLine==0) {
872 fLine = 1;
873 if (prvy==1) {
874 if (yprv<=fEnvelope[x]) fYs = yprv;
875 else fYs = fEnvelope[x];
876 fXs = x;
877 } else {
878 fXs = x;
879 fYs = y;
880 }
881 }
882 if (x!=x2) fEnvelope[x] = y;
883 } else {
884 prvy = 1;
885 yprv = y;
886 if (fLine==1) {
887 fLine = 2;
888 fXe = xold;
889 fYe = yold;
890 }
891 }
892 if (x1==x2) {
893 if (y1!=y2) y += (y1<y2)? +1:-1;
894 if (y!=y2) goto l1;
895 } else {
896 x += krok;
897 fy1 = y2-y1;
898 fx1 = x2-x1;
899 fx = x-x1;
900 fy = fy1*fx/fx1;
901 y = (Int_t)(y1+fy);
902 if (((x<=x2)&&(x1<x2)) || ((x>=x2)&&(x1>x2))) goto l1;
903 }
904 return;
905}
906
907
908////////////////////////////////////////////////////////////////////////////////
909/// Ensures hidden surface removal for Bars, BarsX and BarsY
910/// display modes.
911
913{
914 Int_t x,y,krok,xold=0,yold=0,prvy,xprv,yprv=0;
915 Double_t fx,fy,fx1,fy1;
916 if (x1==x2) {
917 if ((y1>=fEnvelope[x1]) && (y2>=fEnvelope[x1])) {
918 fLine = 0;
919 return;
920 }
921 if ((y1<fEnvelope[x1]) && (y2<fEnvelope[x1])) {
922 fLine = 2;
923 fXs = x1;
924 fYs = y1;
925 fXe = x2;
926 fYe = y2;
927 if (y1<y2) fEnvelope[x1] = y1;
928 else fEnvelope[x1] = y2;
929 return;
930 }
931 if (y1<y2) {
932 fLine = 2;
933 fXs = x1;
934 fYs = y1;
935 fXe = x1;
936 fYe = fEnvelope[x1];
937 fEnvelope[x1] = y1;
938 return;
939 } else {
940 fLine = 2;
941 fXs = x1;
942 fYs = y2;
943 fXe = x1;
944 fYe = fEnvelope[x1];
945 fEnvelope[x1] = y2;
946 return;
947 }
948 }
949 krok = (x1<x2)? 1:-1;
950 fLine = 0;
951 prvy = 0;
952 x = x1;
953 y = y1;
954l1:
955 if (y<=fEnvelope[x]) {
956 xold = x;
957 yold = y;
958 if (fLine==0) {
959 fLine = 1;
960 if (prvy==1) {
961 xprv = x;
962 fXs = xprv;
963 fYs = yprv;
964 } else {
965 fXs = x;
966 fYs = y;
967 }
968 }
969 if (x!=x2) fEnvelope[x] = y;
970 } else {
971 prvy = 1;
972 xprv = x;
973 yprv = y;
974 if (fLine==1) {
975 fLine = 2;
976 fXe = xold;
977 fYe = yold;
978 }
979 }
980 if (x1==x2) {
981 if (y1!=y2) y+=(y1<y2)? +1:-1;
982 if (y!=y2) goto l1;
983 } else {
984 x += krok;
985 fy1 = y2-y1;
986 fx1 = x2-x1;
987 fx = x-x1;
988 fy = fy1*fx/fx1;
989 y = (Int_t)(y1+fy);
990 if (((x<=x2)&&(x1<x2)) || ((x>=x2)&&(x1>x2))) goto l1;
991 }
992 return;
993}
994
995
996////////////////////////////////////////////////////////////////////////////////
997/// Draws channel mark at the screen coordinates x, y. Width of
998/// the mark is w, height is h and the type of the mark is determined by the
999/// parameter type.
1000
1002{
1003 TLine *line=new TLine();
1004 TEllipse *ellipse=new TEllipse();
1006 line->SetLineWidth(1);
1008 ellipse->SetLineColor(fChanmarkColor);
1009 ellipse->SetLineWidth(1);
1010 ellipse->SetLineStyle(kPenStyleSolid);
1011 switch (type) {
1013 ellipse->SetX1(gPad->PixeltoX(x));
1014 ellipse->SetY1(gPad->PixeltoY(y)+1);
1015 ellipse->SetR1(gPad->PixeltoX(w/2));
1016 ellipse->SetR2(gPad->PixeltoY(h/2));
1017 ellipse->SetPhimin(0);
1018 ellipse->SetPhimax(360);
1019 ellipse->SetTheta(0);
1020 ellipse->Paint("");
1021 break;
1023 line->PaintLine(gPad->PixeltoX(x-w/2),gPad->PixeltoY(y)+1,
1024 gPad->PixeltoX(x+w/2),gPad->PixeltoY(y)+1);
1025 line->PaintLine(gPad->PixeltoX(x) ,gPad->PixeltoY(y-h/2)+1,
1026 gPad->PixeltoX(x) ,gPad->PixeltoY(y+h/2+1)+1);
1027 break;
1029 line->PaintLine(gPad->PixeltoX(x-w/2) ,gPad->PixeltoY(y)+1,
1030 gPad->PixeltoX(x+w/2+1),gPad->PixeltoY(y)+1);
1031 line->PaintLine(gPad->PixeltoX(x) ,gPad->PixeltoY(y-h/2)+1,
1032 gPad->PixeltoX(x) ,gPad->PixeltoY(y+h/2+1)+1);
1033 line->PaintLine(gPad->PixeltoX(x-w/2) ,gPad->PixeltoY(y-h/2)+1,
1034 gPad->PixeltoX(x+w/2+1),gPad->PixeltoY(y+h/2+1)+1);
1035 line->PaintLine(gPad->PixeltoX(x-w/2) ,gPad->PixeltoY(y+h/2)+1,
1036 gPad->PixeltoX(x+w/2+1),gPad->PixeltoY(y-h/2-1)+1);
1037 break;
1039 line->PaintLine(gPad->PixeltoX(x-w/2),gPad->PixeltoY(y-h/2)+1,
1040 gPad->PixeltoX(x-w/2),gPad->PixeltoY(y+h/2)+1);
1041 line->PaintLine(gPad->PixeltoX(x-w/2),gPad->PixeltoY(y+h/2)+1,
1042 gPad->PixeltoX(x+w/2),gPad->PixeltoY(y+h/2)+1);
1043 line->PaintLine(gPad->PixeltoX(x+w/2),gPad->PixeltoY(y+h/2)+1,
1044 gPad->PixeltoX(x+w/2),gPad->PixeltoY(y-h/2)+1);
1045 line->PaintLine(gPad->PixeltoX(x+w/2),gPad->PixeltoY(y-h/2)+1,
1046 gPad->PixeltoX(x-w/2),gPad->PixeltoY(y-h/2)+1);
1047 break;
1049 line->PaintLine(gPad->PixeltoX(x-w/2) ,gPad->PixeltoY(y-h/2)+1,
1050 gPad->PixeltoX(x+w/2+1),gPad->PixeltoY(y+h/2+1)+1);
1051 line->PaintLine(gPad->PixeltoX(x-w/2) ,gPad->PixeltoY(y+h/2)+1,
1052 gPad->PixeltoX(x+w/2+1),gPad->PixeltoY(y-h/2-1)+1);
1053 break;
1055 line->PaintLine(gPad->PixeltoX(x) ,gPad->PixeltoY(y-h/2)+1,
1056 gPad->PixeltoX(x-w/2),gPad->PixeltoY(y)+1);
1057 line->PaintLine(gPad->PixeltoX(x-w/2),gPad->PixeltoY(y)+1,
1058 gPad->PixeltoX(x) ,gPad->PixeltoY(y+h/2)+1);
1059 line->PaintLine(gPad->PixeltoX(x) ,gPad->PixeltoY(y+h/2)+1,
1060 gPad->PixeltoX(x+w/2),gPad->PixeltoY(y)+1);
1061 line->PaintLine(gPad->PixeltoX(x+w/2),gPad->PixeltoY(y)+1,
1062 gPad->PixeltoX(x) ,gPad->PixeltoY(y-h/2)+1);
1063 break;
1065 line->PaintLine(gPad->PixeltoX(x) ,gPad->PixeltoY(y-h/2)+1,
1066 gPad->PixeltoX(x-w/2),gPad->PixeltoY(y+h/2)+1);
1067 line->PaintLine(gPad->PixeltoX(x-w/2),gPad->PixeltoY(y+h/2)+1,
1068 gPad->PixeltoX(x+w/2),gPad->PixeltoY(y+h/2)+1);
1069 line->PaintLine(gPad->PixeltoX(x+w/2),gPad->PixeltoY(y+h/2)+1,
1070 gPad->PixeltoX(x) ,gPad->PixeltoY(y-h/2)+1);
1071 break;
1072 }
1073 delete line;
1074 delete ellipse;
1075 return;
1076}
1077
1078
1079////////////////////////////////////////////////////////////////////////////////
1080/// Calculates screen coordinates of the line given by two
1081/// nodes for contours display mode. The line is given by two points
1082/// xr, yr, xs, ys. Finally it draws the line.
1083
1085 Double_t ys, TLine *line)
1086{
1087 Int_t krok,xi,yi,xj,yj,a,b,as,bs,pr,ae,be;
1088 Double_t fx,fy,fx1,fy1;
1089 xi = (Int_t)(fTxx*(xr-fXmin)/fKx+fTxy*(yr-fYmin)/fKy+fVx);
1090 xj = (Int_t)(fTxx*(xs-fXmin)/fKx+fTxy*(ys-fYmin)/fKy+fVx);
1091 yi = (Int_t)(fTyx*(xr-fXmin)/fKx+fTyy*(yr-fYmin)/fKy+fTyz*fZ+fVy);
1092 yj = (Int_t)(fTyx*(xs-fXmin)/fKx+fTyy*(ys-fYmin)/fKy+fTyz*fZ+fVy);
1093 as = xi;
1094 bs = yi;
1095 ae = xj;
1096 be = yj;
1097 a = xi;
1098 b = yi;
1099 pr = 0;
1100 krok = (xi<xj)? 1:-1;
1101l1:
1102 if (b<=fEnvelope[a]) {
1103 fEnvelopeContour[a] = b;
1104 if (pr==0) {
1105 pr = 1;
1106 as = a;
1107 bs = b;
1108 }
1109 } else {
1110 if (pr==1) {
1111 pr = 2;
1112 ae = a;
1113 be = b;
1114 }
1115 }
1116 if (xi==xj) {
1117 if (yi!=yj) b += (yi<yj)? +1:-1;
1118 if (b!=yj) goto l1;
1119 } else {
1120 a += krok;
1121 fy1 = yj-yi;
1122 fx1 = xj-xi;
1123 fx = a-xi;
1124 fy = fy1*fx/fx1;
1125 b = (Int_t)(yi+fy);
1126 if (a!=xj) goto l1;
1127 }
1128 if (pr!=0) {
1129 if (pr==1) {
1130 ae = xj;
1131 be = yj;
1132 }
1133 line->PaintLine(gPad->PixeltoX(as),gPad->PixeltoY(bs)+1,
1134 gPad->PixeltoX(ae),gPad->PixeltoY(be)+1);
1135 }
1136 return;
1137}
1138
1139
1140////////////////////////////////////////////////////////////////////////////////
1141/// Copies envelope vector, which ensures hidden surface removal for the
1142/// contours display mode.
1143
1145 Double_t ys)
1146{
1147 Int_t xi,xj,a;
1148 xi = (Int_t)(fTxx*(xr-fXmin)/fKx+fTxy*(yr-fYmin)/fKy+fVx);
1149 xj = (Int_t)(fTxx*(xs-fXmin)/fKx+fTxy*(ys-fYmin)/fKy+fVx);
1150 if (xi<xj) {
1151 for (a=xi;a<=xj;a++) {
1155 }
1156 } else if (xj<xi) {
1157 for (a=xj;a<=xi;a++) {
1161 }
1162 }
1163 return;
1164}
1165
1166
1167////////////////////////////////////////////////////////////////////////////////
1168/// Paints histogram according to preset parameters.
1169/// ### Visualization
1170/// #### Goal: to present 2-dimensional spectra in suitable visual form
1171/// This package has several display mode groups and display modes, which can be
1172/// employed for the presentation of 2-dimensional histograms
1173/// #### Display modes groups:
1174///
1175/// - `kModeGroupSimple` - it covers simple display modes using one
1176/// color only
1177/// - `kModeGroupLight` - in this group the shading is carried out
1178/// according to the position of the fictive
1179/// light source
1180/// - `kModeGroupHeight` - in this group the shading is carried out
1181/// according to the channel contents
1182/// - `kModeGroupLightHeight` - combination of two previous shading
1183/// algorithms. One can control the weight
1184/// between both algorithms.
1185///
1186/// #### Display modes:
1187///
1188/// - `kDisplayModePoints, `
1189/// - `kDisplayModeGrid, `
1190/// - `kDisplayModeContours,`
1191/// - `kDisplayModeBars,`
1192/// - `kDisplayModeLinesX,`
1193/// - `kDisplayModeLinesY,`
1194/// - `kDisplayModeBarsX,`
1195/// - `kDisplayModeBarsY,`
1196/// - `kDisplayModeNeedles,`
1197/// - `kDisplayModeSurface,`
1198/// - `kDisplayModeTriangles.`
1199///
1200/// one can combine the above given modes groups and display modes. The meaningful
1201/// combinations (denoted by x) are given in the next table.
1202///
1203/// | | Simple | Light | Height | Light-Height |
1204/// |-----------|--------|-------|--------|--------------|
1205/// | Points | X | X | X | X |
1206/// | Grid | X | X | X | X |
1207/// | Contours | X | - | X | - |
1208/// | Bars | X | - | X | - |
1209/// | LinesX | X | X | X | X |
1210/// | LinesY | X | X | X | X |
1211/// | BarsX | X | - | X | - |
1212/// | BarsY | X | - | X | - |
1213/// | Needles | X | - | - | - |
1214/// | Surface | - | X | X | X |
1215/// | Triangles | X | X | X | X |
1216///
1217/// #### Function: void TSpectrum2Painter::SetDisplayMode (Int_t modeGroup, Int_t displayMode)
1218///
1219/// This function controls the display mode group and display mode of the
1220/// histogram drawing. To illustrate the possible effects of the various display
1221/// modes we introduce a set of examples. Default values:
1222///
1223/// - `modeGroup = kModeGroupLightHeight `
1224/// - `displayMode = kDisplayModeSurface `
1225///
1226/// \image html spectrumpainter001.jpg
1227///
1228/// Simple modes group, display mode = points, 256 x 256 channels.
1229/// \image html spectrumpainter002.jpg
1230///
1231/// Simple modes group, display mode = grid, 64 x 64 channels.
1232/// \image html spectrumpainter003.jpg
1233///
1234/// Simple modes group, display mode = contours, 64 x 64 channels.
1235/// \image html spectrumpainter004.jpg
1236///
1237/// Simple modes group, display mode = bars, 64 x 64 channels.
1238/// \image html spectrumpainter005.jpg
1239///
1240/// Simple modes group, display mode = linesX, 64 x 64 channels.
1241/// \image html spectrumpainter006.jpg
1242///
1243/// Simple modes group, display mode = linesY, 64 x 64 channels.
1244/// \image html spectrumpainter007.jpg
1245///
1246/// Simple modes group, display mode = barsX, 64 x 64 channels.
1247/// \image html spectrumpainter008.jpg
1248///
1249/// Simple modes group, display mode = barsY, 64 x 64 channels.
1250/// \image html spectrumpainter009.jpg
1251///
1252/// Simple modes group, display mode = needles, 64 x 64 channels.
1253/// \image html spectrumpainter010.jpg
1254///
1255/// Simple modes group, display mode = triangles, 64 x 64 channels.
1256/// \image html spectrumpainter011.jpg
1257///
1258/// Light modes group, display mode = points, 256 x 256 channels.
1259/// \image html spectrumpainter012.jpg
1260///
1261/// Light modes group, display mode = grid, 256 x 256 channels.
1262/// \image html spectrumpainter013.jpg
1263///
1264/// Light modes group, display mode = surface, 64 x 64 channels.
1265/// \image html spectrumpainter014.jpg
1266///
1267/// Light modes group, display mode = triangles, 64 x 64 channels.
1268/// \image html spectrumpainter015.jpg
1269///
1270/// Height modes group, display mode = points, 256 x 256 channels.
1271/// \image html spectrumpainter016.jpg
1272///
1273/// Height modes group, display mode = grid, 256 x 256 channels.
1274/// \image html spectrumpainter017.jpg
1275///
1276/// Height modes group, display mode = contours, 64 x 64 channels.
1277/// \image html spectrumpainter018.jpg
1278///
1279/// Height modes group, display mode = bars, 64 x 64 channels.
1280/// \image html spectrumpainter019.jpg
1281///
1282/// Height modes group, display mode = surface, 64 x 64 channels.
1283/// \image html spectrumpainter020.jpg
1284///
1285/// Height modes group, display mode = triangles, 64 x 64 channels.
1286/// \image html spectrumpainter021.jpg
1287///
1288/// Light - height modes group, display mode = surface, 64 x 64 channels. The weight
1289/// between both shading algorithms is set to 0.5. One can observe the influence of
1290/// both shadings.
1291///
1292/// #### Function: TSpectrum2Painter::SetPenAttr(Int_t color,Int_t style,Int_t width)
1293///
1294/// Using this function one can change pen color, pen style and pen width.
1295/// Possible pen styles are:
1296///
1297/// - ` kPenStyleSolid,`
1298/// - ` kPenStyleDash,`
1299/// - ` kPenStyleDot,`
1300/// - ` kPenStyleDashDot.`
1301///
1302/// Default values:
1303///
1304/// - ` color = kBlack`
1305/// - ` style = kPenStyleSolid`
1306/// - ` width = 1`
1307///
1308/// \image html spectrumpainter022.jpg
1309///
1310/// Simple modes group, display mode = linesX, 64 x 64 channels. Pen width = 3.
1311///
1312/// #### Function: TSpectrum2Painter::SetNodes(Int_t nodesx,Int_t nodesy)
1313///
1314/// Sometimes the displayed region is rather large. When displaying all channels
1315/// pictures become very dense and complicated. It is very difficult to understand
1316/// overall shape of the data. Therefore in the package we have implemented the
1317/// possibility to change the density of displayed channels. Only channels
1318/// coinciding with given nodes are displayed. In the next figure we introduce the
1319/// example of the above presented spectrum with number of nodes set to 64x64.
1320///
1321/// Default values:
1322///
1323/// - ` nodesx = Xmax-Xmin+1`
1324/// - ` nodesy = Ymax-Ymin+1`
1325///
1326/// \image html spectrumpainter023.jpg
1327///
1328/// Simple modes group, display mode = grid, 256 x 256 channels.
1329/// Number of nodes is 64x64.
1330///
1331/// #### Function: void TSpectrum2Painter::SetAngles (Int_t alpha,Int_t beta, Int_t view)
1332///
1333/// One can change the angles of the position of 3-d space and to rotate the
1334/// space. Alpha parameter defines the angle between bottom horizontal screen line
1335/// and the displayed space on the right side of the picture and beta on the left
1336/// side, respectively. One can rotate the 3-d space around vertical axis going
1337/// through the center of it employing the view parameter.
1338/// Allowed values are 0, 90, 180 and 270 degrees respectively.
1339///
1340/// Default values:
1341///
1342/// - ` alpha = 20`
1343/// - ` beta = 60`
1344/// - ` view = 0`
1345///
1346/// \image html spectrumpainter024.jpg
1347///
1348/// Light modes group, display mode = surface, 256 x 256 channels. Angles are
1349/// set as follows: alpha=40, beta=30, view=0.
1350/// \image html spectrumpainter025.jpg
1351///
1352/// Light modes group, display mode = surface, 256 x 256 channels. Angles are
1353/// set as follows: alpha=30, beta=30, view=90.
1354///
1355/// #### Function: TSpectrum2Painter::SetZScale(Int_t scale)
1356///
1357/// One can change the scale of z-axis. Possible values are:
1358///
1359/// - ` kZScaleLinear`
1360/// - ` kZScaleLog`
1361/// - ` kZScaleSqrt`
1362///
1363/// Default value is:
1364///
1365/// - ` scale = kZScaleLinear`
1366///
1367/// \image html spectrumpainter026.jpg
1368///
1369/// Height modes group, display mode = surface, 64 x 64 channels, log scale.
1370///
1371/// #### Function: TSpectrum2Painter::SetColorIncrements(Double_t r,Double_t g,Double_t b);
1372///
1373/// For sophisticated shading (in kModeGroupLight, kModeGroupHeight
1374/// and kModeGroupLightHeight display modes groups) the color palette starts
1375/// from the basic pen color (see SetPenAttr function). There is a predefined number
1376/// of color levels (256). Color in every level is calculated by adding the
1377/// increments of the r, g, b components to the previous level. Using this function
1378/// one can change the color increments between two neighbouring color levels. The
1379/// function does not apply for kModeGroupSimple display modes group.
1380/// Default values: r=1, g=1, b=1;
1381/// \image html spectrumpainter027.jpg
1382///
1383/// Light modes group, display mode = surface, 64 x 64 channels, color increments
1384/// r=1, g=2, b=3.
1385/// \image html spectrumpainter028.jpg
1386///
1387/// Light modes group, display mode = surface, 64 x 64 channels, color
1388/// increments r=4, g=2, b=1.
1389///
1390/// #### Function: TSpectrum2Painter::SetColorAlgorithm(Int_t colorAlgorithm)
1391///
1392/// To define the colors one can employ one of the following color algorithms
1393/// (rgb, cmy, cie, yiq, hvs models [1], [2]). When the level of a component
1394/// achieves the limit value one can choose either smooth transition (by decreasing
1395/// the limit value) or sharp - modulo transition (continuing with 0 value). This
1396/// makes possible to realize various visual effects. One can choose from the
1397/// following set of the algorithms:
1398///
1399/// - ` kColorAlgRgbSmooth `
1400/// - ` kColorAlgRgbModulo `
1401/// - ` kColorAlgCmySmooth `
1402/// - ` kColorAlgCmyModulo `
1403/// - ` kColorAlgCieSmooth `
1404/// - ` kColorAlgCieModulo `
1405/// - ` kColorAlgYiqSmooth `
1406/// - ` kColorAlgYiqModulo `
1407/// - ` kColorAlgHvsSmooth `
1408/// - ` kColorAlgHvsModulo `
1409///
1410/// The function does not apply for kModeGroupSimple display modes group.
1411/// Default value is:
1412///
1413/// - ` colorAlgorithm = kColorAlgRgbSmooth`
1414///
1415/// \image html spectrumpainter029.jpg
1416///
1417/// Light modes group, display mode = surface, 64 x 64 channels, color algorithm
1418/// is cmy smooth.
1419/// \image html spectrumpainter030.jpg
1420///
1421/// Light modes group, display mode = surface, 64 x 64 channels, color algorithm
1422/// is hvs smooth.
1423/// \image html spectrumpainter031.jpg
1424///
1425/// Light modes group, display mode = surface, 64 x 64 channels, color algorithm
1426/// is yiq smooth.
1427/// \image html spectrumpainter032.jpg
1428///
1429/// Light modes group, display mode = surface, 64 x 64 channels, color algorithm
1430/// is rgb modulo.
1431/// \image html spectrumpainter033.jpg
1432///
1433/// Height modes group, display mode = surface, 256 x 256 channels, color
1434/// algorithm is rgb modulo, increments r=5, g=5, b=5, angles alpha=0, beta=90,
1435/// view=0.
1436///
1437/// #### Function: TSpectrum2Painter::SetLightPosition(Int_t x, Int_t y, Int_t z)
1438///
1439/// In kModeGroupLight and kModeGroupLightHeight display modes
1440/// groups the color palette is calculated according to the fictive light source
1441/// position in 3-d space. Using this function one can change the position of the
1442/// source and thus to achieve various graphical effects. The function does not
1443/// apply for kModeGroupSimple and kModeGroupHeight display modes
1444/// groups. Default values are: x=1000, y=1000, z=1000.
1445/// \image html spectrumpainter034.jpg
1446///
1447/// Light modes group, display mode = surface, 64 x 64 channels. Position of the
1448/// light source was set to x=0, y=1000, z=1000.
1449///
1450/// #### Function: TSpectrum2Painter::SetShading(Int_t shading,Int_t shadow)
1451///
1452/// Surface of the picture is composed of triangles. If desired the edges of the
1453/// neighbouring triangles can be smoothed (shaded). If desired the display of the
1454/// shadow can be painted as well. The function does not apply for
1455/// kModeGroupSimple display modes group.
1456///
1457/// Possible values for shading are:
1458///
1459/// - ` kNotShaded`
1460/// - ` kShaded.`
1461///
1462/// Possible values for shadow are:
1463///
1464/// - ` kShadowsNotPainted`
1465/// - ` kShadowsPainted`
1466///
1467/// Default values:
1468///
1469/// - ` shading = kShaded`
1470/// - ` shadow = kShadowsNotPainted`
1471///
1472/// \image html spectrumpainter035.jpg
1473///
1474/// Light modes group, display mode = surface, 64 x 64 channels, not shaded.
1475/// \image html spectrumpainter036.jpg
1476///
1477/// Light modes group, display mode = surface, 64 x 64 channels, shaded, with
1478/// shadow.
1479///
1480/// #### Function: TSpectrum2Painter::SetBezier(Int_t bezier)
1481///
1482/// For kModeGroupSimple display modes group and for kDisplayModeGrid,
1483/// kDisplayModeLinesX >and kDisplayModeLinesY display modes one
1484/// can smooth data using Bezier smoothing algorithm. The function does not apply
1485/// for other display modes groups and display modes. Possible values are:
1486///
1487/// - ` kNoBezierInterpol`
1488/// - ` kBezierInterpol`
1489///
1490/// Default value is:
1491///
1492/// - ` bezier = kNoBezierInterpol.`
1493///
1494/// \image html spectrumpainter005.jpg
1495///
1496/// Simple modes group, display mode = linesX, 64 x 64 channels with Bezier
1497/// smoothing.
1498///
1499/// #### Function: TSpectrum2Painter::SetContourWidth(Int_t width)
1500///
1501/// This function applies only for kDisplayModeContours display mode.
1502/// One can change the width between horizontal slices and thus their density.
1503/// Default value: width=50.
1504/// \image html spectrumpainter037.jpg
1505///
1506/// Simple modes group, display mode = contours, 64 x 64 channels. Width between
1507/// slices was set to 30.
1508///
1509/// #### Function: TSpectrum2Painter::SetLightHeightWeight(Double_t weight)
1510///
1511/// For kModeGroupLightHeight display modes group one can change the
1512/// weight between both shading algorithm. The function does not apply for other
1513/// display modes groups. Default value is: weight=0.5.
1514/// \image html spectrumpainter038.jpg
1515///
1516/// Light - height modes group, display mode = surface, 64 x 64 channels.
1517/// The weight between both shading algorithms is set to 0.7.
1518///
1519/// #### Function: TSpectrum2Painter::SetChanMarks(Int_t enable,Int_t color,Int_t width,Int_t height,Int_t style)
1520/// In addition to the surface drawn using any above given algorithm one can display
1521/// channel marks. One can control the color as well as the width, height
1522/// (in pixels) and the style of the marks. The parameter enable can be set to:
1523///
1524/// - `kChannelMarksNotDrawn`
1525/// - `kChannelMarksDrawn.`
1526///
1527/// The possible styles can be chosen from the set:
1528///
1529/// - ` kChannelMarksStyleDot`
1530/// - ` kChannelMarksStyleCross`
1531/// - ` kChannelMarksStyleStar`
1532/// - ` kChannelMarksStyleRectangle`
1533/// - ` kChannelMarksStyleX`
1534/// - ` kChannelMarksStyleDiamond`
1535/// - ` kChannelMarksStyleTriangle.`
1536///
1537/// \image html spectrumpainter039.jpg
1538///
1539/// Light modes group, display mode = surface, 64 x 64 channels,
1540/// with marks (red circles).</p>
1541///
1542/// #### Function: TSpectrum2Painter::SetChanGrid(Int_t enable,Int_t color)
1543///
1544/// In addition to the surface drawn using any above given algorithm one can
1545/// display grid using the color parameter. The parameter enable can be set to:
1546///
1547/// - ` kChannelGridNotDrawn`
1548/// - ` kChannelGridDrawn.`
1549///
1550/// \image html spectrumpainter040.jpg
1551///
1552/// Height modes group, display mode = surface, 64 x 64 channels, with blue grid.
1553/// \image html spectrumpainter041.jpg
1554///
1555/// Height modes group, display mode = surface, 64 x 64 channels, with marks
1556/// (red circles) and blue grid.
1557/// #### References:
1558///
1559/// [1] Morhac M., Kliman J., Matouoek V., Turzo I.,
1560/// Sophisticated visualization algorithms for analysis of multidimensional
1561/// experimental nuclear data, Acta Physica Slovaca 54 (2004) 385.
1562///
1563/// [2] D. Hearn, M. P. Baker: Computer Graphics, Prentice Hall International,
1564/// Inc. 1994.
1565/// #### Script:
1566///
1567/// Example to draw source spectrum (class TSpectrum2Painter).
1568/// To execute this example, do:
1569/// ~~~
1570/// root > .x VisA.C
1571/// ~~~
1572/// ~~~ {.cpp}
1573/// #include "TSpectrum2Painter.h"
1574///
1575/// void VisA() {
1576/// TFile *f = new TFile("TSpectrum2.root");
1577/// TH2F *graph=(TH2F*) f->Get("graph2;1");
1578/// TCanvas *Graph2 = new TCanvas("Graph2","Illustration of 2D graphics",10,10,1000,700);
1579/// graph->Draw("SPEC");
1580/// }
1581/// ~~~
1582
1584{
1585
1586
1587 Int_t turni,turnj,w1,w2,x,y;
1588 Int_t q1=0,q2=0,qv=0,smer=0,flag=0,i=0,j=0,x1=0,y1=0,x2=0,y2=0,x3=0,y3=0,x4=0,y4=0,uhl=0,xp1=0,yp1=0,xp2=0,yp2=0;
1589 Int_t ix5,iy5,x6,y6,x7,y7,y8,x1d,y1d,x2d=0,y2d=0;
1590 Int_t i1=0,i2=0,i3=0,i4=0,j1=0,j2=0,j3=0,j4=0;
1591 Int_t s1=0,s2=0,s3=0,s4=0,t1=0,t2=0,t3=0,t4=0;
1592 Double_t dx1,dx2,dx3,dx4,dy1,dy2,dy3,dy4,z1,z2,z3,z4,zl,zh;
1593 Double_t xa,xb=0,ya,yb=0,x5=0,y5=0;
1594 Double_t da=0,db=0,dc=0,dd=0,xtaz,ytaz,ztaz,v,shad_noise;
1595 Int_t iv=0,ekv,stvor,sx1,sx2,sx3,sx4,sx5,sy1,sy2,sy3,sy4,sy5;
1596 Double_t pom1,pom2,sdx1,sdy1,sdx2=0,sdy2,sdx3,sdy3,sdy4,spriz;
1597 Int_t sr1=0,sr2=0,sr3=0,sr4=0,sr5=0,sr6=0,sr7=0,sr8=0;
1598 Int_t tr1=0,tr2=0,tr3=0,tr4=0,tr5=0,tr6=0,tr7=0,tr8=0;
1599 Int_t il,iv1=0,iv2=0,iv3=0,iv4=0;
1600 Double_t v1=0,v2=0,v3=0,v4=0,dxr1,dxr2,dyr1,dyr2,zr1,zr2,bezf;
1601 Double_t dcount_reg,z1l,z2l,z3l,z4l,sdx2p,sdy2p,dap,dbp,dcp,ddp;
1602 Int_t sx1p,sy1p,sx3p,uip=0;
1603 Double_t bezx1,bezy1,bezx2,bezy2;
1604 Double_t p000x,p000y,p100x,p100y,p010x,p010y,p110x,p110y;
1605 Double_t p001x,p001y,p101x,p101y,p011x,p011y,p111x,p111y;
1606 Int_t ibezx1=0,ibezy1=0,ibezx2,ibezy2;
1607 unsigned ui1,ui2,ui3;
1608 Double_t fi,alfa,beta,x3max,y3max,mul,movx,movy;
1609 Double_t xmin,xmax,ymin,ymax,zmin,zmax,mx,my,mz;
1610 Double_t mxx,mxy,myx,myy,myz,px,py,kx,ky;
1611 Double_t bxl,bxh,byl,byh,xd,yd,a,b,rotx,roty;
1612 TLine *line = new TLine();
1613 TBox *box = new TBox();
1614 TColor *pen_col;
1615 pen_col = (TColor*)(gROOT->GetListOfColors()->At(fPenColor));
1616 ui1 = (Int_t)(256*pen_col->GetRed());
1617 ui2 = (Int_t)(256*pen_col->GetGreen());
1618 ui3 = (Int_t)(256*pen_col->GetBlue());
1619
1621 printf("The canvas size exceed the maximum X screen resolution.\n");
1622 printf("Use the option bf() to increase the buffer size (it should be greater than %d).\n",fBx2);
1623 return;
1624 }
1625
1626 for (i=fBx1;i<fBx2;i++) {
1627 fEnvelope[i] = fBy2;
1629 }
1630
1631// gPad->Range(0, 0, 1 ,1);
1632
1633 // Set the histogram's parameters.
1634 fBx1 = gPad->XtoPixel(0.1);
1635 fBx2 = gPad->XtoPixel(0.99);
1636 fBy1 = gPad->YtoPixel(0.99);
1637 fBy2 = gPad->YtoPixel(0.05);
1638 fXmin = fH2->GetXaxis()->GetFirst();
1639 fXmax = fH2->GetXaxis()->GetLast();
1640 fYmin = fH2->GetYaxis()->GetFirst();
1641 fYmax = fH2->GetYaxis()->GetLast();
1642 fZmax = fH2->GetMaximum();
1643 fZmin = fH2->GetMinimum();
1644
1645 // Calculation of display parameters.
1646 xmin = fXmin;
1647 xmax = fXmax;
1648 ymin = fYmin;
1649 ymax = fYmax;
1650 zmin = fZmin;
1651 zmax = fZmax;
1652 xd = (xmax-xmin)/2;
1653 yd = (ymax-ymin)/2;
1654 a = (xmax+xmin)/2;
1655 b = (ymax+ymin)/2;
1656 fi = (fViewAngle*3.1415927)/180;
1657 alfa = (fAlpha*3.1415927)/180;
1658 beta = (fBeta*3.1415927)/180;
1659 rotx = (-1)*a*cos(fi)+b*sin(fi)+xd*TMath::Abs(cos(fi))+yd*TMath::Abs(sin(fi));
1660 roty = (-1)*a*sin(fi)-b*cos(fi)+xd*TMath::Abs(sin(fi))+yd*TMath::Abs(cos(fi));
1661 x3max = (xmax-xmin)*TMath::Abs(cos(fi))+(ymax-ymin)*TMath::Abs(sin(fi));
1662 y3max = (xmax-xmin)*TMath::Abs(sin(fi))+(ymax-ymin)*TMath::Abs(cos(fi));
1663 bxl = fBx1;
1664 bxh = fBx2;
1665 byl = fBy1;
1666 byh = fBy2;
1667 mx = (bxh-bxl)/(x3max*(cos(alfa)+cos(beta)));
1668 my = (bxh-bxl)/(y3max*(cos(alfa)+cos(beta)));
1669 mul = (byh-byl)/(bxh-bxl);
1670 movx = bxl+my*cos(alfa)*y3max;
1671 mxx = mx*cos(beta)*cos(fi)-my*cos(alfa)*sin(fi);
1672 mxy = (-1)*mx*cos(beta)*sin(fi)-my*cos(alfa)*cos(fi);
1673 myx = mul*(mx*sin(beta)*cos(fi)+my*sin(alfa)*sin(fi));
1674 myy = mul*((-1)*mx*sin(beta)*sin(fi)+my*sin(alfa)*cos(fi));
1675 px = rotx*mx*cos(beta)-roty*my*cos(alfa)+movx;
1676 kx = (xmax-xmin)/(fNodesx-1);
1677 ky = (ymax-ymin)/(fNodesy-1);
1678 fKx = kx;
1679 fKy = ky;
1680 fMxx = mxx;
1681 fMxy = mxy;
1682 fMyx = myx;
1683 fMyy = myy;
1684 fTxx = mxx*kx;
1685 fTxy = mxy*ky;
1686 fTyx = myx*kx;
1687 fTyy = myy*ky;
1688 fVx = mxx*xmin+mxy*ymin+px;
1689 if (fZscale==kZScaleLinear) {
1690 mz = (bxh-bxl)*(cos(alfa)+cos(beta)-sin(alfa)-sin(beta));
1691 mz = mz/((zmax-zmin)*(cos(alfa)+cos(beta)));
1692 movy = byl+mul*mz*zmax;
1693 myz = (-1)*mz*mul;
1694 py = mul*(rotx*mx*sin(beta)+roty*my*sin(alfa))+movy;
1695 fTyz = myz;
1696 fVy = myx*xmin+myy*ymin+py;
1697 fNuSli = (zmax-zmin)/(Double_t)fContWidth;
1698 } else if (fZscale==kZScaleLog) {
1699 if (zmin>=1) zmin = log(zmin);
1700 else zmin = 0;
1701 if (zmax>=1) zmax = log(zmax);
1702 else zmax = 0;
1703 if ((zmax-zmin)<0.000001) zmax = zmin+0.000001;
1704 mz = (bxh-bxl)*(cos(alfa)+cos(beta)-sin(alfa)-sin(beta));
1705 mz = mz/((zmax-zmin)*(cos(alfa)+cos(beta)));
1706 movy = byl+mul*mz*zmax;
1707 myz = (-1)*mz*mul;
1708 py = mul*(rotx*mx*sin(beta)+roty*my*sin(alfa))+movy;
1709 fTyz = myz;
1710 fVy = myx*xmin+myy*ymin+py;
1711 fNuSli = (zmax-zmin)/(Double_t)fContWidth;
1712 } else if (fZscale==kZScaleSqrt) {
1713 if (zmin>=1) zmin = sqrt(zmin);
1714 else zmin = 0;
1715 if (zmax>=1) zmax = sqrt(zmax);
1716 else zmax = 0;
1717 if ((zmax-zmin)<0.000001) zmax = zmin+0.000001;
1718 mz = (bxh-bxl)*(cos(alfa)+cos(beta)-sin(alfa)-sin(beta));
1719 mz = mz/((zmax-zmin)*(cos(alfa)+cos(beta)));
1720 movy = byl+mul*mz*zmax;
1721 myz = (-1)*mz*mul;
1722 py = mul*(rotx*mx*sin(beta)+roty*my*sin(alfa))+movy;
1723 fTyz = myz;
1724 fVy = myx*xmin+myy*ymin+py;
1725 fNuSli = (zmax-zmin)/(Double_t)fContWidth;
1726 }
1727
1728 // End of calculations of display parameters.
1729 dcount_reg=fContWidth;
1730 switch (fZscale) {
1731 case kZScaleLog:
1732 dcount_reg=log(dcount_reg);
1733 break;
1734 case kZScaleSqrt:
1735 dcount_reg=sqrt(dcount_reg);
1736 break;
1737 }
1738 shad_noise = fZmax;
1739 shad_noise /= 100.;
1740 w1 = fNodesx-1;
1741 w2 = fNodesy-1;
1742
1743 // Drawing axis in backplanes.
1744 Transform(0,0,-1);
1745 p000x = gPad->PixeltoX(fXt);
1746 p000y = gPad->PixeltoY(fYt)+1;
1747 Transform(w1,0,-1);
1748 p100x = gPad->PixeltoX(fXt);
1749 p100y = gPad->PixeltoY(fYt)+1;
1750 Transform(0,w2,-1);
1751 p010x = gPad->PixeltoX(fXt);
1752 p010y = gPad->PixeltoY(fYt)+1;
1753 Transform(w1,w2,-1);
1754 p110x = gPad->PixeltoX(fXt);
1755 p110y = gPad->PixeltoY(fYt)+1;
1757 Transform(0,0,-2);
1758 p001x = gPad->PixeltoX(fXt);
1759 p001y = gPad->PixeltoY(fYt)+1;
1760 Transform(w1,0,-2);
1761 p101x = gPad->PixeltoX(fXt);
1762 p101y = gPad->PixeltoY(fYt)+1;
1763 Transform(0,w2,-2);
1764 p011x = gPad->PixeltoX(fXt);
1765 p011y = gPad->PixeltoY(fYt)+1;
1766 Transform(w1,w2,-2);
1767 p111x = gPad->PixeltoX(fXt);
1768 p111y = gPad->PixeltoY(fYt)+1;
1769 Double_t bmin, bmax, binLow, binHigh, binWidth;
1770 Double_t axisLevel, gridDist, gridY1, gridY2;
1771 Int_t ndivx = 0, ndivy, ndivz, nbins;
1772 TGaxis *axis = new TGaxis();
1773 TGaxis *xaxis = new TGaxis();
1774 TGaxis *yaxis = new TGaxis();
1775 TGaxis *zaxis = new TGaxis();
1777 if (fViewAngle==0) {
1778 axis->PaintAxis(p000x, p000y, p100x, p100y, bmin, bmax, ndivx, "");
1779 axis->PaintAxis(p000x, p000y, p010x, p010y, bmin, bmax, ndivx, "");
1780 if(fAlpha+fBeta<90)
1781 axis->PaintAxis(p000x, p000y, p001x, p001y, bmin, bmax, ndivx, "");
1782 if(fAlpha+fBeta<90)
1783 axis->PaintAxis(p100x, p100y, p101x, p101y, bmin, bmax, ndivx, "");
1784 axis->PaintAxis(p101x, p101y, p001x, p001y, bmin, bmax, ndivx, "");
1785 axis->PaintAxis(p001x, p001y, p011x, p011y, bmin, bmax, ndivx, "");
1786 if (fZscale==kZScaleLinear) {
1787 bmin = fZmin;
1788 bmax = fZmax;
1789 ndivz = 10;
1790 THLimitsFinder::Optimize(bmin, bmax, ndivz, binLow, binHigh,
1791 nbins, binWidth, " ");
1792 for (i = 0; i < nbins + 1; i++) {
1793 axisLevel = binLow+i*binWidth;
1794 gridDist = (axisLevel-bmin)*(p001y-p000y)/(bmax-bmin);
1795 gridY1 = p000y + gridDist, gridY2 = p100y + gridDist;
1796 line->PaintLine(p000x,gridY1,p100x,gridY2);
1797 gridY2 = p010y + gridDist;
1798 line->PaintLine(p000x,gridY1,p010x,gridY2);
1799 }
1800 }
1801 } else if (fViewAngle==90) {
1802 axis->PaintAxis(p010x, p010y, p000x, p000y, bmin, bmax, ndivx, "");
1803 axis->PaintAxis(p010x, p010y, p110x, p110y, bmin, bmax, ndivx, "");
1804 if(fAlpha+fBeta<90)
1805 axis->PaintAxis(p010x, p010y, p011x, p011y, bmin, bmax, ndivx, "");
1806 if(fAlpha+fBeta<90)
1807 axis->PaintAxis(p000x, p000y, p001x, p001y, bmin, bmax, ndivx, "");
1808 axis->PaintAxis(p001x, p001y, p011x, p011y, bmin, bmax, ndivx, "");
1809 axis->PaintAxis(p011x, p011y, p111x, p111y, bmin, bmax, ndivx, "");
1810 if (fZscale==kZScaleLinear) {
1811 bmin = fZmin;
1812 bmax = fZmax;
1813 ndivz = 10;
1814 THLimitsFinder::Optimize(bmin, bmax, ndivz, binLow, binHigh,
1815 nbins, binWidth, " ");
1816 for (i = 0; i < nbins + 1; i++) {
1817 axisLevel = binLow+i*binWidth;
1818 gridDist = (axisLevel-bmin)*(p011y-p010y)/(bmax-bmin);
1819 gridY1 = p010y + gridDist, gridY2 = p000y + gridDist;
1820 line->PaintLine(p010x,gridY1,p000x,gridY2);
1821 gridY2 = p110y + gridDist;
1822 line->PaintLine(p010x,gridY1,p110x,gridY2);
1823 }
1824 }
1825 } else if (fViewAngle==180) {
1826 axis->PaintAxis(p110x, p110y, p010x, p010y, bmin, bmax, ndivx, "");
1827 axis->PaintAxis(p110x, p110y, p100x, p100y, bmin, bmax, ndivx, "");
1828 if(fAlpha+fBeta<90)
1829 axis->PaintAxis(p110x, p110y, p111x, p111y, bmin, bmax, ndivx, "");
1830 if(fAlpha+fBeta<90)
1831 axis->PaintAxis(p010x, p010y, p011x, p011y, bmin, bmax, ndivx, "");
1832 axis->PaintAxis(p011x, p011y, p111x, p111y, bmin, bmax, ndivx, "");
1833 axis->PaintAxis(p111x, p111y, p101x, p101y, bmin, bmax, ndivx, "");
1834 if (fZscale==kZScaleLinear) {
1835 bmin = fZmin;
1836 bmax = fZmax;
1837 ndivz = 10;
1838 THLimitsFinder::Optimize(bmin, bmax, ndivz, binLow, binHigh,
1839 nbins, binWidth, " ");
1840 for (i = 0; i < nbins + 1; i++) {
1841 axisLevel = binLow+i*binWidth;
1842 gridDist = (axisLevel-bmin)*(p111y-p110y)/(bmax-bmin);
1843 gridY1 = p110y + gridDist, gridY2 = p010y + gridDist;
1844 line->PaintLine(p110x,gridY1,p010x,gridY2);
1845 gridY2 = p100y + gridDist;
1846 line->PaintLine(p110x,gridY1,p100x,gridY2);
1847 }
1848 }
1849 } else if (fViewAngle==270) {
1850 axis->PaintAxis(p100x, p100y, p110x, p110y, bmin, bmax, ndivx, "");
1851 axis->PaintAxis(p100x, p100y, p000x, p000y, bmin, bmax, ndivx, "");
1852 if(fAlpha+fBeta<90)
1853 axis->PaintAxis(p100x, p100y, p101x, p101y, bmin, bmax, ndivx, "");
1854 if(fAlpha+fBeta<90)
1855 axis->PaintAxis(p110x, p110y, p111x, p111y, bmin, bmax, ndivx, "");
1856 axis->PaintAxis(p111x, p111y, p101x, p101y, bmin, bmax, ndivx, "");
1857 axis->PaintAxis(p101x, p101y, p001x, p001y, bmin, bmax, ndivx, "");
1858 if (fZscale==kZScaleLinear) {
1859 bmin = fZmin;
1860 bmax = fZmax;
1861 ndivz = 10;
1862 THLimitsFinder::Optimize(bmin, bmax, ndivz, binLow, binHigh,
1863 nbins, binWidth, " ");
1864 for (i = 0; i < nbins + 1; i++) {
1865 axisLevel = binLow+i*binWidth;
1866 gridDist = (axisLevel-bmin)*(p101y-p100y)/(bmax-bmin);
1867 gridY1 = p100y + gridDist, gridY2 = p110y + gridDist;
1868 line->PaintLine(p100x,gridY1,p110x,gridY2);
1869 gridY2 = p000y + gridDist;
1870 line->PaintLine(p100x,gridY1,p000x,gridY2);
1871 }
1872 }
1873 }
1874
1875 // End.
1876 line->ResetAttLine("");
1880 turni = 0;
1881 turnj = 0;
1882 Transform(w1,0,0);
1883 x1 = fXt;
1884 Transform(0,0,0);
1885 x2 = fXt;
1886 Transform(0,w2,0);
1887 x3 = fXt;
1888 if (x2>=x1) turnj = 1;
1889 if (x3>=x2) turni = 1;
1890 q1 = 1;
1891 q2 = 0;
1892 qv = 1;
1893 do {
1894 uhl = 0;
1895 smer = 0;
1896 flag = 0;
1897l2:
1898 if (turni==1) {
1899 i = q1;
1900 } else {
1901 i = w1-q1;
1902 }
1903 if (turnj==1) {
1904 j = q2;
1905 } else {
1906 j = w2-q2;
1907 }
1908 Transform(i,j,0);
1909 x1 = fXt;
1910 y1 = fYt;
1911 Transform(i,j,-1);
1912 x1d = fXt;
1913 y1d = fYt;
1914 do {
1915 if (flag==0) {
1916 flag = 1;
1917 if (smer==0) q1 -= 1;
1918 else q2 -= 1;
1919 } else {
1920 flag = 0;
1921 if (smer==0) q2 += 1;
1922 else q1 += 1;
1923 }
1924 if (turni==1) {
1925 i = q1;
1926 } else {
1927 i = w1-q1;
1928 }
1929 if (turnj==1) {
1930 j = q2;
1931 } else {
1932 j = w2-q2;
1933 }
1934 Transform(i,j,0);
1935 x2 = fXt;
1936 y2 = fYt;
1937 if (flag==1) {
1938 x = x1;
1939 y = y1;
1940 x1 = x2;
1941 y1 = y2;
1942 x2 = x;
1943 y2 = y;
1944 }
1945 switch (fDisplayMode) {
1946 case kDisplayModePoints:
1948 Envelope(x1,y1,x2,y2);
1949 if (y1<=fEnvelope[x1]) {
1950 line->PaintLine(gPad->PixeltoX(x1) ,gPad->PixeltoY(y1)+1,
1951 gPad->PixeltoX(x1+1),gPad->PixeltoY(y1)+1);
1952 }
1953 if (y2<=fEnvelope[x2]) {
1954 line->PaintLine(gPad->PixeltoX(x2) ,gPad->PixeltoY(y2)+1,
1955 gPad->PixeltoX(x2+1),gPad->PixeltoY(y2)+1);
1956 }
1957 } else {
1958 if ((q1!=q2||smer!=0) && flag==1) {
1959 s1 = q1+1;
1960 t1 = q2;
1961 s2 = q1;
1962 t2 = q2;
1963 s3 = q1;
1964 t3 = q2+1;
1965 s4 = q1+1;
1966 t4 = q2+1;
1967 if (fShading==kShaded) {
1968 sr1 = s1;
1969 tr1 = (Int_t)TMath::Max(t1-1,0);
1970 sr2 = s2;
1971 tr2 = (Int_t)TMath::Max(t2-1,0);
1972 sr3 = (Int_t)TMath::Max(s2-1,0);
1973 tr3 = t2;
1974 sr4 = (Int_t)TMath::Max(s3-1,0);
1975 tr4 = t3;
1976 sr5 = s3;
1977 tr5 = t3+1;
1978 sr6 = s4;
1979 tr6 = t4+1;
1980 sr7 = s4+1;
1981 tr7 = t4;
1982 sr8 = s1+1;
1983 tr8 = t1;
1984 }
1985 if (turni==1) {
1986 i1 = s1;
1987 i2 = s2;
1988 i3 = s3;
1989 i4 = s4;
1990 } else {
1991 i1 = (Int_t)TMath::Max(w1-s1,0);
1992 i2 = (Int_t)TMath::Max(w1-s2,0);
1993 i3 = (Int_t)TMath::Max(w1-s3,0);
1994 i4 = (Int_t)TMath::Max(w1-s4,0);
1995 if (fShading==kShaded) {
1996 sr1 = (Int_t)TMath::Max(w1-sr1,0);
1997 sr2 = (Int_t)TMath::Max(w1-sr2,0);
1998 sr3 = (Int_t)TMath::Max(w1-sr3,0);
1999 sr4 = (Int_t)TMath::Max(w1-sr4,0);
2000 sr5 = (Int_t)TMath::Max(w1-sr5,0);
2001 sr6 = (Int_t)TMath::Max(w1-sr6,0);
2002 sr7 = (Int_t)TMath::Max(w1-sr7,0);
2003 sr8 = (Int_t)TMath::Max(w1-sr8,0);
2004 }
2005 }
2006 if (turnj==1) {
2007 j1 = t1;
2008 j2 = t2;
2009 j3 = t3;
2010 j4 = t4;
2011 } else {
2012 j1 = (Int_t)TMath::Max(w2-t1,0);
2013 j2 = (Int_t)TMath::Max(w2-t2,0);
2014 j3 = (Int_t)TMath::Max(w2-t3,0);
2015 j4 = (Int_t)TMath::Max(w2-t4,0);
2016 if (fShading==kShaded) {
2017 tr1 = (Int_t)TMath::Max(w2-tr1,0);
2018 tr2 = (Int_t)TMath::Max(w2-tr2,0);
2019 tr3 = (Int_t)TMath::Max(w2-tr3,0);
2020 tr4 = (Int_t)TMath::Max(w2-tr4,0);
2021 tr5 = (Int_t)TMath::Max(w2-tr5,0);
2022 tr6 = (Int_t)TMath::Max(w2-tr6,0);
2023 tr7 = (Int_t)TMath::Max(w2-tr7,0);
2024 tr8 = (Int_t)TMath::Max(w2-tr8,0);
2025 }
2026 }
2027 Transform(i1,j1,0);
2028 x1 = fXt;
2029 y1 = fYt;
2030 dx1 = fDxspline;
2031 dy1 = fDyspline;
2032 z1 = fZ;
2033 Transform(i2,j2,0);
2034 x2 = fXt;
2035 y2 = fYt;
2036 dx2 = fDxspline;
2037 dy2 = fDyspline;
2038 z2 = fZ;
2039 Transform(i3,j3,0);
2040 x3 = fXt;
2041 y3 = fYt;
2042 dx3 = fDxspline;
2043 dy3 = fDyspline;
2044 z3 = fZ;
2045 Transform(i4,j4,0);
2046 x4 = fXt;
2047 y4 = fYt;
2048 dx4 = fDxspline;
2049 dy4 = fDyspline;
2050 z4 = fZ;
2051 Envelope(x1,y1,x2,y2);
2052 Envelope(x2,y2,x3,y3);
2053 xtaz = (dx1+dx2+dx4)/3;
2054 ytaz = (dy1+dy2+dy4)/3;
2055 ztaz = (z1+z2+z4)/3;
2056 v = ColorCalculation(dx1,dy1,z1,dx2,dy2,z2,dx4,dy4,z4);
2057 if (fShading==kShaded) {
2059 if (sr1<0||sr1>w1||tr1<0||tr1>w2) Transform(sr1,tr1,-1);
2060 else Transform(sr1,tr1,0);
2061 dxr1 = fDxspline;
2062 dyr1 = fDyspline;
2063 zr1 = fZ;
2064 if (sr8<0||sr8>w1||tr8<0||tr8>w2) Transform(sr8,tr8,-1);
2065 else Transform(sr8,tr8,0);
2066 dxr2 = fDxspline;
2067 dyr2 = fDyspline;
2068 zr2 = fZ;
2069 v = v+ColorCalculation(dxr1,dyr1,zr1,dx2,dy2,z2,dx1,dy1,z1);
2070 v = v+ColorCalculation(dxr1,dyr1,zr1,dx1,dy1,z1,dxr2,dyr2,zr2);
2071 v = v+ColorCalculation(dxr2,dyr2,zr2,dx1,dy1,z1,dx4,dy4,z4);
2072 v1 = v/4;
2073 if (sr3<0||sr3>w1||tr3<0||tr3>w2) Transform(sr3,tr3,-1);
2074 else Transform(sr3,tr3,0);
2075 dxr1 = fDxspline;
2076 dyr1 = fDyspline;
2077 zr1 = fZ;
2078 if (sr2<0||sr2>w1||tr2<0||tr2>w2) Transform(sr2,tr2,-1);
2079 else Transform(sr2,tr2,0);
2080 dxr2 = fDxspline;
2081 dyr2 = fDyspline;
2082 zr2 = fZ;
2083 v = ColorCalculation(dx1,dy1,z1,dx2,dy2,z2,dx3,dy3,z3);
2084 v = v+ColorCalculation(dx2,dy2,z2,dxr1,dyr1,zr1,dx3,dy3,z3);
2085 v = v+ColorCalculation(dx2,dy2,z2,dxr2,dyr2,zr2,dxr1,dyr1,zr1);
2086 v = v+ColorCalculation(dxr2,dyr2,zr2,dx2,dy2,z2,dx1,dy1,z1);
2087 v2 = v/4;
2088 if (sr5<0||sr5>w1||tr5<0||tr5>w2) Transform(sr5,tr5,-1);
2089 else Transform(sr5,tr5,0);
2090 dxr1 = fDxspline;
2091 dyr1 = fDyspline;
2092 zr1 = fZ;
2093 if (sr4<0||sr4>w1||tr4<0||tr4>w2) Transform(sr4,tr4,-1);
2094 else Transform(sr4,tr4,0);
2095 dxr2 = fDxspline;
2096 dyr2 = fDyspline;
2097 zr2 = fZ;
2098 v = ColorCalculation(dx2,dy2,z2,dx3,dy3,z3,dx4,dy4,z4);
2099 v = v+ColorCalculation(dx4,dy4,z4,dx3,dy3,z3,dxr1,dyr1,zr1);
2100 v = v+ColorCalculation(dx3,dy3,z3,dxr2,dyr2,zr2,dxr1,dyr1,zr1);
2101 v = v+ColorCalculation(dx2,dy2,z2,dxr2,dyr2,zr2,dx3,dy3,z3);
2102 v3 = v/4;
2103 if (sr7<0||sr7>w1||tr7<0||tr7>w2) Transform(sr7,tr7,-1);
2104 else Transform(sr7,tr7,0);
2105 dxr1 = fDxspline;
2106 dyr1 = fDyspline;
2107 zr1 = fZ;
2108 if (sr6<0||sr6>w1||tr6<0||tr6>w2) Transform(sr6,tr6,-1);
2109 else Transform(sr6,tr6,0);
2110 dxr2 = fDxspline;
2111 dyr2 = fDyspline;
2112 zr2 = fZ;
2113 v = ColorCalculation(dx1,dy1,z1,dx3,dy3,z3,dx4,dy4,z4);
2114 v = v+ColorCalculation(dx4,dy4,z4,dx3,dy3,z3,dxr2,dyr2,zr2);
2115 v = v+ColorCalculation(dx4,dy4,z4,dxr2,dyr2,zr2,dxr1,dyr1,zr1);
2116 v = v+ColorCalculation(dx1,dy1,z1,dx4,dy4,z4,dxr1,dyr1,zr1);
2117 v4 = v/4;
2118 } else {
2119 spriz = ShadowColorCalculation(xtaz,ytaz,ztaz,shad_noise);
2120 v = v+spriz;
2121 v = v/2;
2122 if (sr1<0||sr1>w1||tr1<0||tr1>w2) Transform(sr1,tr1,-1);
2123 else Transform(sr1,tr1,0);
2124 dxr1 = fDxspline;
2125 dyr1 = fDyspline;
2126 zr1 = fZ;
2127 if (sr8<0||sr8>w1||tr8<0||tr8>w2) Transform(sr8,tr8,-1);
2128 else Transform(sr8,tr8,0);
2129 dxr2 = fDxspline;
2130 dyr2 = fDyspline;
2131 zr2 = fZ;
2132 da = (dxr1+dx2+dx1)/3;
2133 db = (dyr1+dy2+dy1)/3;
2134 dc = (zr1+z2+z1)/3;
2135 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2136 v = v+(ColorCalculation(dxr1,dyr1,zr1,dx2,dy2,z2,dx1,dy1,z1)+spriz)/2;
2137 da = (dxr1+dxr2+dx1)/3;
2138 db = (dyr1+dyr2+dy1)/3;
2139 dc = (zr1+zr2+z1)/3;
2140 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2141 v = v+(ColorCalculation(dxr1,dyr1,zr1,dx1,dy1,z1,dxr2,dyr2,zr2)+spriz)/2;
2142 da = (dxr2+dx1+dx4)/3;
2143 db = (dyr2+dy1+dy4)/3;
2144 dc = (zr2+z1+z4)/3;
2145 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2146 v = v+(ColorCalculation(dxr2,dyr2,zr2,dx1,dy1,z1,dx4,dy4,z4)+spriz)/2;
2147 v1 = v/4;
2148 if (sr3<0||sr3>w1||tr3<0||tr3>w2) Transform(sr3,tr3,-1);
2149 else Transform(sr3,tr3,0);
2150 dxr1 = fDxspline;
2151 dyr1 = fDyspline;
2152 zr1 = fZ;
2153 if (sr2<0||sr2>w1||tr2<0||tr2>w2) Transform(sr2,tr2,-1);
2154 else Transform(sr2,tr2,0);
2155 dxr2 = fDxspline;
2156 dyr2 = fDyspline;
2157 zr2 = fZ;
2158 da = (dx1+dx2+dx3)/3;
2159 db = (dy1+dy2+dy3)/3;
2160 dc = (z1+z2+z3)/3;
2161 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2162 v = (ColorCalculation(dx1,dy1,z1,dx2,dy2,z2,dx3,dy3,z3)+spriz)/2;
2163 da = (dx2+dxr1+dx3)/3;
2164 db = (dy2+dyr1+dy3)/3;
2165 dc = (z2+zr1+z3)/3;
2166 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2167 v = v+(ColorCalculation(dx2,dy2,z2,dxr1,dyr1,zr1,dx3,dy3,z3)+spriz)/2;
2168 da = (dx2+dxr2+dxr1)/3;
2169 db = (dy2+dyr2+dyr1)/3;
2170 dc = (z2+zr2+zr1)/3;
2171 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2172 v = v+(ColorCalculation(dx2,dy2,z2,dxr2,dyr2,zr2,dxr1,dyr1,zr1)+spriz)/2;
2173 da = (dxr2+dx2+dx1)/3;
2174 db = (dyr2+dy2+dy1)/3;
2175 dc = (zr2+z2+z1)/3;
2176 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2177 v = v+(ColorCalculation(dxr2,dyr2,zr2,dx2,dy2,z2,dx1,dy1,z1)+spriz)/2;
2178 v2 = v/4;
2179 if (sr5<0||sr5>w1||tr5<0||tr5>w2) Transform(sr5,tr5,-1);
2180 else Transform(sr5,tr5,0);
2181 dxr1 = fDxspline;
2182 dyr1 = fDyspline;
2183 zr1 = fZ;
2184 if (sr4<0||sr4>w1||tr4<0||tr4>w2) Transform(sr4,tr4,-1);
2185 else Transform(sr4,tr4,0);
2186 dxr2 = fDxspline;
2187 dyr2 = fDyspline;
2188 zr2 = fZ;
2189 da = (dx2+dx3+dx4)/3;
2190 db = (dy2+dy3+dy4)/3;
2191 dc = (z2+z3+z4)/3;
2192 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2193 v = (ColorCalculation(dx2,dy2,z2,dx3,dy3,z3,dx4,dy4,z4)+spriz)/2;
2194 da = (dx4+dx3+dxr1)/3;
2195 db = (dy4+dy3+dyr1)/3;
2196 dc = (z4+z3+zr1)/3;
2197 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2198 v = v+(ColorCalculation(dx4,dy4,z4,dx3,dy3,z3,dxr1,dyr1,zr1)+spriz)/2;
2199 da = (dx3+dxr2+dxr1)/3;
2200 db = (dy3+dyr2+dyr1)/3;
2201 dc = (z3+zr2+zr1)/3;
2202 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2203 v = v+(ColorCalculation(dx3,dy3,z3,dxr2,dyr2,zr2,dxr1,dyr1,zr1)+spriz)/2;
2204 da = (dx2+dxr2+dx3)/3;
2205 db = (dy2+dyr2+dy3)/3;
2206 dc = (z2+zr2+z3)/3;
2207 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2208 v = v+(ColorCalculation(dx2,dy2,z2,dxr2,dyr2,zr2,dx3,dy3,z3)+spriz)/2;
2209 v3 = v/4;
2210 if (sr7<0||sr7>w1||tr7<0||tr7>w2) Transform(sr7,tr7,-1);
2211 else Transform(sr7,tr7,0);
2212 dxr1 = fDxspline;
2213 dyr1 = fDyspline;
2214 zr1 = fZ;
2215 if (sr6<0||sr6>w1||tr6<0||tr6>w2) Transform(sr6,tr6,-1);
2216 else Transform(sr6,tr6,0);
2217 dxr2 = fDxspline;
2218 dyr2 = fDyspline;
2219 zr2 = fZ;
2220 da = (dx1+dx3+dx4)/3;
2221 db = (dy1+dy3+dy4)/3;
2222 dc = (z1+z3+z4)/3;
2223 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2224 v = (ColorCalculation(dx1,dy1,z1,dx3,dy3,z3,dx4,dy4,z4)+spriz)/2;
2225 da = (dx4+dx3+dxr2)/3;
2226 db = (dy4+dy3+dyr2)/3;
2227 dc = (z4+z3+zr2)/3;
2228 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2229 v = v+(ColorCalculation(dx4,dy4,z4,dx3,dy3,z3,dxr2,dyr2,zr2)+spriz)/2;
2230 da = (dx4+dxr2+dxr1)/3;
2231 db = (dy4+dyr2+dyr1)/3;
2232 dc = (z4+zr2+zr1)/3;
2233 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2234 v = v+(ColorCalculation(dx4,dy4,z4,dxr2,dyr2,zr2,dxr1,dyr1,zr1)+spriz)/2;
2235 da = (dx1+dx4+dxr1)/3;
2236 db = (dy1+dy4+dyr1)/3;
2237 dc = (z1+z4+zr1)/3;
2238 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2239 v = v+(ColorCalculation(dx1,dy1,z1,dx4,dy4,z4,dxr1,dyr1,zr1)+spriz)/2;
2240 v4 = v/4;
2241 }
2242 }
2243 spriz = 0;
2245 if (fShading==kNotShaded) {
2246 v = v*fLevels+0.5;
2247 iv = fLevels-(Int_t)v;
2248 } else {
2249 v1 = v1*fLevels;
2250 iv1 = fLevels-(Int_t)v1;
2251 v2 = v2*fLevels;
2252 iv2 = fLevels-(Int_t)v2;
2253 v4 = v4*fLevels;
2254 iv4 = fLevels-(Int_t)v4;
2255 }
2256 } else {
2257 spriz = ShadowColorCalculation(xtaz,ytaz,ztaz,shad_noise);
2258 if (fShading==kNotShaded) {
2259 v = v*fLevels/2.0;
2260 iv = fLevels-(Int_t)(v+0.5);
2261 } else {
2262 v1 = v1*fLevels;
2263 iv1 = fLevels-(Int_t)v1;
2264 v2 = v2*fLevels;
2265 iv2 = fLevels-(Int_t)v2;
2266 v4 = v4*fLevels;
2267 iv4 = fLevels-(Int_t)v4;
2268 }
2269 }
2270 if (fShading==kNotShaded) {
2271 ColorModel(iv,ui1,ui2,ui3);
2273 if (fEnvelope[x1]>=y1) {
2274 line->PaintLine(gPad->PixeltoX(x1),gPad->PixeltoY(y1)+1,gPad->PixeltoX(x1+1),gPad->PixeltoY(y1)+1);
2275 fEnvelope[x1] = y1;
2276 }
2277 if (fEnvelope[x2]>=y2) {
2278 line->PaintLine(gPad->PixeltoX(x2),gPad->PixeltoY(y2)+1,gPad->PixeltoX(x2+1),gPad->PixeltoY(y2)+1);
2279 fEnvelope[x2] = y2;
2280 }
2281 if (fEnvelope[x4]>=y4) {
2282 line->PaintLine(gPad->PixeltoX(x4),gPad->PixeltoY(y4)+1,gPad->PixeltoX(x4+1),gPad->PixeltoY(y4)+1);
2283 fEnvelope[x4] = y4;
2284 }
2285 } else {
2286 if (fEnvelope[x1]>=y1) {
2287 iv = iv1;
2288 ColorModel(iv,ui1,ui2,ui3);
2290 line->PaintLine(gPad->PixeltoX(x1),gPad->PixeltoY(y1)+1,gPad->PixeltoX(x1+1),gPad->PixeltoY(y1)+1);
2291 fEnvelope[x1] = y1;
2292 }
2293 if (fEnvelope[x2]>=y2) {
2294 iv = iv2;
2295 ColorModel(iv,ui1,ui2,ui3);
2297 line->PaintLine(gPad->PixeltoX(x2),gPad->PixeltoY(y2)+1,gPad->PixeltoX(x2+1),gPad->PixeltoY(y2)+1);
2298 fEnvelope[x2]=y2;
2299 }
2300 if (fEnvelope[x4]>=y4) {
2301 iv = iv4;
2302 ColorModel(iv,ui1,ui2,ui3);
2304 line->PaintLine(gPad->PixeltoX(x4),gPad->PixeltoY(y4)+1,gPad->PixeltoX(x4+1),gPad->PixeltoY(y4)+1);
2305 fEnvelope[x4] = y4;
2306 }
2307 }
2308 xtaz = (dx3+dx2+dx4)/3;
2309 ytaz = (dy3+dy2+dy4)/3;
2310 ztaz = (z3+z2+z4)/3;
2311 if (fShading==kNotShaded) v=ColorCalculation(dx2,dy2,z2,dx3,dy3,z3,dx4,dy4,z4);
2312 spriz = 0;
2314 if (fShading==kNotShaded) {
2315 v = v*fLevels;
2316 iv = fLevels-(Int_t)v;
2317 } else {
2318 v3 = v3*fLevels;
2319 iv3 = fLevels-(Int_t)v3;
2320 }
2321 } else {
2322 spriz = ShadowColorCalculation(xtaz,ytaz,ztaz,shad_noise);
2323 if (fShading==kNotShaded) {
2324 v = v*fLevels/2;
2325 iv = fLevels-(Int_t)v;
2326 iv = (Int_t)(iv-fLevels*spriz/2);
2327 } else {
2328 v3 = v3*fLevels;
2329 iv3 = fLevels-(Int_t)v3;
2330 }
2331 }
2332 if (fShading==kNotShaded) {
2333 ColorModel(iv,ui1,ui2,ui3);
2334 line->ResetAttLine("");
2336 if (fEnvelope[x3]>=y3) {
2337 line->PaintLine(gPad->PixeltoX(x3),gPad->PixeltoY(y3)+1,gPad->PixeltoX(x3+1),gPad->PixeltoY(y3)+1);
2338 fEnvelope[x3] = y3;
2339 }
2340 } else {
2341 if (fEnvelope[x3]>=y3) {
2342 iv = iv3;
2343 ColorModel(iv,ui1,ui2,ui3);
2344 line->ResetAttLine("");
2346 line->PaintLine(gPad->PixeltoX(x3),gPad->PixeltoY(y3)+1,gPad->PixeltoX(x3+1),gPad->PixeltoY(y3)+1);
2347 fEnvelope[x3]=y3;
2348 }
2349 }
2350 }
2351 }
2352 break;
2353 case kDisplayModeGrid:
2356 Envelope(x1,y1,x2,y2);
2357 if (fLine!=0) {
2358 if (fLine==1) {
2359 fXe = x2;
2360 fYe = y2;
2361 }
2362 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
2363 }
2364 } else {
2365 if ((q1!=q2||smer!=0)&&flag==1) {
2366 s1 = q1+1;
2367 t1 = q2;
2368 s2 = q1;
2369 t2 = q2;
2370 s3 = q1;
2371 t3 = q2+1;
2372 s4 = q1+1;
2373 t4 = q2+1;
2374 if (fShading==kShaded) {
2375 sr1 = s1;
2376 tr1 = (Int_t)TMath::Max(t1-1,0);
2377 sr2 = s2;
2378 tr2 = (Int_t)TMath::Max(t2-1,0);
2379 sr3 = (Int_t)TMath::Max(s2-1,0);
2380 tr3 = t2;
2381 sr4 = (Int_t)TMath::Max(s3-1,0);
2382 tr4 = t3;
2383 sr5 = s3;
2384 tr5 = t3+1;
2385 sr6 = s4;
2386 tr6 = t4+1;
2387 sr7 = s4+1;
2388 tr7 = t4;
2389 sr8 = s1+1;
2390 tr8 = t1;
2391 }
2392 if (turni==1) {
2393 i1 = s1;
2394 i2 = s2;
2395 i3 = s3;
2396 i4 = s4;
2397 } else {
2398 i1 = (Int_t)TMath::Max(w1-s1,0);
2399 i2 = (Int_t)TMath::Max(w1-s2,0);
2400 i3 = (Int_t)TMath::Max(w1-s3,0);
2401 i4 = (Int_t)TMath::Max(w1-s4,0);
2402 if (fShading==kShaded) {
2403 sr1 = (Int_t)TMath::Max(w1-sr1,0);
2404 sr2 = (Int_t)TMath::Max(w1-sr2,0);
2405 sr3 = (Int_t)TMath::Max(w1-sr3,0);
2406 sr4 = (Int_t)TMath::Max(w1-sr4,0);
2407 sr5 = (Int_t)TMath::Max(w1-sr5,0);
2408 sr6 = (Int_t)TMath::Max(w1-sr6,0);
2409 sr7 = (Int_t)TMath::Max(w1-sr7,0);
2410 sr8 = (Int_t)TMath::Max(w1-sr8,0);
2411 }
2412 }
2413 if (turnj==1) {
2414 j1 = t1;
2415 j2 = t2;
2416 j3 = t3;
2417 j4 = t4;
2418 } else {
2419 j1 = (Int_t)TMath::Max(w2-t1,0);
2420 j2 = (Int_t)TMath::Max(w2-t2,0);
2421 j3 = (Int_t)TMath::Max(w2-t3,0);
2422 j4 = (Int_t)TMath::Max(w2-t4,0);
2423 if (fShading==kShaded) {
2424 tr1 = (Int_t)TMath::Max(w2-tr1,0);
2425 tr2 = (Int_t)TMath::Max(w2-tr2,0);
2426 tr3 = (Int_t)TMath::Max(w2-tr3,0);
2427 tr4 = (Int_t)TMath::Max(w2-tr4,0);
2428 tr5 = (Int_t)TMath::Max(w2-tr5,0);
2429 tr6 = (Int_t)TMath::Max(w2-tr6,0);
2430 tr7 = (Int_t)TMath::Max(w2-tr7,0);
2431 tr8 = (Int_t)TMath::Max(w2-tr8,0);
2432 }
2433 }
2434 Transform(i1,j1,0);
2435 x1 = fXt;
2436 y1 = fYt;
2437 dx1 = fDxspline;
2438 dy1 = fDyspline;
2439 z1 = fZ;
2440 Transform(i2,j2,0);
2441 x2 = fXt;
2442 y2 = fYt;
2443 dx2 = fDxspline;
2444 dy2 = fDyspline;
2445 z2 = fZ;
2446 Transform(i3,j3,0);
2447 x3 = fXt;
2448 y3 = fYt;
2449 dx3 = fDxspline;
2450 dy3 = fDyspline;
2451 z3 = fZ;
2452 Transform(i4,j4,0);
2453 x4 = fXt;
2454 y4 = fYt;
2455 dx4 = fDxspline;
2456 dy4 = fDyspline;
2457 z4 = fZ;
2458 Envelope(x1,y1,x2,y2);
2459 Envelope(x2,y2,x3,y3);
2460 xtaz = (dx1+dx2+dx4)/3;
2461 ytaz = (dy1+dy2+dy4)/3;
2462 ztaz = (z1+z2+z4)/3;
2463 v = ColorCalculation(dx1,dy1,z1,dx2,dy2,z2,dx4,dy4,z4);
2464 if (fShading==kShaded) {
2466 if (sr1<0||sr1>w1||tr1<0||tr1>w2) Transform(sr1,tr1,-1);
2467 else Transform(sr1,tr1,0);
2468 dxr1 = fDxspline;
2469 dyr1 = fDyspline;
2470 zr1 = fZ;
2471 if (sr8<0||sr8>w1||tr8<0||tr8>w2) Transform(sr8,tr8,-1);
2472 else Transform(sr8,tr8,0);
2473 dxr2 = fDxspline;
2474 dyr2 = fDyspline;
2475 zr2 = fZ;
2476 v = v+ColorCalculation(dxr1,dyr1,zr1,dx2,dy2,z2,dx1,dy1,z1);
2477 v = v+ColorCalculation(dxr1,dyr1,zr1,dx1,dy1,z1,dxr2,dyr2,zr2);
2478 v = v+ColorCalculation(dxr2,dyr2,zr2,dx1,dy1,z1,dx4,dy4,z4);
2479 v1 = v/4;
2480 if (sr3<0||sr3>w1||tr3<0||tr3>w2) Transform(sr3,tr3,-1);
2481 else Transform(sr3,tr3,0);
2482 dxr1 = fDxspline;
2483 dyr1 = fDyspline;
2484 zr1 = fZ;
2485 if (sr2<0||sr2>w1||tr2<0||tr2>w2) Transform(sr2,tr2,-1);
2486 else Transform(sr2,tr2,0);
2487 dxr2 = fDxspline;
2488 dyr2 = fDyspline;
2489 zr2 = fZ;
2490 v = ColorCalculation(dx1,dy1,z1,dx2,dy2,z2,dx3,dy3,z3);
2491 v = v+ColorCalculation(dx2,dy2,z2,dxr1,dyr1,zr1,dx3,dy3,z3);
2492 v = v+ColorCalculation(dx2,dy2,z2,dxr2,dyr2,zr2,dxr1,dyr1,zr1);
2493 v = v+ColorCalculation(dxr2,dyr2,zr2,dx2,dy2,z2,dx1,dy1,z1);
2494 v2 = v/4;
2495 if (sr5<0||sr5>w1||tr5<0||tr5>w2) Transform(sr5,tr5,-1);
2496 else Transform(sr5,tr5,0);
2497 dxr1 = fDxspline;
2498 dyr1 = fDyspline;
2499 zr1 = fZ;
2500 if (sr4<0||sr4>w1||tr4<0||tr4>w2) Transform(sr4,tr4,-1);
2501 else Transform(sr4,tr4,0);
2502 dxr2 = fDxspline;
2503 dyr2 = fDyspline;
2504 zr2 = fZ;
2505 v = ColorCalculation(dx2,dy2,z2,dx3,dy3,z3,dx4,dy4,z4);
2506 v = v+ColorCalculation(dx4,dy4,z4,dx3,dy3,z3,dxr1,dyr1,zr1);
2507 v = v+ColorCalculation(dx3,dy3,z3,dxr2,dyr2,zr2,dxr1,dyr1,zr1);
2508 v = v+ColorCalculation(dx2,dy2,z2,dxr2,dyr2,zr2,dx3,dy3,z3);
2509 v3 = v/4;
2510 if (sr7<0||sr7>w1||tr7<0||tr7>w2) Transform(sr7,tr7,-1);
2511 else Transform(sr7,tr7,0);
2512 dxr1 = fDxspline;
2513 dyr1 = fDyspline;
2514 zr1 = fZ;
2515 if (sr6<0||sr6>w1||tr6<0||tr6>w2) Transform(sr6,tr6,-1);
2516 else Transform(sr6,tr6,0);
2517 dxr2 = fDxspline;
2518 dyr2 = fDyspline;
2519 zr2 = fZ;
2520 v = ColorCalculation(dx1,dy1,z1,dx3,dy3,z3,dx4,dy4,z4);
2521 v = v+ColorCalculation(dx4,dy4,z4,dx3,dy3,z3,dxr2,dyr2,zr2);
2522 v = v+ColorCalculation(dx4,dy4,z4,dxr2,dyr2,zr2,dxr1,dyr1,zr1);
2523 v = v+ColorCalculation(dx1,dy1,z1,dx4,dy4,z4,dxr1,dyr1,zr1);
2524 v4 = v/4;
2525 } else {
2526 spriz = ShadowColorCalculation(xtaz,ytaz,ztaz,shad_noise);
2527 v = v+spriz;
2528 v = v/2;
2529 if (sr1<0||sr1>w1||tr1<0||tr1>w2) Transform(sr1,tr1,-1);
2530 else Transform(sr1,tr1,0);
2531 dxr1 = fDxspline;
2532 dyr1 = fDyspline;
2533 zr1 = fZ;
2534 if (sr8<0||sr8>w1||tr8<0||tr8>w2) Transform(sr8,tr8,-1);
2535 else Transform(sr8,tr8,0);
2536 dxr2 = fDxspline;
2537 dyr2 = fDyspline;
2538 zr2 = fZ;
2539 da = (dxr1+dx2+dx1)/3;
2540 db = (dyr1+dy2+dy1)/3;
2541 dc = (zr1+z2+z1)/3;
2542 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2543 v = v+(ColorCalculation(dxr1,dyr1,zr1,dx2,dy2,z2,dx1,dy1,z1)+spriz)/2;
2544 da = (dxr1+dxr2+dx1)/3;
2545 db = (dyr1+dyr2+dy1)/3;
2546 dc = (zr1+zr2+z1)/3;
2547 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2548 v = v+(ColorCalculation(dxr1,dyr1,zr1,dx1,dy1,z1,dxr2,dyr2,zr2)+spriz)/2;
2549 da = (dxr2+dx1+dx4)/3;
2550 db = (dyr2+dy1+dy4)/3;
2551 dc = (zr2+z1+z4)/3;
2552 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2553 v = v+(ColorCalculation(dxr2,dyr2,zr2,dx1,dy1,z1,dx4,dy4,z4)+spriz)/2;
2554 v1 = v/4;
2555 if (sr3<0||sr3>w1||tr3<0||tr3>w2) Transform(sr3,tr3,-1);
2556 else Transform(sr3,tr3,0);
2557 dxr1 = fDxspline;
2558 dyr1 = fDyspline;
2559 zr1 = fZ;
2560 if (sr2<0||sr2>w1||tr2<0||tr2>w2) Transform(sr2,tr2,-1);
2561 else Transform(sr2,tr2,0);
2562 dxr2 = fDxspline;
2563 dyr2 = fDyspline;
2564 zr2 = fZ;
2565 da = (dx1+dx2+dx3)/3;
2566 db = (dy1+dy2+dy3)/3;
2567 dc = (z1+z2+z3)/3;
2568 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2569 v = (ColorCalculation(dx1,dy1,z1,dx2,dy2,z2,dx3,dy3,z3)+spriz)/2;
2570 da = (dx2+dxr1+dx3)/3;
2571 db = (dy2+dyr1+dy3)/3;
2572 dc = (z2+zr1+z3)/3;
2573 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2574 v = v+(ColorCalculation(dx2,dy2,z2,dxr1,dyr1,zr1,dx3,dy3,z3)+spriz)/2;
2575 da = (dx2+dxr2+dxr1)/3;
2576 db = (dy2+dyr2+dyr1)/3;
2577 dc = (z2+zr2+zr1)/3;
2578 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2579 v = v+(ColorCalculation(dx2,dy2,z2,dxr2,dyr2,zr2,dxr1,dyr1,zr1)+spriz)/2;
2580 da = (dxr2+dx2+dx1)/3;
2581 db = (dyr2+dy2+dy1)/3;
2582 dc = (zr2+z2+z1)/3;
2583 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2584 v = v+(ColorCalculation(dxr2,dyr2,zr2,dx2,dy2,z2,dx1,dy1,z1)+spriz)/2;
2585 v2 = v/4;
2586 if (sr5<0||sr5>w1||tr5<0||tr5>w2) Transform(sr5,tr5,-1);
2587 else Transform(sr5,tr5,0);
2588 dxr1 = fDxspline;
2589 dyr1 = fDyspline;
2590 zr1 = fZ;
2591 if (sr4<0||sr4>w1||tr4<0||tr4>w2) Transform(sr4,tr4,-1);
2592 else Transform(sr4,tr4,0);
2593 dxr2 = fDxspline;
2594 dyr2 = fDyspline;
2595 zr2 = fZ;
2596 da = (dx2+dx3+dx4)/3;
2597 db = (dy2+dy3+dy4)/3;
2598 dc = (z2+z3+z4)/3;
2599 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2600 v = (ColorCalculation(dx2,dy2,z2,dx3,dy3,z3,dx4,dy4,z4)+spriz)/2;
2601 da = (dx4+dx3+dxr1)/3;
2602 db = (dy4+dy3+dyr1)/3;
2603 dc = (z4+z3+zr1)/3;
2604 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2605 v = v+(ColorCalculation(dx4,dy4,z4,dx3,dy3,z3,dxr1,dyr1,zr1)+spriz)/2;
2606 da = (dx3+dxr2+dxr1)/3;
2607 db = (dy3+dyr2+dyr1)/3;
2608 dc = (z3+zr2+zr1)/3;
2609 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2610 v = v+(ColorCalculation(dx3,dy3,z3,dxr2,dyr2,zr2,dxr1,dyr1,zr1)+spriz)/2;
2611 da = (dx2+dxr2+dx3)/3;
2612 db = (dy2+dyr2+dy3)/3;
2613 dc = (z2+zr2+z3)/3;
2614 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2615 v = v+(ColorCalculation(dx2,dy2,z2,dxr2,dyr2,zr2,dx3,dy3,z3)+spriz)/2;
2616 v3 = v/4;
2617 if (sr7<0||sr7>w1||tr7<0||tr7>w2) Transform(sr7,tr7,-1);
2618 else Transform(sr7,tr7,0);
2619 dxr1 = fDxspline;
2620 dyr1 = fDyspline;
2621 zr1 = fZ;
2622 if (sr6<0||sr6>w1||tr6<0||tr6>w2) Transform(sr6,tr6,-1);
2623 else Transform(sr6,tr6,0);
2624 dxr2 = fDxspline;
2625 dyr2 = fDyspline;
2626 zr2 = fZ;
2627 da = (dx1+dx3+dx4)/3;
2628 db = (dy1+dy3+dy4)/3;
2629 dc = (z1+z3+z4)/3;
2630 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2631 v = (ColorCalculation(dx1,dy1,z1,dx3,dy3,z3,dx4,dy4,z4)+spriz)/2;
2632 da = (dx4+dx3+dxr2)/3;
2633 db = (dy4+dy3+dyr2)/3;
2634 dc = (z4+z3+zr2)/3;
2635 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2636 v = v+(ColorCalculation(dx4,dy4,z4,dx3,dy3,z3,dxr2,dyr2,zr2)+spriz)/2;
2637 da = (dx4+dxr2+dxr1)/3;
2638 db = (dy4+dyr2+dyr1)/3;
2639 dc = (z4+zr2+zr1)/3;
2640 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2641 v = v+(ColorCalculation(dx4,dy4,z4,dxr2,dyr2,zr2,dxr1,dyr1,zr1)+spriz)/2;
2642 da = (dx1+dx4+dxr1)/3;
2643 db = (dy1+dy4+dyr1)/3;
2644 dc = (z1+z4+zr1)/3;
2645 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
2646 v = v+(ColorCalculation(dx1,dy1,z1,dx4,dy4,z4,dxr1,dyr1,zr1)+spriz)/2;
2647 v4 = v/4;
2648 }
2649 }
2650 spriz = 0;
2652 if (fShading==kNotShaded) {
2653 v = v*fLevels+0.5;
2654 iv = fLevels-(Int_t)v;
2655 } else {
2656 v1 = v1*fLevels;
2657 iv1 = fLevels-(Int_t)v1;
2658 v2 = v2*fLevels;
2659 iv2 = fLevels-(Int_t)v2;
2660 v4 = v4*fLevels;
2661 iv4 = fLevels-(Int_t)v4;
2662 }
2663 } else {
2664 spriz = ShadowColorCalculation(xtaz,ytaz,ztaz,shad_noise);
2665 if (fShading==kNotShaded) {
2666 v = v*fLevels/2.0;
2667 iv = fLevels-(Int_t)(v+0.5);
2668 } else {
2669 v1 = v1*fLevels;
2670 iv1 = fLevels-(Int_t)v1;
2671 v2 = v2*fLevels;
2672 iv2 = fLevels-(Int_t)v2;
2673 v4 = v4*fLevels;
2674 iv4 = fLevels-(Int_t)v4;
2675 }
2676 }
2677 if (fShading==kNotShaded) {
2678 ColorModel(iv,ui1,ui2,ui3);
2680 } else {
2681 dx1 = x1;
2682 dy1 = y1;
2683 dx2 = x2;
2684 dy2 = y2;
2685 dx3 = x4;
2686 dy3 = y4;
2687 z1 = iv1;
2688 z2 = iv2;
2689 z3 = iv4;
2690 da = (dy2-dy1)*(z3-z1)-(dy3-dy1)*(z2-z1);
2691 db = (z2-z1)*(dx3-dx1)-(z3-z1)*(dx2-dx1);
2692 dc = (dx2-dx1)*(dy3-dy1)-(dx3-dx1)*(dy2-dy1);
2693 dd = -da*dx1-db*dy1-dc*z1;
2694 }
2695 sx1 = x1;
2696 sy1 = y1;
2697 sx2 = x2;
2698 sy2 = y2;
2699 if (sx2<sx1) {
2700 sx4 = sx1;
2701 sy4 = sy1;
2702 sx1 = sx2;
2703 sy1 = sy2;
2704 sx2 = sx4;
2705 sy2 = sy4;
2706 }
2707 sdx1 = 0;
2708 pom1 = sy2-sy1;
2709 pom2 = sx2-sx1;
2710 if (pom2!=0) sdx1 = pom1/pom2;
2711 pom1 = sy1;
2712 pom2 = sx1;
2713 sdy1 = pom1-sdx1*pom2;
2714 for (sx4=sx1,sx5=sx1,sy5=sy1;sx4<=sx2;sx4++) {
2715 pom1 = sx4;
2716 sdy4 = sdx1*pom1+sdy1;
2717 sy4 = (Int_t)(sdy4);
2718 if (sy4<=fEnvelope[sx4]) {
2719 fEnvelope[sx4] = sy4;
2720 if (fShading==kNotShaded) {
2721 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx5),gPad->PixeltoY(sy5)+1);
2722 } else {
2723 dx1 = sx4;
2724 dy1 = sy4;
2725 if (dc!=0) v = (-dd-da*dx1-db*dy1)/dc;
2726 else v = (iv1+iv2+iv4)/3;
2727 iv = (Int_t)v;
2728 ColorModel(iv,ui1,ui2,ui3);
2730 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx5),gPad->PixeltoY(sy5)+1);
2731 }
2732 sy5 = sy4;
2733 } else {
2734 sy4 = fEnvelope[sx4];
2735 if (fShading==kNotShaded&&sy5<=fEnvelope[sx5]) {
2736 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx5),gPad->PixeltoY(sy5)+1);
2737 } else if (sy5<=fEnvelope[sx5]) {
2738 dx1 = sx4;
2739 dy1 = sy4;
2740 if (dc!=0) v = (-dd-da*dx1-db*dy1)/dc;
2741 else v = (iv1+iv2+iv4)/3;
2742 iv = (Int_t)v;
2743 ColorModel(iv,ui1,ui2,ui3);
2745 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx5),gPad->PixeltoY(sy5)+1);
2746 }
2747 sy5 = fEnvelope[sx4];
2748 }
2749 sx5 = sx4;
2750 }
2751 sx1 = x1;
2752 sy1 = y1;
2753 sx3 = x4;
2754 sy3 = y4;
2755 if (sx3<sx1) {
2756 sx4 = sx1;
2757 sy4 = sy1;
2758 sx1 = sx3;
2759 sy1 = sy3;
2760 sx3 = sx4;
2761 sy3 = sy4;
2762 }
2763 pom1 = sy3-sy1;
2764 pom2 = sx3-sx1;
2765 if (pom2!=0) sdx2 = pom1/pom2;
2766 pom1 = sy1;
2767 pom2 = sx1;
2768 sdy2 = pom1-sdx2*pom2;
2769 sx1p = sx1;
2770 sy1p = sy1;
2771 sx3p = sx3;
2772 sdx2p = sdx2;
2773 sdy2p = sdy2;
2774 dap = da;
2775 dbp = db;
2776 dcp = dc;
2777 ddp = dd;
2778 uip = fNewColorIndex;
2779 xtaz = (dx3+dx2+dx4)/3;
2780 ytaz = (dy3+dy2+dy4)/3;
2781 ztaz = (z3+z2+z4)/3;
2782 if (fShading==kNotShaded) v = ColorCalculation(dx2,dy2,z2,dx3,dy3,z3,dx4,dy4,z4);
2783 spriz = 0;
2785 if (fShading==kNotShaded) {
2786 v = v*fLevels;
2787 iv = fLevels-(Int_t)v;
2788 } else {
2789 v3 = v3*fLevels;
2790 iv3 = fLevels-(Int_t)v3;
2791 }
2792 } else {
2793 spriz = ShadowColorCalculation(xtaz,ytaz,ztaz,shad_noise);
2794 if (fShading==kNotShaded) {
2795 v = v*fLevels/2;
2796 iv = fLevels-(Int_t)v;
2797 iv = (Int_t)(iv-fLevels*spriz/2);
2798 } else {
2799 v3 = v3*fLevels;
2800 iv3 = fLevels-(Int_t)v3;
2801 }
2802 }
2803 if (fShading==kNotShaded) {
2804 ColorModel(iv,ui1,ui2,ui3);
2806 } else {
2807 dx1 = x2;
2808 dy1 = y2;
2809 dx2 = x3;
2810 dy2 = y3;
2811 dx3 = x4;
2812 dy3 = y4;
2813 z1 = iv2;
2814 z2 = iv3;
2815 z3 = iv4;
2816 da = (dy2-dy1)*(z3-z1)-(dy3-dy1)*(z2-z1);
2817 db = (z2-z1)*(dx3-dx1)-(z3-z1)*(dx2-dx1);
2818 dc = (dx2-dx1)*(dy3-dy1)-(dx3-dx1)*(dy2-dy1);
2819 dd = -da*dx1-db*dy1-dc*z1;
2820 }
2821 sx1 = x2;
2822 sy1 = y2;
2823 sx2 = x3;
2824 sy2 = y3;
2825 if (sx2<sx1) {
2826 sx4 = sx1;
2827 sy4 = sy1;
2828 sx1 = sx2;
2829 sy1 = sy2;
2830 sx2 = sx4;
2831 sy2 = sy4;
2832 }
2833 pom1 = sy2-sy1;
2834 pom2 = sx2-sx1;
2835 sdx1 = 0;
2836 if (pom2!=0) sdx1 = pom1/pom2;
2837 pom1 = sy1;
2838 pom2 = sx1;
2839 sdy1 = pom1-sdx1*pom2;
2840 for (sx4=sx1,sx5=sx1,sy5=sy1;sx4<=sx2;sx4++) {
2841 pom1 = sx4;
2842 sdy4 = sdx1*pom1+sdy1;
2843 sy4 = (Int_t)sdy4;
2844 if (sy4<=fEnvelope[sx4]) {
2845 fEnvelope[sx4] = sy4;
2846 if (fShading==kNotShaded) {
2847 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx5),gPad->PixeltoY(sy5)+1);
2848 } else {
2849 dx1 = sx4;
2850 dy1 = sy4;
2851 if (dc!=0) v = (-dd-da*dx1-db*dy1)/dc;
2852 else v = (iv1+iv2+iv4)/3;
2853 iv = (Int_t)v;
2854 ColorModel(iv,ui1,ui2,ui3);
2856 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx5),gPad->PixeltoY(sy5)+1);
2857 }
2858 sy5 = sy4;
2859 } else {
2860 sy4 = fEnvelope[sx4];
2861 if (fShading==kNotShaded&&sy5<=fEnvelope[sx5]) {
2862 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx5),gPad->PixeltoY(sy5)+1);
2863 } else if (sy5<=fEnvelope[sx5]) {
2864 dx1 = sx4;
2865 dy1 = sy4;
2866 if (dc!=0) v = (-dd-da*dx1-db*dy1)/dc;
2867 else v = (iv1+iv2+iv4)/3;
2868 iv = (Int_t)v;
2869 ColorModel(iv,ui1,ui2,ui3);
2871 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx5),gPad->PixeltoY(sy5)+1);
2872 }
2873 sy5 = fEnvelope[sx4];
2874 }
2875 sx5 = sx4;
2876 }
2877 for (sx4=sx1p,sx5=sx1p,sy5=sy1p;sx4<=sx3p;sx4++) {
2878 pom1 = sx4;
2879 sdy4 = sdx2p*pom1+sdy2p;
2880 sy4 = (Int_t)sdy4;
2881 if (sy4<=fEnvelope[sx4]) {
2882 fEnvelope[sx4]=sy4;
2883 if (fShading==kNotShaded) {
2884 line->SetLineColor(uip);
2885 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx5),gPad->PixeltoY(sy5)+1);
2886 } else {
2887 dx1 = sx4;
2888 dy1 = sy4;
2889 if (dcp!=0) v = (-ddp-dap*dx1-dbp*dy1)/dcp;
2890 else v = (iv1+iv2+iv4)/3;
2891 iv = (Int_t)v;
2892 ColorModel(iv,ui1,ui2,ui3);
2894 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx5),gPad->PixeltoY(sy5)+1);
2895 }
2896 sy5 = sy4;
2897 } else {
2898 sy4 = fEnvelope[sx4];
2899 if (fShading==kNotShaded&&sy5<=fEnvelope[sx5]) {
2900 line->SetLineColor(uip);
2901 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx5),gPad->PixeltoY(sy5)+1);
2902 } else if (sy5<=fEnvelope[sx5]) {
2903 dx1 = sx4;
2904 dy1 = sy4;
2905 if (dcp!=0) v = (-ddp-dap*dx1-dbp*dy1)/dcp;
2906 else v = (iv1+iv2+iv4)/3;
2907 iv = (Int_t)v;
2908 ColorModel(iv,ui1,ui2,ui3);
2910 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx5),gPad->PixeltoY(sy5)+1);
2911 }
2912 sy5 = fEnvelope[sx4];
2913 }
2914 sx5 = sx4;
2915 }
2916 sx2 = x3;
2917 sy2 = y3;
2918 sx3 = x4;
2919 sy3 = y4;
2920 if (sx3<sx2) {
2921 sx4 = sx2;
2922 sy4 = sy2;
2923 sx2 = sx3;
2924 sy2 = sy3;
2925 sx3 = sx4;
2926 sy3 = sy4;
2927 }
2928 sdx2 = 0;
2929 pom1 = sy3-sy2;
2930 pom2 = sx3-sx2;
2931 if (pom2!=0) sdx2 = pom1/pom2;
2932 pom1 = sy2;
2933 pom2 = sx2;
2934 sdy2 = pom1-sdx2*pom2;
2935 for (sx4=sx2,sx5=sx2,sy5=sy2;sx4<=sx3;sx4++) {
2936 pom1 = sx4;
2937 sdy4 = sdx2*pom1+sdy2;
2938 sy4 = (Int_t)sdy4;
2939 if (sy4<=fEnvelope[sx4]) {
2940 fEnvelope[sx4] = sy4;
2941 if (fShading==kNotShaded) {
2942 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx5),gPad->PixeltoY(sy5)+1);
2943 } else {
2944 dx1 = sx4;
2945 dy1 = sy4;
2946 if (dc!=0) v = (-dd-da*dx1-db*dy1)/dc;
2947 else v =(iv1+iv2+iv4)/3;
2948 iv = (Int_t)v;
2949 ColorModel(iv,ui1,ui2,ui3);
2951 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx5),gPad->PixeltoY(sy5)+1);
2952 }
2953 sy5 = sy4;
2954 } else {
2955 sy4 = fEnvelope[sx4];
2956 if (fShading==kNotShaded&&sy5<=fEnvelope[sx5]) {
2957 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx5),gPad->PixeltoY(sy5)+1);
2958 } else if (sy5<=fEnvelope[sx5]) {
2959 dx1 = sx4;
2960 dy1 = sy4;
2961 if (dc!=0) v = (-dd-da*dx1-db*dy1)/dc;
2962 else v =(iv1+iv2+iv4)/3;
2963 iv = (Int_t)v;
2964 ColorModel(iv,ui1,ui2,ui3);
2966 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx5),gPad->PixeltoY(sy5)+1);
2967 }
2968 sy5 = fEnvelope[sx4];
2969 }
2970 sx5 = sx4;
2971 }
2972 }
2973 }
2974 } else {
2975 if (((flag==0)&&(smer==0))||((flag!=0)&&(smer!=0))) {
2976 s1 = q1;
2977 t1 = (Int_t)TMath::Max(q2-1,0);
2978 s2 = q1;
2979 t2 = (Int_t)TMath::Min(q2+2,w2);
2980 } else if (((flag!=0)&&(smer==0))||((flag==0)&&(smer!=0))) {
2981 s1 = (Int_t)TMath::Max(q1-1,0);
2982 t1 = q2;
2983 s2 = (Int_t)TMath::Min(q1+2,w1);
2984 t2 = q2;
2985 }
2986 if (turni==1) {
2987 i1 = s1;
2988 i2 = s2;
2989 } else {
2990 i1 = w1-s1;
2991 i2 = w1-s2;
2992 }
2993 if (turnj==1) {
2994 j1 = t1;
2995 j2 = t2;
2996 } else {
2997 j1 = w2-t1;
2998 j2 = w2-t2;
2999 }
3000 Transform(i1,j1,0);
3001 x3 = fXt;
3002 y3 = fYt;
3003 Transform(i2,j2,0);
3004 x4 = fXt;
3005 y4 = fYt;
3006 bezx1 = x1+(x2-x1)/3;
3007 bezx2 = x1+2*(x2-x1)/3;
3008 bezy1 = y1+(y2-y3)/6;
3009 bezy2 = y2-(y4-y1)/6;
3010 if (x1<=x2) {
3011 if (bezx1<=x1) {
3012 bezx1 = x1;
3013 bezy1 = y1;
3014 }
3015 if (bezx1>=x2) {
3016 bezx1 = x2;
3017 bezy1 = y2;
3018 }
3019 if (bezx2<=x1) {
3020 bezx2 = x1;
3021 bezy2 = y1;
3022 }
3023 if (bezx2>=x2) {
3024 bezx2 = x2;
3025 bezy2 = y2;
3026 }
3027 fBzX[0] = x1;
3028 fBzY[0] = y1;
3029 fBzX[1] = (Int_t)bezx1;
3030 fBzY[1] = (Int_t)bezy1;
3031 fBzX[2] = (Int_t)bezx2;
3032 fBzY[2] = (Int_t)bezy2;
3033 fBzX[3] = x2;
3034 fBzY[3] = y2;
3035 for (bezf=0;bezf<1.01;bezf+=0.1) {
3036 BezierSmoothing(bezf);
3037 if (bezf==0) {
3038 ibezx1 = (Int_t)(fGbezx+0.5);
3039 ibezy1 = (Int_t)(fGbezy+0.5);
3040 } else {
3041 ibezx2 = ibezx1;
3042 ibezy2 = ibezy1;
3043 ibezx1 = (Int_t)(fGbezx+0.5);
3044 ibezy1 = (Int_t)(fGbezy+0.5);
3045 Envelope(ibezx2,ibezy2,ibezx1,ibezy1);
3046 if (fLine!=0) {
3047 if (fLine==1) {
3048 fXe = ibezx1;
3049 fYe = ibezy1;
3050 }
3051 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
3052 }
3053 }
3054 }
3055 } else if (x1>x2) {
3056 if (bezx1>=x1) {
3057 bezx1 = x1;
3058 bezy1 = y1;
3059 }
3060 if (bezx1<=x2) {
3061 bezx1 = x2;
3062 bezy1 = y2;
3063 }
3064 if (bezx2>=x1) {
3065 bezx2 = x1;
3066 bezy2 = y1;
3067 }
3068 if (bezx2<=x2) {
3069 bezx2 = x2;
3070 bezy2 = y2;
3071 }
3072 fBzX[0] = x1;
3073 fBzY[0] = y1;
3074 fBzX[1] = (Int_t)bezx1;
3075 fBzY[1] = (Int_t)bezy1;
3076 fBzX[2] = (Int_t)bezx2;
3077 fBzY[2] = (Int_t)bezy2;
3078 fBzX[3] = x2;
3079 fBzY[3] = y2;
3080 for (bezf=0;bezf<1.01;bezf+=0.1) {
3081 BezierSmoothing(bezf);
3082 if (bezf==0) {
3083 ibezx1 = (Int_t)(fGbezx+0.5);
3084 ibezy1 = (Int_t)(fGbezy+0.5);
3085 } else {
3086 ibezx2 = ibezx1;
3087 ibezy2 = ibezy1;
3088 ibezx1 = (Int_t)(fGbezx+0.5);
3089 ibezy1 = (Int_t)(fGbezy+0.5);
3090 Envelope(ibezx1,ibezy1,ibezx2,ibezy2);
3091 if (fLine!=0) {
3092 if (fLine==1) {
3093 fXe = ibezx2;
3094 fYe = ibezy2;
3095 }
3096 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
3097 }
3098 }
3099 }
3100 }
3101 }
3102 break;
3104 if ((q1!=q2||smer!=0)&&flag==1) {
3105 s1 = q1+1;
3106 t1 = q2;
3107 s2 = q1;
3108 t2 = q2;
3109 s3 = q1;
3110 t3 = q2+1;
3111 s4 = q1+1;
3112 t4 = q2+1;
3113 if (turni==1) {
3114 i1 = (Int_t)TMath::Min(w1,s1);
3115 i2 = (Int_t)TMath::Min(w1,s2);
3116 i3 = (Int_t)TMath::Min(w1,s3);
3117 i4 = (Int_t)TMath::Min(w1,s4);
3118 } else {
3119 i1 = (Int_t)TMath::Max(w1-s1,0);
3120 i2 = (Int_t)TMath::Max(w1-s2,0);
3121 i3 = (Int_t)TMath::Max(w1-s3,0);
3122 i4 = (Int_t)TMath::Max(w1-s4,0);
3123 }
3124 if (turnj==1) {
3125 j1 = (Int_t)TMath::Min(w2,t1);
3126 j2 = (Int_t)TMath::Min(w2,t2);
3127 j3 = (Int_t)TMath::Min(w2,t3);
3128 j4 = (Int_t)TMath::Min(w2,t4);
3129 } else {
3130 j1 = (Int_t)TMath::Max(w2-t1,0);
3131 j2 = (Int_t)TMath::Max(w2-t2,0);
3132 j3 = (Int_t)TMath::Max(w2-t3,0);
3133 j4 = (Int_t)TMath::Max(w2-t4,0);
3134 }
3135 Transform(i1,j1,0);
3136 dx1 = fDxspline;
3137 dy1 = fDyspline;
3138 z1 = fZ;
3139 z1l = fZeq;
3140 Transform(i2,j2,0);
3141 dx2 = fDxspline;
3142 dy2 = fDyspline;
3143 z2 = fZ;
3144 z2l = fZeq;
3145 Transform(i3,j3,0);
3146 dx3 = fDxspline;
3147 dy3 = fDyspline;
3148 z3 = fZ;
3149 z3l = fZeq;
3150 Transform(i4,j4,0);
3151 dx4 = fDxspline;
3152 dy4 = fDyspline;
3153 z4 = fZ;
3154 z4l = fZeq;
3155 zh = (Double_t)TMath::Max(z1,z2);
3156 zh = (Double_t)TMath::Max(zh,z3);
3157 zh = (Double_t)TMath::Max(zh,z4);
3158 zl = (Double_t)TMath::Min(z1l,z2l);
3159 zl = (Double_t)TMath::Min(zl,z3l);
3160 zl = (Double_t)TMath::Min(zl,z4l);
3161 i1 = (Int_t)(zl/dcount_reg+1);
3162 if (z1!=z2||z2!=z3||z3!=z4) {
3163 do {
3164 fZ = i1*dcount_reg;
3165 switch (fZscale) {
3166 case kZScaleLog:
3167 if (fZ>=1.0) fZ = log(fZ);
3168 else fZ = 0;
3169 break;
3170 case kZScaleSqrt:
3171 if (fZ>0) fZ = sqrt(fZ);
3172 else fZ = 0;
3173 break;
3174 }
3176 v = ColorCalculation(dx1,dy1,fZ,dx2,dy2,fZ,dx4,dy4,fZ);
3177 v = v*fLevels+0.5;
3178 iv = fLevels-(Int_t)v;
3179 ColorModel(iv,ui1,ui2,ui3);
3181 }
3182 if (fZ>zh) goto eqend;
3183 i1 += 1;
3184 ekv = 0;
3185 stvor = 0;
3186 if ((z2<=fZ&&fZ<z1)||(z2<fZ&&fZ<=z1)) {
3187 xb = (fZ-z2)*(dx1-dx2)/(z1-z2)+dx2;
3188 goto ekvi1;
3189 }
3190 if ((z1<=fZ&&fZ<z2)||(z1<fZ&&fZ<=z2)) {
3191 xb = (fZ-z1)*(dx2-dx1)/(z2-z1)+dx1;
3192 goto ekvi1;
3193 }
3194 if (z2==fZ&&fZ==z1) {
3195 xb = dx2;
3196ekvi1:
3197 yb = dy2;
3198 ekv = 1;
3199 x5 = xb;
3200 y5 = yb;
3201 stvor += 1;
3202 }
3203 if ((z1<=fZ&&fZ<z4)||(z1<fZ&&fZ<=z4)) {
3204 ya = (fZ-z1)*(dy4-dy1)/(z4-z1)+dy1;
3205 goto ekvi2;
3206 }
3207 if ((z4<=fZ&&fZ<z1)||(z4<fZ&&fZ<=z1)) {
3208 ya = (fZ-z4)*(dy1-dy4)/(z1-z4)+dy4;
3209 goto ekvi2;
3210 }
3211 if (z4==fZ&&fZ==z1) {
3212 ya = dy1;
3213ekvi2:
3214 xa = dx1;
3215 if (ekv==1) {
3216 Slice(xa,ya,xb,yb,line);
3217 stvor += 1;
3218 }
3219 xb = xa;
3220 yb = ya;
3221 ekv = 1;
3222 }
3223 if ((z3<=fZ&&fZ<z4)||(z3<fZ&&fZ<=z4)) {
3224 xa = (fZ-z3)*(dx4-dx3)/(z4-z3)+dx3;
3225 goto ekvi3;
3226 }
3227 if ((z4<=fZ&&fZ<z3)||(z4<fZ&&fZ<=z3)) {
3228 xa = (fZ-z4)*(dx3-dx4)/(z3-z4)+dx4;
3229 goto ekvi3;
3230 }
3231 if (z4==fZ&&fZ==z3) {
3232 xa = dx4;
3233ekvi3:
3234 ya = dy4;
3235 if (ekv==1) {
3236 Slice(xa,ya,xb,yb,line);
3237 stvor += 1;
3238 }
3239 xb = xa;
3240 yb = ya;
3241 ekv = 1;
3242 }
3243 if ((z2<=fZ&&fZ<z3)||(z2<fZ&&fZ<=z3)) {
3244 ya = (fZ-z2)*(dy3-dy2)/(z3-z2)+dy2;
3245 goto ekvi4;
3246 }
3247 if ((z3<=fZ&&fZ<z2)||(z3<fZ&&fZ<=z2)) {
3248 ya = (fZ-z3)*(dy2-dy3)/(z2-z3)+dy3;
3249 goto ekvi4;
3250 }
3251 if (z3==fZ&&fZ==z2) {
3252 ya = dy3;
3253ekvi4:
3254 xa = dx3;
3255 if (ekv==1) {
3256 Slice(xa,ya,xb,yb,line);
3257 stvor += 1;
3258 }
3259 if (stvor==4) Slice(xa,ya,x5,y5,line);
3260 }
3261 } while (fZ<=zh);
3262eqend:
3263 CopyEnvelope(dx1,dx3,dy1,dy3);
3264 }
3265 }
3266 break;
3267 case kDisplayModeBars:
3268 case kDisplayModeBarsX:
3269 case kDisplayModeBarsY:
3270 if ((q1!=q2||smer!=0)&&flag==1) {
3271 s1 = q1+1;
3272 t1 = q2;
3273 s2 = q1;
3274 t2 = q2;
3275 s3 = q1;
3276 t3 = q2+1;
3277 s4 = q1+1;
3278 t4 = q2+1;
3279 }
3280 if (turni==1) {
3282 if (s1<=w1&&s2<=w1&&s3<=w1&&s4<=w1) {
3283 i1 = s1;
3284 i2 = s2;
3285 i3 = s3;
3286 i4 = s4;
3287 }
3288 } else {
3289 i1 = (Int_t)TMath::Min(w1,s1);
3290 i2 = (Int_t)TMath::Min(w1,s2);
3291 i3 = (Int_t)TMath::Min(w1,s3);
3292 i4 = (Int_t)TMath::Min(w1,s4);
3293 }
3294 } else {
3296 if (s1<=w1&&s2<=w1&&s3<=w1&&s4<=w1) {
3297 i1 = w1-s1;
3298 i2 = w1-s2;
3299 i3 = w1-s3;
3300 i4 = w1-s4;
3301 }
3302 } else {
3303 i1 = (Int_t)TMath::Max(w1-s1,0);
3304 i2 = (Int_t)TMath::Max(w1-s2,0);
3305 i3 = (Int_t)TMath::Max(w1-s3,0);
3306 i4 = (Int_t)TMath::Max(w1-s4,0);
3307 }
3308 }
3309 if (turnj==1) {
3311 if (t1<=w2&&t2<=w2&&t3<=w2&&t4<=w2) {
3312 j1 = t1;
3313 j2 = t2;
3314 j3 = t3;
3315 j4 = t4;
3316 }
3317 } else {
3318 j1 = (Int_t)TMath::Min(w2,t1);
3319 j2 = (Int_t)TMath::Min(w2,t2);
3320 j3 = (Int_t)TMath::Min(w2,t3);
3321 j4 = (Int_t)TMath::Min(w2,t4);
3322 }
3323 } else {
3325 if (t1<=w2&&t2<=w2&&t3<=w2&&t4<=w2) {
3326 j1 = w2-t1;
3327 j2 = w2-t2;
3328 j3 = w2-t3;
3329 j4 = w2-t4;
3330 }
3331 } else {
3332 j1 = (Int_t)TMath::Max(w2-t1,0);
3333 j2 = (Int_t)TMath::Max(w2-t2,0);
3334 j3 = (Int_t)TMath::Max(w2-t3,0);
3335 j4 = (Int_t)TMath::Max(w2-t4,0);
3336 }
3337 }
3338 Transform(i1,j1,0);
3339 x1 = fXt;
3340 dx1 = fDxspline;
3341 dy1 = fDyspline;
3342 z1 = fZ;
3343 Transform(i2,j2,0);
3344 x2 = fXt;
3345 dx2 = fDxspline;
3346 dy2 = fDyspline;
3347 z2 = fZ;
3348 Transform(i3,j3,0);
3349 x3 = fXt;
3350 dx3 = fDxspline;
3351 dy3 = fDyspline;
3352 z3 = fZ;
3353 Transform(i4,j4,0);
3354 x4 = fXt;
3355 y4 = fYt;
3356 dx4 = fDxspline;
3357 dy4 = fDyspline;
3358 z4 = fZ;
3359 Transform(i1,j1,-1);
3360 ix5 = fXt;
3361 iy5 = fYt;
3362 Transform(i2,j2,-1);
3363 x6 = fXt;
3364 y6 = fYt;
3365 Transform(i3,j3,-1);
3366 x7 = fXt;
3367 y7 = fYt;
3368 Transform(i4,j4,-1);
3369 y8 = fYt;
3370 y1 = iy5+(y4-y8);
3371 y2 = y6+(y4-y8);
3372 y3 = y7+(y4-y8);
3373 if ((fDisplayMode==kDisplayModeBars)&&(q1!=q2||smer!=0)&&(flag==1)) {
3374 EnvelopeBars(ix5,iy5,x6,y6);
3375 if (fLine!=0) {
3376 if (fLine==1) {
3377 fXe = x6;
3378 fYe = y6;
3379 }
3380 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
3381 }
3382 EnvelopeBars(x6,y6,x7,y7);
3383 if (fLine!=0) {
3384 if (fLine==1) {
3385 fXe = x7;
3386 fYe = y7;
3387 }
3388 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
3389 }
3390 EnvelopeBars(ix5,iy5,x1,y1);
3391 if (fLine!=0) {
3392 if (fLine==1) {
3393 fXe = x1;
3394 fYe = y1;
3395 }
3396 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
3397 }
3398 EnvelopeBars(x6,y6,x2,y2);
3399 if (fLine!=0) {
3400 if (fLine==1) {
3401 fXe = x2;
3402 fYe = y2;
3403 }
3404 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
3405 }
3406 EnvelopeBars(x7,y7,x3,y3);
3407 if (fLine!=0) {
3408 if (fLine==1) {
3409 fXe = x3;
3410 fYe = y3;
3411 }
3412 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
3413 }
3415 v = ColorCalculation(dx1,dy1,z4,dx2,dy2,z4,dx4,dy4,z4);
3416 v = v*fLevels+0.5;
3417 iv = fLevels-(Int_t)v;
3418 uip = fNewColorIndex;
3419 ColorModel(iv,ui1,ui2,ui3);
3421 sx1 = x1;
3422 sy1 = y1;
3423 sx2 = x2;
3424 sy2 = y2;
3425 sx3 = x4;
3426 sy3 = y4;
3427 if (sx2<sx1) {
3428 sx4 = sx1;
3429 sy4 = sy1;
3430 sx1 = sx2;
3431 sy1 = sy2;
3432 sx2 = sx4;
3433 sy2 = sy4;
3434 }
3435 if (sx3<sx1) {
3436 sx4 = sx1;
3437 sy4 = sy1;
3438 sx1 = sx3;
3439 sy1 = sy3;
3440 sx3 = sx4;
3441 sy3 = sy4;
3442 }
3443 if (sy2<sy3) {
3444 sx4 = sx2;
3445 sy4 = sy2;
3446 sx2 = sx3;
3447 sy2 = sy3;
3448 sx3 = sx4;
3449 sy3 = sy4;
3450 }
3451 sdx1 = 0;
3452 sdx2 = 0;
3453 sdx3 = 0;
3454 pom1 = sy2-sy1;
3455 pom2 = sx2-sx1;
3456 if (pom2!=0) sdx1 = pom1/pom2;
3457 pom1 = sy1;
3458 pom2 = sx1;
3459 sdy1 = pom1-sdx1*pom2;
3460 pom1 = sy3-sy1;
3461 pom2 = sx3-sx1;
3462 if (pom2!=0) sdx2 = pom1/pom2;
3463 pom1 = sy1;
3464 pom2 = sx1;
3465 sdy2 = pom1-sdx2*pom2;
3466 pom1 = sy3-sy2;
3467 pom2 = sx3-sx2;
3468 if (pom2!=0) sdx3 = pom1/pom2;
3469 pom1 = sy2;
3470 pom2 = sx2;
3471 sdy3 = pom1-sdx3*pom2;
3472 if (sx2<sx3) {
3473 if (sx1!=sx2) {
3474 for (sx4=sx1;sx4<=sx2;sx4++) {
3475 pom1 = sx4;
3476 sdy4 = sdx1*pom1+sdy1;
3477 sy4 = (Int_t)sdy4;
3478 if (sx3!=sx1) {
3479 sdy4 = sdx2*pom1+sdy2;
3480 sy5 = (Int_t)sdy4;
3481 y5 = fEnvelope[sx4];
3482 if (sy4<sy5) {
3483 pom1 = sy4;
3484 sy4 = sy5;
3485 sy5 = (Int_t)pom1;
3486 }
3487 if ((sy4<=y5)||(sy5<y5)) {
3488 sy4 = (Int_t)TMath::Min(sy4,(Int_t)y5);
3489 sy5 = (Int_t)TMath::Min(sy5,(Int_t)y5);
3490 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx4),gPad->PixeltoY(sy5)+1);
3491 }
3492 }
3493 }
3494 }
3495 if (sx2!=sx3) {
3496 for (sx4=sx2;sx4<=sx3;sx4++) {
3497 pom1 = sx4;
3498 sdy4 = sdx3*pom1+sdy3;
3499 sy4 = (Int_t)sdy4;
3500 if (sx3!=sx1) {
3501 sdy4 = sdx2*pom1+sdy2;
3502 sy5 = (Int_t)sdy4;
3503 y5 = fEnvelope[sx4];
3504 if (sy4<sy5) {
3505 pom1 = sy4;
3506 sy4 = sy5;
3507 sy5 = (Int_t)pom1;
3508 }
3509 if ((sy4<=y5)||(sy5<y5)) {
3510 sy4 = (Int_t)TMath::Min(sy4,(Int_t)y5);
3511 sy5 = (Int_t)TMath::Min(sy5,(Int_t)y5);
3512 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx4),gPad->PixeltoY(sy5)+1);
3513 }
3514 }
3515 }
3516 }
3517 } else {
3518 if (sx3!=sx1) {
3519 for (sx4=sx1;sx4<=sx3;sx4++) {
3520 pom1 = sx4;
3521 sdy4 = sdx2*pom1+sdy2;
3522 sy4 = (Int_t)sdy4;
3523 if (sx2!=sx1) {
3524 sdy4 = sdx1*pom1+sdy1;
3525 sy5 = (Int_t)sdy4;
3526 y5 = fEnvelope[sx4];
3527 if (sy4<sy5) {
3528 pom1 = sy4;
3529 sy4 = sy5;
3530 sy5 = (Int_t)pom1;
3531 }
3532 if ((sy4<=y5)||(sy5<y5)) {
3533 sy4 = (Int_t)TMath::Min(sy4,(Int_t)y5);
3534 sy5 = (Int_t)TMath::Min(sy5,(Int_t)y5);
3535 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx4),gPad->PixeltoY(sy5)+1);
3536 }
3537 }
3538 }
3539 }
3540 if (sx2!=sx3) {
3541 for (sx4=sx3;sx4<=sx2;sx4++) {
3542 pom1 = sx4;
3543 sdy4 = sdx3*pom1+sdy3;
3544 sy4 = (Int_t)sdy4;
3545 if (sx2!=sx1) {
3546 sdy4 = sdx1*pom1+sdy1;
3547 sy5 = (Int_t)sdy4;
3548 y5 = fEnvelope[sx4];
3549 if (sy4<sy5) {
3550 pom1 = sy4;
3551 sy4 = sy5;
3552 sy5 = (Int_t)pom1;
3553 }
3554 if ((sy4<=y5)||(sy5<y5)) {
3555 sy4 = (Int_t)TMath::Min(sy4,(Int_t)y5);
3556 sy5 = (Int_t)TMath::Min(sy5,(Int_t)y5);
3557 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx4),gPad->PixeltoY(sy5)+1);
3558 }
3559 }
3560 }
3561 }
3562 }
3563 sx1 = x2;
3564 sy1 = y2;
3565 sx2 = x3;
3566 sy2 = y3;
3567 sx3 = x4;
3568 sy3 = y4;
3569 if (sx2<sx1) {
3570 sx4 = sx1;
3571 sy4 = sy1;
3572 sx1 = sx2;
3573 sy1 = sy2;
3574 sx2 = sx4;
3575 sy2 = sy4;
3576 }
3577 if (sx3<sx1) {
3578 sx4 = sx1;
3579 sy4 = sy1;
3580 sx1 = sx3;
3581 sy1 = sy3;
3582 sx3 = sx4;
3583 sy3 = sy4;
3584 }
3585 if (sy2<sy3) {
3586 sx4 = sx2;
3587 sy4 = sy2;
3588 sx2 = sx3;
3589 sy2 = sy3;
3590 sx3 = sx4;
3591 sy3 = sy4;
3592 }
3593 sdx1 = 0;
3594 sdx2 = 0;
3595 sdx3 = 0;
3596 pom1 = sy2-sy1;
3597 pom2 = sx2-sx1;
3598 if (pom2!=0) sdx1 = pom1/pom2;
3599 pom1 = sy1;
3600 pom2 = sx1;
3601 sdy1 = pom1-sdx1*pom2;
3602 pom1 = sy3-sy1;
3603 pom2 = sx3-sx1;
3604 if (pom2!=0) sdx2 = pom1/pom2;
3605 pom1 = sy1;
3606 pom2 = sx1;
3607 sdy2 = pom1-sdx2*pom2;
3608 pom1 = sy3-sy2;
3609 pom2 = sx3-sx2;
3610 if (pom2!=0) sdx3 = pom1/pom2;
3611 pom1 = sy2;
3612 pom2 = sx2;
3613 sdy3 = pom1-sdx3*pom2;
3614 if (sx2<sx3) {
3615 if (sx1!=sx2) {
3616 for (sx4=sx1;sx4<=sx2;sx4++) {
3617 pom1 = sx4;
3618 sdy4 = sdx1*pom1+sdy1;
3619 sy4 = (Int_t)sdy4;
3620 if (sx3!=sx1) {
3621 sdy4 = sdx2*pom1+sdy2;
3622 sy5 = (Int_t)sdy4;
3623 y5 = fEnvelope[sx4];
3624 if (sy4<sy5) {
3625 pom1 = sy4;
3626 sy4 = sy5;
3627 sy5 = (Int_t)pom1;
3628 }
3629 if ((sy4<=y5)||(sy5<y5)) {
3630 sy4 = (Int_t)TMath::Min(sy4,(Int_t)y5);
3631 sy5 = (Int_t)TMath::Min(sy5,(Int_t)y5);
3632 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx4),gPad->PixeltoY(sy5)+1);
3633 }
3634 }
3635 }
3636 }
3637 if (sx2!=sx3) {
3638 for (sx4=sx2;sx4<=sx3;sx4++) {
3639 pom1 = sx4;
3640 sdy4 = sdx3*pom1+sdy3;
3641 sy4 = (Int_t)sdy4;
3642 if (sx3!=sx1) {
3643 sdy4 = sdx2*pom1+sdy2;
3644 sy5 = (Int_t)sdy4;
3645 y5 = fEnvelope[sx4];
3646 if (sy4<sy5) {
3647 pom1 = sy4;
3648 sy4 = sy5;
3649 sy5 = (Int_t)pom1;
3650 }
3651 if ((sy4<=y5)||(sy5<y5)) {
3652 sy4 = (Int_t)TMath::Min(sy4,(Int_t)y5);
3653 sy5 = (Int_t)TMath::Min(sy5,(Int_t)y5);
3654 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx4),gPad->PixeltoY(sy5)+1);
3655 }
3656 }
3657 }
3658 }
3659 } else {
3660 if (sx3!=sx1) {
3661 for (sx4=sx1;sx4<=sx3;sx4++) {
3662 pom1 = sx4;
3663 sdy4 = sdx2*pom1+sdy2;
3664 sy4 = (Int_t)sdy4;
3665 if (sx2!=sx1) {
3666 sdy4 = sdx1*pom1+sdy1;
3667 sy5 = (Int_t)sdy4;
3668 y5 = fEnvelope[sx4];
3669 if (sy4<sy5) {
3670 pom1 = sy4;
3671 sy4 = sy5;
3672 sy5 = (Int_t)pom1;
3673 }
3674 if ((sy4<=y5)||(sy5<y5)) {
3675 sy4 = (Int_t)TMath::Min(sy4,(Int_t)y5);
3676 sy5 = (Int_t)TMath::Min(sy5,(Int_t)y5);
3677 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx4),gPad->PixeltoY(sy5)+1);
3678 }
3679 }
3680 }
3681 }
3682 if (sx2!=sx3) {
3683 for (sx4=sx3;sx4<=sx2;sx4++) {
3684 pom1 = sx4;
3685 sdy4 = sdx3*pom1+sdy3;
3686 sy4 = (Int_t)sdy4;
3687 if (sx2!=sx1) {
3688 sdy4 = sdx1*pom1+sdy1;
3689 sy5 = (Int_t)sdy4;
3690 y5 = fEnvelope[sx4];
3691 if (sy4<sy5) {
3692 pom1 = sy4;
3693 sy4 = sy5;
3694 sy5 = (Int_t)pom1;
3695 }
3696 if ((sy4<=y5)||(sy5<y5)) {
3697 sy4 = (Int_t)TMath::Min(sy4,(Int_t)y5);
3698 sy5 = (Int_t)TMath::Min(sy5,(Int_t)y5);
3699 line->PaintLine(gPad->PixeltoX(sx4),gPad->PixeltoY(sy4)+1,gPad->PixeltoX(sx4),gPad->PixeltoY(sy5)+1);
3700 }
3701 }
3702 }
3703 }
3704 }
3705 line->SetLineColor(uip);
3706 }
3708 if (fLine!=0) {
3709 if (fLine==1) {
3710 fXe = x2;
3711 fYe = y2;
3712 }
3713 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
3714 }
3715 EnvelopeBars(x2,y2,x3,y3);
3716 if (fLine!=0) {
3717 if (fLine==1) {
3718 fXe = x3;
3719 fYe = y3;
3720 }
3721 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
3722 }
3723 EnvelopeBars(x1,y1,x4,y4);
3724 if (fLine!=0) {
3725 if (fLine==1) {
3726 fXe = x4;
3727 fYe = y4;
3728 }
3729 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
3730 }
3731 EnvelopeBars(x4,y4,x3,y3);
3732 if (fLine!=0) {
3733 if (fLine==1) {
3734 fXe = x3;
3735 fYe = y3;
3736 }
3737 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
3738 }
3739 } else if ((fDisplayMode==kDisplayModeBarsY)&&(((flag!=0)&&(smer==0))||((flag==0)&&(smer!=0)))) {
3740 EnvelopeBars(ix5,iy5,x6,y6);
3741 if (fLine!=0) {
3742 if (fLine==1) {
3743 fXe = x6;
3744 fYe = y6;
3745 }
3746 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
3747 }
3748 EnvelopeBars(x1,y1,ix5,iy5);
3749 if (fLine!=0) {
3750 if (fLine==1) {
3751 fXe = ix5;
3752 fYe = iy5;
3753 }
3754 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
3755 }
3756 EnvelopeBars(x2,y2,x6,y6);
3757 if (fLine!=0) {
3758 if (fLine==1) {
3759 fXe = x6;
3760 fYe = y6;
3761 }
3762 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
3763 }
3765 v = ColorCalculation(dx1,dy1,z4,dx2,dy2,z4,dx4,dy4,z4);
3766 v = v*fLevels+0.5;
3767 iv = fLevels-(Int_t)v;
3768 uip = fNewColorIndex;
3769 ColorModel(iv,ui1,ui2,ui3);
3771 }
3773 if (fLine!=0) {
3774 if (fLine==1) {
3775 fXe = x2;
3776 fYe = y2;
3777 }
3778 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
3779 }
3781 line->SetLineColor(uip);
3782 }
3783 } else if ((fDisplayMode==kDisplayModeBarsX)&&(((flag==0)&&(smer==0))||((flag!=0)&&(smer!=0)))) {
3784 EnvelopeBars(x7,y7,x6,y6);
3785 if (fLine!=0) {
3786 if (fLine==1) {
3787 fXe = x6;
3788 fYe = y6;
3789 }
3790 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
3791 }
3792 EnvelopeBars(x2,y2,x6,y6);
3793 if (fLine!=0) {
3794 if (fLine==1) {
3795 fXe = x6;
3796 fYe = y6;
3797 }
3798 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
3799 }
3800 EnvelopeBars(x3,y3,x7,y7);
3801 if (fLine!=0) {
3802 if (fLine==1) {
3803 fXe = x7;
3804 fYe = y7;
3805 }
3806 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
3807 }
3809 v = ColorCalculation(dx1,dy1,z4,dx2,dy2,z4,dx4,dy4,z4);
3810 v = v*fLevels+0.5;
3811 iv = fLevels-(Int_t)v;
3812 uip = fNewColorIndex;
3813 ColorModel(iv,ui1,ui2,ui3);
3815 }
3816 EnvelopeBars(x3,y3,x2,y2);
3817 if (fLine!=0) {
3818 if (fLine==1) {
3819 fXe = x2;
3820 fYe = y2;
3821 }
3822 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
3823 }
3825 line->SetLineColor(uip);
3826 }
3827 }
3828 break;
3829 case kDisplayModeLinesX:
3831 if (((flag==0)&&(smer==0))||((flag!=0)&&(smer!=0))) {
3833 Envelope(x1,y1,x2,y2);
3834 if (fLine!=0) {
3835 if (fLine==1) {
3836 fXe = x2;
3837 fYe = y2;
3838 }
3839 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
3840 }
3841 } else {
3842 s1 = q1;
3843 t1 = (Int_t)TMath::Max(q2-1,0);
3844 s2 = q1;
3845 t2 = (Int_t)TMath::Min(q2+2,w2);
3846 if (turni==1) {
3847 i1 = s1;
3848 i2 = s2;
3849 } else {
3850 i1 = (Int_t)TMath::Max(w1-s1,0);
3851 i2 = (Int_t)TMath::Max(w1-s2,0);
3852 }
3853 if (turnj==1) {
3854 j1 = t1;
3855 j2 = t2;
3856 } else {
3857 j1 = (Int_t)TMath::Max(w2-t1,0);
3858 j2 = (Int_t)TMath::Max(w2-t2,0);
3859 }
3860 Transform(i1,j1,0);
3861 x3 = fXt;
3862 y3 = fYt;
3863 Transform(i2,j2,0);
3864 x4 = fXt;
3865 y4 = fYt;
3866 bezx1 = x1+(x2-x1)/3;
3867 bezx2 = x1+2*(x2-x1)/3;
3868 bezy1 = y1+(y2-y3)/6;
3869 bezy2 = y2-(y4-y1)/6;
3870 if (x1<=x2) {
3871 if (bezx1<=x1) {
3872 bezx1 = x1;
3873 bezy1 = y1;
3874 }
3875 if (bezx1>=x2) {
3876 bezx1 = x2;
3877 bezy1 = y2;
3878 }
3879 if (bezx2<=x1) {
3880 bezx2 = x1;
3881 bezy2 = y1;
3882 }
3883 if (bezx2>=x2) {
3884 bezx2 = x2;
3885 bezy2 = y2;
3886 }
3887 fBzX[0] = x1;
3888 fBzY[0] = y1;
3889 fBzX[1] = (Int_t)bezx1;
3890 fBzY[1] = (Int_t)bezy1;
3891 fBzX[2] = (Int_t)bezx2;
3892 fBzY[2] = (Int_t)bezy2;
3893 fBzX[3] = x2;
3894 fBzY[3] = y2;
3895 for (bezf=0;bezf<1.01;bezf+=0.1) {
3896 BezierSmoothing(bezf);
3897 if (bezf==0) {
3898 ibezx1 = (Int_t)(fGbezx+0.5);
3899 ibezy1 = (Int_t)(fGbezy+0.5);
3900 } else {
3901 ibezx2 = ibezx1;
3902 ibezy2 = ibezy1;
3903 ibezx1 = (Int_t)(fGbezx+0.5);
3904 ibezy1 = (Int_t)(fGbezy+0.5);
3905 Envelope(ibezx2,ibezy2,ibezx1,ibezy1);
3906 if (fLine!=0) {
3907 if (fLine==1) {
3908 fXe = ibezx1;
3909 fYe = ibezy1;
3910 }
3911 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
3912 }
3913 }
3914 }
3915 } else if (x1>x2) {
3916 if (bezx1>=x1) {
3917 bezx1 = x1;
3918 bezy1 = y1;
3919 }
3920 if (bezx1<=x2) {
3921 bezx1 = x2;
3922 bezy1 = y2;
3923 }
3924 if (bezx2>=x1) {
3925 bezx2 = x1;
3926 bezy2 = y1;
3927 }
3928 if (bezx2<=x2) {
3929 bezx2 = x2;
3930 bezy2 = y2;
3931 }
3932 fBzX[0] = x1;
3933 fBzY[0] = y1;
3934 fBzX[1] = (Int_t)bezx1;
3935 fBzY[1] = (Int_t)bezy1;
3936 fBzX[2] = (Int_t)bezx2;
3937 fBzY[2] = (Int_t)bezy2;
3938 fBzX[3] = x2;
3939 fBzY[3] = y2;
3940 for (bezf=0;bezf<1.01;bezf+=0.1) {
3941 BezierSmoothing(bezf);
3942 if (bezf==0) {
3943 ibezx1 = (Int_t)(fGbezx+0.5);
3944 ibezy1 = (Int_t)(fGbezy+0.5);
3945 } else {
3946 ibezx2 = ibezx1;
3947 ibezy2 = ibezy1;
3948 ibezx1 = (Int_t)(fGbezx+0.5);
3949 ibezy1 = (Int_t)(fGbezy+0.5);
3950 Envelope(ibezx1,ibezy1,ibezx2,ibezy2);
3951 if (fLine!=0) {
3952 if (fLine==1) {
3953 fXe = ibezx2;
3954 fYe = ibezy2;
3955 }
3956 line->PaintLine(gPad->PixeltoX(fXs),gPad->PixeltoY(fYs)+1,gPad->PixeltoX(fXe),gPad->PixeltoY(fYe)+1);
3957 }
3958 }
3959 }
3960 }
3961 }
3962 }
3963 } else {
3964 if ((q1!=q2||smer!=0)&&flag==1) {
3965 s1 = q1+1;
3966 t1 = q2;
3967 s2 = q1;
3968 t2 = q2;
3969 s3 = q1;
3970 t3 = q2+1;
3971 s4 = q1+1;
3972 t4 = q2+1;
3973 if (fShading==kShaded) {
3974 sr1 = s1;
3975 tr1 = (Int_t)TMath::Max(t1-1,0);
3976 sr2 = s2;
3977 tr2 = (Int_t)TMath::Max(t2-1,0);
3978 sr3 = (Int_t)TMath::Max(s2-1,0);
3979 tr3 = t2;
3980 sr4 = (Int_t)TMath::Max(s3-1,0);
3981 tr4 = t3;
3982 sr5 = s3;
3983 tr5 = t3+1;
3984 sr6 = s4;
3985 tr6 = t4+1;
3986 sr7 = s4+1;
3987 tr7 = t4;
3988 sr8 = s1+1;
3989 tr8 = t1;
3990 }
3991 if (turni==1) {
3992 i1 = s1;
3993 i2 = s2;
3994 i3 = s3;
3995 i4 = s4;
3996 } else {
3997 i1 = (Int_t)TMath::Max(w1-s1,0);
3998 i2 = (Int_t)TMath::Max(w1-s2,0);
3999 i3 = (Int_t)TMath::Max(w1-s3,0);
4000 i4 = (Int_t)TMath::Max(w1-s4,0);
4001 if (fShading==kShaded) {
4002 sr1 = (Int_t)TMath::Max(w1-sr1,0);
4003 sr2 = (Int_t)TMath::Max(w1-sr2,0);
4004 sr3 = (Int_t)TMath::Max(w1-sr3,0);
4005 sr4 = (Int_t)TMath::Max(w1-sr4,0);
4006 sr5 = (Int_t)TMath::Max(w1-sr5,0);
4007 sr6 = (Int_t)TMath::Max(w1-sr6,0);
4008 sr7 = (Int_t)TMath::Max(w1-sr7,0);
4009 sr8 = (Int_t)TMath::Max(w1-sr8,0);
4010 }
4011 }
4012 if (turnj==1) {
4013 j1 = t1;
4014 j2 = t2;
4015 j3 = t3;
4016 j4 = t4;
4017 } else {
4018 j1 = (Int_t)TMath::Max(w2-t1,0);
4019 j2 = (Int_t)TMath::Max(w2-t2,0);
4020 j3 = (Int_t)TMath::Max(w2-t3,0);
4021 j4 = (Int_t)TMath::Max(w2-t4,0);
4022 if (fShading==kShaded) {
4023 tr1 = (Int_t)TMath::Max(w2-tr1,0);
4024 tr2 = (Int_t)TMath::Max(w2-tr2,0);
4025 tr3 = (Int_t)TMath::Max(w2-tr3,0);
4026 tr4 = (Int_t)TMath::Max(w2-tr4,0);
4027 tr5 = (Int_t)TMath::Max(w2-tr5,0);
4028 tr6 = (Int_t)TMath::Max(w2-tr6,0);
4029 tr7 = (Int_t)TMath::Max(w2-tr7,0);
4030 tr8 = (Int_t)TMath::Max(w2-tr8,0);
4031 }
4032 }
4033 Transform(i1,j1,0);
4034 x1 = fXt;
4035 y1 = fYt;
4036 dx1 = fDxspline;
4037 dy1 = fDyspline;
4038 z1 = fZ;
4039 Transform(i2,j2,0);
4040 x2 = fXt;
4041 y2 = fYt;
4042 dx2 = fDxspline;
4043 dy2 = fDyspline;
4044 z2 = fZ;
4045 Transform(i3,j3,0);
4046 x3 = fXt;
4047 y3 = fYt;
4048 dx3 = fDxspline;
4049 dy3 = fDyspline;
4050 z3 = fZ;
4051 Transform(i4,j4,0);
4052 x4 = fXt;
4053 y4 = fYt;
4054 dx4 = fDxspline;
4055 dy4 = fDyspline;
4056 z4 = fZ;
4057 Envelope(x1,y1,x2,y2);
4058 Envelope(x2,y2,x3,y3);
4059 xtaz = (dx1+dx2+dx4)/3;
4060 ytaz = (dy1+dy2+dy4)/3;
4061 ztaz = (z1+z2+z4)/3;
4062 v = ColorCalculation(dx1,dy1,z1,dx2,dy2,z2,dx4,dy4,z4);
4063 if (fShading==kShaded) {
4065 if (sr1<0||sr1>w1||tr1<0||tr1>w2) Transform(sr1,tr1,-1);
4066 else Transform(sr1,tr1,0);
4067 dxr1 = fDxspline;
4068 dyr1 = fDyspline;
4069 zr1 = fZ;
4070 if (sr8<0||sr8>w1||tr8<0||tr8>w2) Transform(sr8,tr8,-1);
4071 else Transform(sr8,tr8,0);
4072 dxr2 = fDxspline;
4073 dyr2 = fDyspline;
4074 zr2 = fZ;
4075 v = v+ColorCalculation(dxr1,dyr1,zr1,dx2,dy2,z2,dx1,dy1,z1);
4076 v = v+ColorCalculation(dxr1,dyr1,zr1,dx1,dy1,z1,dxr2,dyr2,zr2);
4077 v = v+ColorCalculation(dxr2,dyr2,zr2,dx1,dy1,z1,dx4,dy4,z4);
4078 v1 = v/4;
4079 if (sr3<0||sr3>w1||tr3<0||tr3>w2) Transform(sr3,tr3,-1);
4080 else Transform(sr3,tr3,0);
4081 dxr1 = fDxspline;
4082 dyr1 = fDyspline;
4083 zr1 = fZ;
4084 if (sr2<0||sr2>w1||tr2<0||tr2>w2) Transform(sr2,tr2,-1);
4085 else Transform(sr2,tr2,0);
4086 dxr2 = fDxspline;
4087 dyr2 = fDyspline;
4088 zr2 = fZ;
4089 v = ColorCalculation(dx1,dy1,z1,dx2,dy2,z2,dx3,dy3,z3);
4090 v = v+ColorCalculation(dx2,dy2,z2,dxr1,dyr1,zr1,dx3,dy3,z3);
4091 v = v+ColorCalculation(dx2,dy2,z2,dxr2,dyr2,zr2,dxr1,dyr1,zr1);
4092 v = v+ColorCalculation(dxr2,dyr2,zr2,dx2,dy2,z2,dx1,dy1,z1);
4093 v2 = v/4;
4094 if (sr5<0||sr5>w1||tr5<0||tr5>w2) Transform(sr5,tr5,-1);
4095 else Transform(sr5,tr5,0);
4096 dxr1 = fDxspline;
4097 dyr1 = fDyspline;
4098 zr1 = fZ;
4099 if (sr4<0||sr4>w1||tr4<0||tr4>w2) Transform(sr4,tr4,-1);
4100 else Transform(sr4,tr4,0);
4101 dxr2 = fDxspline;
4102 dyr2 = fDyspline;
4103 zr2 = fZ;
4104 v = ColorCalculation(dx2,dy2,z2,dx3,dy3,z3,dx4,dy4,z4);
4105 v = v+ColorCalculation(dx4,dy4,z4,dx3,dy3,z3,dxr1,dyr1,zr1);
4106 v = v+ColorCalculation(dx3,dy3,z3,dxr2,dyr2,zr2,dxr1,dyr1,zr1);
4107 v = v+ColorCalculation(dx2,dy2,z2,dxr2,dyr2,zr2,dx3,dy3,z3);
4108 v3 = v/4;
4109 if (sr7<0||sr7>w1||tr7<0||tr7>w2) Transform(sr7,tr7,-1);
4110 else Transform(sr7,tr7,0);
4111 dxr1 = fDxspline;
4112 dyr1 = fDyspline;
4113 zr1 = fZ;
4114 if (sr6<0||sr6>w1||tr6<0||tr6>w2) Transform(sr6,tr6,-1);
4115 else Transform(sr6,tr6,0);
4116 dxr2 = fDxspline;
4117 dyr2 = fDyspline;
4118 zr2 = fZ;
4119 v = ColorCalculation(dx1,dy1,z1,dx3,dy3,z3,dx4,dy4,z4);
4120 v = v+ColorCalculation(dx4,dy4,z4,dx3,dy3,z3,dxr2,dyr2,zr2);
4121 v = v+ColorCalculation(dx4,dy4,z4,dxr2,dyr2,zr2,dxr1,dyr1,zr1);
4122 v = v+ColorCalculation(dx1,dy1,z1,dx4,dy4,z4,dxr1,dyr1,zr1);
4123 v4 = v/4;
4124 } else {
4125 spriz = ShadowColorCalculation(xtaz,ytaz,ztaz,shad_noise);
4126 v = v+spriz;
4127 v = v/2;
4128 if (sr1<0||sr1>w1||tr1<0||tr1>w2) Transform(sr1,tr1,-1);
4129 else Transform(sr1,tr1,0);
4130 dxr1 = fDxspline;
4131 dyr1 = fDyspline;
4132 zr1 = fZ;
4133 if (sr8<0||sr8>w1||tr8<0||tr8>w2) Transform(sr8,tr8,-1);
4134 else Transform(sr8,tr8,0);
4135 dxr2 = fDxspline;
4136 dyr2 = fDyspline;
4137 zr2 = fZ;
4138 da = (dxr1+dx2+dx1)/3;
4139 db = (dyr1+dy2+dy1)/3;
4140 dc = (zr1+z2+z1)/3;
4141 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
4142 v = v+(ColorCalculation(dxr1,dyr1,zr1,dx2,dy2,z2,dx1,dy1,z1)+spriz)/2;
4143 da = (dxr1+dxr2+dx1)/3;
4144 db = (dyr1+dyr2+dy1)/3;
4145 dc = (zr1+zr2+z1)/3;
4146 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
4147 v = v+(ColorCalculation(dxr1,dyr1,zr1,dx1,dy1,z1,dxr2,dyr2,zr2)+spriz)/2;
4148 da = (dxr2+dx1+dx4)/3;
4149 db = (dyr2+dy1+dy4)/3;
4150 dc = (zr2+z1+z4)/3;
4151 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
4152 v = v+(ColorCalculation(dxr2,dyr2,zr2,dx1,dy1,z1,dx4,dy4,z4)+spriz)/2;
4153 v1 = v/4;
4154 if (sr3<0||sr3>w1||tr3<0||tr3>w2) Transform(sr3,tr3,-1);
4155 else Transform(sr3,tr3,0);
4156 dxr1 = fDxspline;
4157 dyr1 = fDyspline;
4158 zr1 = fZ;
4159 if (sr2<0||sr2>w1||tr2<0||tr2>w2) Transform(sr2,tr2,-1);
4160 else Transform(sr2,tr2,0);
4161 dxr2 = fDxspline;
4162 dyr2 = fDyspline;
4163 zr2 = fZ;
4164 da = (dx1+dx2+dx3)/3;
4165 db = (dy1+dy2+dy3)/3;
4166 dc = (z1+z2+z3)/3;
4167 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
4168 v = (ColorCalculation(dx1,dy1,z1,dx2,dy2,z2,dx3,dy3,z3)+spriz)/2;
4169 da = (dx2+dxr1+dx3)/3;
4170 db = (dy2+dyr1+dy3)/3;
4171 dc = (z2+zr1+z3)/3;
4172 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
4173 v = v+(ColorCalculation(dx2,dy2,z2,dxr1,dyr1,zr1,dx3,dy3,z3)+spriz)/2;
4174 da = (dx2+dxr2+dxr1)/3;
4175 db = (dy2+dyr2+dyr1)/3;
4176 dc = (z2+zr2+zr1)/3;
4177 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
4178 v = v+(ColorCalculation(dx2,dy2,z2,dxr2,dyr2,zr2,dxr1,dyr1,zr1)+spriz)/2;
4179 da = (dxr2+dx2+dx1)/3;
4180 db = (dyr2+dy2+dy1)/3;
4181 dc = (zr2+z2+z1)/3;
4182 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
4183 v = v+(ColorCalculation(dxr2,dyr2,zr2,dx2,dy2,z2,dx1,dy1,z1)+spriz)/2;
4184 v2 = v/4;
4185 if (sr5<0||sr5>w1||tr5<0||tr5>w2) Transform(sr5,tr5,-1);
4186 else Transform(sr5,tr5,0);
4187 dxr1 = fDxspline;
4188 dyr1 = fDyspline;
4189 zr1 = fZ;
4190 if (sr4<0||sr4>w1||tr4<0||tr4>w2) Transform(sr4,tr4,-1);
4191 else Transform(sr4,tr4,0);
4192 dxr2 = fDxspline;
4193 dyr2 = fDyspline;
4194 zr2 = fZ;
4195 da = (dx2+dx3+dx4)/3;
4196 db = (dy2+dy3+dy4)/3;
4197 dc = (z2+z3+z4)/3;
4198 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
4199 v = (ColorCalculation(dx2,dy2,z2,dx3,dy3,z3,dx4,dy4,z4)+spriz)/2;
4200 da = (dx4+dx3+dxr1)/3;
4201 db = (dy4+dy3+dyr1)/3;
4202 dc = (z4+z3+zr1)/3;
4203 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
4204 v = v+(ColorCalculation(dx4,dy4,z4,dx3,dy3,z3,dxr1,dyr1,zr1)+spriz)/2;
4205 da = (dx3+dxr2+dxr1)/3;
4206 db = (dy3+dyr2+dyr1)/3;
4207 dc = (z3+zr2+zr1)/3;
4208 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
4209 v = v+(ColorCalculation(dx3,dy3,z3,dxr2,dyr2,zr2,dxr1,dyr1,zr1)+spriz)/2;
4210 da = (dx2+dxr2+dx3)/3;
4211 db = (dy2+dyr2+dy3)/3;
4212 dc = (z2+zr2+z3)/3;
4213 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
4214 v = v+(ColorCalculation(dx2,dy2,z2,dxr2,dyr2,zr2,dx3,dy3,z3)+spriz)/2;
4215 v3 = v/4;
4216 if (sr7<0||sr7>w1||tr7<0||tr7>w2) Transform(sr7,tr7,-1);
4217 else Transform(sr7,tr7,0);
4218 dxr1 = fDxspline;
4219 dyr1 = fDyspline;
4220 zr1 = fZ;
4221 if (sr6<0||sr6>w1||tr6<0||tr6>w2) Transform(sr6,tr6,-1);
4222 else Transform(sr6,tr6,0);
4223 dxr2 = fDxspline;
4224 dyr2 = fDyspline;
4225 zr2 = fZ;
4226 da = (dx1+dx3+dx4)/3;
4227 db = (dy1+dy3+dy4)/3;
4228 dc = (z1+z3+z4)/3;
4229 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
4230 v = (ColorCalculation(dx1,dy1,z1,dx3,dy3,z3,dx4,dy4,z4)+spriz)/2;
4231 da = (dx4+dx3+dxr2)/3;
4232 db = (dy4+dy3+dyr2)/3;
4233 dc = (z4+z3+zr2)/3;
4234 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
4235 v = v+(ColorCalculation(dx4,dy4,z4,dx3,dy3,z3,dxr2,dyr2,zr2)+spriz)/2;
4236 da = (dx4+dxr2+dxr1)/3;
4237 db = (dy4+dyr2+dyr1)/3;
4238 dc = (z4+zr2+zr1)/3;
4239 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
4240 v = v+(ColorCalculation(dx4,dy4,z4,dxr2,dyr2,zr2,dxr1,dyr1,zr1)+spriz)/2;
4241 da = (dx1+dx4+dxr1)/3;
4242 db = (dy1+dy4+dyr1)/3;
4243 dc = (z1+z4+zr1)/3;
4244 spriz = ShadowColorCalculation(da,db,dc,shad_noise);
4245 v = v+(ColorCalculation(dx1,dy1,z1,dx4,dy4,z4,dxr1,dyr1,zr1)+spriz)/2;
4246