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