]> git.uio.no Git - u/mrichter/AliRoot.git/blame - EVE/alice-macros/trd_detectors.C
Error codes for histogram booking taken with getters, no longer hard coded (D. Elia)
[u/mrichter/AliRoot.git] / EVE / alice-macros / trd_detectors.C
CommitLineData
6983e87a 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//
66c3c96b 14#ifdef __CINT__
15class TEvePointSet;
16class TEveElement;
17#else
18#include <TEveManager.h>
19#include <TEvePointSet.h>
20#include <EveBase/AliEveEventManager.h>
9fb2a3ed 21#include <EveDet/AliEveTRDModuleImp.h>
66c3c96b 22
23#include "AliRunLoader.h"
24#include "AliCluster.h"
9fb2a3ed 25#include "TRD/AliTRDcluster.h"
26#include "TRD/AliTRDgeometry.h"
27#include "TRD/AliTRDdigitsManager.h"
66c3c96b 28#endif
29
794d4a46 30TEveElementList* trd_detectors(Int_t sector = -1, TEveElement *cont = 0)
66c3c96b 31{
32 // Link data containers
9fb2a3ed 33 AliEveEventManager::AssertGeometry();
66c3c96b 34
6983e87a 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;
66c3c96b 45
6983e87a 46 // Link TRD containers
47 TObjArray *clusters = 0x0;
48 rl->LoadRecPoints("TRD");
49 TTree *tR = rl->GetTreeR("TRD", kFALSE);
50 tR->SetBranchAddress("TRDcluster", &clusters);
66c3c96b 51
6983e87a 52 rl->LoadDigits("TRD");
53 TTree *tD = rl->GetTreeD("TRD", kFALSE);
54 AliTRDdigitsManager dm; dm.ReadDigits(tD);
66c3c96b 55
6983e87a 56 for(Int_t i=0; i<tR->GetEntries(); i++) {
57 if (!tR->GetEvent(i)) continue;
d4d1b414 58 if(!clusters->GetEntries()) continue;
59 Int_t icl=0; AliTRDcluster *c = 0x0;
6983e87a 60 while(!(c = (AliTRDcluster*)clusters->UncheckedAt(icl++))) {;}
d4d1b414 61 if(!c) continue;
62
63 Int_t idet, ism, istk, ipla;
6983e87a 64 idet = c->GetDetector();
65 ism = geo->GetSector(idet);
794d4a46 66 istk = geo->GetStack(idet);
67 ipla = geo->GetLayer(idet);
68 if(sector>=0 && ism != sector) continue;
9fb2a3ed 69 if(!(sm = (AliEveTRDNode*)list->FindChild(Form("SM%03d", ism)))){
6983e87a 70 list->AddElement(sm = new AliEveTRDNode("SM", ism));
71 sm->SetElementTitle(Form("Supermodule %2d", ism));
66c3c96b 72 }
9fb2a3ed 73 if(!(stk=(AliEveTRDNode*)sm->FindChild(Form("Stack%03d", istk)))){
6983e87a 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);
66c3c96b 81
9b7da201 82 //clusters->Clear();
66c3c96b 83 }
6983e87a 84
9fb2a3ed 85 rl->UnloadDigits("TRD");
86 rl->UnloadRecPoints("TRD");
87
6983e87a 88 gEve->AddElement(list, cont);
66c3c96b 89 gEve->Redraw3D();
90
6983e87a 91 return list;
66c3c96b 92}