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