]> git.uio.no Git - u/mrichter/AliRoot.git/blame - STEER/CopyReference.C
Merge remote-tracking branch 'remotes/origin/master' into TPCdev
[u/mrichter/AliRoot.git] / STEER / CopyReference.C
CommitLineData
d6859b2d 1//
2// This macro translate original Tracks Reference tree
3// into new tree sorted according to track label.
4//
5// TClonesArray connected to a branch, is filled
6// with Track References from one track only.
7//
8// Algorithm:
9//
10// 1) Open files and initialise trees
11// 2) Create and initialise Look Up Table
12// 3) Loop over entries in References tree
13// - store entry number in a LUT
14// - store first and last array index of references for one track
15// 4) Check LUT's consistancy.
16// 5) Fill up new tree with references ordered acording to label.
17//
18//
19// Sylwester Radomski (S.Radomski@gsi.de)
20// GSI, Jan 20, 2003
21//
22
23void CopyReference () {
24
25 // create output file
26 const char *outFileName = "trackRef.root";
27 TFile *outFile = new TFile(outFileName, "RECREATE");
28
29 // connect to gAlice object in input file
30 const char *inFileName = "galice.root";
31 TFile *inFile = new TFile(inFileName, "READ");
32
33 gAlice = (AliRun*)inFile->Get("gAlice");
34 gAlice->GetEvent(0);
35
36 // declare detectors to translate
37 const Int_t nDet = 4;
38 const char *detectors[nDet] = {"ITS","TPC", "TRD", "TOF"};
39
40 // Create LUT's
41
42 Int_t nPart = gAlice->GetNtrack();
43 Int_t *lutTrack = new Int_t[nPart * nDet];
44 Int_t *lutStart = new Int_t[nPart * nDet];
45 Int_t *lutStop = new Int_t[nPart * nDet];
46
47 // clean LUT
48
49 for(Int_t i=0; i<nPart; i++)
50 for (Int_t j=0; j<nDet; j++)
51 lutTrack[i*nDet+j] = lutStart[i*nDet+j] = lutStop[i*nDet+j] = -1;
52
53
54 // Declare tree, branches and CloneArrays
55
56 TTree *refTree = (TTree*)inFile->Get("TreeTR0");
57 TBranch *refBranch[nDet];
58 TClonesArray *refArray[nDet];
59
60 // Connect branches with arrays
61
62 for (Int_t j=0; j<nDet; j++) {
63 refBranch[j] = refTree->GetBranch(detectors[j]);
64 refArray[j] = new TClonesArray("AliTrackReference", 100);
65 refBranch[j]->SetAddress(&(refArray+j));
66 }
67
68 Int_t nTracks = refBranch[0]->GetEntries();
69
70 cout << "N Particles\t" << nPart << endl
71 << "N Tracks\t" << nTracks << endl;
72
73 cout << "Filling Look Up Tables ..." << endl;
74
75 for (Int_t i=0; i<nTracks; i++) {
76
77 cout << i << "\r";
78
79 for (Int_t j=0; j<nDet; j++) {
80
81 Int_t lab, start;
82 AliTrackReference *ref = 0;
83 refBranch[j]->GetEvent(i);
84 Int_t nObj = refArray[j]->GetEntries();
85
86 if (!nObj) continue;
87
88 lab = ((AliTrackReference*)refArray[j]->At(0))->GetTrack();
89 start = 0;
90
91 for (Int_t e=0; e<nObj-1; e++) {
92
93 ref = (AliTrackReference*)refArray[j]->At(e+1);
94
95 if (ref->GetTrack() != lab) {
96 lutTrack[lab*nDet + j] = i;
97 lutStart[lab*nDet + j] = start;
98 lutStop[lab*nDet + j] = e+1;
99
100 start = e+1;
101 lab = ref->GetTrack();
102 }
103 }
104
105 lutTrack[lab*nDet + j] = i;
106 lutStart[lab*nDet + j] = start;
107 lutStop[lab*nDet + j] = nObj;
108 }
109 }
110
111 cout << endl << "done" << endl;
112 cout << "Checking consistancy of LUTs ..." << endl;
113
114 // check consistancy
115
116 for(Int_t i=0; i<nPart; i++) {
117
118 cout << i << "\r";
119 Int_t ctrlTrack = -1;
120
121 for(Int_t j=0; j<nDet; j++) {
122
123 if (ctrlTrack == -1 && lutTrack[i*nDet+j] != -1)
124 ctrlTrack = lutTrack[i*nDet+j];
125
126 if (lutTrack[i*nDet+j] != -1 && lutTrack[i*nDet+j] != ctrlTrack)
127 cout << "Error: " << i << " " << j << endl;
128 }
129 }
130
131 cout << endl << "done" << endl;
132 cout << "Writing to a new Tree ..." << endl;
133
134 // Create a Tree
135
136 outFile->cd();
137 TTree *outTree = new TTree("TreeTR0_Sorted", "Sorted Tracks References");
138 TBranch *outBranch[nDet];
139 TClonesArray *outArray[nDet];
140
141 // Create Branches
142
143 for(Int_t j=0; j<nDet; j++) {
144 outArray[j] = new TClonesArray("AliTrackReference", 200);
145 outBranch[j] = outTree->Branch(detectors[j], &(outArray+j), 64000, 3);
146 }
147
148 // Loop over particles
149 // Fill Tree if entries exists
150
151 for(Int_t i=0; i<nPart; i++) {
152
153 cout << i << "\r";
154
155 for(Int_t j=0; j<nDet; j++) {
156
157 outArray[j]->Clear();
158
159 if (lutTrack[i*nDet+j] != -1) {
160
161 refBranch[j]->GetEvent(lutTrack[i*nDet+j]);
162
163 // rewrite objects in clonearrays
164 for(Int_t k = lutStart[i*nDet+j], en = 0; k<lutStop[i*nDet+j]; k++ ) {
165 AliTrackReference *ref = (AliTrackReference*)refArray[j]->At(k);
166 new( (*(outArray[j]))[en++] ) AliTrackReference(*ref);
167 }
168 }
169 }
170
171 outTree->Fill();
172 }
173
174 cout << endl << "done" << endl;
175
176 // close and clean
177
178 outTree->Write();
179 inFile->Close();
180 outFile->Close();
181
182 //if (refTree) delete refTree;
183 //if (outTree) delete outTree;
184 //if (inFile) delete inFile;
185 //if (outFile) delete outFile;
186}