ROOT logo
/**************************************************************************
 * Copyright(c) 1998-2008, ALICE Experiment at CERN, all rights reserved. *
 * See http://aliceinfo.cern.ch/Offline/AliRoot/License.html for          *
 * full copyright notice.                                                 *
 **************************************************************************/

#if !defined(__CINT__) || defined(__MAKECINT__)
#include <TClonesArray.h>
#include <TBranch.h>
#include <TGeoMatrix.h>
#include <TTree.h>
#include <TStyle.h>
#include <TEveElement.h>
#include <TEveFrameBox.h>
#include <TEveManager.h>
#include <TEvePointSet.h>
#include <TEveRGBAPalette.h>
#include <TEveTrans.h>
#include <TEveQuadSet.h>

#include <AliHMPIDDigit.h>
#include <AliHMPIDv3.h>
#include <AliCluster3D.h>
#include <AliRunLoader.h>
#include <AliEveEventManager.h>
#endif

void hmpid_digits()
{
  const Char_t *name[]={ "HMPID0", "HMPID1", "HMPID2", "HMPID3",
			 "HMPID4", "HMPID5", "HMPID6" };

  AliRunLoader* rl = AliEveEventManager::AssertRunLoader();
  rl->LoadDigits("HMPID");

  TTree *dTree = rl->GetTreeD("HMPID", kFALSE);
  if (!dTree) return;

  TEveElementList* list = new TEveElementList("HMPID Digits");
  gEve->AddElement(list);

  gStyle->SetPalette(1, 0);

  TEveRGBAPalette *pal = new TEveRGBAPalette(0, 3000);
  pal->SetMax(1000);
  TEveFrameBox    *box = new TEveFrameBox();
  box->SetAAQuadXY(0, 0, 0, 144, 121);
  box->SetFrameColor(kGray);

  TClonesArray* digits = new TClonesArray("AliHMPIDDigit");
  for (Int_t iCh = 0; iCh < 7; ++iCh)
  {
    TBranch *br = dTree->GetBranch(name[iCh]);
    br->SetAddress(&digits);
    br->GetEntry(0);

    TEveQuadSet* q = new TEveQuadSet(Form("Chamber %d", iCh));
    q->SetOwnIds(kTRUE);
    q->SetPalette(pal);
    q->SetFrame(box);
    q->SetAntiFlick(kTRUE);
    q->SetPickable(kTRUE);

    q->Reset(TEveQuadSet::kQT_RectangleXYFixedDimZ, kFALSE, 64);
    q->SetDefCoord(0);
    q->SetDefHeight(0.84f);
    q->SetDefWidth(0.8f);

    for(Int_t iDig = 0; iDig < digits->GetEntriesFast(); ++iDig)
    {
      AliHMPIDDigit *pDig = (AliHMPIDDigit*) digits->At(iDig);

      q->AddQuad(pDig->PadChX()*0.8f,  pDig->PadChY()*0.84f);
      q->QuadValue(TMath::Nint(pDig->Q()));
      q->QuadId(new AliHMPIDDigit(*pDig));
    }

    q->RefitPlex();

    TGeoHMatrix mat;
    AliHMPIDv3::IdealPosition(iCh, &mat);
    q->RefMainTrans().SetFrom(mat);
    q->RefMainTrans().Move3LF(-0.5f*144, -0.5f*121, 0);

    list->AddElement(q);
  }

  delete digits;
  rl->UnloadDigits("HMPID");

  gEve->Redraw3D();
}
 hmpid_digits.C:1
 hmpid_digits.C:2
 hmpid_digits.C:3
 hmpid_digits.C:4
 hmpid_digits.C:5
 hmpid_digits.C:6
 hmpid_digits.C:7
 hmpid_digits.C:8
 hmpid_digits.C:9
 hmpid_digits.C:10
 hmpid_digits.C:11
 hmpid_digits.C:12
 hmpid_digits.C:13
 hmpid_digits.C:14
 hmpid_digits.C:15
 hmpid_digits.C:16
 hmpid_digits.C:17
 hmpid_digits.C:18
 hmpid_digits.C:19
 hmpid_digits.C:20
 hmpid_digits.C:21
 hmpid_digits.C:22
 hmpid_digits.C:23
 hmpid_digits.C:24
 hmpid_digits.C:25
 hmpid_digits.C:26
 hmpid_digits.C:27
 hmpid_digits.C:28
 hmpid_digits.C:29
 hmpid_digits.C:30
 hmpid_digits.C:31
 hmpid_digits.C:32
 hmpid_digits.C:33
 hmpid_digits.C:34
 hmpid_digits.C:35
 hmpid_digits.C:36
 hmpid_digits.C:37
 hmpid_digits.C:38
 hmpid_digits.C:39
 hmpid_digits.C:40
 hmpid_digits.C:41
 hmpid_digits.C:42
 hmpid_digits.C:43
 hmpid_digits.C:44
 hmpid_digits.C:45
 hmpid_digits.C:46
 hmpid_digits.C:47
 hmpid_digits.C:48
 hmpid_digits.C:49
 hmpid_digits.C:50
 hmpid_digits.C:51
 hmpid_digits.C:52
 hmpid_digits.C:53
 hmpid_digits.C:54
 hmpid_digits.C:55
 hmpid_digits.C:56
 hmpid_digits.C:57
 hmpid_digits.C:58
 hmpid_digits.C:59
 hmpid_digits.C:60
 hmpid_digits.C:61
 hmpid_digits.C:62
 hmpid_digits.C:63
 hmpid_digits.C:64
 hmpid_digits.C:65
 hmpid_digits.C:66
 hmpid_digits.C:67
 hmpid_digits.C:68
 hmpid_digits.C:69
 hmpid_digits.C:70
 hmpid_digits.C:71
 hmpid_digits.C:72
 hmpid_digits.C:73
 hmpid_digits.C:74
 hmpid_digits.C:75
 hmpid_digits.C:76
 hmpid_digits.C:77
 hmpid_digits.C:78
 hmpid_digits.C:79
 hmpid_digits.C:80
 hmpid_digits.C:81
 hmpid_digits.C:82
 hmpid_digits.C:83
 hmpid_digits.C:84
 hmpid_digits.C:85
 hmpid_digits.C:86
 hmpid_digits.C:87
 hmpid_digits.C:88
 hmpid_digits.C:89
 hmpid_digits.C:90
 hmpid_digits.C:91
 hmpid_digits.C:92
 hmpid_digits.C:93