]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/alice-macros/hmpid_digits.C
doxy: code cleanup: comments and clarifications
[u/mrichter/AliRoot.git] / EVE / alice-macros / hmpid_digits.C
1 /**************************************************************************
2  * Copyright(c) 1998-2008, ALICE Experiment at CERN, all rights reserved. *
3  * See http://aliceinfo.cern.ch/Offline/AliRoot/License.html for          *
4  * full copyright notice.                                                 *
5  **************************************************************************/
6
7 #if !defined(__CINT__) || defined(__MAKECINT__)
8 #include <TClonesArray.h>
9 #include <TBranch.h>
10 #include <TGeoMatrix.h>
11 #include <TTree.h>
12 #include <TStyle.h>
13 #include <TEveElement.h>
14 #include <TEveFrameBox.h>
15 #include <TEveManager.h>
16 #include <TEvePointSet.h>
17 #include <TEveRGBAPalette.h>
18 #include <TEveTrans.h>
19 #include <TEveQuadSet.h>
20
21 #include <AliHMPIDDigit.h>
22 #include <AliHMPIDv3.h>
23 #include <AliCluster3D.h>
24 #include <AliRunLoader.h>
25 #include <AliEveEventManager.h>
26 #endif
27
28 void hmpid_digits()
29 {
30   const Char_t *name[]={ "HMPID0", "HMPID1", "HMPID2", "HMPID3",
31                          "HMPID4", "HMPID5", "HMPID6" };
32
33   AliRunLoader* rl = AliEveEventManager::AssertRunLoader();
34   rl->LoadDigits("HMPID");
35
36   TTree *dTree = rl->GetTreeD("HMPID", kFALSE);
37   if (!dTree) return;
38
39   TEveElementList* list = new TEveElementList("HMPID Digits");
40   gEve->AddElement(list);
41
42   gStyle->SetPalette(1, 0);
43
44   TEveRGBAPalette *pal = new TEveRGBAPalette(0, 3000);
45   pal->SetMax(1000);
46   TEveFrameBox    *box = new TEveFrameBox();
47   box->SetAAQuadXY(0, 0, 0, 144, 121);
48   box->SetFrameColor(kGray);
49
50   TClonesArray* digits = new TClonesArray("AliHMPIDDigit");
51   for (Int_t iCh = 0; iCh < 7; ++iCh)
52   {
53     TBranch *br = dTree->GetBranch(name[iCh]);
54     br->SetAddress(&digits);
55     br->GetEntry(0);
56
57     TEveQuadSet* q = new TEveQuadSet(Form("Chamber %d", iCh));
58     q->SetOwnIds(kTRUE);
59     q->SetPalette(pal);
60     q->SetFrame(box);
61     q->SetAntiFlick(kTRUE);
62     q->SetPickable(kTRUE);
63
64     q->Reset(TEveQuadSet::kQT_RectangleXYFixedDimZ, kFALSE, 64);
65     q->SetDefCoord(0);
66     q->SetDefHeight(0.84f);
67     q->SetDefWidth(0.8f);
68
69     for(Int_t iDig = 0; iDig < digits->GetEntriesFast(); ++iDig)
70     {
71       AliHMPIDDigit *pDig = (AliHMPIDDigit*) digits->At(iDig);
72
73       q->AddQuad(pDig->PadChX()*0.8f,  pDig->PadChY()*0.84f);
74       q->QuadValue(TMath::Nint(pDig->Q()));
75       q->QuadId(new AliHMPIDDigit(*pDig));
76     }
77
78     q->RefitPlex();
79
80     TGeoHMatrix mat;
81     AliHMPIDv3::IdealPosition(iCh, &mat);
82     q->RefMainTrans().SetFrom(mat);
83     q->RefMainTrans().Move3LF(-0.5f*144, -0.5f*121, 0);
84
85     list->AddElement(q);
86   }
87
88   delete digits;
89   rl->UnloadDigits("HMPID");
90
91   gEve->Redraw3D();
92 }