Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TEveTriangleSet.cxx
Go to the documentation of this file.
1// @(#)root/eve:$Id$
2// Authors: Matevz Tadel & Alja Mrak-Tadel: 2006, 2007
3
4/*************************************************************************
5 * Copyright (C) 1995-2007, 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 "TEveTriangleSet.h"
13#include "TEveRGBAPalette.h"
14#include "TEveManager.h"
15
16#include "TMath.h"
17#include "TVector3.h"
18#include "TRandom3.h"
19
20
21/** \class TEveTriangleSet
22\ingroup TEve
23Made from a list of vertices and a list of triangles (triplets of
24vertex indices).
25
26If input is composed from triangles with direct vertex coordinates
27one should consider finding all occurrences of the same vertex
28and specifying it only once.
29*/
30
32
33////////////////////////////////////////////////////////////////////////////////
34/// Constructor.
35
37 TEveElementList("TEveTriangleSet", "", kTRUE),
38 fNVerts (nv), fVerts(nullptr),
39 fNTrings (nt), fTrings(nullptr), fTringNorms(nullptr), fTringCols(nullptr)
40{
42
43 fVerts = new Float_t[3*fNVerts];
44 fTrings = new Int_t [3*fNTrings];
45 fTringNorms = (norms) ? new Float_t[3*fNTrings] : nullptr;
46 fTringCols = (cols) ? new UChar_t[3*fNTrings] : nullptr;
47}
48
49////////////////////////////////////////////////////////////////////////////////
50/// Destructor.
51
53{
54 delete [] fVerts;
55 delete [] fTrings;
56 delete [] fTringNorms;
57 delete [] fTringCols;
58}
59
60////////////////////////////////////////////////////////////////////////////////
61/// Generate triangle normals via cross product of triangle edges.
62
64{
65 if (fTringNorms == nullptr) fTringNorms = new Float_t[3*fNTrings];
66
67 TVector3 e1, e2, n;
68 Float_t *norm = fTringNorms;
69 Int_t *tring = fTrings;
70 for(Int_t t=0; t<fNTrings; ++t, norm+=3, tring+=3)
71 {
72 Float_t* v0 = Vertex(tring[0]);
73 Float_t* v1 = Vertex(tring[1]);
74 Float_t* v2 = Vertex(tring[2]);
75 e1.SetXYZ(v1[0]-v0[0], v1[1]-v0[1], v1[2]-v0[2]);
76 e2.SetXYZ(v2[0]-v0[0], v2[1]-v0[1], v2[2]-v0[2]);
77 n = e1.Cross(e2);
78 n.SetMag(1);
79 n.GetXYZ(norm);
80 }
81}
82
83////////////////////////////////////////////////////////////////////////////////
84/// Assign random colors to all triangles.
85
87{
88 if (fTringCols == nullptr) fTringCols = new UChar_t[3*fNTrings];
89
90 TRandom r;
91 r.SetSeed(0);
92 UChar_t *col = fTringCols;
93 for(Int_t t=0; t<fNTrings; ++t, col+=3)
94 {
95 col[0] = (UChar_t) r.Uniform(60, 255);
96 col[1] = (UChar_t) r.Uniform(60, 255);
97 col[2] = (UChar_t) r.Uniform(60, 255);
98 }
99}
100
101////////////////////////////////////////////////////////////////////////////////
102/// Generate triangle colors by the z-component of the normal.
103/// Current palette is taken from gStyle.
104
106 Bool_t interp, Bool_t wrap)
107{
108 if (fTringCols == nullptr) fTringCols = new UChar_t[3*fNTrings];
109 if (fTringNorms == nullptr) GenerateTriangleNormals();
110
111 TEveRGBAPalette pal(min, max, interp, wrap);
112 UChar_t *col = fTringCols;
113 Float_t *norm = fTringNorms;
114 for(Int_t t=0; t<fNTrings; ++t, col+=3, norm+=3)
115 {
116 Int_t v = TMath::Nint(fac * norm[2]);
117 pal.ColorFromValue(v, col, kFALSE);
118 }
119 gEve->Redraw3D();
120}
121
122////////////////////////////////////////////////////////////////////////////////
123/// Compute bounding box.
124/// Virtual from TAttBBox.
125
127{
128 if (fNVerts <= 0) {
129 BBoxZero();
130 return;
131 }
132
133 BBoxInit();
134 Float_t* v = fVerts;
135 for (Int_t i=0; i<fNVerts; ++i, v += 3)
137}
138
139////////////////////////////////////////////////////////////////////////////////
140/// Paint this object. Only direct rendering is supported.
141
143{
144 PaintStandard(this);
145}
146
147////////////////////////////////////////////////////////////////////////////////
148/// Read a simple ascii input file describing vertices and triangles.
149
151{
152 static const TEveException kEH("TEveTriangleSet::ReadTrivialFile ");
153
154 FILE* f = fopen(file, "r");
155 if (f == nullptr) {
156 ::Error(kEH, "file '%s' not found.", file);
157 return nullptr;
158 }
159
160 Int_t nv, nt;
161 if (fscanf(f, "%d %d", &nv, &nt) != 2) {
162 fclose(f);
163 throw kEH + "Reading nv, nt failed.";
164 }
165
166 if (nv < 0 || nt < 0) {
167 fclose(f);
168 throw kEH + "Negative number of vertices / triangles specified.";
169 }
170
171 TEveTriangleSet* ts = new TEveTriangleSet(nv, nt);
172
173 Float_t *vtx = ts->Vertex(0);
174 for (Int_t i=0; i<nv; ++i, vtx+=3) {
175 if (fscanf(f, "%f %f %f", &vtx[0], &vtx[1], &vtx[2]) != 3) {
176 fclose(f);
177 throw kEH + TString::Format("Reading vertex data %d failed.", i);
178 }
179 }
180
181 Int_t *tngl = ts->Triangle(0);
182 for (Int_t i=0; i<nt; ++i, tngl+=3) {
183 if (fscanf(f, "%d %d %d", &tngl[0], &tngl[1], &tngl[2]) != 3) {
184 fclose(f);
185 throw kEH + TString::Format("Reading triangle data %d failed.", i);
186 }
187 }
188
189 fclose(f);
190
191 return ts;
192}
#define f(i)
Definition RSha256.hxx:104
unsigned char UChar_t
Definition RtypesCore.h:38
float Float_t
Definition RtypesCore.h:57
constexpr Bool_t kFALSE
Definition RtypesCore.h:94
constexpr Bool_t kTRUE
Definition RtypesCore.h:93
const char Option_t
Definition RtypesCore.h:66
#define ClassImp(name)
Definition Rtypes.h:382
R__EXTERN TEveManager * gEve
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
void BBoxCheckPoint(Float_t x, Float_t y, Float_t z)
Definition TAttBBox.h:69
void BBoxZero(Float_t epsilon=0, Float_t x=0, Float_t y=0, Float_t z=0)
Create cube of volume (2*epsilon)^3 at (x,y,z).
Definition TAttBBox.cxx:42
void BBoxInit(Float_t infinity=1e6)
Dynamic Float_t[6] X(min,max), Y(min,max), Z(min,max)
Definition TAttBBox.cxx:29
A list of TEveElements.
virtual void InitMainTrans(Bool_t can_edit=kTRUE)
Initialize the main transformation to identity matrix.
virtual void PaintStandard(TObject *id)
Paint object – a generic implementation for EVE elements.
Exception class thrown by TEve classes and macros.
Definition TEveUtil.h:102
void Redraw3D(Bool_t resetCameras=kFALSE, Bool_t dropLogicals=kFALSE)
A generic, speed-optimised mapping from value to RGBA color supporting different wrapping and range t...
const UChar_t * ColorFromValue(Int_t val) const
Made from a list of vertices and a list of triangles (triplets of vertex indices).
Float_t * fTringNorms
void ComputeBBox() override
Compute bounding box.
void GenerateZNormalColors(Float_t fac=20, Int_t min=-20, Int_t max=20, Bool_t interp=kFALSE, Bool_t wrap=kFALSE)
Generate triangle colors by the z-component of the normal.
void GenerateRandomColors()
Assign random colors to all triangles.
Float_t * Vertex(Int_t i)
~TEveTriangleSet() override
Destructor.
static TEveTriangleSet * ReadTrivialFile(const char *file)
Read a simple ascii input file describing vertices and triangles.
TEveTriangleSet(const TEveTriangleSet &)
void GenerateTriangleNormals()
Generate triangle normals via cross product of triangle edges.
void Paint(Option_t *option="") override
Paint this object. Only direct rendering is supported.
Int_t * Triangle(Int_t i)
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition TObject.cxx:993
This is the base class for the ROOT Random number generators.
Definition TRandom.h:27
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2378
void SetXYZ(Double_t x, Double_t y, Double_t z)
Definition TVector3.h:227
void SetMag(Double_t)
Definition TVector3.h:376
TVector3 Cross(const TVector3 &) const
Definition TVector3.h:335
const Int_t n
Definition legend1.C:16
Int_t Nint(T x)
Round to nearest integer. Rounds half integers to the nearest even integer.
Definition TMath.h:693