swh:1:snp:af87cd67498ef4fe47c76ed3e7caffe5b61facaf
Raw File
Tip revision: be74382292a29ef72a13d39a6eef7316572e63d3 authored by Fons Rademakers on 24 November 2011, 15:43:23 UTC
tag patch release v5-28-00h.
Tip revision: be74382
quadset.C
// @(#)root/eve:$Id$
// Author: Matevz Tadel

// Demonstates usage of 2D digit class TEveQuadSet.

TEveQuadSet* quadset(Float_t x=0, Float_t y=0, Float_t z=0,
                       Int_t num=100, Bool_t register=kTRUE)
{
   TEveManager::Create();

   TRandom r(0);
   gStyle->SetPalette(1, 0);

   TEveRGBAPalette *pal = new TEveRGBAPalette(0, 130);
   TEveFrameBox    *box = new TEveFrameBox();
   box->SetAAQuadXY(-10, -10, 0, 20, 20);
   box->SetFrameColor(kGray);

   TEveQuadSet* q = new TEveQuadSet("RectangleXY");
   q->SetOwnIds(kTRUE);
   q->SetPalette(pal);
   q->SetFrame(box);
   q->Reset(TEveQuadSet::kQT_RectangleXY, kFALSE, 32);
   for (Int_t i=0; i<num; ++i)
   {
      q->AddQuad(r.Uniform(-10, 9), r.Uniform(-10, 9), 0,
                 r.Uniform(0.2, 1), r.Uniform(0.2, 1));
      q->QuadValue(r.Uniform(0, 130));
      q->QuadId(new TNamed(Form("QuadIdx %d", i),
                           "TNamed assigned to a quad as an indentifier."));
   }
   q->RefitPlex();

   TEveTrans& t = q->RefMainTrans();
   t.RotateLF(1, 3, 0.5*TMath::Pi());
   t.SetPos(x, y, z);

   TGLViewer* v = gEve->GetDefaultGLViewer();
   v->SetCurrentCamera(TGLViewer::kCameraOrthoZOY);
   TGLCameraOverlay* co = v->GetCameraOverlay();
   co->SetShowOrthographic(kTRUE);
   co->SetOrthographicMode(TGLCameraOverlay::kGridFront);

   // Uncomment these two lines to get internal highlight / selection.
   // q->SetPickable(1);
   // q->SetAlwaysSecSelect(1);

   if (register)
   {
      gEve->AddElement(q);
      gEve->Redraw3D(kTRUE);
   }

   Info("quadset", "use alt-left-mouse to select individual digits.");

   return q;
}

TEveQuadSet* quadset_emc(Float_t x=0, Float_t y=0, Float_t z=0,
                           Int_t num=100)
{
   TEveManager::Create();

   TRandom r(0);
   gStyle->SetPalette(1, 0);

   TEveQuadSet* q = new TEveQuadSet("EMC Supermodule");
   q->SetOwnIds(kTRUE);
   q->Reset(TEveQuadSet::kQT_RectangleXZFixedDimY, kFALSE, 32);
   q->SetDefWidth(8);
   q->SetDefHeight(8);

   for (Int_t i=0; i<num; ++i)
   {
      q->AddQuad(r.Uniform(-100, 100), r.Uniform(-100, 100));
      q->QuadValue(r.Uniform(0, 130));
      q->QuadId(new TNamed(Form("Cell %d", i), "Dong!"));
   }
   q->RefitPlex();

   TEveTrans& t = q->RefMainTrans();
   t.SetPos(x, y, z);

   gEve->AddElement(q);
   gEve->Redraw3D();

   return q;
}

TEveQuadSet* quadset_circ()
{
   TEveManager::Create();

   TRandom rnd(0);
   gStyle->SetPalette(1, 0);

   Float_t R = 10, dW = 1, dH = .5;

   TEveFrameBox *box = new TEveFrameBox();
   {
      Float_t  frame[3*36];
      Float_t *p = frame;
      for (Int_t i = 0; i < 36; ++i, p += 3) {
         p[0] = 11 * TMath::Cos(TMath::TwoPi()*i/36);
         p[1] = 11 * TMath::Sin(TMath::TwoPi()*i/36);
         p[2] = 0;
      }
      box->SetQuadByPoints(frame, 36);
   }
   box->SetFrameColor(kGray);

   TEveQuadSet* q = new TEveQuadSet("Pepe");
   q->SetFrame(box);
   q->Reset(TEveQuadSet::kQT_HexagonXY, kFALSE, 32);

   for (Float_t r = R; r > 2; r *= 0.8)
   {
      Int_t maxI = 2.0*TMath::Pi()*r / 2;
      for (Int_t i = 0; i < maxI; ++i)
      {
         Float_t x = r * TMath::Cos(TMath::TwoPi()*i/maxI);
         Float_t y = r * TMath::Sin(TMath::TwoPi()*i/maxI);
         q->AddHexagon(x, y, rnd.Uniform(-1, 1), rnd.Uniform(0.2, 1));
         q->QuadValue(rnd.Uniform(0, 130));
      }
   }
   q->RefitPlex();

   TEveTrans& t = q->RefMainTrans();
   t.RotateLF(1, 3, 0.5*TMath::Pi());

   gEve->AddElement(q);
   gEve->Redraw3D();

   return q;
}

