]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/alice-macros/hmpid_raw.C
Modifications in CMake* files to export all required headers.
[u/mrichter/AliRoot.git] / EVE / alice-macros / hmpid_raw.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 <TGeoMatrix.h>
9 #include <TStyle.h>
10 #include <TEveElement.h>
11 #include <TEveFrameBox.h>
12 #include <TEveManager.h>
13 #include <TEveRGBAPalette.h>
14 #include <TEveTrans.h>
15 #include <TEveQuadSet.h>
16
17 #include <AliHMPIDDigit.h>
18 #include <AliHMPIDv3.h>
19 #include <AliHMPIDRawStream.h>
20 #include <AliRawReader.h>
21 #include <AliEveEventManager.h>
22 #endif
23
24 void hmpid_raw()
25 {
26   const Char_t *name[] = { "HMPID0", "HMPID1", "HMPID2", "HMPID3",
27                            "HMPID4", "HMPID5", "HMPID6" };
28
29   AliRawReader *rawReader = AliEveEventManager::AssertRawReader();
30   AliHMPIDRawStream stream(rawReader);    
31
32   TEveElementList* list = new TEveElementList("HMPID Raw");
33   gEve->AddElement(list);
34
35   gStyle->SetPalette(1, 0);
36
37   TEveRGBAPalette *pal = new TEveRGBAPalette(0, 3000);
38   pal->SetMax(1000);
39   TEveFrameBox    *box = new TEveFrameBox();
40   box->SetAAQuadXY(0, 0, 0, 144, 121);
41   box->SetFrameColor(kGray);
42
43   TEveQuadSet* ms[7];
44   for (Int_t iCh = 0; iCh < 7; ++iCh)
45   {
46     ms[iCh] = new TEveQuadSet(Form("Chamber %d", iCh));
47
48     TEveQuadSet* q = ms[iCh];
49     q->SetOwnIds(kTRUE);
50     q->SetPalette(pal);
51     q->SetFrame(box);
52     q->SetAntiFlick(kTRUE);
53     q->SetPickable(kTRUE);
54
55     q->Reset(TEveQuadSet::kQT_RectangleXYFixedDimZ, kFALSE, 64);
56     q->SetDefCoord(0);
57     q->SetDefHeight(0.84f);
58     q->SetDefWidth(0.8f);
59   }
60
61   while (stream.Next())
62   {
63     Int_t ch = AliHMPIDParam::DDL2C(stream.GetDDLNumber());
64     TEveQuadSet* q = ms[ch];
65
66     for (Int_t iPad = 0; iPad < stream.GetNPads(); ++iPad)
67     {
68       AliHMPIDDigit dig(stream.GetPadArray()[iPad],stream.GetChargeArray()[iPad]);
69
70       q->AddQuad(dig.PadChX()*0.8f,  dig.PadChY()*0.84f);
71       q->QuadValue(TMath::Nint(dig.Q()));
72       q->QuadId(new AliHMPIDDigit(dig));
73     }
74   }
75
76   for (Int_t iCh = 0; iCh < 7; ++iCh)
77   {
78     TEveQuadSet* q = ms[iCh];
79
80     q->RefitPlex();
81
82     TGeoHMatrix mat;
83     AliHMPIDv3::IdealPosition(iCh, &mat);
84     q->RefMainTrans().SetFrom(mat);
85     q->RefMainTrans().Move3LF(-0.5*144, -0.5*121, 0);
86
87     list->AddElement(q);
88   }
89
90   gEve->Redraw3D();
91 }