]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/alice-macros/trd_detectors.C
Protection against missing digits (the case in mainstream reconstruction). Otherwise...
[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   AliEveTRDNode *sm = 0x0, *stk = 0x0; 
41   AliEveTRDChamber *chm=0x0;
42
43   // Link TRD containers
44   TObjArray *clusters = 0x0;
45   rl->LoadRecPoints("TRD");
46   TTree *tR = rl->GetTreeR("TRD", kFALSE);
47   if (!tR) return 0;
48   tR->SetBranchAddress("TRDcluster", &clusters);
49
50   rl->LoadDigits("TRD");
51   TTree *tD = rl->GetTreeD("TRD", kFALSE);
52   AliTRDdigitsManager dm;
53   if (tD) dm.ReadDigits(tD);
54
55   AliTRDgeometry *geo = new AliTRDgeometry();
56   
57   for(Int_t i=0; i<tR->GetEntries(); i++) {
58     if (!tR->GetEvent(i)) continue;
59     if(!clusters->GetEntries()) continue;
60     Int_t icl=0; AliTRDcluster *c = 0x0;
61     while(!(c = (AliTRDcluster*)clusters->UncheckedAt(icl++))) {;}
62     if(!c) continue;
63
64     Int_t idet, ism, istk, ipla; 
65     idet = c->GetDetector();
66     ism  = geo->GetSector(idet);
67     istk = geo->GetStack(idet);
68     ipla = geo->GetLayer(idet);
69     if(sector>=0 && ism != sector) continue;
70     if(!(sm = (AliEveTRDNode*)list->FindChild(Form("SM%03d", ism)))){ 
71       list->AddElement(sm = new AliEveTRDNode("SM", ism));
72       sm->SetElementTitle(Form("Supermodule %2d", ism));
73     }
74     if(!(stk=(AliEveTRDNode*)sm->FindChild(Form("Stack%03d", istk)))){
75       sm->AddElement(stk = new AliEveTRDNode("Stack", istk));
76       stk->SetElementTitle(Form("SM %2d Stack %1d", ism, istk));
77     }
78     stk->AddElement(chm = new AliEveTRDChamber(idet));
79     chm->SetGeometry(geo);
80     chm->LoadClusters(clusters);
81     if (tD) chm->LoadDigits(&dm);
82
83     //clusters->Clear();
84   }
85
86   rl->UnloadDigits("TRD");
87   rl->UnloadRecPoints("TRD");
88
89   gEve->AddElement(list, cont);
90   gEve->Redraw3D();
91
92   return list;
93 }