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