1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 //-------------------------------------------------------------------------
17 // Implementation of the ITS tracker class
18 // It reads AliITSclusterV2 clusters and creates AliITStrackV2 tracks
19 // and fills with them the ESD
20 // Origin: Iouri Belikov, CERN, Jouri.Belikov@cern.ch
21 // dEdx analysis by: Boris Batyunya, JINR, Boris.Batiounia@cern.ch
22 //-------------------------------------------------------------------------
28 #include "AliITSgeom.h"
29 #include "AliITSRecPoint.h"
30 #include "AliTPCtrack.h"
32 #include "AliITSclusterV2.h"
33 #include "AliITStrackerV2.h"
35 ClassImp(AliITStrackerV2)
37 AliITStrackerV2::AliITSlayer AliITStrackerV2::fgLayers[kMaxLayer]; // ITS layers
39 AliITStrackerV2::AliITStrackerV2(const AliITSgeom *geom) : AliTracker() {
40 //--------------------------------------------------------------------
41 //This is the AliITStrackerV2 constructor
42 //--------------------------------------------------------------------
43 AliITSgeom *g=(AliITSgeom*)geom;
47 for (i=1; i<kMaxLayer+1; i++) {
48 Int_t nlad=g->GetNladders(i);
49 Int_t ndet=g->GetNdetectors(i);
51 g->GetTrans(i,1,1,x,y,z);
52 Double_t r=TMath::Sqrt(x*x + y*y);
53 Double_t poff=TMath::ATan2(y,x);
56 g->GetTrans(i,1,2,x,y,z);
57 r += TMath::Sqrt(x*x + y*y);
58 g->GetTrans(i,2,1,x,y,z);
59 r += TMath::Sqrt(x*x + y*y);
60 g->GetTrans(i,2,2,x,y,z);
61 r += TMath::Sqrt(x*x + y*y);
64 new (fgLayers+i-1) AliITSlayer(r,poff,zoff,nlad,ndet);
66 for (Int_t j=1; j<nlad+1; j++) {
67 for (Int_t k=1; k<ndet+1; k++) { //Fill this layer with detectors
68 Float_t x,y,zshift; g->GetTrans(i,j,k,x,y,zshift);
69 Double_t rot[9]; g->GetRotMatrix(i,j,k,rot);
71 Double_t phi=TMath::ATan2(rot[1],rot[0])+TMath::Pi();
73 if (i==1) phi+=TMath::Pi();
74 Double_t cp=TMath::Cos(phi), sp=TMath::Sin(phi);
77 AliITSdetector &det=fgLayers[i-1].GetDetector((j-1)*ndet + k-1);
78 new(&det) AliITSdetector(r,phi);
87 fConstraint[0]=1; fConstraint[1]=0;
89 Double_t xyz[]={kXV,kYV,kZV}, ers[]={kSigmaXV,kSigmaYV,kSigmaZV};
92 for (Int_t i=0; i<kMaxLayer; i++) fLayersNotToSkip[i]=kLayersNotToSkip[i];
93 fLastLayerToTrackTo=kLastLayerToTrackTo;
97 void AliITStrackerV2::SetLayersNotToSkip(Int_t *l) {
98 //--------------------------------------------------------------------
99 //This function set masks of the layers which must be not skipped
100 //--------------------------------------------------------------------
101 for (Int_t i=0; i<kMaxLayer; i++) fLayersNotToSkip[i]=l[i];
104 Int_t AliITStrackerV2::LoadClusters(TTree *cTree) {
105 //--------------------------------------------------------------------
106 //This function loads ITS clusters
107 //--------------------------------------------------------------------
108 TBranch *branch=cTree->GetBranch("Clusters");
110 Error("LoadClusters"," can't get the branch !\n");
114 TClonesArray dummy("AliITSclusterV2",10000), *clusters=&dummy;
115 branch->SetAddress(&clusters);
118 for (Int_t i=0; i<kMaxLayer; i++) {
119 Int_t ndet=fgLayers[i].GetNdetectors();
120 Int_t jmax = j + fgLayers[i].GetNladders()*ndet;
121 for (; j<jmax; j++) {
122 if (!cTree->GetEvent(j)) continue;
123 Int_t ncl=clusters->GetEntriesFast();
125 AliITSclusterV2 *c=(AliITSclusterV2*)clusters->UncheckedAt(ncl);
126 fgLayers[i].InsertCluster(new AliITSclusterV2(*c));
130 fgLayers[i].ResetRoad(); //road defined by the cluster density
136 void AliITStrackerV2::UnloadClusters() {
137 //--------------------------------------------------------------------
138 //This function unloads ITS clusters
139 //--------------------------------------------------------------------
140 for (Int_t i=0; i<kMaxLayer; i++) fgLayers[i].ResetClusters();
143 static Int_t CorrectForDeadZoneMaterial(AliITStrackV2 *t) {
144 //--------------------------------------------------------------------
145 // Correction for the material between the TPC and the ITS
146 // (should it belong to the TPC code ?)
147 //--------------------------------------------------------------------
148 Double_t riw=80., diw=0.0053, x0iw=30; // TPC inner wall ?
149 Double_t rcd=61., dcd=0.0053, x0cd=30; // TPC "central drum" ?
150 Double_t yr=12.8, dr=0.03; // rods ?
151 Double_t zm=0.2, dm=0.40; // membrane
152 //Double_t rr=52., dr=0.19, x0r=24., yyr=7.77; //rails
153 Double_t rs=50., ds=0.001; // something belonging to the ITS (screen ?)
155 if (t->GetX() > riw) {
156 if (!t->PropagateTo(riw,diw,x0iw)) return 1;
157 if (TMath::Abs(t->GetY())>yr) t->CorrectForMaterial(dr);
158 if (TMath::Abs(t->GetZ())<zm) t->CorrectForMaterial(dm);
159 if (!t->PropagateTo(rcd,dcd,x0cd)) return 1;
160 //Double_t x,y,z; t->GetGlobalXYZat(rr,x,y,z);
161 //if (TMath::Abs(y)<yyr) t->PropagateTo(rr,dr,x0r);
162 if (!t->PropagateTo(rs,ds)) return 1;
163 } else if (t->GetX() < rs) {
164 if (!t->PropagateTo(rs,-ds)) return 1;
165 //Double_t x,y,z; t->GetGlobalXYZat(rr,x,y,z);
166 //if (TMath::Abs(y)<yyr) t->PropagateTo(rr,-dr,x0r);
167 if (!t->PropagateTo(rcd,-dcd,x0cd)) return 1;
168 if (!t->PropagateTo(riw,-diw,x0iw)) return 1;
170 ::Error("CorrectForDeadZoneMaterial","track is already in the dead zone !");
177 Int_t AliITStrackerV2::Clusters2Tracks(AliESD *event) {
178 //--------------------------------------------------------------------
179 // This functions reconstructs ITS tracks
180 // The clusters must be already loaded !
181 //--------------------------------------------------------------------
182 TObjArray itsTracks(15000);
184 {/* Read ESD tracks */
185 Int_t nentr=event->GetNumberOfTracks();
186 Info("Clusters2Tracks", "Number of ESD tracks: %d\n", nentr);
188 AliESDtrack *esd=event->GetTrack(nentr);
190 if (esd->GetStatus() != AliESDtrack::kTPCin) continue;
194 t=new AliITStrackV2(*esd);
195 } catch (const Char_t *msg) {
196 Warning("Clusters2Tracks",msg);
200 if (TMath::Abs(t->GetD())>4) continue;
202 if (CorrectForDeadZoneMaterial(t)!=0) {
203 Warning("Clusters2Tracks",
204 "failed to correct for the material in the dead zone !\n");
208 itsTracks.AddLast(t);
210 } /* End Read ESD tracks */
213 Int_t nentr=itsTracks.GetEntriesFast();
216 for (fPass=0; fPass<2; fPass++) {
217 Int_t &constraint=fConstraint[fPass]; if (constraint<0) continue;
218 for (Int_t i=0; i<nentr; i++) {
219 AliITStrackV2 *t=(AliITStrackV2*)itsTracks.UncheckedAt(i);
220 if (t==0) continue; //this track has been already tracked
221 Int_t tpcLabel=t->GetLabel(); //save the TPC track label
223 ResetTrackToFollow(*t);
226 for (FollowProlongation(); fI<kMaxLayer; fI++) {
227 while (TakeNextProlongation()) FollowProlongation();
230 if (fBestTrack.GetNumberOfClusters() == 0) continue;
232 if (fConstraint[fPass]) {
233 ResetTrackToFollow(*t);
234 if (!RefitAt(3.7, &fTrackToFollow, &fBestTrack)) continue;
238 fBestTrack.SetLabel(tpcLabel);
239 fBestTrack.CookdEdx();
240 CookLabel(&fBestTrack,0.); //For comparison only
241 fBestTrack.UpdateESDtrack(AliESDtrack::kITSin);
242 UseClusters(&fBestTrack);
243 delete itsTracks.RemoveAt(i);
250 Info("Clusters2Tracks","Number of prolonged tracks: %d\n",ntrk);
255 Int_t AliITStrackerV2::Clusters2Tracks(TTree *tpcTree, TTree *itsTree) {
256 //--------------------------------------------------------------------
257 // This functions reconstructs ITS tracks
258 // The clusters must be already loaded !
259 //--------------------------------------------------------------------
260 Int_t nentr=0; TObjArray itsTracks(15000);
262 {/* Read TPC tracks */
263 AliTPCtrack *itrack=new AliTPCtrack;
264 TBranch *branch=tpcTree->GetBranch("tracks");
266 Error("Clusters2Tracks","Can't get the branch !");
269 tpcTree->SetBranchAddress("tracks",&itrack);
270 nentr=(Int_t)tpcTree->GetEntries();
272 Info("Clusters2Tracks","Number of TPC tracks: %d\n",nentr);
274 for (Int_t i=0; i<nentr; i++) {
275 tpcTree->GetEvent(i);
278 t=new AliITStrackV2(*itrack);
279 } catch (const Char_t *msg) {
280 Warning("Clusters2Tracks",msg);
284 if (TMath::Abs(t->GetD())>4) continue;
286 if (CorrectForDeadZoneMaterial(t)!=0) {
287 Warning("Clusters2Tracks",
288 "failed to correct for the material in the dead zone !\n");
292 itsTracks.AddLast(t);
297 nentr=itsTracks.GetEntriesFast();
300 AliITStrackV2 *otrack=&fBestTrack;
301 TBranch *branch=itsTree->GetBranch("tracks");
302 if (!branch) itsTree->Branch("tracks","AliITStrackV2",&otrack,32000,3);
303 else branch->SetAddress(&otrack);
305 for (fPass=0; fPass<2; fPass++) {
306 Int_t &constraint=fConstraint[fPass]; if (constraint<0) continue;
307 for (Int_t i=0; i<nentr; i++) {
308 AliITStrackV2 *t=(AliITStrackV2*)itsTracks.UncheckedAt(i);
309 if (t==0) continue; //this track has been already tracked
310 Int_t tpcLabel=t->GetLabel(); //save the TPC track label
312 ResetTrackToFollow(*t);
315 for (FollowProlongation(); fI<kMaxLayer; fI++) {
316 while (TakeNextProlongation()) FollowProlongation();
319 if (fBestTrack.GetNumberOfClusters() == 0) continue;
321 if (fConstraint[fPass]) {
322 ResetTrackToFollow(*t);
323 if (!RefitAt(3.7, &fTrackToFollow, &fBestTrack)) continue;
327 fBestTrack.SetLabel(tpcLabel);
328 fBestTrack.CookdEdx();
329 CookLabel(&fBestTrack,0.); //For comparison only
331 UseClusters(&fBestTrack);
332 delete itsTracks.RemoveAt(i);
336 nentr=(Int_t)itsTree->GetEntries();
337 Info("Clusters2Tracks","Number of prolonged tracks: %d\n",nentr);
344 Int_t AliITStrackerV2::PropagateBack(AliESD *event) {
345 //--------------------------------------------------------------------
346 // This functions propagates reconstructed ITS tracks back
347 // The clusters must be loaded !
348 //--------------------------------------------------------------------
349 Int_t nentr=event->GetNumberOfTracks();
350 Info("PropagateBack", "Number of ESD tracks: %d\n", nentr);
353 for (Int_t i=0; i<nentr; i++) {
354 AliESDtrack *esd=event->GetTrack(i);
356 if (esd->GetStatus()!=(AliESDtrack::kTPCin|AliESDtrack::kITSin)) continue;
360 t=new AliITStrackV2(*esd);
361 } catch (const Char_t *msg) {
362 Warning("PropagateBack",msg);
367 ResetTrackToFollow(*t);
369 // propagete to vertex [SR, GSI 17.02.2003]
370 fTrackToFollow.PropagateTo(3.,0.0028,65.19);
371 fTrackToFollow.PropagateToVertex();
373 // Start Time measurement [SR, GSI 17.02.2003]
374 fTrackToFollow.StartTimeIntegral();
375 fTrackToFollow.PropagateTo(3.,-0.0028,65.19);
377 fTrackToFollow.ResetCovariance(); fTrackToFollow.ResetClusters();
378 if (RefitAt(49.,&fTrackToFollow,t)) {
379 if (CorrectForDeadZoneMaterial(&fTrackToFollow)!=0) {
380 Warning("PropagateBack",
381 "failed to correct for the material in the dead zone !\n");
385 fTrackToFollow.SetLabel(t->GetLabel());
386 fTrackToFollow.CookdEdx();
387 CookLabel(&fTrackToFollow,0.); //For comparison only
388 fTrackToFollow.UpdateESDtrack(AliESDtrack::kITSout);
389 UseClusters(&fTrackToFollow);
395 Info("PropagateBack","Number of back propagated ITS tracks: %d\n",ntrk);
400 Int_t AliITStrackerV2::PropagateBack(TTree *inp, TTree *out) {
401 //--------------------------------------------------------------------
402 //This functions propagates reconstructed ITS tracks back
403 //--------------------------------------------------------------------
404 Error("PropagateBack","This method is not converted to NewIO yet (Args: %x , %x )\n",inp,out);
407 TFile *in=(TFile*)inp;
408 TDirectory *savedir=gDirectory;
410 if (LoadClusters()!=0) return 1;
413 Error("PropagateBack","file with ITS tracks is not open !\n");
417 if (!out->IsOpen()) {
418 Error("PropagateBack","file for back propagated ITS tracks is not open !\n");
425 sprintf(tname,"TreeT_ITS_%d",GetEventNumber());
426 TTree *itsTree=(TTree*)in->Get(tname);
428 Error("PropagateBack","can't get a tree with ITS tracks !\n");
431 AliITStrackV2 *itrack=new AliITStrackV2;
432 itsTree->SetBranchAddress("tracks",&itrack);
436 sprintf(tname,"TreeT_ITSb_%d",GetEventNumber());
437 TTree backTree(tname,"Tree with back propagated ITS tracks");
438 AliTPCtrack *otrack=0;
439 backTree.Branch("tracks","AliTPCtrack",&otrack,32000,2);
441 Int_t nentr=(Int_t)itsTree->GetEntries();
443 for (i=0; i<nentr; i++) {
444 itsTree->GetEvent(i);
445 Int_t itsLabel=itrack->GetLabel(); //save the ITS track label
446 ResetTrackToFollow(*itrack);
448 // propagete to vertex [SR, GSI 17.02.2003]
449 fTrackToFollow.PropagateTo(3.,0.0028,65.19);
450 fTrackToFollow.PropagateToVertex();
452 // Start Time measurement [SR, GSI 17.02.2003]
453 fTrackToFollow.StartTimeIntegral();
454 fTrackToFollow.PropagateTo(3.,-0.0028,65.19);
457 fTrackToFollow.ResetCovariance(); fTrackToFollow.ResetClusters();
458 if (!RefitAt(49.,&fTrackToFollow,itrack)) continue;
460 if (CorrectForDeadZoneMaterial(&fTrackToFollow)!=0) {
461 Warning("PropagateBack",
462 "failed to correct for the material in the dead zone !\n");
466 fTrackToFollow.SetLabel(itsLabel);
467 otrack=new AliTPCtrack(fTrackToFollow,fTrackToFollow.GetAlpha());
468 backTree.Fill(); delete otrack;
469 UseClusters(&fTrackToFollow);
471 i=(Int_t)backTree.GetEntries();
474 Info("PropagateBack","Number of ITS tracks: %d\n",nentr);
475 Info("PropagateBack","Number of back propagated ITS tracks: %d\n",i);
478 delete itsTree; //Thanks to Mariana Bondila
488 Int_t AliITStrackerV2::RefitInward(AliESD *event) {
489 //--------------------------------------------------------------------
490 // This functions refits ITS tracks using the
491 // "inward propagated" TPC tracks
492 // The clusters must be loaded !
493 //--------------------------------------------------------------------
494 Int_t nentr=event->GetNumberOfTracks();
495 Info("RefitInward", "Number of ESD tracks: %d\n", nentr);
498 for (Int_t i=0; i<nentr; i++) {
499 AliESDtrack *esd=event->GetTrack(i);
501 ULong_t flags = AliESDtrack::kITSin | AliESDtrack::kTPCrefit;
503 if ( (esd->GetStatus() & flags) != flags ) continue;
504 if ( esd->GetStatus() & AliESDtrack::kITSrefit) continue;
508 t=new AliITStrackV2(*esd);
509 } catch (const Char_t *msg) {
510 Warning("RefitInward",msg);
515 if (CorrectForDeadZoneMaterial(t)!=0) {
516 Warning("RefitInward",
517 "failed to correct for the material in the dead zone !\n");
522 ResetTrackToFollow(*t);
523 fTrackToFollow.ResetClusters();
526 if (RefitAt(3.7, &fTrackToFollow, t)) {
527 fTrackToFollow.SetLabel(t->GetLabel());
528 fTrackToFollow.CookdEdx();
529 CookLabel(&fTrackToFollow,0.); //For comparison only
531 fTrackToFollow.PropagateTo(3.,0.0028,65.19); //The beam pipe
532 fTrackToFollow.PropagateToVertex();
534 fTrackToFollow.UpdateESDtrack(AliESDtrack::kITSrefit);
535 UseClusters(&fTrackToFollow);
541 Info("RefitInward","Number of refitted tracks: %d\n",ntrk);
546 Int_t AliITStrackerV2::RefitInward(TTree *in, TTree *out) {
547 //--------------------------------------------------------------------
548 // This functions refits ITS tracks using the
549 // "inward propagated" TPC tracks
550 //--------------------------------------------------------------------
551 Error("RefitInward","This method is not converted to NewIO yet (Args: %x , %x )\n",in,out);
554 TFile *in=(TFile*)inp;
555 TDirectory *savedir=gDirectory;
557 if (LoadClusters()!=0) return 1;
559 if (!inSeeds->IsOpen()) {
560 Error("RefitInward","file with inward TPC tracks is not open !\n");
565 Error("RefitInward","file with ITS tracks is not open !\n");
569 if (!out->IsOpen()) {
570 Error("RefitInward","file for inward ITS tracks is not open !\n");
576 //LUT used for the track matching (S.Radomski's idea)
577 const Int_t nLab = 400000; // limit on the number of track labels
579 for(i=0; i<nLab; i++) lut[i] = -1;
583 TObjArray itsTracks(15000);
584 {// Read the ITS tracks
585 sprintf(tname,"TreeT_ITS_%d",GetEventNumber());
586 TTree *itsTree=(TTree*)in->Get(tname);
588 Error("RefitInward","can't get a tree with ITS tracks !\n");
591 AliITStrackV2 *itrack=new AliITStrackV2;
592 itsTree->SetBranchAddress("tracks",&itrack);
593 Int_t nits=(Int_t)itsTree->GetEntries();
595 Info("RefitInward","Number of ITS tracks: %d\n",nits);
597 for (Int_t i=0; i<nits; i++) {
598 itsTree->GetEvent(i);
599 Int_t lab=TMath::Abs(itrack->GetLabel());
601 if (lut[lab]>=0) Warning("RefitInward","double track %d\n",lab);
604 Warning("RefitInward","Too big ITS track label: %d\n!",lab);
607 itsTracks.AddLast(new AliITStrackV2(*itrack));
613 //Get the input seeds tree
615 sprintf(tname,"tracksTPC_%d",GetEventNumber());
616 TTree *tpcTree=(TTree*)inSeeds->Get(tname);
618 Error("RefitInward","can't get a tree with TPC tracks !\n");
621 AliTPCtrack *itrack=new AliTPCtrack;
622 tpcTree->SetBranchAddress("tracks",&itrack);
623 Int_t ntpc=(Int_t)tpcTree->GetEntries();
625 Info("RefitInward","Number of TPC tracks: %d\n",ntpc);
627 //Create the output tree
629 sprintf(tname,"TreeT_ITSinward_%d",GetEventNumber());
630 TTree outTree(tname,"Tree with inward refitted ITS tracks");
631 AliITStrackV2 *otrack=0;
632 outTree.Branch("tracks","AliITStrackV2",&otrack,32000,0);
634 for (i=0; i<ntpc; i++) {
635 tpcTree->GetEvent(i);
637 otrack=new AliITStrackV2(*itrack);
638 } catch (const Char_t *msg) {
639 Warning("RefitInward",msg);
643 //check if this track was reconstructed in the ITS
644 Int_t lab=TMath::Abs(otrack->GetLabel());
646 Warning("RefitInward","Too big TPC track label: %d\n!",lab);
650 if (idx<0) continue; //no prolongation in the ITS for this track
652 if (CorrectForDeadZoneMaterial(otrack)!=0) {
653 Warning("RefitInward",
654 "failed to correct for the material in the dead zone !\n");
660 AliITStrackV2 *ctrack=(AliITStrackV2*)itsTracks.UncheckedAt(idx);
661 if (!RefitAt(3.7, otrack, ctrack)) continue;
663 otrack->SetLabel(itrack->GetLabel()); //For comparison only
665 CookLabel(otrack,0.); //For comparison only
669 i=(Int_t)outTree.GetEntries();
670 Info("RefitInward","Number of inward refitted ITS tracks: %d\n",i);
684 AliCluster *AliITStrackerV2::GetCluster(Int_t index) const {
685 //--------------------------------------------------------------------
686 // Return pointer to a given cluster
687 //--------------------------------------------------------------------
688 Int_t l=(index & 0xf0000000) >> 28;
689 Int_t c=(index & 0x0fffffff) >> 00;
690 return fgLayers[l].GetCluster(c);
694 void AliITStrackerV2::FollowProlongation() {
695 //--------------------------------------------------------------------
696 //This function finds a track prolongation
697 //--------------------------------------------------------------------
698 while (fI>fLastLayerToTrackTo) {
701 AliITSlayer &layer=fgLayers[i];
702 AliITStrackV2 &track=fTracks[i];
704 Double_t r=layer.GetR();
707 Double_t rs=0.5*(fgLayers[i+1].GetR() + r);
708 Double_t d=0.0034, x0=38.6;
709 if (i==1) {rs=9.; d=0.0097; x0=42;}
710 if (!fTrackToFollow.PropagateTo(rs,d,x0)) {
711 //Warning("FollowProlongation","propagation failed !\n");
718 if (!fTrackToFollow.GetGlobalXYZat(r,x,y,z)) {
719 //Warning("FollowProlongation","failed to estimate track !\n");
722 Double_t phi=TMath::ATan2(y,x);
724 Int_t idet=layer.FindDetectorIndex(phi,z);
726 //Warning("FollowProlongation","failed to find a detector !\n");
730 //propagate to the intersection
731 const AliITSdetector &det=layer.GetDetector(idet);
733 if (!fTrackToFollow.Propagate(phi,det.GetR())) {
734 //Warning("FollowProlongation","propagation failed !\n");
737 fTrackToFollow.SetDetectorIndex(idet);
739 //Select possible prolongations and store the current track estimation
740 track.~AliITStrackV2(); new(&track) AliITStrackV2(fTrackToFollow);
741 Double_t dz=7*TMath::Sqrt(track.GetSigmaZ2() + kSigmaZ2[i]);
742 Double_t dy=7*TMath::Sqrt(track.GetSigmaY2() + kSigmaY2[i]);
743 Double_t road=layer.GetRoad();
744 if (dz*dy>road*road) {
745 Double_t dd=TMath::Sqrt(dz*dy), scz=dz/dd, scy=dy/dd;
746 dz=road*scz; dy=road*scy;
749 //Double_t dz=4*TMath::Sqrt(track.GetSigmaZ2() + kSigmaZ2[i]);
750 if (dz < 0.5*TMath::Abs(track.GetTgl())) dz=0.5*TMath::Abs(track.GetTgl());
752 //Warning("FollowProlongation","too broad road in Z !\n");
756 if (TMath::Abs(fTrackToFollow.GetZ()-GetZ()) > r+dz) return;
758 //Double_t dy=4*TMath::Sqrt(track.GetSigmaY2() + kSigmaY2[i]);
759 if (dy < 0.5*TMath::Abs(track.GetSnp())) dy=0.5*TMath::Abs(track.GetSnp());
761 //Warning("FollowProlongation","too broad road in Y !\n");
765 Double_t zmin=track.GetZ() - dz;
766 Double_t zmax=track.GetZ() + dz;
767 Double_t ymin=track.GetY() + r*phi - dy;
768 Double_t ymax=track.GetY() + r*phi + dy;
769 layer.SelectClusters(zmin,zmax,ymin,ymax);
772 //take another prolongation
773 if (!TakeNextProlongation())
774 if (fLayersNotToSkip[fI]) return;
778 //deal with the best track
779 Int_t ncl=fTrackToFollow.GetNumberOfClusters();
780 Int_t nclb=fBestTrack.GetNumberOfClusters();
783 Double_t chi2=fTrackToFollow.GetChi2();
784 if (chi2/ncl < kChi2PerCluster) {
785 if (ncl > nclb || chi2 < fBestTrack.GetChi2()) {
793 Int_t AliITStrackerV2::TakeNextProlongation() {
794 //--------------------------------------------------------------------
795 // This function takes another track prolongation
797 // dEdx analysis by: Boris Batyunya, JINR, Boris.Batiounia@cern.ch
798 //--------------------------------------------------------------------
799 AliITSlayer &layer=fgLayers[fI];
800 ResetTrackToFollow(fTracks[fI]);
802 Double_t dz=7*TMath::Sqrt(fTrackToFollow.GetSigmaZ2() + kSigmaZ2[fI]);
803 Double_t dy=7*TMath::Sqrt(fTrackToFollow.GetSigmaY2() + kSigmaY2[fI]);
804 Double_t road=layer.GetRoad();
805 if (dz*dy>road*road) {
806 Double_t dd=TMath::Sqrt(dz*dy), scz=dz/dd, scy=dy/dd;
807 dz=road*scz; dy=road*scy;
810 const AliITSclusterV2 *c=0; Int_t ci=-1;
811 Double_t chi2=12345.;
812 while ((c=layer.GetNextCluster(ci))!=0) {
813 Int_t idet=c->GetDetectorIndex();
815 if (fTrackToFollow.GetDetectorIndex()!=idet) {
816 const AliITSdetector &det=layer.GetDetector(idet);
817 ResetTrackToFollow(fTracks[fI]);
818 if (!fTrackToFollow.Propagate(det.GetPhi(),det.GetR())) {
819 //Warning("TakeNextProlongation","propagation failed !\n");
822 fTrackToFollow.SetDetectorIndex(idet);
823 if (TMath::Abs(fTrackToFollow.GetZ()-GetZ())>layer.GetR()+dz) continue;
826 if (TMath::Abs(fTrackToFollow.GetZ() - c->GetZ()) > dz) continue;
827 if (TMath::Abs(fTrackToFollow.GetY() - c->GetY()) > dy) continue;
829 chi2=fTrackToFollow.GetPredictedChi2(c); if (chi2<kMaxChi2) break;
832 if (chi2>=kMaxChi2) return 0;
835 if (!fTrackToFollow.Update(c,chi2,(fI<<28)+ci)) {
836 //Warning("TakeNextProlongation","filtering failed !\n");
840 if (fTrackToFollow.GetNumberOfClusters()>1)
841 if (TMath::Abs(fTrackToFollow.GetD())>4) return 0;
844 SetSampledEdx(c->GetQ(),fTrackToFollow.GetNumberOfClusters()-1); //b.b.
848 Double_t d=layer.GetThickness(fTrackToFollow.GetY(),fTrackToFollow.GetZ(),x0);
849 fTrackToFollow.CorrectForMaterial(d,x0);
852 if (fConstraint[fPass]) {
853 Double_t d=GetEffectiveThickness(0,0); //Think of this !!!!
854 Double_t xyz[]={GetX(),GetY(),GetZ()};
855 Double_t ers[]={GetSigmaX(),GetSigmaY(),GetSigmaZ()};
856 fTrackToFollow.Improve(d,xyz,ers);
863 AliITStrackerV2::AliITSlayer::AliITSlayer() {
864 //--------------------------------------------------------------------
865 //default AliITSlayer constructor
866 //--------------------------------------------------------------------
871 AliITStrackerV2::AliITSlayer::
872 AliITSlayer(Double_t r,Double_t p,Double_t z,Int_t nl,Int_t nd) {
873 //--------------------------------------------------------------------
874 //main AliITSlayer constructor
875 //--------------------------------------------------------------------
876 fR=r; fPhiOffset=p; fZOffset=z;
877 fNladders=nl; fNdetectors=nd;
878 fDetectors=new AliITSdetector[fNladders*fNdetectors];
883 fRoad=2*fR*TMath::Sqrt(3.14/1.);//assuming that there's only one cluster
886 AliITStrackerV2::AliITSlayer::~AliITSlayer() {
887 //--------------------------------------------------------------------
888 // AliITSlayer destructor
889 //--------------------------------------------------------------------
891 for (Int_t i=0; i<fN; i++) delete fClusters[i];
894 void AliITStrackerV2::AliITSlayer::ResetClusters() {
895 //--------------------------------------------------------------------
896 // This function removes loaded clusters
897 //--------------------------------------------------------------------
898 for (Int_t i=0; i<fN; i++) delete fClusters[i];
903 void AliITStrackerV2::AliITSlayer::ResetRoad() {
904 //--------------------------------------------------------------------
905 // This function calculates the road defined by the cluster density
906 //--------------------------------------------------------------------
908 for (Int_t i=0; i<fN; i++) {
909 if (TMath::Abs(fClusters[i]->GetZ())<fR) n++;
911 if (n>1) fRoad=2*fR*TMath::Sqrt(3.14/n);
914 Int_t AliITStrackerV2::AliITSlayer::InsertCluster(AliITSclusterV2 *c) {
915 //--------------------------------------------------------------------
916 //This function adds a cluster to this layer
917 //--------------------------------------------------------------------
918 if (fN==kMaxClusterPerLayer) {
919 ::Error("InsertCluster","Too many clusters !\n");
923 if (fN==0) {fClusters[fN++]=c; return 0;}
924 Int_t i=FindClusterIndex(c->GetZ());
925 memmove(fClusters+i+1 ,fClusters+i,(fN-i)*sizeof(AliITSclusterV2*));
926 fClusters[i]=c; fN++;
931 Int_t AliITStrackerV2::AliITSlayer::FindClusterIndex(Double_t z) const {
932 //--------------------------------------------------------------------
933 // This function returns the index of the nearest cluster
934 //--------------------------------------------------------------------
936 if (z <= fClusters[0]->GetZ()) return 0;
937 if (z > fClusters[fN-1]->GetZ()) return fN;
938 Int_t b=0, e=fN-1, m=(b+e)/2;
939 for (; b<e; m=(b+e)/2) {
940 if (z > fClusters[m]->GetZ()) b=m+1;
946 void AliITStrackerV2::AliITSlayer::
947 SelectClusters(Double_t zmin,Double_t zmax,Double_t ymin, Double_t ymax) {
948 //--------------------------------------------------------------------
949 // This function sets the "window"
950 //--------------------------------------------------------------------
951 fI=FindClusterIndex(zmin); fZmax=zmax;
952 Double_t circle=2*TMath::Pi()*fR;
953 if (ymax>circle) { ymax-=circle; ymin-=circle; }
954 fYmin=ymin; fYmax=ymax;
957 const AliITSclusterV2 *AliITStrackerV2::AliITSlayer::GetNextCluster(Int_t &ci){
958 //--------------------------------------------------------------------
959 // This function returns clusters within the "window"
960 //--------------------------------------------------------------------
961 const AliITSclusterV2 *cluster=0;
962 for (Int_t i=fI; i<fN; i++) {
963 const AliITSclusterV2 *c=fClusters[i];
964 if (c->GetZ() > fZmax) break;
965 if (c->IsUsed()) continue;
966 const AliITSdetector &det=GetDetector(c->GetDetectorIndex());
967 Double_t y=fR*det.GetPhi() + c->GetY();
969 if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi();
970 if (y>1.*fR*TMath::Pi() && fYmax<y) y -= 2*fR*TMath::Pi();
972 if (y<fYmin) continue;
973 if (y>fYmax) continue;
982 Int_t AliITStrackerV2::AliITSlayer::
983 FindDetectorIndex(Double_t phi, Double_t z) const {
984 //--------------------------------------------------------------------
985 //This function finds the detector crossed by the track
986 //--------------------------------------------------------------------
987 Double_t dphi=-(phi-fPhiOffset);
988 if (dphi < 0) dphi += 2*TMath::Pi();
989 else if (dphi >= 2*TMath::Pi()) dphi -= 2*TMath::Pi();
990 Int_t np=Int_t(dphi*fNladders*0.5/TMath::Pi()+0.5);
991 if (np>=fNladders) np-=fNladders;
992 if (np<0) np+=fNladders;
994 Double_t dz=fZOffset-z;
995 Int_t nz=Int_t(dz*(fNdetectors-1)*0.5/fZOffset+0.5);
996 if (nz>=fNdetectors) return -1;
999 return np*fNdetectors + nz;
1003 AliITStrackerV2::AliITSlayer::GetThickness(Double_t y,Double_t z,Double_t &x0)
1005 //--------------------------------------------------------------------
1006 //This function returns the layer thickness at this point (units X0)
1007 //--------------------------------------------------------------------
1011 if (43<fR&&fR<45) { //SSD2
1014 if (TMath::Abs(y-0.00)>3.40) d+=dd;
1015 if (TMath::Abs(y-1.90)<0.45) {d+=(0.013-0.0034);}
1016 if (TMath::Abs(y+1.90)<0.45) {d+=(0.013-0.0034);}
1017 for (Int_t i=0; i<12; i++) {
1018 if (TMath::Abs(z-3.9*(i+0.5))<0.15) {
1019 if (TMath::Abs(y-0.00)>3.40) d+=dd;
1023 if (TMath::Abs(z+3.9*(i+0.5))<0.15) {
1024 if (TMath::Abs(y-0.00)>3.40) d+=dd;
1028 if (TMath::Abs(z-3.4-3.9*i)<0.50) {d+=(0.016-0.0034); break;}
1029 if (TMath::Abs(z+0.5+3.9*i)<0.50) {d+=(0.016-0.0034); break;}
1032 if (37<fR&&fR<41) { //SSD1
1035 if (TMath::Abs(y-0.00)>3.40) d+=dd;
1036 if (TMath::Abs(y-1.90)<0.45) {d+=(0.013-0.0034);}
1037 if (TMath::Abs(y+1.90)<0.45) {d+=(0.013-0.0034);}
1038 for (Int_t i=0; i<11; i++) {
1039 if (TMath::Abs(z-3.9*i)<0.15) {
1040 if (TMath::Abs(y-0.00)>3.40) d+=dd;
1044 if (TMath::Abs(z+3.9*i)<0.15) {
1045 if (TMath::Abs(y-0.00)>3.40) d+=dd;
1049 if (TMath::Abs(z-1.85-3.9*i)<0.50) {d+=(0.016-0.0034); break;}
1050 if (TMath::Abs(z+2.05+3.9*i)<0.50) {d+=(0.016-0.0034); break;}
1053 if (13<fR&&fR<26) { //SDD
1056 if (TMath::Abs(y-0.00)>3.30) d+=dd;
1058 if (TMath::Abs(y-1.80)<0.55) {
1060 for (Int_t j=0; j<20; j++) {
1061 if (TMath::Abs(z+0.7+1.47*j)<0.12) {d+=0.08; x0=9.; break;}
1062 if (TMath::Abs(z-0.7-1.47*j)<0.12) {d+=0.08; x0=9.; break;}
1065 if (TMath::Abs(y+1.80)<0.55) {
1067 for (Int_t j=0; j<20; j++) {
1068 if (TMath::Abs(z-0.7-1.47*j)<0.12) {d+=0.08; x0=9.; break;}
1069 if (TMath::Abs(z+0.7+1.47*j)<0.12) {d+=0.08; x0=9.; break;}
1073 for (Int_t i=0; i<4; i++) {
1074 if (TMath::Abs(z-7.3*i)<0.60) {
1076 if (TMath::Abs(y-0.00)>3.30) d+=dd;
1079 if (TMath::Abs(z+7.3*i)<0.60) {
1081 if (TMath::Abs(y-0.00)>3.30) d+=dd;
1086 if (6<fR&&fR<8) { //SPD2
1087 Double_t dd=0.0063; x0=21.5;
1089 if (TMath::Abs(y-3.08)>0.5) d+=dd;
1090 //if (TMath::Abs(y-3.08)>0.45) d+=dd;
1091 if (TMath::Abs(y-3.03)<0.10) {d+=0.014;}
1093 if (3<fR&&fR<5) { //SPD1
1094 Double_t dd=0.0063; x0=21.5;
1096 if (TMath::Abs(y+0.21)>0.6) d+=dd;
1097 //if (TMath::Abs(y+0.21)>0.45) d+=dd;
1098 if (TMath::Abs(y+0.10)<0.10) {d+=0.014;}
1104 Double_t AliITStrackerV2::GetEffectiveThickness(Double_t y,Double_t z) const
1106 //--------------------------------------------------------------------
1107 //Returns the thickness between the current layer and the vertex (units X0)
1108 //--------------------------------------------------------------------
1109 Double_t d=0.0028*3*3; //beam pipe
1112 Double_t xn=fgLayers[fI].GetR();
1113 for (Int_t i=0; i<fI; i++) {
1114 Double_t xi=fgLayers[i].GetR();
1115 d+=fgLayers[i].GetThickness(y,z,x0)*xi*xi;
1124 Double_t xi=0.5*(fgLayers[3].GetR()+fgLayers[4].GetR());
1131 Int_t AliITStrackerV2::AliITSlayer::InRoad() const {
1132 //--------------------------------------------------------------------
1133 // This function returns number of clusters within the "window"
1134 //--------------------------------------------------------------------
1136 for (Int_t i=fI; i<fN; i++) {
1137 const AliITSclusterV2 *c=fClusters[i];
1138 if (c->GetZ() > fZmax) break;
1139 if (c->IsUsed()) continue;
1140 const AliITSdetector &det=GetDetector(c->GetDetectorIndex());
1141 Double_t y=fR*det.GetPhi() + c->GetY();
1143 if (y>2.*fR*TMath::Pi()) y -= 2*fR*TMath::Pi();
1144 if (y>1.*fR*TMath::Pi() && fYmax<y) y -= 2*fR*TMath::Pi();
1146 if (y<fYmin) continue;
1147 if (y>fYmax) continue;
1154 AliITStrackerV2::RefitAt(Double_t xx,AliITStrackV2 *t,const AliITStrackV2 *c) {
1155 //--------------------------------------------------------------------
1156 // This function refits the track "t" at the position "x" using
1157 // the clusters from "c"
1158 //--------------------------------------------------------------------
1159 Int_t index[kMaxLayer];
1161 for (k=0; k<kMaxLayer; k++) index[k]=-1;
1162 Int_t nc=c->GetNumberOfClusters();
1163 for (k=0; k<nc; k++) {
1164 Int_t idx=c->GetClusterIndex(k),nl=(idx&0xf0000000)>>28;
1168 Int_t from, to, step;
1169 if (xx > t->GetX()) {
1170 from=0; to=kMaxLayer;
1173 from=kMaxLayer-1; to=-1;
1177 for (Int_t i=from; i != to; i += step) {
1178 AliITSlayer &layer=fgLayers[i];
1179 Double_t r=layer.GetR();
1182 Double_t hI=i-0.5*step;
1183 if (TMath::Abs(hI-1.5)<0.01 || TMath::Abs(hI-3.5)<0.01) {
1184 Double_t rs=0.5*(fgLayers[i-step].GetR() + r);
1185 Double_t d=0.0034, x0=38.6;
1186 if (TMath::Abs(hI-1.5)<0.01) {rs=9.; d=0.0097; x0=42;}
1187 if (!t->PropagateTo(rs,-step*d,x0)) {
1193 // remember old position [SR, GSI 18.02.2003]
1194 Double_t oldX=0., oldY=0., oldZ=0.;
1195 if (t->IsStartedTimeIntegral() && step==1) {
1196 t->GetGlobalXYZat(t->GetX(),oldX,oldY,oldZ);
1201 if (!t->GetGlobalXYZat(r,x,y,z)) {
1204 Double_t phi=TMath::ATan2(y,x);
1205 Int_t idet=layer.FindDetectorIndex(phi,z);
1209 const AliITSdetector &det=layer.GetDetector(idet);
1211 if (!t->Propagate(phi,det.GetR())) {
1214 t->SetDetectorIndex(idet);
1216 const AliITSclusterV2 *cl=0;
1217 Double_t maxchi2=kMaxChi2;
1221 const AliITSclusterV2 *c=(AliITSclusterV2 *)GetCluster(idx);
1222 if (idet != c->GetDetectorIndex()) {
1223 idet=c->GetDetectorIndex();
1224 const AliITSdetector &det=layer.GetDetector(idet);
1225 if (!t->Propagate(det.GetPhi(),det.GetR())) {
1228 t->SetDetectorIndex(idet);
1230 Double_t chi2=t->GetPredictedChi2(c);
1240 if (t->GetNumberOfClusters()>2) {
1241 Double_t dz=4*TMath::Sqrt(t->GetSigmaZ2()+kSigmaZ2[i]);
1242 Double_t dy=4*TMath::Sqrt(t->GetSigmaY2()+kSigmaY2[i]);
1243 Double_t zmin=t->GetZ() - dz;
1244 Double_t zmax=t->GetZ() + dz;
1245 Double_t ymin=t->GetY() + phi*r - dy;
1246 Double_t ymax=t->GetY() + phi*r + dy;
1247 layer.SelectClusters(zmin,zmax,ymin,ymax);
1249 const AliITSclusterV2 *c=0; Int_t ci=-1;
1250 while ((c=layer.GetNextCluster(ci))!=0) {
1251 if (idet != c->GetDetectorIndex()) continue;
1252 Double_t chi2=t->GetPredictedChi2(c);
1253 if (chi2<maxchi2) { cl=c; maxchi2=chi2; idx=ci; }
1258 if (!t->Update(cl,maxchi2,idx)) {
1261 t->SetSampledEdx(cl->GetQ(),t->GetNumberOfClusters()-1);
1266 Double_t d=layer.GetThickness(t->GetY(),t->GetZ(),x0);
1267 t->CorrectForMaterial(-step*d,x0);
1270 // track time update [SR, GSI 17.02.2003]
1271 if (t->IsStartedTimeIntegral() && step==1) {
1272 Double_t newX, newY, newZ;
1273 t->GetGlobalXYZat(t->GetX(),newX,newY,newZ);
1274 Double_t dL2 = (oldX-newX)*(oldX-newX) + (oldY-newY)*(oldY-newY) +
1275 (oldZ-newZ)*(oldZ-newZ);
1276 t->AddTimeStep(TMath::Sqrt(dL2));
1282 if (!t->PropagateTo(xx,0.,0.)) return kFALSE;
1286 void AliITStrackerV2::UseClusters(const AliKalmanTrack *t, Int_t from) const {
1287 //--------------------------------------------------------------------
1288 // This function marks clusters assigned to the track
1289 //--------------------------------------------------------------------
1290 AliTracker::UseClusters(t,from);
1292 AliITSclusterV2 *c=(AliITSclusterV2 *)GetCluster(t->GetClusterIndex(0));
1293 //if (c->GetQ()>2) c->Use();
1294 if (c->GetSigmaZ2()>0.1) c->Use();
1295 c=(AliITSclusterV2 *)GetCluster(t->GetClusterIndex(1));
1296 //if (c->GetQ()>2) c->Use();
1297 if (c->GetSigmaZ2()>0.1) c->Use();