]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TRD/AliTRDdisplayTracks.C
Introduce parameter class
[u/mrichter/AliRoot.git] / TRD / AliTRDdisplayTracks.C
CommitLineData
fa148e6c 1void AliTRDdisplayTracks(Int_t track)
2{
3
4 c1 = new TCanvas( "RecPoints", "RecPoints Display", 10, 10, 710, 740);
5 c1->SetFillColor(1);
6 TView *v=new TView(1);
7
8 v->SetRange(-350.,-350.,-400.,350.,350.,400.); // full
9 // v->SetRange(0.,0.,0.,350.,350.,400.); // top right
10 // v->SetRange(-350.,0.,0.,0.,350.,400.); // top left
11 // v->SetRange(0.,-350.,0.,350.,0.,400.); // bottom right
12 // v->SetRange(-350.,-350.,0.,0.,0.,400.); // bottom left
13 // v->Side();
14 v->Top();
15 cerr<<"psi = "<<v->GetPsi()<<endl;
16
17 // Dynamically link some shared libs
18 if (gClassTable->GetID("AliRun") < 0) {
19 gROOT->LoadMacro("loadlibs.C");
20 loadlibs();
21 cout << "Loaded shared libraries" << endl;
22 }
23
24
25// Load tracks
26
27 TFile *tf=TFile::Open("AliTRDtracks.root");
28 if (!tf->IsOpen()) {cerr<<"Can't open AliTRDtracks.root !\n"; return 3;}
29 TObjArray tarray(2000);
30 TTree *tracktree=(TTree*)tf->Get("TreeT");
31 TBranch *tbranch=tracktree->GetBranch("tracks");
32 Int_t nentr=tracktree->GetEntries();
33 cerr<<"Found "<<nentr<<" entries in the track tree"<<endl;
34 for (Int_t i=0; i<nentr; i++) {
35 AliTRDtrack *iotrack=new AliTRDtrack;
36 tbranch->SetAddress(&iotrack);
37 tracktree->GetEvent(i);
38 tarray.AddLast(iotrack);
39 cerr<<"got track "<<i<<": index ="<<iotrack->GetLabel()<<endl;
40 }
41 tf->Close();
42
43
44 // Load clusters
45 Char_t *alifile = "AliTRDclusters.root";
46 Int_t nEvent = 0;
47 TObjArray rparray(2000);
48 TObjArray *RecPointsArray = &rparray;
49 AliTRDtracker *Tracker = new AliTRDtracker("dummy","dummy");
50 Tracker->ReadClusters(RecPointsArray,alifile,-1);
51 Int_t nRecPoints = RecPointsArray->GetEntriesFast();
52 cerr<<"Found "<<nRecPoints<<" rec. points"<<endl;
53
54
55 // Connect the AliRoot file containing Geometry, Kine, Hits, and Digits
56 alifile = "galice.root";
57 TFile *gafl = (TFile*) gROOT->GetListOfFiles()->FindObject(alifile);
58 if (!gafl) {
59 cout << "Open the ALIROOT-file " << alifile << endl;
60 gafl = new TFile(alifile);
61 }
62 else {
63 cout << alifile << " is already open" << endl;
64 }
65
66 // Get AliRun object from file or create it if not on file
67 gAlice = (AliRun*) gafl->Get("gAlice");
68 if (gAlice)
69 cout << "AliRun object found on file" << endl;
70 else
71 gAlice = new AliRun("gAlice","Alice test program");
72
73
74 AliTRDv1 *TRD = (AliTRDv1*) gAlice->GetDetector("TRD");
75 AliTRDgeometry *fGeom = TRD->GetGeometry();
76
77 Int_t i,j,index,det,sector, ti[3];
78 Double_t x,y,z, cs,sn,tmp;
79 Float_t global[3], local[3];
80
81 // Display all TRD RecPoints
82 TPolyMarker3D *pm = new TPolyMarker3D(nRecPoints);
83
84 for (Int_t i = 0; i < nRecPoints; i++) {
85 printf("\r point %d out of %d",i,nRecPoints);
86 AliTRDrecPoint *rp = (AliTRDrecPoint *) RecPointsArray->UncheckedAt(i);
87
88 local[0]=rp->GetLocalRow();
89 local[1]=rp->GetLocalCol();
90 local[2]=rp->GetLocalTime();
91 det=rp->GetDetector();
92
93 for (j = 0; j < 3; j++) { ti[j] = rp->GetTrackIndex(j); }
94
95 if((track < 0) ||
96 ((ti[0]==track)||(ti[1]==track)||(ti[2]==track))) {
97
98 if(fGeom->Local2Global(det,local,global)) {
99 x=global[0]; y=global[1]; z=global[2];
100 pm->SetPoint(i,x,y,z);
101 }
102 }
103 }
104
105 pm->SetMarkerSize(1); pm->SetMarkerColor(10); pm->SetMarkerStyle(1);
106 pm->Draw();
107
108
109 Int_t ntracks = tarray.GetEntriesFast();
110
111 for (i = 0; i < ntracks; i++) {
112 AliTRDtrack *pt = (AliTRDtrack *) tarray.UncheckedAt(i), &t=*pt;
113 Int_t nclusters = t.GetNclusters();
114 cerr<<"in track "<<i<<" found "<<nclusters<<" clusters"<<endl;
115
116 TPolyMarker3D *pm = new TPolyMarker3D(nclusters);
117 for(j = 0; j < nclusters; j++) {
118 index = t.GetClusterIndex(j);
119 AliTRDrecPoint *rp = (AliTRDrecPoint *) RecPointsArray->UncheckedAt(index);
120 local[0]=rp->GetLocalRow();
121 local[1]=rp->GetLocalCol();
122 local[2]=rp->GetLocalTime();
123 det=rp->GetDetector();
124
125 if(fGeom->Local2Global(det,local,global)) {
126 x=global[0]; y=global[1]; z=global[2];
127 pm->SetPoint(j,x,y,z);
128 }
129 }
130 pm->SetMarkerSize(1); pm->SetMarkerColor(i%6+3); pm->SetMarkerStyle(1);
131 pm->Draw();
132 }
133
134 TGeometry *geom=(TGeometry*)gafl->Get("AliceGeom");
135 geom->Draw("same");
136
137 c1->Modified();
138
139 c1->Update();
140
141 gafl->Close();
142
143}