]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/alice-macros/trd_detectors.C
0742933e3a130f9d5354b43fd11bc11a3ddfeac6
[u/mrichter/AliRoot.git] / EVE / alice-macros / trd_detectors.C
1 //
2 // How to steer the basic TRD data containers from a macro.
3 // 
4 // The loading and looping over events is performed by the
5 // AliEve mechanism. For a stand alone TRD loop control check 
6 // the macro "trd_loader.C"
7 // 
8 // Usage:
9 // .x trd_detectors.C(sector)
10 // 
11 // Author:
12 // Alex Bercuci (A.Bercuci@gsi.de)
13 // 
14 #ifdef __CINT__
15 class TEvePointSet;
16 class TEveElement;
17 #else
18 #include <TEveManager.h>
19 #include <TEvePointSet.h>
20 #include <EveBase/AliEveEventManager.h>
21 #include <EveDet/AliEveTRDModuleImp.h>
22
23 #include "AliRunLoader.h"
24 #include "AliCluster.h"
25 #include "TRD/AliTRDcluster.h"
26 #include "TRD/AliTRDgeometry.h"
27 #include "TRD/AliTRDdigitsManager.h"
28 #endif
29
30 TEveElementList* trd_detectors(Int_t sector = -1, TEveElement *cont = 0)
31 {
32   // Link data containers
33   AliEveEventManager::AssertGeometry();
34
35   AliRunLoader *rl = AliEveEventManager::AssertRunLoader();
36   
37   // define EVE containers
38   TEveElementList *list = new TEveElementList("TRD Detectors");
39   
40   AliTRDgeometry *geo = new AliTRDgeometry();
41   //geo->CreateClusterMatrixArray();
42   
43   AliEveTRDNode *sm = 0x0, *stk = 0x0; 
44   AliEveTRDChamber *chm=0x0;
45
46   // Link TRD containers
47   TObjArray *clusters = 0x0;
48   rl->LoadRecPoints("TRD");
49   TTree *tR = rl->GetTreeR("TRD", kFALSE);
50   tR->SetBranchAddress("TRDcluster", &clusters);
51
52   rl->LoadDigits("TRD");
53   TTree *tD = rl->GetTreeD("TRD", kFALSE);
54   AliTRDdigitsManager dm; dm.ReadDigits(tD);
55
56   for(Int_t i=0; i<tR->GetEntries(); i++) {
57     if (!tR->GetEvent(i)) continue;
58     if(!clusters->GetEntries()) continue;
59     Int_t icl=0; AliTRDcluster *c = 0x0;
60     while(!(c = (AliTRDcluster*)clusters->UncheckedAt(icl++))) {;}
61     if(!c) continue;
62
63     Int_t idet, ism, istk, ipla; 
64     idet = c->GetDetector();
65     ism  = geo->GetSector(idet);
66     istk = geo->GetStack(idet);
67     ipla = geo->GetLayer(idet);
68     if(sector>=0 && ism != sector) continue;
69     if(!(sm = (AliEveTRDNode*)list->FindChild(Form("SM%03d", ism)))){ 
70       list->AddElement(sm = new AliEveTRDNode("SM", ism));
71       sm->SetElementTitle(Form("Supermodule %2d", ism));
72     }
73     if(!(stk=(AliEveTRDNode*)sm->FindChild(Form("Stack%03d", istk)))){
74       sm->AddElement(stk = new AliEveTRDNode("Stack", istk));
75       stk->SetElementTitle(Form("SM %2d Stack %1d", ism, istk));
76     }
77     stk->AddElement(chm = new AliEveTRDChamber(idet));
78     chm->SetGeometry(geo);
79     chm->LoadClusters(clusters);
80     chm->LoadDigits(&dm);
81
82     //clusters->Clear();
83   }
84
85   rl->UnloadDigits("TRD");
86   rl->UnloadRecPoints("TRD");
87
88   gEve->AddElement(list, cont);
89   gEve->Redraw3D();
90
91   return list;
92 }