Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TGraph2DPainter.cxx
Go to the documentation of this file.
1// @(#)root/histpainter:$Id: TGraph2DPainter.cxx,v 1.00
2// Author: Olivier Couet
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#include "TGraph2DPainter.h"
13#include "TMath.h"
14#include "TList.h"
15#include "TGraph.h"
16#include "TGraph2D.h"
17#include "TGraphDelaunay.h"
18#include "TGraphDelaunay2D.h"
19#include "TPolyLine.h"
20#include "TPolyMarker.h"
21#include "TVirtualPad.h"
22#include "TView.h"
23#include "THLimitsFinder.h"
24#include "TStyle.h"
25#include "Hoption.h"
26#include "TH1.h"
27#include <vector>
28
31
32
33
34////////////////////////////////////////////////////////////////////////////////
35/*! \class TGraph2DPainter
36 \ingroup Histpainter
37 \brief The TGraphDelaunay painting class.
38
39TGraph2DPainter paints a TGraphDelaunay using triangles or clouds of points.
40
41See documentation of TGraph2D, TGraph2DErrors and TGraph2DAsymmErrors to get the list of
42drawing options for these classes.
43
44*/
45
46
47////////////////////////////////////////////////////////////////////////////////
48/// TGraph2DPainter default constructor
49
51{
52 fX = nullptr;
53 fY = nullptr;
54 fZ = nullptr;
55 fEXlow = nullptr;
56 fEXhigh = nullptr;
57 fEYlow = nullptr;
58 fEYhigh = nullptr;
59 fEZlow = nullptr;
60 fEZhigh = nullptr;
61 fXN = nullptr;
62 fYN = nullptr;
63 fPTried = nullptr;
64 fNTried = nullptr;
65 fMTried = nullptr;
66 fGraph2D = nullptr;
67 fDelaunay = nullptr;
68 fDelaunay2D = nullptr;
69 fXmin = 0.;
70 fXmax = 0.;
71 fYmin = 0.;
72 fYmax = 0.;
73 fZmin = 0.;
74 fZmax = 0.;
75 fXNmin = 0;
76 fXNmax = 0;
77 fYNmin = 0;
78 fYNmax = 0;
79 fNdt = 0;
80 fNpoints = 0;
81}
82
83
84////////////////////////////////////////////////////////////////////////////////
85/// TGraph2DPainter constructor using the old Delaunay algorithm
86
88{
89 fDelaunay = gd;
90 fDelaunay2D = nullptr;
93 fNdt = 0;
94 fXN = nullptr;
95 fYN = nullptr;
96 fXNmin = 0;
97 fXNmax = 0;
98 fYNmin = 0;
99 fYNmax = 0;
100 fPTried = nullptr;
101 fNTried = nullptr;
102 fMTried = nullptr;
103 fXmin = 0.;
104 fXmax = 0.;
105 fYmin = 0.;
106 fYmax = 0.;
107 fZmin = 0.;
108 fZmax = 0.;
109}
110
111////////////////////////////////////////////////////////////////////////////////
112/// TGraph2DPainter constructor using the new Delaunay algorithm
113
115{
116 fDelaunay = nullptr;
117 fDelaunay2D = gd;
120 fNdt = 0;
121 fXN = nullptr;
122 fYN = nullptr;
123 fXNmin = 0;
124 fXNmax = 0;
125 fYNmin = 0;
126 fYNmax = 0;
127 fPTried = nullptr;
128 fNTried = nullptr;
129 fMTried = nullptr;
130 fXmin = 0.;
131 fXmax = 0.;
132 fYmin = 0.;
133 fYmax = 0.;
134 fZmin = 0.;
135 fZmax = 0.;
136}
137
138
139////////////////////////////////////////////////////////////////////////////////
140/// TGraph2DPainter destructor.
141
145
146
147////////////////////////////////////////////////////////////////////////////////
148/// Protected method to get all TGraph2D properties
149
151{
153 fX = fGraph2D->GetX();
154 fY = fGraph2D->GetY();
155 fZ = fGraph2D->GetZ();
157 else fEXlow = fGraph2D->GetEX();
159 else fEXhigh = fGraph2D->GetEX();
161 else fEYlow = fGraph2D->GetEY();
163 else fEYhigh = fGraph2D->GetEY();
165 else fEZlow = fGraph2D->GetEZ();
167 else fEZhigh = fGraph2D->GetEZ();
168}
169
170////////////////////////////////////////////////////////////////////////////////
171/// Find triangles in fDelaunay and initialise the TGraph2DPainter values
172/// needed to paint triangles or find contours.
173
198
199
200////////////////////////////////////////////////////////////////////////////////
201/// Returns the X and Y graphs building a contour. A contour level may
202/// consist in several parts not connected to each other. This function
203/// finds them and returns them in a graphs' list.
204
206{
207 // Exit if the contour is outside the Z range.
210 if (Hoption.Logz) {
211 if (zmin > 0) {
212 zmin = TMath::Log10(zmin);
213 zmax = TMath::Log10(zmax);
214 } else {
215 return nullptr;
216 }
217 }
219 Error("GetContourList", "Contour level (%g) outside the Z scope [%g,%g]",
220 contour,zmin,zmax);
221 return nullptr;
222 }
223
224 if (!fNdt) FindTriangles();
225
226 TGraph *graph = nullptr; // current graph
227 Int_t npg = 0; // number of points in the current graph
228
229 // Find all the segments making the contour
230
232 Int_t p0, p1, p2;
233 Double_t x0, y0, z0;
234 Double_t x1, y1, z1;
235 Double_t x2, y2, z2;
236 Int_t t[3],i,it,i0,i1,i2;
237
238 // Allocate space to store the segments. They cannot be more than the
239 // number of triangles.
241 std::vector<Double_t> xs0(fNdt);
242 std::vector<Double_t> ys0(fNdt);
243 std::vector<Double_t> xs1(fNdt);
244 std::vector<Double_t> ys1(fNdt);
245 for (i=0;i<fNdt;i++) {
246 xs0[i] = 0.;
247 ys0[i] = 0.;
248 xs1[i] = 0.;
249 ys1[i] = 0.;
250 }
251 Int_t nbSeg = 0;
252
253 // Loop over all the triangles in order to find all the line segments
254 // making the contour.
255
256 // old implementation
257 if (fDelaunay) {
258 for(it=0; it<fNdt; it++) {
259 t[0] = fPTried[it];
260 t[1] = fNTried[it];
261 t[2] = fMTried[it];
262 p0 = t[0]-1;
263 p1 = t[1]-1;
264 p2 = t[2]-1;
265 x0 = fX[p0]; x2 = fX[p0];
266 y0 = fY[p0]; y2 = fY[p0];
267 z0 = fZ[p0]; z2 = fZ[p0];
268
269 // Order along Z axis the points (xi,yi,zi) where "i" belongs to {0,1,2}
270 // After this z0 < z1 < z2
271 i0 = i1 = i2 = 0;
272 if (fZ[p1]<=z0) {z0=fZ[p1]; x0=fX[p1]; y0=fY[p1]; i0=1;}
273 if (fZ[p1]>z2) {z2=fZ[p1]; x2=fX[p1]; y2=fY[p1]; i2=1;}
274 if (fZ[p2]<=z0) {z0=fZ[p2]; x0=fX[p2]; y0=fY[p2]; i0=2;}
275 if (fZ[p2]>z2) {z2=fZ[p2]; x2=fX[p2]; y2=fY[p2]; i2=2;}
276 if (i0==0 && i2==0) {
277 Error("GetContourList", "wrong vertices ordering");
278 return nullptr;
279 } else {
280 i1 = 3-i2-i0;
281 }
282 x1 = fX[t[i1]-1];
283 y1 = fY[t[i1]-1];
284 z1 = fZ[t[i1]-1];
285
286 if (Hoption.Logz) {
287 z0 = TMath::Log10(z0);
288 z1 = TMath::Log10(z1);
289 z2 = TMath::Log10(z2);
290 }
291
292 if(contour >= z0 && contour <=z2) {
293 r20 = (contour-z0)/(z2-z0);
294 xs0c = r20*(x2-x0)+x0;
295 ys0c = r20*(y2-y0)+y0;
296 if(contour >= z1 && contour <=z2) {
297 r21 = (contour-z1)/(z2-z1);
298 xs1c = r21*(x2-x1)+x1;
299 ys1c = r21*(y2-y1)+y1;
300 } else {
301 r10 = (contour-z0)/(z1-z0);
302 xs1c = r10*(x1-x0)+x0;
303 ys1c = r10*(y1-y0)+y0;
304 }
305 // do not take the segments equal to a point
306 if(xs0c != xs1c || ys0c != ys1c) {
307 nbSeg++;
308 xs0[nbSeg-1] = xs0c;
309 ys0[nbSeg-1] = ys0c;
310 xs1[nbSeg-1] = xs1c;
311 ys1[nbSeg-1] = ys1c;
312 }
313 }
314 }
315 }
316 else if (fDelaunay2D) {
317
318 Int_t p[3];
319 for(const auto & face : *fDelaunay2D) {
320 p[0] = face.idx[0];
321 p[1] = face.idx[1];
322 p[2] = face.idx[2];
323 x0 = fX[p[0]]; x2 = fX[p[0]];
324 y0 = fY[p[0]]; y2 = fY[p[0]];
325 z0 = fZ[p[0]]; z2 = fZ[p[0]];
326
327 // Order along Z axis the points (xi,yi,zi) where "i" belongs to {0,1,2}
328 // After this z0 < z1 < z2
329 i0 = i1 = i2 = 0;
330 if (fZ[p[1]]<=z0) {z0=fZ[p[1]]; x0=fX[p[1]]; y0=fY[p[1]]; i0=1;}
331 if (fZ[p[1]]>z2) {z2=fZ[p[1]]; x2=fX[p[1]]; y2=fY[p[1]]; i2=1;}
332 if (fZ[p[2]]<=z0) {z0=fZ[p[2]]; x0=fX[p[2]]; y0=fY[p[2]]; i0=2;}
333 if (fZ[p[2]]>z2) {z2=fZ[p[2]]; x2=fX[p[2]]; y2=fY[p[2]]; i2=2;}
334 if (i0==0 && i2==0) {
335 Error("GetContourList", "wrong vertices ordering");
336 return nullptr;
337 } else {
338 i1 = 3-i2-i0;
339 }
340 x1 = fX[p[i1]];
341 y1 = fY[p[i1]];
342 z1 = fZ[p[i1]];
343
344 if (Hoption.Logz) {
345 z0 = TMath::Log10(z0);
346 z1 = TMath::Log10(z1);
347 z2 = TMath::Log10(z2);
348 }
349
350 if(contour >= z0 && contour <=z2) {
351 r20 = (contour-z0)/(z2-z0);
352 xs0c = r20*(x2-x0)+x0;
353 ys0c = r20*(y2-y0)+y0;
354 if(contour >= z1 && contour <=z2) {
355 r21 = (contour-z1)/(z2-z1);
356 xs1c = r21*(x2-x1)+x1;
357 ys1c = r21*(y2-y1)+y1;
358 } else {
359 r10 = (contour-z0)/(z1-z0);
360 xs1c = r10*(x1-x0)+x0;
361 ys1c = r10*(y1-y0)+y0;
362 }
363 // do not take the segments equal to a point
364 if(xs0c != xs1c || ys0c != ys1c) {
365 nbSeg++;
366 xs0[nbSeg-1] = xs0c;
367 ys0[nbSeg-1] = ys0c;
368 xs1[nbSeg-1] = xs1c;
369 ys1[nbSeg-1] = ys1c;
370 }
371 }
372 }
373 }
374
375 TList *list = new TList(); // list holding all the graphs
376
377 std::vector<Bool_t> segUsed(fNdt);
378 for(i=0; i<fNdt; i++) segUsed[i] = kFALSE;
379
380 // Find all the graphs making the contour. There is two kind of graphs,
381 // either they are "opened" or they are "closed"
382
383 // Find the opened graphs
384 Double_t xc=0, yc=0, xnc=0, ync=0;
386 Bool_t s0, s1;
387 Int_t is, js;
388 for (is=0; is<nbSeg; is++) {
389 if (segUsed[is]) continue;
390 s0 = s1 = kFALSE;
391
392 // Find to which segment is is connected. It can be connected
393 // via 0, 1 or 2 vertices.
394 for (js=0; js<nbSeg; js++) {
395 if (is==js) continue;
396 if (xs0[is]==xs0[js] && ys0[is]==ys0[js]) s0 = kTRUE;
397 if (xs0[is]==xs1[js] && ys0[is]==ys1[js]) s0 = kTRUE;
398 if (xs1[is]==xs0[js] && ys1[is]==ys0[js]) s1 = kTRUE;
399 if (xs1[is]==xs1[js] && ys1[is]==ys1[js]) s1 = kTRUE;
400 }
401
402 // Segment is is alone, not connected. It is stored in the
403 // list and the next segment is examined.
404 if (!s0 && !s1) {
405 graph = new TGraph();
406 graph->SetPoint(npg,xs0[is],ys0[is]); npg++;
407 graph->SetPoint(npg,xs1[is],ys1[is]); npg++;
408 segUsed[is] = kTRUE;
409 list->Add(graph); npg = 0;
410 continue;
411 }
412
413 // Segment is is connected via 1 vertex only and can be considered
414 // as the starting point of an opened contour.
415 if (!s0 || !s1) {
416 // Find all the segments connected to segment is
417 graph = new TGraph();
418 if (s0) {xc = xs0[is]; yc = ys0[is]; xnc = xs1[is]; ync = ys1[is];}
419 if (s1) {xc = xs1[is]; yc = ys1[is]; xnc = xs0[is]; ync = ys0[is];}
420 graph->SetPoint(npg,xnc,ync); npg++;
421 segUsed[is] = kTRUE;
422 js = 0;
423L01:
424 findNew = kFALSE;
425 if (js < nbSeg && segUsed[js]) {
426 js++;
427 goto L01;
428 } else if (xc==xs0[js] && yc==ys0[js]) {
429 xc = xs1[js];
430 yc = ys1[js];
431 findNew = kTRUE;
432 } else if (xc==xs1[js] && yc==ys1[js]) {
433 xc = xs0[js];
434 yc = ys0[js];
435 findNew = kTRUE;
436 }
437 if (findNew) {
438 segUsed[js] = kTRUE;
439 graph->SetPoint(npg,xc,yc); npg++;
440 js = 0;
441 goto L01;
442 }
443 js++;
444 if (js<nbSeg) goto L01;
445 list->Add(graph); npg = 0;
446 }
447 }
448
449
450 // Find the closed graphs. At this point all the remaining graphs
451 // are closed. Any segment can be used to start the search.
452 for (is=0; is<nbSeg; is++) {
453 if (segUsed[is]) continue;
454
455 // Find all the segments connected to segment is
456 graph = new TGraph();
457 segUsed[is] = kTRUE;
458 xc = xs0[is];
459 yc = ys0[is];
460 js = 0;
461 graph->SetPoint(npg,xc,yc); npg++;
462L02:
463 findNew = kFALSE;
464 if (js < nbSeg && segUsed[js]) {
465 js++;
466 goto L02;
467 } else if (xc==xs0[js] && yc==ys0[js]) {
468 xc = xs1[js];
469 yc = ys1[js];
470 findNew = kTRUE;
471 } else if (xc==xs1[js] && yc==ys1[js]) {
472 xc = xs0[js];
473 yc = ys0[js];
474 findNew = kTRUE;
475 }
476 if (findNew) {
477 segUsed[js] = kTRUE;
478 graph->SetPoint(npg,xc,yc); npg++;
479 js = 0;
480 goto L02;
481 }
482 js++;
483 if (js<nbSeg) goto L02;
484 // Close the contour
485 graph->SetPoint(npg,xs0[is],ys0[is]); npg++;
486 list->Add(graph); npg = 0;
487 }
488
489 return list;
490}
491
492
493////////////////////////////////////////////////////////////////////////////////
494/// Paint a TGraphDelaunay according to the value of "option":
495///
496/// | Option | Description |
497/// |----------|---------------------------------------------------------------|
498/// | "TRI" | The Delaunay triangles are drawn using filled area. An hidden surface drawing technique is used. The surface is painted with the current fill area color. The edges of each triangles are painted with the current line color. |
499/// | "TRIW" | The Delaunay triangles are drawn as wire frame. |
500/// | "TRI1" | The Delaunay triangles are painted with color levels. The edges of each triangles are painted with the current line color. |
501/// | "TRI2" | the Delaunay triangles are painted with color levels. |
502/// | "P" | Draw a marker at each vertex. |
503/// | "P0" | Draw a circle at each vertex. Each circle background is white. |
504/// | "PCOL" | Draw a marker at each vertex. The color of each marker is defined according to its Z position. |
505/// | "CONT" | Draw contours. |
506/// | "LINE" | Draw a 3D polyline. |
507
509{
510 TString opt = option;
511 opt.ToLower();
512 Bool_t triangles = opt.Contains("tri") ||
513 opt.Contains("tri1") ||
514 opt.Contains("tri2");
515 if (opt.Contains("tri0")) triangles = kFALSE;
516
517 Bool_t markers = opt.Contains("p") && !triangles;
518 Bool_t contour = opt.Contains("cont");
519 Bool_t line = opt.Contains("line");
520 Bool_t err = opt.Contains("err");
521
523
524 fGraph2D->TAttLine::Modify();
525 fGraph2D->TAttFill::Modify();
526 fGraph2D->TAttMarker::Modify();
527
528 // Compute minimums and maximums
530 Int_t first = xaxis->GetFirst();
531 fXmin = xaxis->GetBinLowEdge(first);
532 if (Hoption.Logx && fXmin <= 0) fXmin = xaxis->GetBinUpEdge(xaxis->FindFixBin(0.01*xaxis->GetBinWidth(first)));
533 fXmax = xaxis->GetBinUpEdge(xaxis->GetLast());
535 first = yaxis->GetFirst();
536 fYmin = yaxis->GetBinLowEdge(first);
537 if (Hoption.Logy && fYmin <= 0) fYmin = yaxis->GetBinUpEdge(yaxis->FindFixBin(0.01*yaxis->GetBinWidth(first)));
538 fYmax = yaxis->GetBinUpEdge(yaxis->GetLast());
541 if (Hoption.Logz && fZmin <= 0) fZmin = TMath::Min((Double_t)1, (Double_t)0.001*fGraph2D->GetZmax());
542
544 if (contour) PaintContour(option);
546 if (err) PaintErrors(option);
548}
549
550
551////////////////////////////////////////////////////////////////////////////////
552/// Paints the 2D graph as a contour plot. Delaunay triangles are used
553/// to compute the contours.
554
556{
557 // Initialize the levels on the Z axis
558 Int_t ncolors = gStyle->GetNumberOfColors();
559 Int_t ndiv = gCurrentHist->GetContour();
560 if (ndiv == 0 ) {
561 ndiv = gStyle->GetNumberContours();
563 }
564 Int_t ndivz = TMath::Abs(ndiv);
566
567 if (!fNdt) FindTriangles();
568
569 for (Int_t k=0; k<ndiv; k++) {
572 TIter next(l);
573 while (auto obj = next()) {
574 if(obj->InheritsFrom(TGraph::Class()) ) {
575 TGraph *g = (TGraph*)obj;
577 g->SetLineStyle(fGraph2D->GetLineStyle());
578 Int_t theColor = Int_t((k+0.99)*Float_t(ncolors)/Float_t(ndivz));
579 g->SetLineColor(gStyle->GetColorPalette(theColor));
580 g->Paint("l");
581 }
582 }
583 if (l) { l->Delete(); delete l; }
584 }
585}
586
587
588////////////////////////////////////////////////////////////////////////////////
589/// Paints the 2D graph as error bars
590
592{
593 Double_t temp1[3],temp2[3];
594
595 TView *view = gPad ? gPad->GetView() : nullptr;
596 if (!view) {
597 Error("PaintErrors", "No TView in current pad");
598 return;
599 }
600
601 Int_t it;
602 Double_t xm[2];
603 Double_t ym[2];
604 Double_t err = 0;
605
609 fGraph2D->TAttLine::Modify();
610
611 for (it=0; it<fNpoints; it++) {
612 if(fX[it] < fXmin || fX[it] > fXmax) continue;
613 if(fY[it] < fYmin || fY[it] > fYmax) continue;
614
615 err = 0;
616 if (fEXlow) err = fEXlow[it];
617 temp1[0] = fX[it]-err;
618 temp1[1] = fY[it];
619 temp1[2] = fZ[it];
620 temp1[0] = TMath::Max(temp1[0],fXmin);
621 temp1[1] = TMath::Max(temp1[1],fYmin);
622 temp1[2] = TMath::Max(temp1[2],fZmin);
623 temp1[2] = TMath::Min(temp1[2],fZmax);
624 if (Hoption.Logx) temp1[0] = TMath::Log10(temp1[0]);
625 if (Hoption.Logy) temp1[1] = TMath::Log10(temp1[1]);
626 if (Hoption.Logz) temp1[2] = TMath::Log10(temp1[2]);
627 view->WCtoNDC(temp1, &temp2[0]);
628 xm[0] = temp2[0];
629 ym[0] = temp2[1];
630 err = 0;
631 if (fEXhigh) err = fEXhigh[it];
632 temp1[0] = fX[it]+err;
633 temp1[0] = TMath::Max(temp1[0],fXmin);
634 if (Hoption.Logx) temp1[0] = TMath::Log10(temp1[0]);
635 view->WCtoNDC(temp1, &temp2[0]);
636 xm[1] = temp2[0];
637 ym[1] = temp2[1];
638 gPad->PaintPolyLine(2,xm,ym);
639
640 err = 0;
641 if (fEYlow) err = fEYlow[it];
642 temp1[0] = fX[it];
643 temp1[1] = fY[it]-err;
644 temp1[2] = fZ[it];
645 temp1[0] = TMath::Max(temp1[0],fXmin);
646 temp1[1] = TMath::Max(temp1[1],fYmin);
647 temp1[2] = TMath::Max(temp1[2],fZmin);
648 temp1[2] = TMath::Min(temp1[2],fZmax);
649 if (Hoption.Logx) temp1[0] = TMath::Log10(temp1[0]);
650 if (Hoption.Logy) temp1[1] = TMath::Log10(temp1[1]);
651 if (Hoption.Logz) temp1[2] = TMath::Log10(temp1[2]);
652 view->WCtoNDC(temp1, &temp2[0]);
653 xm[0] = temp2[0];
654 ym[0] = temp2[1];
655 err = 0;
656 if (fEYhigh) err = fEYhigh[it];
657 temp1[1] = fY[it]+err;
658 temp1[1] = TMath::Max(temp1[1],fYmin);
659 if (Hoption.Logy) temp1[1] = TMath::Log10(temp1[1]);
660 view->WCtoNDC(temp1, &temp2[0]);
661 xm[1] = temp2[0];
662 ym[1] = temp2[1];
663 gPad->PaintPolyLine(2,xm,ym);
664
665 err = 0;
666 if (fEZlow) err = fEZlow[it];
667 temp1[0] = fX[it];
668 temp1[1] = fY[it];
669 temp1[2] = fZ[it]-err;
670 temp1[0] = TMath::Max(temp1[0],fXmin);
671 temp1[1] = TMath::Max(temp1[1],fYmin);
672 temp1[2] = TMath::Max(temp1[2],fZmin);
673 temp1[2] = TMath::Min(temp1[2],fZmax);
674 if (Hoption.Logx) temp1[0] = TMath::Log10(temp1[0]);
675 if (Hoption.Logy) temp1[1] = TMath::Log10(temp1[1]);
676 if (Hoption.Logz) temp1[2] = TMath::Log10(temp1[2]);
677 view->WCtoNDC(temp1, &temp2[0]);
678 xm[0] = temp2[0];
679 ym[0] = temp2[1];
680 err = 0;
681 if (fEZhigh) err = fEZhigh[it];
682 temp1[2] = fZ[it]+err;
683 temp1[2] = TMath::Max(temp1[2],fZmin);
684 temp1[2] = TMath::Min(temp1[2],fZmax);
685 if (Hoption.Logz) temp1[2] = TMath::Log10(temp1[2]);
686 view->WCtoNDC(temp1, &temp2[0]);
687 xm[1] = temp2[0];
688 ym[1] = temp2[1];
689 gPad->PaintPolyLine(2,xm,ym);
690 }
691}
692
693
694////////////////////////////////////////////////////////////////////////////////
695/// Paints one triangle.
696///
697/// - nblev = 0 : paint the color levels
698/// - nblev != 0 : paint the grid
699
702{
703 Int_t i, fillColor, ncolors, theColor0, theColor2;
704
705 Int_t p[3];
706 if (fDelaunay) {
707 p[0]=t[0]-1;
708 p[1]=t[1]-1;
709 p[2]=t[2]-1;
710 }
711 else {
712 p[0]=t[0];
713 p[1]=t[1];
714 p[2]=t[2];
715 }
716 Double_t xl[2],yl[2];
717 Double_t zl, r21, r20, r10;
718 Double_t x0 = x[0] , x2 = x[0];
719 Double_t y0 = y[0] , y2 = y[0];
720 Double_t z0 = fZ[p[0]], z2 = fZ[p[0]];
721 Double_t zmin = fGraph2D->GetMinimum();
722 Double_t zmax = fGraph2D->GetMaximum();
723 if (zmin==-1111 && zmax==-1111) {
724 zmin = TMath::Min(fZmin, 0.);
725 if (Hoption.Logz && zmin <= 0) zmin = TMath::Min((Double_t)1, (Double_t)0.001*fGraph2D->GetZmax());
726 zmax = fZmax;
727 }
728
729 // Order along Z axis the points (xi,yi,zi) where "i" belongs to {0,1,2}
730 // After this z0 < z1 < z2
731 Int_t i0=0, i1=0, i2=0;
732 if (fZ[p[1]]<=z0) {z0=fZ[p[1]]; x0=x[1]; y0=y[1]; i0=1;}
733 if (fZ[p[1]]>z2) {z2=fZ[p[1]]; x2=x[1]; y2=y[1]; i2=1;}
734 if (fZ[p[2]]<=z0) {z0=fZ[p[2]]; x0=x[2]; y0=y[2]; i0=2;}
735 if (fZ[p[2]]>z2) {z2=fZ[p[2]]; x2=x[2]; y2=y[2]; i2=2;}
736 i1 = 3-i2-i0;
737 Double_t x1 = x[i1];
738 Double_t y1 = y[i1];
739 Double_t z1 = fZ[p[i1]];
740
741 if (z0>zmax) z0 = zmax;
742 if (z2>zmax) z2 = zmax;
743 if (z0<zmin) z0 = zmin;
744 if (z2<zmin) z2 = zmin;
745 if (z1>zmax) z1 = zmax;
746 if (z1<zmin) z1 = zmin;
747
748 if (Hoption.Logz) {
749 z0 = TMath::Log10(z0);
750 z1 = TMath::Log10(z1);
751 z2 = TMath::Log10(z2);
752 zmin = TMath::Log10(zmin);
753 zmax = TMath::Log10(zmax);
754 }
755
756 // zi = Z values of the stripe number i
757 // zip = Previous zi
758 Double_t zi=0, zip=0;
759
760 if (nblev <= 0) {
761 // Paint the colors levels
762
763 // Compute the color associated to z0 (theColor0) and z2 (theColor2)
764 ncolors = gStyle->GetNumberOfColors();
765 theColor0 = (Int_t)( ((z0-zmin)/(zmax-zmin))*(ncolors-1) );
766 theColor2 = (Int_t)( ((z2-zmin)/(zmax-zmin))*(ncolors-1) );
767
768 // The stripes drawn to fill the triangles may have up to 5 points
769 Double_t xp[5], yp[5];
770
771 // rl = Ratio between z0 and z2 (long)
772 // rs = Ratio between z0 and z1 or z1 and z2 (short)
773 Double_t rl,rs;
774
775 // ci = Color of the stripe number i
776 // npf = number of point needed to draw the current stripe
777 Int_t ci,npf;
778
780
781 // If the z0's color and z2's colors are the same, the whole triangle
782 // can be painted in one go.
783 if(theColor0 == theColor2) {
785 fGraph2D->TAttFill::Modify();
786 gPad->PaintFillArea(3,x,y);
787
788 // The triangle must be painted with several colors
789 } else {
790 for(ci=theColor0; ci<=theColor2; ci++) {
792 fGraph2D->TAttFill::Modify();
793 if (ci==theColor0) {
794 zi = (((ci+1)*(zmax-zmin))/(ncolors-1))+zmin;
795 xp[0] = x0;
796 yp[0] = y0;
797 rl = (zi-z0)/(z2-z0);
798 xp[1] = rl*(x2-x0)+x0;
799 yp[1] = rl*(y2-y0)+y0;
800 if (zi>=z1 || z0==z1) {
801 rs = (zi-z1)/(z2-z1);
802 xp[2] = rs*(x2-x1)+x1;
803 yp[2] = rs*(y2-y1)+y1;
804 xp[3] = x1;
805 yp[3] = y1;
806 npf = 4;
807 } else {
808 rs = (zi-z0)/(z1-z0);
809 xp[2] = rs*(x1-x0)+x0;
810 yp[2] = rs*(y1-y0)+y0;
811 npf = 3;
812 }
813 } else if (ci==theColor2) {
814 xp[0] = xp[1];
815 yp[0] = yp[1];
816 xp[1] = x2;
817 yp[1] = y2;
818 if (zi<z1 || z2==z1) {
819 xp[3] = xp[2];
820 yp[3] = yp[2];
821 xp[2] = x1;
822 yp[2] = y1;
823 npf = 4;
824 } else {
825 npf = 3;
826 }
827 } else {
828 zi = (((ci+1)*(zmax-zmin))/(ncolors-1))+zmin;
829 xp[0] = xp[1];
830 yp[0] = yp[1];
831 rl = (zi-z0)/(z2-z0);
832 xp[1] = rl*(x2-x0)+x0;
833 yp[1] = rl*(y2-y0)+y0;
834 if ( zi>=z1 && zip<=z1) {
835 xp[3] = x1;
836 yp[3] = y1;
837 xp[4] = xp[2];
838 yp[4] = yp[2];
839 npf = 5;
840 } else {
841 xp[3] = xp[2];
842 yp[3] = yp[2];
843 npf = 4;
844 }
845 if (zi<z1) {
846 rs = (zi-z0)/(z1-z0);
847 xp[2] = rs*(x1-x0)+x0;
848 yp[2] = rs*(y1-y0)+y0;
849 } else {
850 rs = (zi-z1)/(z2-z1);
851 xp[2] = rs*(x2-x1)+x1;
852 yp[2] = rs*(y2-y1)+y1;
853 }
854 }
855 zip = zi;
856 // Paint a stripe
857 gPad->PaintFillArea(npf,xp,yp);
858 }
859 }
861 fGraph2D->TAttFill::Modify();
862
863 } else {
864 // Paint the grid levels
866 fGraph2D->TAttLine::Modify();
867 for(i=0; i<nblev; i++){
868 zl=glev[i];
869 if(zl >= z0 && zl <=z2) {
870 r21=(zl-z1)/(z2-z1);
871 r20=(zl-z0)/(z2-z0);
872 r10=(zl-z0)/(z1-z0);
873 xl[0]=r20*(x2-x0)+x0;
874 yl[0]=r20*(y2-y0)+y0;
875 if(zl >= z1 && zl <=z2) {
876 xl[1]=r21*(x2-x1)+x1;
877 yl[1]=r21*(y2-y1)+y1;
878 } else {
879 xl[1]=r10*(x1-x0)+x0;
880 yl[1]=r10*(y1-y0)+y0;
881 }
882 gPad->PaintPolyLine(2,xl,yl);
883 }
884 }
886 fGraph2D->TAttLine::Modify();
887 }
888}
889
890
891////////////////////////////////////////////////////////////////////////////////
892/// Paints the 2D graph as PaintPolyMarker
893
895{
896 Double_t temp1[3],temp2[3];
897
898 TView *view = gPad ? gPad->GetView() : nullptr;
899 if (!view) {
900 Error("PaintPolyMarker", "No TView in current pad");
901 return;
902 }
903
904 TString opt = option;
905 opt.ToLower();
906 Bool_t markers0 = opt.Contains("p0");
907 Bool_t colors = opt.Contains("pcol");
908 Int_t ncolors = gStyle->GetNumberOfColors();
909 Int_t it, theColor;
910
911 // Initialize the levels on the Z axis
912 if (colors) {
913 Int_t ndiv = gCurrentHist->GetContour();
914 if (ndiv == 0 ) {
915 ndiv = gStyle->GetNumberContours();
917 }
919 }
920
921 std::vector<Double_t> xm(fNpoints);
922 std::vector<Double_t> ym(fNpoints);
923 std::vector<Double_t> zm(fNpoints);
926
927 // min and max for colors
930 if (hzmincol==-1111 && hzmaxcol==-1111) {
933 hzmaxcol = fZmax;
934 }
935 if (Hoption.Logz) {
938 }
939
940 Double_t Xeps = (fXmax-fXmin)*0.0001;
941 Double_t Yeps = (fYmax-fYmin)*0.0001;
942 Double_t Zeps = (hzmax-hzmin)*0.0001;
943
944 Int_t npd = 0;
945 for (it=0; it<fNpoints; it++) {
946 xm[it] = 0;
947 ym[it] = 0;
948 if(fXmin - fX[it] > Xeps || fX[it] - fXmax > Xeps) continue;
949 if(fYmin - fY[it] > Yeps || fY[it] - fYmax > Yeps) continue;
950 if(hzmin - fZ[it] > Zeps || fZ[it] - hzmax > Zeps) continue;
951 temp1[0] = fX[it];
952 temp1[1] = fY[it];
953 temp1[2] = fZ[it];
954 temp1[0] = TMath::Max(temp1[0],fXmin);
955 temp1[1] = TMath::Max(temp1[1],fYmin);
956 temp1[2] = TMath::Max(temp1[2],hzmin);
957 temp1[2] = TMath::Min(temp1[2],hzmax);
958 if (Hoption.Logx) temp1[0] = TMath::Log10(temp1[0]);
959 if (Hoption.Logy) temp1[1] = TMath::Log10(temp1[1]);
960 if (Hoption.Logz) temp1[2] = TMath::Log10(temp1[2]);
961 view->WCtoNDC(temp1, &temp2[0]);
962 xm[npd] = temp2[0];
963 ym[npd] = temp2[1];
964 zm[npd] = temp1[2];
965 npd++;
966 }
967 if (markers0) {
968 PaintPolyMarker0(npd,xm.data(),ym.data());
969 } else if (colors) {
971 for (it=0; it<npd; it++) {
972 theColor = (Int_t)( ((zm[it]-hzmincol)/(hzmaxcol-hzmincol))*(ncolors-1) );
974 fGraph2D->TAttMarker::Modify();
975 gPad->PaintPolyMarker(1,xm.data()+it,ym.data()+it);
976 }
978 } else {
982 fGraph2D->TAttMarker::Modify();
983 gPad->PaintPolyMarker(npd,xm.data(),ym.data());
984 }
985}
986
987
988////////////////////////////////////////////////////////////////////////////////
989/// Paints the 2D graph as PaintPolyLine
990
992{
993 Double_t temp1[3],temp2[3];
994
995 TView *view = gPad ? gPad->GetView() : nullptr;
996 if (!view) {
997 Error("PaintPolyLine", "No TView in current pad");
998 return;
999 }
1000
1001 Int_t it;
1002 Double_t Xeps = (fXmax - fXmin) * 0.0001;
1003 Double_t Yeps = (fYmax - fYmin) * 0.0001;
1004
1005 std::vector<Double_t> xm(fNpoints);
1006 std::vector<Double_t> ym(fNpoints);
1007 Int_t npd = 0;
1008
1009 for (it=0; it<fNpoints; it++) {
1010 if (fXmin - fX[it] > Xeps || fX[it] - fXmax > Xeps)
1011 continue;
1012 if (fYmin - fY[it] > Yeps || fY[it] - fYmax > Yeps)
1013 continue;
1014 npd++;
1015 temp1[0] = fX[it];
1016 temp1[1] = fY[it];
1017 temp1[2] = fZ[it];
1018 temp1[0] = TMath::Max(temp1[0],fXmin);
1019 temp1[1] = TMath::Max(temp1[1],fYmin);
1020 temp1[2] = TMath::Max(temp1[2],fZmin);
1021 temp1[2] = TMath::Min(temp1[2],fZmax);
1022 if (Hoption.Logx) temp1[0] = TMath::Log10(temp1[0]);
1023 if (Hoption.Logy) temp1[1] = TMath::Log10(temp1[1]);
1024 if (Hoption.Logz) temp1[2] = TMath::Log10(temp1[2]);
1025 view->WCtoNDC(temp1, &temp2[0]);
1026 xm[npd - 1] = temp2[0];
1027 ym[npd - 1] = temp2[1];
1028 }
1032 fGraph2D->TAttLine::Modify();
1033 gPad->PaintPolyLine(npd,xm.data(),ym.data());
1034}
1035
1036
1037////////////////////////////////////////////////////////////////////////////////
1038/// Paints a circle at each vertex. Each circle background is white.
1039
1041{
1045 for (Int_t i=0; i<n; i++) {
1048 fGraph2D->TAttMarker::Modify();
1049 gPad->PaintPolyMarker(1,&x[i],&y[i]);
1052 fGraph2D->TAttMarker::Modify();
1053 gPad->PaintPolyMarker(1,&x[i],&y[i]);
1054 }
1056}
1057
1058
1059////////////////////////////////////////////////////////////////////////////////
1060/// Paints the 2D graph as triangles
1061
1069
1070
1071////////////////////////////////////////////////////////////////////////////////
1072/// Paints the 2D graph as triangles (old implementation)
1073
1075{
1076 Double_t x[4], y[4], temp1[3],temp2[3];
1077 Int_t it,t[3];
1078 std::vector<Int_t> order;
1079 std::vector<Double_t> dist;
1080
1081 TView *view = gPad ? gPad->GetView() : nullptr;
1082 if (!view) {
1083 Error("PaintTriangles", "No TView in current pad");
1084 return;
1085 }
1086
1087 TString opt = option;
1088 opt.ToLower();
1089 Bool_t tri1 = opt.Contains("tri1");
1090 Bool_t tri2 = opt.Contains("tri2");
1091 Bool_t markers = opt.Contains("p");
1092 Bool_t markers0 = opt.Contains("p0");
1093 Bool_t wire = opt.Contains("w");
1094
1095 // Define the grid levels drawn on the triangles.
1096 // The grid levels are aligned on the Z axis' main tick marks.
1097 Int_t nblev = 0;
1098 std::vector<Double_t> glev;
1099 if (!tri1 && !tri2 && !wire) {
1101 Int_t nbins;
1102 Double_t binLow = 0, binHigh = 0, binWidth = 0;
1103
1104 // Find the main tick marks positions.
1105 Double_t *r0 = view->GetRmin();
1106 Double_t *r1 = view->GetRmax();
1107 if (!r0 || !r1) return;
1108
1109 if (ndivz > 0) {
1111 binLow, binHigh, nbins, binWidth, " ");
1112 } else {
1113 nbins = TMath::Abs(ndivz);
1114 binLow = r0[2];
1115 binHigh = r1[2];
1116 binWidth = (binHigh-binLow)/nbins;
1117 }
1118 // Define the grid levels
1119 nblev = nbins+1;
1120 glev.resize(nblev);
1121 for (Int_t i = 0; i < nblev; ++i) glev[i] = binLow+i*binWidth;
1122 }
1123
1124 // Initialize the levels on the Z axis
1125 if (tri1 || tri2) {
1126 Int_t ndiv = gCurrentHist->GetContour();
1127 if (ndiv == 0 ) {
1128 ndiv = gStyle->GetNumberContours();
1129 gCurrentHist->SetContour(ndiv);
1130 }
1132 }
1133
1134 // For each triangle, compute the distance between the triangle centre
1135 // and the back planes. Then these distances are sorted in order to draw
1136 // the triangles from back to front.
1137 if (!fNdt) FindTriangles();
1138 Double_t cp = TMath::Cos(view->GetLongitude()*TMath::Pi()/180.);
1139 Double_t sp = TMath::Sin(view->GetLongitude()*TMath::Pi()/180.);
1140 order.resize(fNdt);
1141 dist.resize(fNdt);
1142 Double_t xd,yd;
1143 Int_t p, n, m;
1144 Bool_t o = kFALSE;
1145 for (it=0; it<fNdt; it++) {
1146 p = fPTried[it];
1147 n = fNTried[it];
1148 m = fMTried[it];
1149 xd = (fXN[p]+fXN[n]+fXN[m])/3;
1150 yd = (fYN[p]+fYN[n]+fYN[m])/3;
1151 if ((cp >= 0) && (sp >= 0.)) {
1152 dist[it] = -(fXNmax-xd+fYNmax-yd);
1153 } else if ((cp <= 0) && (sp >= 0.)) {
1154 dist[it] = -(fXNmax-xd+yd-fYNmin);
1155 o = kTRUE;
1156 } else if ((cp <= 0) && (sp <= 0.)) {
1157 dist[it] = -(xd-fXNmin+yd-fYNmin);
1158 } else {
1159 dist[it] = -(xd-fXNmin+fYNmax-yd);
1160 o = kTRUE;
1161 }
1162 }
1163 TMath::Sort(fNdt, dist.data(), order.data(), o);
1164
1165 // Draw the triangles and markers if requested
1168 fGraph2D->SetFillStyle(1001);
1169 fGraph2D->TAttFill::Modify();
1171 fGraph2D->TAttLine::Modify();
1172 int lst = fGraph2D->GetLineStyle();
1175 for (it=0; it<fNdt; it++) {
1176 t[0] = fPTried[order[it]];
1177 t[1] = fNTried[order[it]];
1178 t[2] = fMTried[order[it]];
1179 for (Int_t k=0; k<3; k++) {
1180 if(fX[t[k]-1] < fXmin || fX[t[k]-1] > fXmax) goto endloop;
1181 if(fY[t[k]-1] < fYmin || fY[t[k]-1] > fYmax) goto endloop;
1182 if(fZ[t[k]-1] < zmin || fZ[t[k]-1] > zmax) goto endloop;
1183 temp1[0] = fX[t[k]-1];
1184 temp1[1] = fY[t[k]-1];
1185 temp1[2] = fZ[t[k]-1];
1186 temp1[0] = TMath::Max(temp1[0],fXmin);
1187 temp1[1] = TMath::Max(temp1[1],fYmin);
1188 temp1[2] = TMath::Max(temp1[2],fZmin);
1189 temp1[2] = TMath::Min(temp1[2],fZmax);
1190 if (Hoption.Logx) temp1[0] = TMath::Log10(temp1[0]);
1191 if (Hoption.Logy) temp1[1] = TMath::Log10(temp1[1]);
1192 if (Hoption.Logz) temp1[2] = TMath::Log10(temp1[2]);
1193 view->WCtoNDC(temp1, &temp2[0]);
1194 x[k] = temp2[0];
1195 y[k] = temp2[1];
1196 }
1197 x[3] = x[0];
1198 y[3] = y[0];
1199 if (tri1 || tri2) PaintLevels(t,x,y);
1200 if (!tri1 && !tri2 && !wire) {
1201 gPad->PaintFillArea(3,x,y);
1202 PaintLevels(t,x,y,nblev,glev.data());
1203 }
1204 if (!tri2) gPad->PaintPolyLine(4,x,y);
1205 if (markers) {
1206 if (markers0) {
1207 PaintPolyMarker0(3,x,y);
1208 } else {
1212 fGraph2D->TAttMarker::Modify();
1213 gPad->PaintPolyMarker(3,x,y);
1214 }
1215 }
1216endloop:
1217 continue;
1218 }
1221 fGraph2D->TAttLine::Modify();
1222 fGraph2D->TAttFill::Modify();
1223}
1224
1225
1226////////////////////////////////////////////////////////////////////////////////
1227/// Paints the 2D graph as triangles (new implementation)
1228
1230{
1231
1232 Double_t x[4], y[4], temp1[3],temp2[3];
1233 Int_t p[3];
1234
1235 TView *view = gPad ? gPad->GetView() : nullptr;
1236 if (!view) {
1237 Error("PaintTriangles", "No TView in current pad");
1238 return;
1239 }
1240
1241 TString opt = option;
1242 opt.ToLower();
1243 Bool_t tri1 = opt.Contains("tri1");
1244 Bool_t tri2 = opt.Contains("tri2");
1245 Bool_t markers = opt.Contains("p");
1246 Bool_t markers0 = opt.Contains("p0");
1247 Bool_t wire = opt.Contains("w");
1248
1249 // Define the grid levels drawn on the triangles.
1250 // The grid levels are aligned on the Z axis' main tick marks.
1251 Int_t nblev=0;
1252 std::vector<Double_t> glev;
1253 if (!tri1 && !tri2 && !wire) {
1255 Int_t nbins;
1256 Double_t binLow = 0, binHigh = 0, binWidth = 0;
1257
1258 // Find the main tick marks positions.
1259 Double_t *r0 = view->GetRmin();
1260 Double_t *r1 = view->GetRmax();
1261 if (!r0 || !r1) return;
1262
1263 if (ndivz > 0) {
1265 binLow, binHigh, nbins, binWidth, " ");
1266 } else {
1267 nbins = TMath::Abs(ndivz);
1268 binLow = r0[2];
1269 binHigh = r1[2];
1270 binWidth = (binHigh-binLow)/nbins;
1271 }
1272 // Define the grid levels
1273 nblev = nbins+1;
1274 glev.resize(nblev);
1275 for (Int_t i = 0; i < nblev; ++i) glev[i] = binLow+i*binWidth;
1276 }
1277
1278 // Initialize the levels on the Z axis
1279 if (tri1 || tri2) {
1280 Int_t ndiv = gCurrentHist->GetContour();
1281 if (ndiv == 0 ) {
1282 ndiv = gStyle->GetNumberContours();
1283 gCurrentHist->SetContour(ndiv);
1284 }
1286 }
1287
1288 // For each triangle, compute the distance between the triangle centre
1289 // and the back planes. Then these distances are sorted in order to draw
1290 // the triangles from back to front.
1291 if (!fNdt) FindTriangles();
1292 Double_t cp = TMath::Cos(view->GetLongitude()*TMath::Pi()/180.);
1293 Double_t sp = TMath::Sin(view->GetLongitude()*TMath::Pi()/180.);
1294
1295 Bool_t reverse = kFALSE;
1296 std::function<Double_t(Double_t, Double_t)> fDist;
1297
1298 if ((cp >= 0) && (sp >= 0.)) {
1299 fDist = [&](Double_t xd, Double_t yd) -> Double_t { return -(fXNmax-xd+fYNmax-yd);};
1300 } else if ((cp <= 0) && (sp >= 0.)) {
1301 fDist = [&](Double_t xd, Double_t yd) -> Double_t { return -(fXNmax-xd+yd-fYNmin);};
1302 reverse = kTRUE;
1303 } else if ((cp <= 0) && (sp <= 0.)) {
1304 fDist = [&](Double_t xd, Double_t yd) -> Double_t { return -(xd-fXNmin+yd-fYNmin);};
1305 } else {
1306 fDist = [&](Double_t xd, Double_t yd) -> Double_t { return -(xd-fXNmin+fYNmax-yd);};
1307 reverse = kTRUE;
1308 }
1309
1310 typedef std::pair<Double_t, TGraphDelaunay2D::Triangles::const_iterator> DistEntry;
1311 std::vector<DistEntry> dist;
1312 for(auto it = fDelaunay2D->begin(); it != fDelaunay2D->end(); ++it){
1313 auto face = *it;
1314 Double_t xd = (face.x[0] + face.x[1] + face.x[2]) / 3;
1315 Double_t yd = (face.y[0] + face.y[1] + face.y[2]) / 3;
1316
1317 dist.emplace_back(fDist(xd, yd), it);
1318 }
1319
1320 std::sort(dist.begin(), dist.end(),
1321 [&](const DistEntry & a, const DistEntry & b){ return !reverse ? (a.first < b.first) : (b.first < a .first); });
1322
1323 // Draw the triangles and markers if requested
1326 fGraph2D->SetFillStyle(1001);
1327 fGraph2D->TAttFill::Modify();
1329 fGraph2D->TAttLine::Modify();
1330 int lst = fGraph2D->GetLineStyle();
1333 for (const auto & it : dist) {
1334 p[0] = it.second->idx[0];
1335 p[1] = it.second->idx[1];
1336 p[2] = it.second->idx[2];
1337 for (Int_t k=0; k<3; k++) {
1338 if(fX[p[k]] < fXmin || fX[p[k]] > fXmax) goto endloop;
1339 if(fY[p[k]] < fYmin || fY[p[k]] > fYmax) goto endloop;
1340 if(fZ[p[k]] < zmin || fZ[p[k]] > zmax) goto endloop;
1341 temp1[0] = fX[p[k]];
1342 temp1[1] = fY[p[k]];
1343 temp1[2] = fZ[p[k]];
1344 temp1[0] = TMath::Max(temp1[0],fXmin);
1345 temp1[1] = TMath::Max(temp1[1],fYmin);
1346 temp1[2] = TMath::Max(temp1[2],fZmin);
1347 temp1[2] = TMath::Min(temp1[2],fZmax);
1348 if (Hoption.Logx) temp1[0] = TMath::Log10(temp1[0]);
1349 if (Hoption.Logy) temp1[1] = TMath::Log10(temp1[1]);
1350 if (Hoption.Logz) temp1[2] = TMath::Log10(temp1[2]);
1351 view->WCtoNDC(temp1, &temp2[0]);
1352 x[k] = temp2[0];
1353 y[k] = temp2[1];
1354 }
1355 x[3] = x[0];
1356 y[3] = y[0];
1357 if (tri1 || tri2) PaintLevels(p,x,y);
1358 if (!tri1 && !tri2 && !wire) {
1359 gPad->PaintFillArea(3,x,y);
1360 PaintLevels(p,x,y,nblev,glev.data());
1361 }
1362 if (!tri2) gPad->PaintPolyLine(4,x,y);
1363 if (markers) {
1364 if (markers0) {
1365 PaintPolyMarker0(3,x,y);
1366 } else {
1370 fGraph2D->TAttMarker::Modify();
1371 gPad->PaintPolyMarker(3,x,y);
1372 }
1373 }
1374endloop:
1375 continue;
1376 }
1379 fGraph2D->TAttLine::Modify();
1380 fGraph2D->TAttFill::Modify();
1381}
#define R__EXTERN
Definition DllImport.h:26
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define g(i)
Definition RSha256.hxx:105
#define s0(x)
Definition RSha256.hxx:90
#define a(i)
Definition RSha256.hxx:99
#define s1(x)
Definition RSha256.hxx:91
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
float Float_t
Float 4 bytes (float)
Definition RtypesCore.h:71
constexpr Bool_t kFALSE
Definition RtypesCore.h:108
double Double_t
Double 8 bytes.
Definition RtypesCore.h:73
constexpr Bool_t kTRUE
Definition RtypesCore.h:107
const char Option_t
Option string (const char)
Definition RtypesCore.h:80
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t option
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize fs
Option_t Option_t TPoint TPoint const char y1
R__EXTERN TH1 * gCurrentHist
R__EXTERN Hoption_t Hoption
R__EXTERN TStyle * gStyle
Definition TStyle.h:442
#define gPad
Color * colors
Definition X3DBuffer.c:21
virtual Int_t GetNdivisions() const
Definition TAttAxis.h:37
virtual Color_t GetFillColor() const
Return the fill area color.
Definition TAttFill.h:31
virtual Style_t GetFillStyle() const
Return the fill area style.
Definition TAttFill.h:32
virtual void SetFillColor(Color_t fcolor)
Set the fill area color.
Definition TAttFill.h:38
virtual void SetFillStyle(Style_t fstyle)
Set the fill area style.
Definition TAttFill.h:40
virtual Color_t GetLineColor() const
Return the line color.
Definition TAttLine.h:35
virtual void SetLineStyle(Style_t lstyle)
Set the line style.
Definition TAttLine.h:44
virtual Width_t GetLineWidth() const
Return the line width.
Definition TAttLine.h:37
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:45
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:42
virtual Style_t GetLineStyle() const
Return the line style.
Definition TAttLine.h:36
virtual Style_t GetMarkerStyle() const
Return the marker style.
Definition TAttMarker.h:33
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition TAttMarker.h:39
virtual Color_t GetMarkerColor() const
Return the marker color.
Definition TAttMarker.h:32
virtual Size_t GetMarkerSize() const
Return the marker size.
Definition TAttMarker.h:34
virtual void SetMarkerStyle(Style_t mstyle=1)
Set the marker style.
Definition TAttMarker.h:41
virtual void SetMarkerSize(Size_t msize=1)
Set the marker size.
Definition TAttMarker.h:46
Class to manage histogram axis.
Definition TAxis.h:32
void PaintContour(Option_t *option)
Paints the 2D graph as a contour plot.
Int_t * fNTried
! Pointer to fDelaunay->fNTried
Double_t * fYN
! Pointer to fDelaunay->fYN
TGraphDelaunay * fDelaunay
! Pointer to the TGraphDelaunay2D to be painted
void Paint(Option_t *option) override
Paint a TGraphDelaunay according to the value of "option":
Double_t fYmin
! fGraph2D->fHistogram Ymin
Double_t * fZ
! Pointer to fGraph2D->fZ
void PaintLevels(Int_t *v, Double_t *x, Double_t *y, Int_t nblev=0, Double_t *glev=nullptr)
Paints one triangle.
TGraphDelaunay2D * fDelaunay2D
! Pointer to the TGraphDelaunay2D to be painted
TGraph2D * fGraph2D
! Pointer to the TGraph2D in fDelaunay
void PaintErrors(Option_t *option)
Paints the 2D graph as error bars.
Double_t * fXN
! Pointer to fDelaunay->fXN
Int_t * fMTried
! Pointer to fDelaunay->fMTried
void FindTriangles()
Find triangles in fDelaunay and initialise the TGraph2DPainter values needed to paint triangles or fi...
Int_t * fPTried
! Pointer to fDelaunay->fPTried
Double_t fXNmin
! Equal to fDelaunay->fXNmin
Double_t fXmax
! fGraph2D->fHistogram Xmax
Double_t * fEZhigh
! Pointer to fGraph2D->fZEhigh
Int_t fNpoints
! Equal to fGraph2D->fNpoints
Double_t fZmin
! fGraph2D->fHistogram Zmin
Int_t fNdt
! Equal to fDelaunay->fNdt
Double_t * fEYlow
! Pointer to fGraph2D->fYElow
TList * GetContourList(Double_t contour)
Returns the X and Y graphs building a contour.
Double_t * fX
! Pointer to fGraph2D->fX
Double_t fYmax
! fGraph2D->fHistogram Ymax
Double_t * fEXhigh
! Pointer to fGraph2D->fXEhigh
void PaintPolyMarker0(Int_t n, Double_t *x, Double_t *y)
Paints a circle at each vertex. Each circle background is white.
void PaintTriangles(Option_t *option)
Paints the 2D graph as triangles.
void PaintPolyMarker(Option_t *option)
Paints the 2D graph as PaintPolyMarker.
Double_t * fY
! Pointer to fGraph2D->fY
Double_t fYNmin
! Equal to fDelaunay->fYNmin
Double_t * fEXlow
! Pointer to fGraph2D->fXElow
void PaintPolyLine(Option_t *option)
Paints the 2D graph as PaintPolyLine.
Double_t * fEYhigh
! Pointer to fGraph2D->fYEhigh
TGraph2DPainter()
TGraph2DPainter default constructor.
Double_t fZmax
! fGraph2D->fHistogram Zmax
void PaintTriangles_old(Option_t *option)
Paints the 2D graph as triangles (old implementation)
void GetGraph2dProperties()
Protected method to get all TGraph2D properties.
~TGraph2DPainter() override
TGraph2DPainter destructor.
Double_t fYNmax
! Equal to fDelaunay->fYNmax
Double_t * fEZlow
! Pointer to fGraph2D->fZElow
Double_t fXmin
! fGraph2D->fHistogram Xmin
Double_t fXNmax
! Equal to fDelaunay->fXNmax
void PaintTriangles_new(Option_t *option)
Paints the 2D graph as triangles (new implementation)
virtual Double_t * GetEZhigh() const
Definition TGraph2D.h:133
Double_t GetMaximum() const
Definition TGraph2D.h:116
virtual Double_t GetZminE() const
Definition TGraph2D.h:145
Double_t GetMinimum() const
Definition TGraph2D.h:117
Double_t * GetY() const
Definition TGraph2D.h:123
Double_t * GetX() const
Definition TGraph2D.h:122
virtual Double_t GetZmaxE() const
Definition TGraph2D.h:144
virtual Double_t * GetEZ() const
Definition TGraph2D.h:127
Double_t GetZmax() const
Returns the Z maximum.
virtual Double_t * GetEYhigh() const
Definition TGraph2D.h:131
virtual Double_t * GetEY() const
Definition TGraph2D.h:126
Int_t GetN() const
Definition TGraph2D.h:121
virtual Double_t * GetEXhigh() const
Definition TGraph2D.h:129
virtual Double_t * GetEXlow() const
Definition TGraph2D.h:128
virtual Double_t * GetEX() const
Definition TGraph2D.h:125
virtual Double_t * GetEZlow() const
Definition TGraph2D.h:132
virtual Double_t * GetEYlow() const
Definition TGraph2D.h:130
Double_t * GetZ() const
Definition TGraph2D.h:124
TGraphDelaunay2D generates a Delaunay triangulation of a TGraph2D.
Triangles::const_iterator begin() const
Double_t GetXNmax() const
TGraph2D * GetGraph2D() const
Triangles::const_iterator end() const
Double_t GetXNmin() const
Int_t GetNdt() const
Double_t GetYNmax() const
Double_t GetYNmin() const
TGraphDelaunay generates a Delaunay triangulation of a TGraph2D.
Int_t GetNdt() const
Double_t GetYNmax() const
Int_t * GetMTried() const
Double_t GetXNmin() const
TGraph2D * GetGraph2D() const
Double_t GetXNmax() const
void FindAllTriangles()
Attempt to find all the Delaunay triangles of the point set.
Double_t * GetXN() const
Int_t * GetPTried() const
Int_t * GetNTried() const
Double_t GetYNmin() const
Double_t * GetYN() const
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
static TClass * Class()
virtual void SetPoint(Int_t i, Double_t x, Double_t y)
Set x and y values for point number i.
Definition TGraph.cxx:2290
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:109
TAxis * GetZaxis()
Definition TH1.h:574
virtual Double_t GetContourLevelPad(Int_t level) const
Return the value of contour number "level" in Pad coordinates.
Definition TH1.cxx:8468
@ kUserContour
User specified contour levels.
Definition TH1.h:405
TAxis * GetXaxis()
Definition TH1.h:572
virtual Double_t GetMaximum(Double_t maxval=FLT_MAX) const
Return maximum value smaller than maxval of bins in the range, unless the value has been overridden b...
Definition TH1.cxx:8573
TAxis * GetYaxis()
Definition TH1.h:573
virtual Int_t GetContour(Double_t *levels=nullptr)
Return contour values into array levels if pointer levels is non zero.
Definition TH1.cxx:8439
virtual void SetContour(Int_t nlevels, const Double_t *levels=nullptr)
Set the number and values of contour levels.
Definition TH1.cxx:8511
virtual Double_t GetMinimum(Double_t minval=-FLT_MAX) const
Return minimum value larger than minval of bins in the range, unless the value has been overridden by...
Definition TH1.cxx:8663
static void Optimize(Double_t A1, Double_t A2, Int_t nold, Double_t &BinLow, Double_t &BinHigh, Int_t &nbins, Double_t &BWID, Option_t *option="")
Static function to compute reasonable axis limits.
A doubly linked list.
Definition TList.h:38
R__ALWAYS_INLINE Bool_t TestBit(UInt_t f) const
Definition TObject.h:202
virtual void Delete(Option_t *option="")
Delete this object.
Definition TObject.cxx:267
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:1071
Basic string class.
Definition TString.h:138
void ToLower()
Change string to lower-case.
Definition TString.cxx:1189
Bool_t Contains(const char *pat, ECaseCompare cmp=kExact) const
Definition TString.h:640
Int_t GetColorPalette(Int_t i) const
Return color number i in current palette.
Definition TStyle.cxx:1102
Int_t GetNumberOfColors() const
Return number of colors in the color palette.
Definition TStyle.cxx:1176
Int_t GetNumberContours() const
Definition TStyle.h:243
See TView3D.
Definition TView.h:25
virtual Double_t * GetRmax()=0
virtual Double_t * GetRmin()=0
virtual Double_t GetLongitude()=0
virtual void WCtoNDC(const Float_t *pw, Float_t *pn)=0
TLine * line
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
Short_t Max(Short_t a, Short_t b)
Returns the largest of a and b.
Definition TMathBase.h:251
Short_t Min(Short_t a, Short_t b)
Returns the smallest of a and b.
Definition TMathBase.h:199
Double_t Cos(Double_t)
Returns the cosine of an angle of x radians.
Definition TMath.h:605
constexpr Double_t Pi()
Definition TMath.h:40
Double_t Sin(Double_t)
Returns the sine of an angle of x radians.
Definition TMath.h:599
void Sort(Index n, const Element *a, Index *index, Bool_t down=kTRUE)
Sort the n elements of the array a of generic templated type Element.
Definition TMathBase.h:432
Double_t Log10(Double_t x)
Returns the common (base-10) logarithm of x.
Definition TMath.h:773
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:124
Histograms' drawing options structure.
Definition Hoption.h:24
int Logx
log scale in X. Also set by histogram option
Definition Hoption.h:70
int Logz
log scale in Z. Also set by histogram option
Definition Hoption.h:72
int Logy
log scale in Y. Also set by histogram option
Definition Hoption.h:71
TMarker m
Definition textangle.C:8
TLine l
Definition textangle.C:4