TEveQuadSet* quadset_hex(Float_t x=0, Float_t y=0, Float_t z=0,
                           Int_t num=100, Bool_t register=kTRUE)
{
   TEveManager::Create();

   TRandom r(0);
   gStyle->SetPalette(1, 0);

   {
      TEveQuadSet* q = new TEveQuadSet("HexagonXY");
      q->Reset(TEveQuadSet::kQT_HexagonXY, kFALSE, 32);
      for (Int_t i=0; i<num; ++i)
      {
         q->AddHexagon(r.Uniform(-10, 10),
                       r.Uniform(-10, 10),
                       r.Uniform(-10, 10),
                       r.Uniform(0.2, 1));
         q->QuadValue(r.Uniform(0, 120));
      }
      q->RefitPlex();

      TEveTrans& t = q->RefMainTrans();
      t.SetPos(x, y, z);

      if (register)
      {
         gEve->AddElement(q);
         gEve->Redraw3D();
      }
   }

   {
      TEveQuadSet* q = new TEveQuadSet("HexagonYX");
      q->Reset(TEveQuadSet::kQT_HexagonYX, kFALSE, 32);
      for (Int_t i=0; i<num; ++i)
      {
         q->AddHexagon(r.Uniform(-10, 10),
                       r.Uniform(-10, 10),
                       r.Uniform(-10, 10),
                       r.Uniform(0.2, 1));
         q->QuadValue(r.Uniform(0, 120));
      }
      q->RefitPlex();

      TEveTrans& t = q->RefMainTrans();
      t.SetPos(x, y, z);

      if (register)
      {
         gEve->AddElement(q);
         gEve->Redraw3D();
      }
   }

   return q;
}

TEveQuadSet* quadset_hexid(Float_t x=0, Float_t y=0, Float_t z=0,
                             Int_t num=100, Bool_t register=kTRUE)
{
   TEveManager::Create();

   TRandom r(0);
   gStyle->SetPalette(1, 0);

   {
      TEveQuadSet* q = new TEveQuadSet("HexagonXY");
      q->SetOwnIds(kTRUE);
      q->Reset(TEveQuadSet::kQT_HexagonXY, kFALSE, 32);
      for (Int_t i=0; i<num; ++i)
      {
         q->AddHexagon(r.Uniform(-10, 10),
                       r.Uniform(-10, 10),
                       r.Uniform(-10, 10),
                       r.Uniform(0.2, 1));
         q->QuadValue(r.Uniform(0, 120));
         q->QuadId(new TNamed(Form("Quad with idx=%d", i),
                              "This title is not confusing."));
      }
      q->RefitPlex();

      TEveTrans& t = q->RefMainTrans();
      t.SetPos(x, y, z);

      if (register)
      {
         gEve->AddElement(q);
         gEve->Redraw3D();
      }
   }

   // This show another way of getting notified about
   // secondary selection hit. The callback function and the
   // setting of it must be done in compiled code.
   gROOT->ProcessLine(".L quadset_callback.cxx+");
   quadset_set_callback(q);

   return q;
}

void quadset_hierarchy(Int_t n=4)
{
   TEveManager::Create();

   gStyle->SetPalette(1, 0);

   TEveRGBAPalette* pal = new TEveRGBAPalette(20, 100);
   pal->SetLimits(0, 120);

   TEveFrameBox* box = new TEveFrameBox();
   box->SetAABox(-10, -10, -10, 20, 20, 20);
   box->SetFrameColor(33);

   TEveElementList* l = new TEveElementList("Parent/Dir");
   l->SetTitle("Tooltip");
   //  l->SetMainColor(3);
   gEve->AddElement(l);

   for (Int_t i=0; i<n; ++i)
   {
      TEveQuadSet* qs = quadset_hexid(0, 0, 50*i, 50, kFALSE);
      qs->SetPalette(pal);
      qs->SetFrame(box);
      l->AddElement(qs);
   }

   gEve->Redraw3D();
}
back to top