]>
Commit | Line | Data |
---|---|---|
13918578 | 1 | /************************************************************************** |
2 | * Copyright(c) 1998-2003, ALICE Experiment at CERN, All rights reserved. * | |
3 | * * | |
4 | * Author: The ALICE Off-line Project. * | |
5 | * Contributors are mentioned in the code where appropriate. * | |
6 | * * | |
7 | * Permission to use, copy, modify and distribute this software and its * | |
8 | * documentation strictly for non-commercial purposes is hereby granted * | |
9 | * without fee, provided that the above copyright notice appears in all * | |
10 | * copies and that both the copyright notice and this permission notice * | |
11 | * appear in the supporting documentation. The authors make no claims * | |
12 | * about the suitability of this software for any purpose. It is * | |
13 | * provided "as is" without express or implied warranty. * | |
14 | **************************************************************************/ | |
15 | ||
16 | //////////////////////////////////////////////////// | |
17 | // Stand alone tracker class // | |
18 | // Origin: Elisabetta Crescio // | |
19 | // e-mail: crescio@to.infn.it // | |
20 | // tracks are saved as AliITStrackV2 objects // | |
21 | //////////////////////////////////////////////////// | |
22 | ||
23 | //#include <stdlib.h> | |
24 | #include "TArrayI.h" | |
25 | #include <TBranch.h> | |
26 | #include <TMath.h> | |
27 | #include <TObjArray.h> | |
28 | #include <TTree.h> | |
29 | #include "AliRun.h" | |
30 | #include "AliITSclusterTable.h" | |
31 | #include "AliITSclusterV2.h" | |
32 | #include "AliITSgeom.h" | |
33 | #include "AliITSRiemannFit.h" | |
34 | #include "AliITStrackerSA.h" | |
35 | #include "AliITStrackSA.h" | |
36 | #include "AliITSVertexer.h" | |
37 | #include "AliITSVertex.h" | |
38 | ||
39 | ||
40 | ClassImp(AliITStrackerSA) | |
41 | ||
42 | //____________________________________________________________________________ | |
43 | AliITStrackerSA::AliITStrackerSA():AliITStrackerV2(){ | |
44 | // Default constructor | |
45 | Init(); | |
46 | } | |
47 | //____________________________________________________________________________ | |
48 | AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom, AliITSVertex *vert):AliITStrackerV2(geom) | |
49 | { | |
50 | // Standard constructor (Vertex is known and passed to this obj.) | |
51 | Init(); | |
52 | fVert = vert; | |
53 | fGeom = geom; | |
54 | } | |
55 | ||
56 | //______________________________________________________________________ | |
57 | AliITStrackerSA::AliITStrackerSA(const AliITStrackerSA &trkr) : | |
58 | AliITStrackerV2(trkr) { | |
59 | // Copy constructor | |
60 | // Copies are not allowed. The method is protected to avoid misuse. | |
61 | Error("AliITStrackerSA","Copy constructor not allowed\n"); | |
62 | } | |
63 | ||
64 | //______________________________________________________________________ | |
65 | AliITStrackerSA& AliITStrackerSA::operator=(const | |
66 | AliITStrackerSA& /* trkr */){ | |
67 | // Assignment operator | |
68 | // Assignment is not allowed. The method is protected to avoid misuse. | |
69 | Error("= operator","Assignment operator not allowed\n"); | |
70 | return *this; | |
71 | } | |
72 | ||
73 | //____________________________________________________________________________ | |
74 | AliITStrackerSA::AliITStrackerSA(AliITSgeom *geom, AliITSVertexer *vertexer):AliITStrackerV2(geom) | |
75 | { | |
76 | // Standard constructor (Vertex is unknown - vertexer is passed to this obj) | |
77 | Init(); | |
78 | fVertexer = vertexer; | |
79 | fGeom = geom; | |
80 | } | |
81 | ||
82 | //____________________________________________________________________________ | |
83 | AliITStrackerSA::AliITStrackerSA(AliITStrackerSA& tracker):AliITStrackerV2(){ | |
84 | // Copy constructor | |
85 | fPhiEstimate = tracker.fPhiEstimate; | |
86 | for(Int_t i=0;i<2;i++){ | |
87 | fPoint1[i]=tracker.fPoint1[i]; | |
88 | fPoint2[i]=tracker.fPoint2[i]; | |
89 | fPoint3[i]=tracker.fPoint3[i]; | |
90 | fPointc[i]=tracker.fPointc[i]; | |
91 | } | |
92 | fLambdac = tracker.fLambdac; | |
93 | fPhic = tracker.fPhic; | |
94 | fCoef1 = tracker.fCoef1; | |
95 | fCoef2 = tracker.fCoef2; | |
96 | fCoef3 = tracker.fCoef3; | |
97 | fNloop = tracker.fNloop; | |
98 | fPhiWin = tracker.fPhiWin; | |
99 | fLambdaWin = tracker.fLambdaWin; | |
100 | if(tracker.fVertexer && tracker.fVert){ | |
101 | fVert = new AliITSVertex(*tracker.fVert); | |
102 | } | |
103 | else { | |
104 | fVert = tracker.fVert; | |
105 | } | |
106 | fVertexer = tracker.fVertexer; | |
107 | fGeom = tracker.fGeom; | |
108 | fFlagLoad = tracker.fFlagLoad; | |
109 | fTable = tracker.fTable; | |
110 | } | |
111 | ||
112 | //____________________________________________________________________________ | |
113 | AliITStrackerSA::~AliITStrackerSA(){ | |
114 | // destructor | |
115 | // if fVertexer is not null, the AliITSVertex obj. is owned by this class | |
116 | // and is deleted here | |
117 | if(fVertexer){ | |
118 | if(fVert)delete fVert; | |
119 | } | |
120 | fVert = 0; | |
121 | fVertexer = 0; | |
122 | fGeom = 0; | |
123 | fFlagLoad = 0; | |
124 | if(fPhiWin)delete []fPhiWin; | |
125 | if(fLambdaWin)delete []fLambdaWin; | |
126 | fTable =0; | |
127 | } | |
128 | ||
129 | //____________________________________________________________________________ | |
130 | void AliITStrackerSA::Init(){ | |
131 | // Reset all data members | |
132 | fPhiEstimate=0; | |
133 | for(Int_t i=0;i<3;i++){fPoint1[i]=0;fPoint2[i]=0;fPoint3[i]=0;} | |
134 | fLambdac=0; | |
135 | fPhic=0; | |
136 | fCoef1=0; | |
137 | fCoef2=0; | |
138 | fCoef3=0; | |
139 | fPointc[0]=0; | |
140 | fPointc[1]=0; | |
141 | fVert = 0; | |
142 | fVertexer = 0; | |
143 | fGeom = 0; | |
144 | fFlagLoad = 0; | |
145 | SetWindowSizes(); | |
146 | fTable = 0; | |
147 | } | |
148 | //_______________________________________________________________________ | |
149 | void AliITStrackerSA::ResetForFinding(){ | |
150 | // Reset data members used in all loops during track finding | |
151 | fPhiEstimate=0; | |
152 | for(Int_t i=0;i<3;i++){fPoint1[i]=0;fPoint2[i]=0;fPoint3[i]=0;} | |
153 | fLambdac=0; | |
154 | fPhic=0; | |
155 | fCoef1=0; | |
156 | fCoef2=0; | |
157 | fCoef3=0; | |
158 | fPointc[0]=0; | |
159 | fPointc[1]=0; | |
160 | } | |
161 | //____________________________________________________________________________ | |
162 | void AliITStrackerSA::FindTracks(TTree *clusterTree, TTree *out,Int_t evnumber,char *opt){ | |
163 | ||
164 | /************************************************************************** | |
165 | * Options: to use only ITS execute FindTracks * | |
166 | * to define good tracks with 6 points out of 6 in the ITS write * | |
167 | * *opt="6/6". The default is *opt="6/6" * | |
168 | * * | |
169 | * * | |
170 | * Example: to execute function with only the ITS (no combined tracking * | |
171 | * with TPC+ITS) and requiring 5/6 points to define a good track * | |
172 | * use: FindTracks(treein,treeout,evnumber,"5/6") * | |
173 | * to execute combined tracking, before using FindTracks use * | |
174 | * UseFoundTracksV2 * | |
175 | *************************************************************************/ | |
176 | ||
177 | //options | |
178 | TString choice(opt); | |
179 | Bool_t sixpoints= choice.Contains("6/6"); | |
180 | ||
181 | //Fill array with cluster indices for each module | |
182 | if(!fTable){ | |
183 | fTable = new AliITSclusterTable(fGeom,this); | |
184 | fTable->FillArray(clusterTree,evnumber); | |
185 | } | |
186 | //Get primary vertex | |
187 | if(fVertexer){ | |
188 | if(fVert)delete fVert; | |
189 | fVert = fVertexer->FindVertexForCurrentEvent(evnumber); | |
190 | } | |
191 | else { | |
192 | gAlice->GetEvent(evnumber); | |
193 | if(!fVert){ | |
194 | Fatal("FindTracks","Vertex is missing\n"); | |
195 | return; | |
196 | } | |
197 | } | |
198 | Double_t primaryVertex[3]; | |
199 | Double_t errorsprimvert[3]; | |
200 | fVert->GetXYZ(primaryVertex); | |
201 | fVert->GetSigmaXYZ(errorsprimvert); | |
202 | if(errorsprimvert[0]==0 || errorsprimvert[1]==0){ | |
203 | Warning("FindTracks","Set errors on vertex positions x and y at 0.0001"); | |
204 | errorsprimvert[0]=0.0001; | |
205 | errorsprimvert[1]=0.0001; | |
206 | } | |
207 | fVert->PrintStatus(); | |
208 | ||
209 | SetEventNumber(evnumber); | |
210 | if(GetFlagLoadedClusters()==0) LoadClusters(clusterTree); | |
211 | ||
212 | //Fill tree for found tracks | |
213 | AliITStrackV2* outrack=0; | |
214 | TBranch* branch=out->Branch("Tracks","AliITStrackV2",&outrack,32000,0); | |
215 | if (!branch) out->Branch("tracks","AliITStrackV2",&outrack,32000,3); | |
216 | else branch->SetAddress(&outrack); | |
217 | ||
218 | ||
219 | Int_t firstmod[fGeom->GetNlayers()]; | |
220 | for(Int_t i=0;i<fGeom->GetNlayers();i++){ | |
221 | firstmod[i]=fGeom->GetModuleIndex(i+1,1,1); | |
222 | } | |
223 | // firstmod [i] number of the first module in the ITS layer i. | |
224 | ||
225 | AliITSlayer &layer=fgLayers[0]; // first layer | |
226 | Int_t ntrack=0; | |
227 | //loop on the different windows | |
228 | for(Int_t nloop=0;nloop<fNloop;nloop++){ | |
229 | Int_t ncl=layer.GetNumberOfClusters(); | |
230 | while(ncl--){ //loop starting from layer 0 | |
231 | ResetForFinding(); | |
232 | Int_t pflag=0; | |
233 | AliITSclusterV2* cl = layer.GetCluster(ncl); | |
234 | if(cl->IsUsed()==1) continue; | |
235 | Double_t x,y,z; | |
236 | Int_t module = cl->GetDetectorIndex()+firstmod[0]; | |
237 | GetClusterGCoordinates(cl,module,x,y,z); | |
238 | fPhic = TMath::ATan2(y,x); | |
239 | fLambdac = TMath::ATan2(z-primaryVertex[2],TMath::Sqrt((x-primaryVertex[0])*(x-primaryVertex[0])+(y-primaryVertex[1])*(y-primaryVertex[1]))); | |
240 | fPhiEstimate = fPhic; | |
241 | AliITStrackSA* trs = new AliITStrackSA(); | |
242 | fPoint1[0]=primaryVertex[0]; | |
243 | fPoint1[1]=primaryVertex[1]; | |
244 | fPoint2[0]=x; | |
245 | fPoint2[1]=y; | |
246 | ||
247 | Int_t nn[fGeom->GetNlayers()];//counter for clusters on each layer | |
248 | for(Int_t i=0;i<fGeom->GetNlayers();i++){ nn[i]=0;} | |
249 | nn[0] = SearchClusters(0,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],primaryVertex[1],primaryVertex[0],pflag,fTable); | |
250 | nn[1] = SearchClusters(1,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],primaryVertex[1],primaryVertex[0],pflag,fTable); | |
251 | ||
252 | if(nn[1]>0){ | |
253 | pflag=1; | |
254 | fPoint3[0] = fPointc[0]; | |
255 | fPoint3[1] = fPointc[1]; | |
256 | } | |
257 | nn[2] = SearchClusters(2,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],primaryVertex[1],primaryVertex[0],pflag,fTable); | |
258 | if(nn[1]==0 && nn[2]==0) pflag=0; | |
259 | if(nn[2]!=0 && nn[1]!=0){ pflag=1; UpdatePoints();} | |
260 | if(nn[2]!=0 && nn[1]==0){ | |
261 | pflag=1; | |
262 | fPoint3[0]=fPointc[0]; | |
263 | fPoint3[1]=fPointc[1]; | |
264 | } | |
265 | ||
266 | nn[3] = SearchClusters(3,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],primaryVertex[1],primaryVertex[0],pflag,fTable); | |
267 | pflag=1; | |
268 | if(nn[3]!=0) UpdatePoints(); | |
269 | nn[4] = SearchClusters(4,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],primaryVertex[1],primaryVertex[0],pflag,fTable); | |
270 | pflag=1; | |
271 | if(nn[4]!=0) UpdatePoints(); | |
272 | nn[5] = SearchClusters(5,fPhiWin[nloop],fLambdaWin[nloop],trs,primaryVertex[2],primaryVertex[1],primaryVertex[0],pflag,fTable); | |
273 | ||
274 | ||
275 | Int_t layOK=0; | |
276 | Int_t numberofpoints; | |
277 | if(sixpoints) numberofpoints=6; //check of the candidate track | |
278 | else numberofpoints=5; //if track is good (with the required number | |
279 | for(Int_t nnp=0;nnp<fGeom->GetNlayers();nnp++){ //of points) it is written on file | |
280 | if(nn[nnp]!=0) layOK+=1; | |
281 | } | |
282 | if(layOK>=numberofpoints){ | |
283 | AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert,opt); | |
284 | if(tr2==0){ | |
285 | Int_t nct = trs->GetNumberOfClustersSA(); | |
286 | while(nct--){ | |
287 | Int_t index = trs->GetClusterIndexSA(nct); | |
288 | AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index); | |
289 | if(kl->IsUsed()==1) kl->Use(); | |
290 | } | |
291 | continue; | |
292 | } | |
293 | outrack=tr2; | |
294 | out->Fill(); | |
295 | ntrack++; | |
296 | Int_t nct = tr2->GetNumberOfClusters(); | |
297 | ||
298 | while(nct--){ | |
299 | Int_t index = tr2->GetClusterIndex(nct); | |
300 | AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index); | |
301 | if(kl->IsUsed()==0) kl->Use(); | |
302 | ||
303 | } | |
304 | } | |
305 | else{ | |
306 | Int_t nct = trs->GetNumberOfClustersSA(); | |
307 | while(nct--){ | |
308 | Int_t index = trs->GetClusterIndexSA(nct); | |
309 | AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index); | |
310 | if(kl->IsUsed()==1) kl->Use(); | |
311 | } | |
312 | } | |
313 | delete trs; | |
314 | ||
315 | }//end loop on clusters of layer1 | |
316 | ||
317 | }//end loop2 | |
318 | ||
319 | //if 5/6 points are required, second loop starting | |
320 | //from second layer, to find tracks with point of | |
321 | //layer 1 missing | |
322 | ||
323 | if(!sixpoints){ | |
324 | // counter for clusters on each layer | |
325 | Int_t *nn = new Int_t[fGeom->GetNlayers()-1]; | |
326 | for(Int_t nloop=0;nloop<fNloop;nloop++){ | |
327 | AliITSlayer &layer2=fgLayers[1]; //loop on layer 2 | |
328 | Int_t ncl2=layer2.GetNumberOfClusters(); | |
329 | while(ncl2--){ //loop starting from layer 2 | |
330 | ResetForFinding(); | |
331 | Int_t pflag=0; | |
332 | AliITSclusterV2* cl = layer2.GetCluster(ncl2); | |
333 | if(cl->IsUsed()==1) continue; | |
334 | Double_t x,y,z; | |
335 | Int_t module = cl->GetDetectorIndex()+firstmod[1]; | |
336 | GetClusterGCoordinates(cl,module,x,y,z); | |
337 | fPhic = TMath::ATan2(y,x); | |
338 | fLambdac = TMath::ATan2(z-primaryVertex[2],TMath::Sqrt((x-primaryVertex[0])*(x-primaryVertex[0])+(y-primaryVertex[1])*(y-primaryVertex[1]))); | |
339 | ||
340 | fPhiEstimate = fPhic; | |
341 | AliITStrackSA* trs = new AliITStrackSA(); | |
342 | ||
343 | fPoint1[0]=primaryVertex[0]; | |
344 | fPoint1[1]=primaryVertex[1]; | |
345 | fPoint2[0]=x; | |
346 | fPoint2[1]=y; | |
347 | ||
348 | for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++)nn[kk] = 0; | |
349 | for(Int_t kk=0;kk<fGeom->GetNlayers()-1;kk++){ | |
350 | nn[kk] = SearchClusters(kk+1,fPhiWin[nloop],fLambdaWin[nloop], | |
351 | trs,primaryVertex[2],primaryVertex[1],primaryVertex[0],pflag,fTable); | |
352 | if(nn[kk]==0)break; | |
353 | if(kk>0){ | |
354 | UpdatePoints(); | |
355 | pflag = 1; | |
356 | } | |
357 | } | |
358 | Int_t fl=0; | |
359 | for(Int_t nnp=0;nnp<fGeom->GetNlayers()-1;nnp++){ | |
360 | if(nn[nnp]!=0) fl+=1; | |
361 | } | |
362 | if(fl>=5){ // 5/6 | |
363 | AliITStrackV2* tr2 = FitTrack(trs,primaryVertex,errorsprimvert,opt); | |
364 | if(tr2==0){ | |
365 | Int_t nct = trs->GetNumberOfClustersSA(); | |
366 | while(nct--){ | |
367 | Int_t index = trs->GetClusterIndexSA(nct); | |
368 | AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index); | |
369 | if(kl->IsUsed()==1) kl->Use(); | |
370 | } | |
371 | continue; | |
372 | } | |
373 | outrack=tr2; | |
374 | out->Fill(); | |
375 | Int_t nct = tr2->GetNumberOfClusters(); | |
376 | while(nct--){ | |
377 | Int_t index = tr2->GetClusterIndex(nct); | |
378 | AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index); | |
379 | if(kl==0) continue; | |
380 | if(kl->IsUsed()==0) kl->Use(); | |
381 | } | |
382 | } | |
383 | else{ | |
384 | Int_t nct = trs->GetNumberOfClustersSA(); | |
385 | while(nct--){ | |
386 | Int_t index = trs->GetClusterIndexSA(nct); | |
387 | AliITSclusterV2* kl = (AliITSclusterV2*)GetCluster(index); | |
388 | if(kl==0) continue; | |
389 | if(kl->IsUsed()==1) kl->Use(); | |
390 | } | |
391 | } | |
392 | delete trs; | |
393 | }//end loop on clusters of layer2 | |
394 | } | |
395 | delete []nn; | |
396 | } //end opt="5/6" | |
397 | ||
398 | ||
399 | UnloadClusters(); | |
400 | } | |
401 | //_______________________________________________________________________ | |
402 | void AliITStrackerSA::UseFoundTracksV2(Int_t evnum,TTree* treev2, TTree* clustertree){ | |
403 | // Marks as used clusters belonging to tracks found with V2 TPC+ITS tracking | |
404 | if(!fTable){ | |
405 | fTable = new AliITSclusterTable(fGeom,this); | |
406 | fTable->FillArray(clustertree,evnum); | |
407 | } | |
408 | LoadClusters(clustertree); | |
409 | SetFlagLoadedClusters(1); | |
410 | ||
411 | TBranch* bra = (TBranch*)treev2->GetBranch("tracks"); | |
412 | if(!bra) Warning("UseFoundTracksV2","No branch for track tree"); | |
413 | AliITStrackV2* ttrrt = new AliITStrackV2; | |
414 | bra->SetAddress(&ttrrt); | |
415 | ||
416 | for(Int_t nj=0;nj<treev2->GetEntries();nj++){ | |
417 | treev2->GetEvent(nj); | |
418 | Int_t ncl = ttrrt->GetNumberOfClusters(); | |
419 | for(Int_t k=0;k<ncl;k++){ | |
420 | Int_t index = ttrrt->GetClusterIndex(k); | |
421 | AliITSclusterV2* clui = (AliITSclusterV2*)GetCluster(index); | |
422 | if(clui->IsUsed()==0) clui->Use(); | |
423 | } | |
424 | } | |
425 | delete ttrrt; | |
426 | ||
427 | } | |
428 | ||
429 | //________________________________________________________________________ | |
430 | ||
431 | AliITStrackV2* AliITStrackerSA::FitTrack(AliITStrackSA* tr,Double_t *primaryVertex,Double_t *errorsprimvert,char *opt){ | |
432 | //fit of the found track | |
433 | ||
434 | //options | |
435 | TString choice(opt); | |
436 | Bool_t sixpoints= choice.Contains("6/6"); | |
437 | ||
438 | Int_t firstmod[fGeom->GetNlayers()]; | |
439 | for(Int_t i=0;i<fGeom->GetNlayers();i++){ | |
440 | firstmod[i]=fGeom->GetModuleIndex(i+1,1,1); | |
441 | } | |
442 | TObjArray* listoftracks=new TObjArray(0,0); | |
443 | AliITStrackV2* otrack2; | |
444 | Int_t nclusters = tr->GetNumberOfClustersSA(); | |
445 | TObjArray** listlayer = new TObjArray*[fGeom->GetNlayers()]; | |
446 | for(Int_t i=0;i<fGeom->GetNlayers();i++){ | |
447 | listlayer[i] = new TObjArray(0,0); | |
448 | } | |
449 | ||
450 | TArrayI clind0(20); | |
451 | TArrayI clind1(20); | |
452 | TArrayI clind2(20); | |
453 | TArrayI clind3(20); | |
454 | TArrayI clind4(20); | |
455 | TArrayI clind5(20); | |
456 | ||
457 | Int_t nnn[fGeom->GetNlayers()]; | |
458 | for(Int_t i=0;i<fGeom->GetNlayers();i++)nnn[i]=0; | |
459 | ||
460 | for(Int_t ncl=0;ncl<nclusters;ncl++){ | |
461 | Int_t index = tr->GetClusterIndexSA(ncl); | |
462 | AliITSclusterV2* cl = (AliITSclusterV2*)GetCluster(index); | |
463 | if(cl->IsUsed()==1) cl->Use(); | |
464 | Int_t lay = (index & 0xf0000000) >> 28; | |
465 | if(lay==0) { listlayer[0]->AddLast(cl); clind0[nnn[0]]=index;nnn[0]++;} | |
466 | if(lay==1) { listlayer[1]->AddLast(cl); clind1[nnn[1]]=index;nnn[1]++;} | |
467 | if(lay==2) { listlayer[2]->AddLast(cl); clind2[nnn[2]]=index;nnn[2]++;} | |
468 | if(lay==3) { listlayer[3]->AddLast(cl); clind3[nnn[3]]=index;nnn[3]++;} | |
469 | if(lay==4) { listlayer[4]->AddLast(cl); clind4[nnn[4]]=index;nnn[4]++;} | |
470 | if(lay==5) { listlayer[5]->AddLast(cl); clind5[nnn[5]]=index;nnn[5]++;} | |
471 | } | |
472 | ||
473 | ||
474 | Int_t end[fGeom->GetNlayers()]; | |
475 | for(Int_t i=0;i<fGeom->GetNlayers();i++){ | |
476 | if(listlayer[i]->GetEntries()==0) end[i]=1; | |
477 | else end[i]=listlayer[i]->GetEntries(); | |
478 | } | |
479 | ||
480 | for(Int_t l1=0;l1<end[0];l1++){//loop on layer 1 | |
481 | AliITSclusterV2* cl0 = (AliITSclusterV2*)listlayer[0]->At(l1); | |
482 | TVector3** recp = new TVector3*[3]; | |
483 | TVector3** errs = new TVector3*[3]; | |
484 | recp[0] = new TVector3(primaryVertex[0],primaryVertex[1],primaryVertex[2]); | |
485 | errs[0] = new TVector3(errorsprimvert[0],errorsprimvert[1],errorsprimvert[2]); | |
486 | Double_t x1,y1,z1,sx1,sy1,sz1; | |
487 | Double_t x2,y2,z2,sx2,sy2,sz2; | |
488 | AliITSclusterV2* p1=0; | |
489 | AliITSclusterV2* p2=0; | |
490 | ||
491 | for(Int_t l2=0;l2<end[1];l2++){//loop on layer 2 | |
492 | AliITSclusterV2* cl1 = (AliITSclusterV2*)listlayer[1]->At(l2); | |
493 | for(Int_t l3=0;l3<end[2];l3++){ //loop on layer 3 | |
494 | AliITSclusterV2* cl2 = (AliITSclusterV2*)listlayer[2]->At(l3); | |
495 | ||
496 | if(cl0==0 && cl1!=0) { | |
497 | p2 = cl2; | |
498 | p1=cl1; | |
499 | ||
500 | } | |
501 | if(cl0!=0 && cl1==0){ | |
502 | p1=cl0; | |
503 | p2=cl2; | |
504 | } | |
505 | if(cl0!=0 && cl1!=0){ | |
506 | p1=cl0; | |
507 | p2=cl1; | |
508 | } | |
509 | Int_t module1 = p1->GetDetectorIndex()+firstmod[0]; | |
510 | Int_t module2 = p2->GetDetectorIndex()+firstmod[1]; | |
511 | GetClusterGCoordinates(p1,module1,x1,y1,z1); | |
512 | GetClusterGCoordinates(p2,module2,x2,y2,z2); | |
513 | GetClusterGErrors(p1,module1,sx1,sy1,sz1); | |
514 | GetClusterGErrors(p2,module2,sx2,sy2,sz2); | |
515 | Double_t phi1=TMath::ATan2(y1,x1); | |
516 | recp[1] = new TVector3(x1,y1,z1); | |
517 | errs[1] = new TVector3(sx1,sy1,sz1); | |
518 | recp[2] = new TVector3(x2,y2,z2); | |
519 | errs[2] = new TVector3(sx2,sy2,sz2); | |
520 | ||
521 | //fit on the Riemann sphere | |
522 | Float_t seed1,seed2,seed3; | |
523 | AliITSRiemannFit fit; | |
524 | Int_t rf = fit.FitHelix(3,recp,errs,seed1,seed2,seed3); //this gives phi,tgl,curvature to start Kalman Filter | |
525 | if(rf==0) continue; | |
526 | Double_t phi=seed1; | |
527 | Double_t tgl=seed2; | |
528 | ||
529 | if(phi1>0){ | |
530 | if(seed1>-TMath::Pi() && seed1<-0.5*TMath::Pi()){ | |
531 | phi=seed1+1.5*TMath::Pi(); | |
532 | tgl=seed2; | |
533 | } | |
534 | if(seed1>-0.5*TMath::Pi() && seed1<0.5*TMath::Pi()){ | |
535 | phi=seed1+0.5*TMath::Pi(); | |
536 | tgl=(-1)*seed2; | |
537 | } | |
538 | if(seed1>0.5*TMath::Pi() && seed1<TMath::Pi()){ | |
539 | phi=seed1-0.5*TMath::Pi(); | |
540 | tgl=seed2; | |
541 | } | |
542 | } | |
543 | if(phi1<0){ | |
544 | if(seed1>-TMath::Pi() && seed1<-0.5*TMath::Pi()){ | |
545 | phi=seed1+0.5*TMath::Pi(); | |
546 | tgl=(-1)*seed2; | |
547 | } | |
548 | if(seed1>-0.5*TMath::Pi() && seed1<0.5*TMath::Pi()){ | |
549 | phi=seed1-0.5*TMath::Pi(); | |
550 | tgl=seed2; | |
551 | } | |
552 | if(seed1>0.5*TMath::Pi() && seed1<TMath::Pi()){ | |
553 | phi=seed1-1.5*TMath::Pi(); | |
554 | tgl=(-1)*seed2; | |
555 | } | |
556 | } | |
557 | ||
558 | Int_t layer,ladder,detector; | |
559 | fGeom->GetModuleId(module1,layer,ladder,detector); | |
560 | Float_t yclu1 = p1->GetY(); | |
561 | Float_t zclu1 = p1->GetZ(); | |
562 | Double_t cv=Curvature(primaryVertex[0],primaryVertex[1],x1,y1,x2,y2); | |
563 | ||
564 | for(Int_t l4=0;l4<end[3];l4++){ //loop on layer 4 | |
565 | AliITSclusterV2* cl3 = (AliITSclusterV2*)listlayer[3]->At(l4); | |
566 | for(Int_t l5=0;l5<end[4];l5++){ //loop on layer 5 | |
567 | AliITSclusterV2* cl4 = (AliITSclusterV2*)listlayer[4]->At(l5); | |
568 | for(Int_t l6=0;l6<end[5];l6++){ //loop on layer 6 | |
569 | AliITSclusterV2* cl5 = (AliITSclusterV2*)listlayer[5]->At(l6); | |
570 | AliITStrackSA* trac = new AliITStrackSA(layer,ladder,detector,yclu1,zclu1,phi,tgl,cv,1); | |
571 | ||
572 | if(cl5!=0) trac->AddClusterV2(5,(clind5[l6] & 0x0fffffff)>>0); | |
573 | if(cl4!=0) trac->AddClusterV2(4,(clind4[l5] & 0x0fffffff)>>0); | |
574 | if(cl3!=0) trac->AddClusterV2(3,(clind3[l4] & 0x0fffffff)>>0); | |
575 | if(cl2!=0) trac->AddClusterV2(2,(clind2[l3] & 0x0fffffff)>>0); | |
576 | if(cl1!=0) trac->AddClusterV2(1,(clind1[l2] & 0x0fffffff)>>0); | |
577 | if(cl0!=0) trac->AddClusterV2(0,(clind0[l1] & 0x0fffffff)>>0); | |
578 | ||
579 | //fit with Kalman filter using AliITStrackerV2::RefitAt() | |
580 | ||
581 | AliITStrackV2* ot = new AliITStrackV2(*trac); | |
582 | ||
583 | ot->ResetCovariance(); | |
584 | ot->ResetClusters(); | |
585 | ||
586 | if(RefitAt(49.,ot,trac)){ //fit from layer 1 to layer 6 | |
587 | ||
588 | otrack2 = new AliITStrackV2(*ot); | |
589 | otrack2->ResetCovariance(); | |
590 | otrack2->ResetClusters(); | |
591 | //fit from layer 6 to layer 1 | |
592 | if(RefitAt(3.7,otrack2,ot)) listoftracks->AddLast(otrack2); | |
593 | ||
594 | } | |
595 | ||
596 | delete ot; | |
597 | delete trac; | |
598 | }//end loop layer 6 | |
599 | }//end loop layer 5 | |
600 | }//end loop layer 4 | |
601 | ||
602 | for(Int_t i=1;i<3;i++){ | |
603 | delete recp[i]; | |
604 | delete errs[i]; | |
605 | } | |
606 | }//end loop layer 3 | |
607 | }//end loop layer 2 | |
608 | delete recp[0]; | |
609 | delete errs[0]; | |
610 | delete[] recp; | |
611 | delete[] errs; | |
612 | }//end loop layer 1 | |
613 | ||
614 | Int_t dim=listoftracks->GetEntries(); | |
615 | if(dim==0){ | |
616 | delete listoftracks; | |
617 | for(Int_t i=0;i<fGeom->GetNlayers();i++){ | |
618 | delete listlayer[i]; | |
619 | } | |
620 | delete listlayer; | |
621 | return 0; | |
622 | } | |
623 | ||
624 | AliITStrackV2* otrack =(AliITStrackV2*)FindTrackLowChiSquare(listoftracks,dim); | |
625 | ||
626 | if(otrack==0) return 0; | |
627 | Int_t indexc[fGeom->GetNlayers()]; | |
628 | for(Int_t i=0;i<fGeom->GetNlayers();i++) indexc[i]=0; | |
629 | for(Int_t nind=0;nind<otrack->GetNumberOfClusters();nind++){ | |
630 | indexc[nind] = otrack->GetClusterIndex(nind); | |
631 | } | |
632 | AliITSclusterV2* cl0 = (AliITSclusterV2*)GetCluster(indexc[0]); | |
633 | AliITSclusterV2* cl1 = (AliITSclusterV2*)GetCluster(indexc[1]); | |
634 | AliITSclusterV2* cl2 = (AliITSclusterV2*)GetCluster(indexc[2]); | |
635 | AliITSclusterV2* cl3 = (AliITSclusterV2*)GetCluster(indexc[3]); | |
636 | AliITSclusterV2* cl4 = (AliITSclusterV2*)GetCluster(indexc[4]); | |
637 | Int_t labl[3]={-1,-1,-1}; | |
638 | if(otrack->GetNumberOfClusters()==fGeom->GetNlayers()){ | |
639 | AliITSclusterV2* cl5 = (AliITSclusterV2*)GetCluster(indexc[5]); | |
640 | labl[0]=cl5->GetLabel(0); | |
641 | labl[1]=cl5->GetLabel(1); | |
642 | labl[2]=cl5->GetLabel(2); | |
643 | } | |
644 | if(otrack->GetNumberOfClusters()==(fGeom->GetNlayers()-1)){ | |
645 | labl[0]=-1; | |
646 | labl[1]=-1; | |
647 | labl[2]=-1; | |
648 | } | |
649 | Int_t numberofpoints; | |
650 | if(sixpoints) numberofpoints=6; | |
651 | else numberofpoints=5; | |
652 | Int_t label = Label(cl0->GetLabel(0),cl1->GetLabel(0), | |
653 | cl2->GetLabel(0),cl3->GetLabel(0), | |
654 | cl4->GetLabel(0),labl[0], | |
655 | cl0->GetLabel(1),cl1->GetLabel(1), | |
656 | cl2->GetLabel(1),cl3->GetLabel(1), | |
657 | cl4->GetLabel(1),labl[1], | |
658 | cl0->GetLabel(2),cl1->GetLabel(2), | |
659 | cl2->GetLabel(2),cl3->GetLabel(2), | |
660 | cl4->GetLabel(2),labl[2],numberofpoints); | |
661 | ||
662 | otrack->SetLabel(label); | |
663 | delete listoftracks; | |
664 | for(Int_t i=0;i<fGeom->GetNlayers();i++){ | |
665 | delete listlayer[i]; | |
666 | } | |
667 | delete listlayer; | |
668 | return otrack; | |
669 | ||
670 | } | |
671 | ||
672 | //________________________________________________________________ | |
673 | void AliITStrackerSA::UpdatePoints(){ | |
674 | //update of points for the estimation of the curvature | |
675 | ||
676 | fPoint1[0]=fPoint2[0]; | |
677 | fPoint1[1]=fPoint2[1]; | |
678 | fPoint2[0]=fPoint3[0]; | |
679 | fPoint2[1]=fPoint3[1]; | |
680 | fPoint3[0]=fPointc[0]; | |
681 | fPoint3[1]=fPointc[1]; | |
682 | ||
683 | ||
684 | } | |
685 | ||
686 | //_________________________________________________________________ | |
687 | ||
688 | Int_t AliITStrackerSA::SearchClusters(Int_t layer,Double_t phiwindow,Double_t lambdawindow, AliITStrackSA* trs,Double_t zvertex,Double_t yvertex,Double_t xvertex,Int_t pflag,AliITSclusterTable* table){ | |
689 | //function used to to find the clusters associated to the track | |
690 | Int_t nc=0; | |
691 | AliITSlayer &lay = fgLayers[layer]; | |
692 | Int_t firstmod[fGeom->GetNlayers()]; | |
693 | for(Int_t i=0;i<fGeom->GetNlayers();i++){ | |
694 | firstmod[i]=fGeom->GetModuleIndex(i+1,1,1); | |
695 | } | |
696 | if(pflag==1){ | |
697 | ||
698 | Float_t cx1,cx2,cy1,cy2; | |
699 | FindEquation(fPoint1[0],fPoint1[1],fPoint2[0],fPoint2[1],fPoint3[0],fPoint3[1],fCoef1,fCoef2,fCoef3); | |
700 | Int_t fun = FindIntersection(fCoef1,fCoef2,fCoef3,-(lay.GetR()*lay.GetR()),cx1,cy1,cx2,cy2); | |
701 | if(fun==0) return 0; | |
702 | Double_t fi1 =TMath::ATan2(cy1,cx1); | |
703 | Double_t fi2 =TMath::ATan2(cy2,cx2); | |
704 | fPhiEstimate = ChoosePoint(fi1,fi2,fPhic); | |
705 | } | |
706 | ||
707 | Double_t zed = TMath::Tan(fLambdac)*lay.GetR()+zvertex; | |
708 | Double_t zed1 = TMath::Tan(fLambdac+lambdawindow)*lay.GetR()+zvertex; | |
709 | Double_t zed2 = TMath::Tan(fLambdac-lambdawindow)*lay.GetR()+zvertex; | |
710 | ||
711 | Double_t fi = fPhiEstimate; | |
712 | Int_t nmod = lay.FindDetectorIndex(fi,zed)+firstmod[layer]; | |
713 | ||
714 | Int_t nm[8]={0,0,0,0,0,0,0,0}; | |
715 | nm[0] = lay.FindDetectorIndex(fi+phiwindow,zed)+firstmod[layer]; | |
716 | nm[1] = lay.FindDetectorIndex(fi-phiwindow,zed)+firstmod[layer]; | |
717 | nm[2] = lay.FindDetectorIndex(fi,zed1)+firstmod[layer]; | |
718 | nm[3] = lay.FindDetectorIndex(fi,zed2)+firstmod[layer]; | |
719 | nm[4] = lay.FindDetectorIndex(fi+phiwindow,zed1)+firstmod[layer]; | |
720 | nm[5] = lay.FindDetectorIndex(fi-phiwindow,zed1)+firstmod[layer]; | |
721 | nm[6] = lay.FindDetectorIndex(fi+phiwindow,zed2)+firstmod[layer]; | |
722 | nm[7] = lay.FindDetectorIndex(fi-phiwindow,zed2)+firstmod[layer]; | |
723 | ||
724 | ||
725 | if(nmod<0) return 0; | |
726 | Int_t nn=0; | |
727 | TArrayI* array =(TArrayI*)table->GetListOfClusters(nmod); | |
728 | TArrayI* list = new TArrayI(array->GetSize()); | |
729 | for(Int_t i=0;i<array->GetSize();i++){ | |
730 | Int_t in=(Int_t)array->At(i); | |
731 | list->AddAt(in,nn); | |
732 | nn++; | |
733 | } | |
734 | ||
735 | for(Int_t ii=0;ii<8;ii++){ | |
736 | if(nm[ii]!=nmod && nm[ii]>=0){ | |
737 | TArrayI* ar =(TArrayI*)table->GetListOfClusters(nm[ii]); | |
738 | list->Set(list->GetSize()+ar->GetSize()); | |
739 | for(Int_t j=0;j<ar->GetSize();j++){ | |
740 | Int_t in=(Int_t)ar->At(j); | |
741 | list->AddAt(in,nn); | |
742 | nn++; | |
743 | } | |
744 | } | |
745 | } | |
746 | for(Int_t i=0;i<list->GetSize();i++){ | |
747 | Int_t index = (Int_t)list->At(i); | |
748 | AliITSclusterV2* cllay = lay.GetCluster(index); | |
749 | if(cllay==0) continue; | |
750 | if(cllay->IsUsed()==1) continue; | |
751 | ||
752 | Double_t x,y,z; | |
753 | Int_t modlay = cllay->GetDetectorIndex()+firstmod[layer]; | |
754 | GetClusterGCoordinates(cllay,modlay,x,y,z); | |
755 | Double_t phi = TMath::ATan2(y,x); | |
756 | Double_t lambda=TMath::ATan2(z-zvertex,TMath::Sqrt((x-xvertex)*(x-xvertex)+(y-yvertex)*(y-yvertex))); | |
757 | ||
758 | ||
759 | if(TMath::Abs(fLambdac-lambda)<lambdawindow && | |
760 | TMath::Abs(fPhiEstimate-phi)<phiwindow){ | |
761 | nc+=1; | |
762 | fLambdac = lambda; | |
763 | if(trs->GetNumberOfClustersSA()==20){ | |
764 | delete list; | |
765 | return 0; | |
766 | } | |
767 | trs->AddClusterSA(layer,index); | |
768 | cllay->Use(); | |
769 | fPhiEstimate=phi; | |
770 | fPointc[0]=x;fPointc[1]=y; | |
771 | } | |
772 | ||
773 | } | |
774 | delete list; | |
775 | return nc; | |
776 | ||
777 | } | |
778 | ||
779 | //__________________________________________________________________ | |
780 | void AliITStrackerSA::GetClusterGCoordinates(AliITSclusterV2* cluster,Int_t module,Double_t& x, Double_t& y, Double_t& z){ | |
781 | ||
782 | //returns x,y,z of cluster in global coordinates | |
783 | ||
784 | Double_t rot[9]; fGeom->GetRotMatrix(module,rot); | |
785 | Int_t lay,lad,det; fGeom->GetModuleId(module,lay,lad,det); | |
786 | Float_t tx,ty,tz; fGeom->GetTrans(lay,lad,det,tx,ty,tz); | |
787 | ||
788 | Double_t alpha=TMath::ATan2(rot[1],rot[0])+TMath::Pi(); | |
789 | Double_t phi=TMath::Pi()/2+alpha; | |
790 | if (lay==1) phi+=TMath::Pi(); | |
791 | ||
792 | Double_t cp=TMath::Cos(phi), sp=TMath::Sin(phi); | |
793 | Double_t r=tx*cp+ty*sp; | |
794 | ||
795 | x= r*cp - cluster->GetY()*sp; | |
796 | y= r*sp + cluster->GetY()*cp; | |
797 | z=cluster->GetZ(); | |
798 | ||
799 | ||
800 | ||
801 | } | |
802 | //__________________________________________________________________ | |
803 | void AliITStrackerSA::GetClusterGErrors(AliITSclusterV2* cluster,Int_t module,Double_t& sigmax, Double_t& sigmay, Double_t& sigmaz){ | |
804 | ||
805 | //returns x,y,z of cluster in global coordinates | |
806 | ||
807 | Double_t rot[9]; fGeom->GetRotMatrix(module,rot); | |
808 | Int_t lay,lad,det; fGeom->GetModuleId(module,lay,lad,det); | |
809 | ||
810 | Double_t alpha=TMath::ATan2(rot[1],rot[0])+TMath::Pi(); | |
811 | Double_t phi=TMath::Pi()/2+alpha; | |
812 | if (lay==1) phi+=TMath::Pi(); | |
813 | ||
814 | Double_t cp=TMath::Cos(phi), sp=TMath::Sin(phi); | |
815 | ||
816 | ||
817 | sigmax = TMath::Sqrt(sp*sp*cluster->GetSigmaY2()); | |
818 | sigmay = TMath::Sqrt(cp*cp*cluster->GetSigmaY2()); | |
819 | sigmaz = TMath::Sqrt(cluster->GetSigmaZ2()); | |
820 | } | |
821 | ||
822 | //___________________________________________________________________ | |
823 | Int_t AliITStrackerSA::FindEquation(Float_t x1, Float_t y1, Float_t x2, Float_t y2, Float_t x3, Float_t y3,Float_t& a, Float_t& b, Float_t& c){ | |
824 | ||
825 | //given (x,y) of three recpoints (in global coordinates) | |
826 | //returns the parameters a,b,c of circonference x*x + y*y +a*x + b*y +c | |
827 | ||
828 | Float_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1); | |
829 | if(den==0) return 0; | |
830 | a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den; | |
831 | b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1); | |
832 | c = -x1*x1-y1*y1-a*x1-b*y1; | |
833 | return 1; | |
834 | } | |
835 | //__________________________________________________________________________ | |
836 | Int_t AliITStrackerSA::FindIntersection(Float_t a1, Float_t b1, Float_t c1, Float_t c2,Float_t& x1,Float_t& y1, Float_t& x2, Float_t& y2){ | |
837 | ||
838 | //Finds the intersection between the circonference of the track and the circonference centered in (0,0) represented by one layer | |
839 | //c2 is -rlayer*rlayer | |
840 | ||
841 | if(a1==0) return 0; | |
842 | Float_t m = c2-c1; | |
843 | Float_t aA = (b1*b1)/(a1*a1)+1; | |
844 | Float_t bB = (-2*m*b1/(a1*a1)); | |
845 | Float_t cC = c2+(m*m)/(a1*a1); | |
846 | if((bB*bB-4*aA*cC)<0) return 0; | |
847 | ||
848 | y1 = (-bB+TMath::Sqrt(bB*bB-4*aA*cC))/(2*aA); | |
849 | y2 = (-bB-TMath::Sqrt(bB*bB-4*aA*cC))/(2*aA); | |
850 | x1 = (c2-c1-b1*y1)/a1; | |
851 | x2 = (c2-c1-b1*y2)/a1; | |
852 | ||
853 | return 1; | |
854 | } | |
855 | //____________________________________________________________________ | |
856 | Double_t AliITStrackerSA::Curvature(Double_t x1,Double_t y1,Double_t | |
857 | x2,Double_t y2,Double_t x3,Double_t y3){ | |
858 | ||
859 | //calculates the curvature of track | |
860 | Double_t den = (x3-x1)*(y2-y1)-(x2-x1)*(y3-y1); | |
861 | if(den==0) return 0; | |
862 | Double_t a = ((y3-y1)*(x2*x2+y2*y2-x1*x1-y1*y1)-(y2-y1)*(x3*x3+y3*y3-x1*x1-y1*y1))/den; | |
863 | Double_t b = -(x2*x2-x1*x1+y2*y2-y1*y1+a*(x2-x1))/(y2-y1); | |
864 | Double_t c = -x1*x1-y1*y1-a*x1-b*y1; | |
865 | Double_t xc=-a/2.; | |
866 | ||
867 | if((a*a+b*b-4*c)<0) return 0; | |
868 | Double_t rad = TMath::Sqrt(a*a+b*b-4*c)/2.; | |
869 | if(rad==0) return 0; | |
870 | ||
871 | if((x1>0 && y1>0 && x1<xc)) rad*=-1; | |
872 | if((x1<0 && y1>0 && x1<xc)) rad*=-1; | |
873 | // if((x1<0 && y1<0 && x1<xc)) rad*=-1; | |
874 | // if((x1>0 && y1<0 && x1<xc)) rad*=-1; | |
875 | ||
876 | return 1/rad; | |
877 | ||
878 | } | |
879 | //____________________________________________________________________ | |
880 | Double_t AliITStrackerSA::ChoosePoint(Double_t p1, Double_t p2, Double_t pp){ | |
881 | ||
882 | //Returns the point closest to pp | |
883 | ||
884 | Double_t diff1 = p1-pp; | |
885 | Double_t diff2 = p2-pp; | |
886 | ||
887 | if(TMath::Abs(diff1)<TMath::Abs(diff2)) fPhiEstimate=p1; | |
888 | else fPhiEstimate=p2; | |
889 | return fPhiEstimate; | |
890 | ||
891 | } | |
892 | ||
893 | ||
894 | //_________________________________________________________________ | |
895 | AliITStrackV2* AliITStrackerSA::FindTrackLowChiSquare(TObjArray* tracklist, Int_t dim) const { | |
896 | // returns track with lowes chi square | |
897 | if(dim==1){ | |
898 | AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(0); | |
899 | return trk; | |
900 | } | |
901 | if(dim==0) return 0; | |
902 | Double_t chi2[dim]; | |
903 | Int_t index[dim]; | |
904 | for(Int_t i=0;i<dim;i++){ | |
905 | AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(i); | |
906 | chi2[i]=trk->GetChi2(); | |
907 | index[i]=i; | |
908 | } | |
909 | ||
910 | Int_t w=0;Double_t value; | |
911 | Int_t lp; | |
912 | while(w<dim){ | |
913 | for(Int_t j=w+1;j<dim;j++){ | |
914 | if(chi2[w]<chi2[j]){ | |
915 | value=chi2[w]; | |
916 | chi2[w]=chi2[j]; | |
917 | chi2[j]=value; | |
918 | lp=index[w]; | |
919 | index[w]=index[j]; | |
920 | index[j]=lp; | |
921 | } | |
922 | } | |
923 | w++; | |
924 | } | |
925 | ||
926 | AliITStrackV2* trk = (AliITStrackV2*)tracklist->At(index[dim-1]); | |
927 | return trk; | |
928 | ||
929 | } | |
930 | ||
931 | //__________________________________________________________ | |
932 | Int_t AliITStrackerSA::FindLabel(Int_t l1, Int_t l2, Int_t l3, Int_t l4, Int_t l5, Int_t l6){ | |
933 | ||
934 | //function used to determine the track label | |
935 | ||
936 | Int_t lb[6] = {l1,l2,l3,l4,l5,l6}; | |
937 | Int_t aa[6]={1,1,1,1,1,1}; | |
938 | Int_t ff=0; | |
939 | Int_t ll=0; | |
940 | Int_t k=0;Int_t w=0;Int_t num=6; | |
941 | if(lb[5]==-1) num=5; | |
942 | ||
943 | while(k<num){ | |
944 | ||
945 | for(Int_t i=k+1;i<num;i++){ | |
946 | ||
947 | if(lb[k]==lb[i] && aa[k]!=0){ | |
948 | ||
949 | aa[k]+=1; | |
950 | aa[i]=0; | |
951 | } | |
952 | } | |
953 | k++; | |
954 | } | |
955 | ||
956 | while(w<num){ | |
957 | ||
958 | for(Int_t j=0;j<6;j++){ | |
959 | if(aa[w]<aa[j]){ | |
960 | ff=aa[w]; | |
961 | aa[w]=aa[j]; | |
962 | aa[j]=ff; | |
963 | ll=lb[w]; | |
964 | lb[w]=lb[j]; | |
965 | lb[j]=ll; | |
966 | } | |
967 | } | |
968 | w++; | |
969 | } | |
970 | if(num==6) return lb[5]; | |
971 | else return lb[4]; | |
972 | } | |
973 | ||
974 | //_____________________________________________________________________________ | |
975 | Int_t AliITStrackerSA::Label(Int_t gl1, Int_t gl2, Int_t gl3, Int_t gl4, Int_t gl5, Int_t gl6,Int_t gl7, Int_t gl8, Int_t gl9, Int_t gl10,Int_t gl11, | |
976 | Int_t gl12, Int_t gl13, Int_t gl14,Int_t gl15, Int_t gl16, Int_t gl17, Int_t gl18, Int_t numberofpoints){ | |
977 | ||
978 | ||
979 | //function used to assign label to the found track. If track is fake, the label is negative | |
980 | ||
981 | Int_t lb0[6] = {gl1,gl2,gl3,gl4,gl5,gl6}; | |
982 | Int_t lb1[6] = {gl7,gl8,gl9,gl10,gl11,gl12}; | |
983 | Int_t lb2[6] = {gl13,gl14,gl15,gl16,gl17,gl18}; | |
984 | Int_t ll=FindLabel(lb0[0],lb0[1],lb0[2],lb0[3],lb0[4],lb0[5]); | |
985 | Int_t lflag=0;Int_t num=6; | |
986 | if(lb0[5]==-1 && lb1[5]==-1 && lb2[5]==-1) num=5; | |
987 | ||
988 | for(Int_t i=0;i<num;i++){ | |
989 | if(lb0[i]==ll || lb1[i]==ll || lb2[i]==ll) lflag+=1; | |
990 | } | |
991 | ||
992 | if(lflag>=numberofpoints) return ll; | |
993 | else return -ll; | |
994 | ||
995 | ||
996 | } | |
997 | ||
998 | //_____________________________________________________________________________ | |
999 | void AliITStrackerSA::SetWindowSizes(Int_t n, Double_t *phi, Double_t *lam){ | |
1000 | // Set sizes of the phi and lambda windows used for track finding | |
1001 | fNloop = n; | |
1002 | if(phi){ // user defined values | |
1003 | fPhiWin = new Double_t[fNloop]; | |
1004 | fLambdaWin = new Double_t[fNloop]; | |
1005 | for(Int_t k=0;k<fNloop;k++){ | |
1006 | fPhiWin[k]=phi[k]; | |
1007 | fLambdaWin[k]=lam[k]; | |
1008 | } | |
1009 | } | |
1010 | else { // default values | |
1011 | Double_t phid[46] = {0.001,0.0015,0.002,0.0023,0.0025,0.0027,0.003, | |
1012 | 0.0033,0.0035,0.0037,0.004,0.0043,0.0045,0.0047, | |
1013 | 0.005,0.0053,0.0055, | |
1014 | 0.006,0.0063,0.0065,0.007,0.0073,0.0075,0.0077, | |
1015 | 0.008,0.0083,0.0085,0.0087,0.009,0.0095,0.0097, | |
1016 | 0.01,0.0105,0.011,0.0115,0.012,0.0125,0.013,0.0135, | |
1017 | 0.014,0.0145,0.015,0.0155,0.016,0.0165,0.017}; | |
1018 | Double_t lambdad[46] = {0.001,0.0015,0.002,0.0023,0.0025,0.0027,0.003, | |
1019 | 0.0033,0.0035,0.0037,0.004,0.0043,0.0045,0.0047, | |
1020 | 0.005,0.0053,0.0055, | |
1021 | 0.006,0.0063,0.0065,0.007,0.0073,0.0075,0.0077, | |
1022 | 0.008,0.0083,0.0085,0.0087,0.009,0.0095,0.0097, | |
1023 | 0.01,0.015,0.011,0.0115,0.012,0.0125,0.013,0.0135, | |
1024 | 0.014,0.0145,0.015,0.0155,0.016,0.016,0.016}; | |
1025 | if(fNloop!=46){ | |
1026 | fNloop = 46; | |
1027 | Warning("SetWindowSizes","Number of loop forced to the default value %d",fNloop); | |
1028 | } | |
1029 | ||
1030 | fPhiWin = new Double_t[fNloop]; | |
1031 | fLambdaWin = new Double_t[fNloop]; | |
1032 | ||
1033 | for(Int_t k=0;k<fNloop;k++){ | |
1034 | fPhiWin[k]=phid[k]; | |
1035 | fLambdaWin[k]=lambdad[k]; | |
1036 | } | |
1037 | ||
1038 | } | |
1039 | ||
1040 | } | |
1041 | ||
1042 | ||
1043 | ||
1044 | ||
1045 | ||
1046 | ||
1047 | ||
1048 | ||
1049 | ||
1050 | ||
1051 | ||
1052 | ||
1053 | ||
1054 | ||
1055 | ||
1056 | ||
1057 | ||
1058 | ||
1059 | ||
1060 | ||
1061 | ||
1062 | ||
1063 | ||
1064 | ||
1065 | ||
1066 | ||
1067 |