00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012 #include "TEveTriangleSet.h"
00013 #include "TEveRGBAPalette.h"
00014 #include "TEveManager.h"
00015
00016 #include "TMath.h"
00017 #include "TVector3.h"
00018 #include "TRandom3.h"
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031 ClassImp(TEveTriangleSet);
00032
00033
00034 TEveTriangleSet::TEveTriangleSet(Int_t nv, Int_t nt, Bool_t norms, Bool_t cols) :
00035 TEveElementList("TEveTriangleSet", "", kTRUE),
00036 fNVerts (nv), fVerts(0),
00037 fNTrings (nt), fTrings(0), fTringNorms(0), fTringCols(0)
00038 {
00039
00040
00041 InitMainTrans();
00042
00043 fVerts = new Float_t[3*fNVerts];
00044 fTrings = new Int_t [3*fNTrings];
00045 fTringNorms = (norms) ? new Float_t[3*fNTrings] : 0;
00046 fTringCols = (cols) ? new UChar_t[3*fNTrings] : 0;
00047 }
00048
00049
00050 TEveTriangleSet::~TEveTriangleSet()
00051 {
00052
00053
00054 delete [] fVerts;
00055 delete [] fTrings;
00056 delete [] fTringNorms;
00057 delete [] fTringCols;
00058 }
00059
00060
00061
00062
00063 void TEveTriangleSet::GenerateTriangleNormals()
00064 {
00065
00066
00067 if (fTringNorms == 0) fTringNorms = new Float_t[3*fNTrings];
00068
00069 TVector3 e1, e2, n;
00070 Float_t *norm = fTringNorms;
00071 Int_t *tring = fTrings;
00072 for(Int_t t=0; t<fNTrings; ++t, norm+=3, tring+=3)
00073 {
00074 Float_t* v0 = Vertex(tring[0]);
00075 Float_t* v1 = Vertex(tring[1]);
00076 Float_t* v2 = Vertex(tring[2]);
00077 e1.SetXYZ(v1[0]-v0[0], v1[1]-v0[1], v1[2]-v0[2]);
00078 e2.SetXYZ(v2[0]-v0[0], v2[1]-v0[1], v2[2]-v0[2]);
00079 n = e1.Cross(e2);
00080 n.SetMag(1);
00081 n.GetXYZ(norm);
00082 }
00083 }
00084
00085
00086 void TEveTriangleSet::GenerateRandomColors()
00087 {
00088
00089
00090 if (fTringCols == 0) fTringCols = new UChar_t[3*fNTrings];
00091
00092 TRandom r;
00093 r.SetSeed(0);
00094 UChar_t *col = fTringCols;
00095 for(Int_t t=0; t<fNTrings; ++t, col+=3)
00096 {
00097 col[0] = (UChar_t) r.Uniform(60, 255);
00098 col[1] = (UChar_t) r.Uniform(60, 255);
00099 col[2] = (UChar_t) r.Uniform(60, 255);
00100 }
00101 }
00102
00103
00104 void TEveTriangleSet::GenerateZNormalColors(Float_t fac, Int_t min, Int_t max,
00105 Bool_t interp, Bool_t wrap)
00106 {
00107
00108
00109
00110 if (fTringCols == 0) fTringCols = new UChar_t[3*fNTrings];
00111 if (fTringNorms == 0) GenerateTriangleNormals();
00112
00113 TEveRGBAPalette pal(min, max, interp, wrap);
00114 UChar_t *col = fTringCols;
00115 Float_t *norm = fTringNorms;
00116 for(Int_t t=0; t<fNTrings; ++t, col+=3, norm+=3)
00117 {
00118 Int_t v = TMath::Nint(fac * norm[2]);
00119 pal.ColorFromValue(v, col, kFALSE);
00120 }
00121 gEve->Redraw3D();
00122 }
00123
00124
00125
00126
00127 void TEveTriangleSet::ComputeBBox()
00128 {
00129
00130
00131
00132 if (fNVerts <= 0) {
00133 BBoxZero();
00134 return;
00135 }
00136
00137 BBoxInit();
00138 Float_t* v = fVerts;
00139 for (Int_t i=0; i<fNVerts; ++i, v += 3)
00140 BBoxCheckPoint(v);
00141 }
00142
00143
00144 void TEveTriangleSet::Paint(Option_t*)
00145 {
00146
00147
00148 PaintStandard(this);
00149 }
00150
00151
00152
00153
00154 TEveTriangleSet* TEveTriangleSet::ReadTrivialFile(const char* file)
00155 {
00156
00157
00158 static const TEveException kEH("TEveTriangleSet::ReadTrivialFile ");
00159
00160 FILE* f = fopen(file, "r");
00161 if (f == 0) {
00162 ::Error(kEH, "file '%s' not found.", file);
00163 return 0;
00164 }
00165
00166 Int_t nv, nt;
00167 if (fscanf(f, "%d %d", &nv, &nt) != 2)
00168 throw kEH + "Reading nv, nt failed.";
00169 if (nv < 0 || nt < 0)
00170 throw kEH + "Negative number of vertices / triangles specified.";
00171
00172
00173 TEveTriangleSet* ts = new TEveTriangleSet(nv, nt);
00174
00175 Float_t *vtx = ts->Vertex(0);
00176 for (Int_t i=0; i<nv; ++i, vtx+=3) {
00177 if (fscanf(f, "%f %f %f", &vtx[0], &vtx[1], &vtx[2]) != 3)
00178 throw kEH + TString::Format("Reading vertex data %d failed.", i);
00179 }
00180
00181 Int_t *tngl = ts->Triangle(0);
00182 for (Int_t i=0; i<nt; ++i, tngl+=3) {
00183 if (fscanf(f, "%d %d %d", &tngl[0], &tngl[1], &tngl[2]) != 3)
00184 throw kEH + TString::Format("Reading triangle data %d failed.", i);
00185 }
00186
00187 fclose(f);
00188
00189 return ts;
00190 }