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