]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/alice-macros/trd_tracklets.C
changes for Vertex and Tracks classes
[u/mrichter/AliRoot.git] / EVE / alice-macros / trd_tracklets.C
1
2 #if !defined(__CINT__) || defined(__MAKECINT__)
3 #include <TClonesArray.h>
4 #include <TEveLine.h>
5 #include <TEveManager.h>
6 #include <TEveElement.h>
7
8 #include <AliRunLoader.h>
9 #include <AliLoader.h>
10 #include <AliDataLoader.h>
11 #include <AliTreeLoader.h>
12 #include <AliTRDarrayADC.h>
13 #include <AliTRDtrackletWord.h>
14 #include <AliTRDtrackletMCM.h>
15 #include <AliEveEventManager.h>
16 #include <AliEveTRDData.h>
17 #endif
18
19 TEveElementList *trd_tracklets()
20 {
21   AliRunLoader* rl =  AliEveEventManager::AssertRunLoader();
22   AliLoader *loader = rl ? rl->GetLoader("TRDLoader") : 0x0;
23
24   TTree *trklTree = 0x0;
25
26   AliDataLoader *dl = loader ? loader->GetDataLoader("tracklets") : 0x0;
27   if (!dl) {
28     printf("No tracklet loader\n");
29     return 0x0;
30   }
31
32   gEve->DisableRedraw();
33
34   // ----- simulated tracklets -----
35   dl->Load();
36   trklTree = dl->Tree();
37   
38   if (trklTree) {
39     TBranch *trklBranch = 0x0;
40     if ((trklBranch = trklTree->GetBranch("mcmtrklbranch"))) {
41       AliTRDtrackletMCM *trkl = 0x0; 
42       trklBranch->SetAddress(&trkl);
43       
44       TEveElementList* listOfTracklets = new TEveElementList("TRD tracklets (sim)");
45       gEve->AddElement(listOfTracklets);
46       
47       for (Int_t i = 0; i < trklBranch->GetEntries(); i++) {
48         trklBranch->GetEntry(i);
49         if (!trkl)
50           continue;
51         gEve->AddElement(new AliEveTRDTrackletOnline(trkl), listOfTracklets);
52       }
53     }
54   }
55
56   // raw tracklets
57   AliTreeLoader *tl = (AliTreeLoader*) dl->GetBaseLoader("tracklets-raw");
58   if (tl) {
59     tl->Load();
60     trklTree = tl->Tree();
61   }
62   else 
63     trklTree = 0x0;
64   //  trklTree = tl ? tl->Load(), tl->Tree : 0x0;
65
66   if (trklTree) {
67     TEveElementList* listOfTracklets = new TEveElementList("TRD tracklets (raw)");
68     gEve->AddElement(listOfTracklets);
69     
70     Int_t hc; 
71     TClonesArray *ar = 0x0;
72     trklTree->SetBranchAddress("hc", &hc);
73     trklTree->SetBranchAddress("trkl", &ar);
74
75     for (Int_t iEntry = 0; iEntry < trklTree->GetEntries(); iEntry++) {
76       trklTree->GetEntry(iEntry);
77       //      printf("%i tracklets in HC %i\n", ar->GetEntriesFast(), hc);
78       for (Int_t iTracklet = 0; iTracklet < ar->GetEntriesFast(); iTracklet++) {
79         AliTRDtrackletWord *trklWord = (AliTRDtrackletWord*) (*ar)[iTracklet];
80         AliEveTRDTrackletOnline *evetrkl = new AliEveTRDTrackletOnline(new AliTRDtrackletWord(trklWord->GetTrackletWord(), hc));
81         gEve->AddElement(evetrkl, listOfTracklets);
82       }
83     }
84   }
85
86   gEve->EnableRedraw();
87   gEve->Redraw3D();
88
89   return 0x0;
90 }
91