#include "TEveTriangleSet.h"
#include "TEveRGBAPalette.h"
#include "TEveManager.h"
#include "TMath.h"
#include "TVector3.h"
#include "TRandom3.h"
#include "TVirtualPad.h"
#include "TVirtualViewer3D.h"
#include "TBuffer3D.h"
#include "TBuffer3DTypes.h"
ClassImp(TEveTriangleSet)
TEveTriangleSet::TEveTriangleSet(Int_t nv, Int_t nt, Bool_t norms, Bool_t cols) :
   TEveElement(fColor),
   TNamed("TEveTriangleSet", 0),
   fNVerts  (nv), fVerts(0),
   fNTrings (nt), fTrings(0), fTringNorms(0), fTringCols(0),
   fColor   (2),  fTransp(0),
   fHMTrans ()
{
   
   fVerts  = new Float_t[3*fNVerts];
   fTrings = new Int_t  [3*fNTrings];
   fTringNorms = (norms) ? new Float_t[3*fNTrings] : 0;
   fTringCols  = (cols)  ? new UChar_t[3*fNTrings] : 0;
}
TEveTriangleSet::~TEveTriangleSet()
{
   
   delete [] fVerts;
   delete [] fTrings;
   delete [] fTringNorms;
   delete [] fTringCols;
}
void TEveTriangleSet::GenerateTriangleNormals()
{
   
   if (fTringNorms == 0)  fTringNorms = new Float_t[3*fNTrings];
   TVector3 e1, e2, n;
   Float_t *norm = fTringNorms;
   Int_t   *tring  = fTrings;
   for(Int_t t=0; t<fNTrings; ++t, norm+=3, tring+=3)
   {
      Float_t* v0 = Vertex(tring[0]);
      Float_t* v1 = Vertex(tring[1]);
      Float_t* v2 = Vertex(tring[2]);
      e1.SetXYZ(v1[0]-v0[0], v1[1]-v0[1], v1[2]-v0[2]);
      e2.SetXYZ(v2[0]-v0[0], v2[1]-v0[1], v2[2]-v0[2]);
      n = e1.Cross(e2);
      n.SetMag(1);
      n.GetXYZ(norm);
   }
}
void TEveTriangleSet::GenerateRandomColors()
{
   
   if (fTringCols == 0)  fTringCols = new UChar_t[3*fNTrings];
   TRandom r;
   r.SetSeed(0);
   UChar_t *col = fTringCols;
   for(Int_t t=0; t<fNTrings; ++t, col+=3)
   {
      col[0] = (UChar_t) r.Uniform(60, 255);
      col[1] = (UChar_t) r.Uniform(60, 255);
      col[2] = (UChar_t) r.Uniform(60, 255);
   }
}
void TEveTriangleSet::GenerateZNormalColors(Float_t fac, Int_t min, Int_t max,
                                            Bool_t interp, Bool_t wrap)
{
   
   
   if (fTringCols  == 0)  fTringCols = new UChar_t[3*fNTrings];
   if (fTringNorms == 0)  GenerateTriangleNormals();
   TEveRGBAPalette pal(min, max, interp, wrap);
   UChar_t *col = fTringCols;
   Float_t *norm = fTringNorms;
   for(Int_t t=0; t<fNTrings; ++t, col+=3, norm+=3)
   {
      Int_t v = TMath::Nint(fac * norm[2]);
      pal.ColorFromValue(v, col, kFALSE);
   }
   gEve->Redraw3D();
}
void TEveTriangleSet::ComputeBBox()
{
   
   
   if (fNVerts <= 0) {
      BBoxZero();
      return;
   }
   BBoxInit();
   Float_t* v = fVerts;
   for (Int_t i=0; i<fNVerts; ++i, v += 3)
      BBoxCheckPoint(v);
}
void TEveTriangleSet::Paint(Option_t* )
{
   
   TBuffer3D buffer(TBuffer3DTypes::kGeneric);
   
   buffer.fID           = this;
   buffer.fColor        = fColor;
   buffer.fTransparency = fTransp;
   fHMTrans.SetBuffer3D(buffer);
   buffer.SetSectionsValid(TBuffer3D::kCore);
   
   Int_t reqSections = gPad->GetViewer3D()->AddObject(buffer);
   if (reqSections == TBuffer3D::kNone) {
      return;
   }
   Error("TEveTriangleSet::Paint", "only direct OpenGL rendering supported.");
}
#include <stdio.h>
TEveTriangleSet* TEveTriangleSet::ReadTrivialFile(const char* file)
{
   
   FILE* f = fopen(file, "r");
   if (f == 0) {
      ::Error("TEveTriangleSet::ReadTrivialFile", Form("file '%s' not found.", file));
      return 0;
   }
   Int_t nv, nt;
   fscanf(f, "%d %d", &nv, &nt);
   TEveTriangleSet* ts = new TEveTriangleSet(nv, nt);
   Float_t *vtx = ts->Vertex(0);
   for (Int_t i=0; i<nv; ++i, vtx+=3) {
      fscanf(f, "%f %f %f", &vtx[0], &vtx[1], &vtx[2]);
   }
   Int_t *tngl = ts->Triangle(0);
   for (Int_t i=0; i<nt; ++i, tngl+=3) {
      fscanf(f, "%d %d %d", &tngl[0], &tngl[1], &tngl[2]);
   }
   fclose(f);
   return ts;
}
Last update: Thu Jan 17 08:49:42 2008
This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.