Logo ROOT   6.10/09
Reference Guide
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
23 Made from a list of vertices and a list of triangles (triplets of
24 vertex indices).
25 
26 If input is composed from triangles with direct vertex coordinates
27 one should consider finding all occurrences of the same vertex
28 and specifying it only once.
29 */
30 
32 
33 ////////////////////////////////////////////////////////////////////////////////
34 /// Constructor.
35 
37  TEveElementList("TEveTriangleSet", "", kTRUE),
38  fNVerts (nv), fVerts(0),
39  fNTrings (nt), fTrings(0), fTringNorms(0), fTringCols(0)
40 {
41  InitMainTrans();
42 
43  fVerts = new Float_t[3*fNVerts];
44  fTrings = new Int_t [3*fNTrings];
45  fTringNorms = (norms) ? new Float_t[3*fNTrings] : 0;
46  fTringCols = (cols) ? new UChar_t[3*fNTrings] : 0;
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 == 0) fTringNorms = new Float_t[3*fNTrings];
66 
67  TVector3 e1, e2, n;
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 == 0) 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 == 0) fTringCols = new UChar_t[3*fNTrings];
110 
111  TEveRGBAPalette pal(min, max, interp, wrap);
112  UChar_t *col = fTringCols;
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)
136  BBoxCheckPoint(v);
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 == 0) {
156  ::Error(kEH, "file '%s' not found.", file);
157  return 0;
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 }
A generic, speed-optimised mapping from value to RGBA color supporting different wrapping and range t...
~TEveTriangleSet()
Destructor.
Float_t * fTringNorms
float Float_t
Definition: RtypesCore.h:53
const char Option_t
Definition: RtypesCore.h:62
Int_t * Triangle(Int_t i)
void GenerateTriangleNormals()
Generate triangle normals via cross product of triangle edges.
virtual void PaintStandard(TObject *id)
Paint object – a generic implementation for EVE elements.
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
void Redraw3D(Bool_t resetCameras=kFALSE, Bool_t dropLogicals=kFALSE)
Definition: TEveManager.h:168
A list of TEveElements.
Definition: TEveElement.h:459
virtual void Paint(Option_t *option="")
Paint this object. Only direct rendering is supported.
void BBoxCheckPoint(Float_t x, Float_t y, Float_t z)
Definition: TAttBBox.h:58
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:2345
This is the base class for the ROOT Random number generators.
Definition: TRandom.h:27
virtual void SetSeed(ULong_t seed=0)
Set the random generator seed.
Definition: TRandom.cxx:568
void SetXYZ(Double_t x, Double_t y, Double_t z)
Definition: TVector3.h:227
void GetXYZ(Double_t *carray) const
Definition: TVector3.h:233
const UChar_t * ColorFromValue(Int_t val) const
SEXP wrap(const TString &s)
Definition: RExports.h:67
TVector3 is a general three vector class, which can be used for the description of different vectors ...
Definition: TVector3.h:22
TVector3 Cross(const TVector3 &) const
Definition: TVector3.h:335
virtual void InitMainTrans(Bool_t can_edit=kTRUE)
Initialize the main transformation to identity matrix.
TRandom2 r(17)
SVector< double, 2 > v
Definition: Dict.h:5
R__EXTERN TEveManager * gEve
Definition: TEveManager.h:243
virtual void ComputeBBox()
Compute bounding box.
virtual void Error(const char *method, const char *msgfmt,...) const
Issue error message.
Definition: TObject.cxx:873
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
const Bool_t kFALSE
Definition: RtypesCore.h:92
#define ClassImp(name)
Definition: Rtypes.h:336
double f(double x)
Made from a list of vertices and a list of triangles (triplets of vertex indices).
virtual Double_t Uniform(Double_t x1=1)
Returns a uniform deviate on the interval (0, x1).
Definition: TRandom.cxx:606
Float_t * Vertex(Int_t i)
static TEveTriangleSet * ReadTrivialFile(const char *file)
Read a simple ascii input file describing vertices and triangles.
Definition: file.py:1
void SetMag(Double_t)
Definition: TVector3.h:376
Exception class thrown by TEve classes and macros.
Definition: TEveUtil.h:102
void GenerateRandomColors()
Assign random colors to all triangles.
unsigned char UChar_t
Definition: RtypesCore.h:34
TEveTriangleSet(const TEveTriangleSet &)
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.
const Bool_t kTRUE
Definition: RtypesCore.h:91
double norm(double *x, double *p)
Definition: unuranDistr.cxx:40
Int_t Nint(T x)
Definition: TMath.h:607
const Int_t n
Definition: legend1.C:16
void BBoxInit(Float_t infinity=1e6)
Dynamic Float_t[6] X(min,max), Y(min,max), Z(min,max)
Definition: TAttBBox.cxx:29
UChar_t * fTringCols