]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/alice-macros/trd_digits.C
SetSeed dummy implementations
[u/mrichter/AliRoot.git] / EVE / alice-macros / trd_digits.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
15 #if !defined(__CINT__) || defined(__MAKECINT__)
16 #include <TGeoManager.h>
17 #include <TEveManager.h>
18 #include <TEveElement.h>
19 #include <TEvePointSet.h>
20
21 #include <AliCluster.h>
22 #include <AliRunLoader.h>
23 #include <AliTRDarrayADC.h>
24 #include <AliTRDcluster.h>
25 #include <AliTRDgeometry.h>
26 #include <AliTRDdigitsManager.h>
27 #include <AliEveEventManager.h>
28 #include <AliEveTRDModuleImp.h>
29 #else
30 class TEvePointSet;
31 class TEveElement;
32 #endif
33
34 TEveElementList* trd_digits(Int_t sector = -1, TEveElement *cont = 0)
35 {
36   // Link data containers
37   if(!gGeoManager) AliEveEventManager::AssertGeometry();
38
39   AliRunLoader *rl = AliEveEventManager::AssertRunLoader();
40   
41   // define EVE containers
42   TEveElementList *list = new TEveElementList("TRD Digits");
43   
44   AliEveTRDNode *sm(NULL), *stk(NULL); 
45   AliEveTRDChamber *chm(NULL);
46   // Link TRD digits
47   rl->LoadDigits("TRD");
48   TTree *tD = rl->GetTreeD("TRD", kFALSE);
49   if(!tD){ 
50     Error("trd_digits", "Missing digits tree");
51     return NULL;
52   }
53   AliTRDdigitsManager dm;
54   dm.ReadDigits(tD);
55
56   AliTRDgeometry *geo = new AliTRDgeometry();
57   Int_t sBegin=sector<0?0:sector,
58         sEnd  =sector<0?(AliTRDgeometry::kNsector):sector+1;
59   Int_t jdet(0);
60   for(Int_t isec=sBegin; isec<sEnd; isec++) {
61     sm = new AliEveTRDNode("Sector", isec);
62     for(Int_t istk(0); istk<AliTRDgeometry::kNstack; istk++) {
63       stk = new AliEveTRDNode("Stack", istk);
64       stk->SetTitle(Form("Index %d", isec*AliTRDgeometry::kNstack+istk));
65       for(Int_t ily(0); ily<AliTRDgeometry::kNlayer; ily++) {
66         Int_t det=AliTRDgeometry::GetDetector(ily, istk, isec);
67         if(!(dm.GetDigits(det)->GetDim())) continue;
68         chm=new AliEveTRDChamber(det);
69         chm->SetGeometry(geo);
70         chm->LoadDigits(&dm);
71         stk->AddElement(chm);
72         jdet++;
73       }
74       if(!stk->HasChildren()){
75         delete stk;
76         continue;
77       }
78       sm->AddElement(stk);
79     }
80     if(!sm->HasChildren()){
81       delete sm;
82       continue;
83     }
84     list->AddElement(sm);
85   }
86   rl->UnloadDigits("TRD");
87   gEve->AddElement(list, cont);
88   gEve->Redraw3D();
89
90   Info("trd_digits", "TRD chambers with data for current selection %d.", jdet);
91   return list;
92 }