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 HLT ITS tracker class
18 // It reads AliITSclusterV2 clusters and HLT ESD tracks and creates
19 // AliITStrackV2 tracks. For details, see also RunHLTITS.C macro.
20 // Origin: Cvetan Cheshkov, CERN, Cvetan.Cheshkov@cern.ch
21 //-------------------------------------------------------------------------
24 #include "AliHLTITStrack.h"
25 #include "AliHLTITStracker.h"
27 ClassImp(AliHLTITStracker)
29 static Int_t CorrectForDeadZoneMaterial(AliITStrackV2 *t) {
30 //--------------------------------------------------------------------
31 // Correction for the material between the TPC and the ITS
32 // (should it belong to the TPC code ?)
33 //--------------------------------------------------------------------
34 Double_t riw=80., diw=0.0053, x0iw=30; // TPC inner wall ?
35 Double_t rcd=61., dcd=0.0053, x0cd=30; // TPC "central drum" ?
36 Double_t yr=12.8, dr=0.03; // rods ?
37 Double_t zm=0.2, dm=0.40; // membrane
38 //Double_t rr=52., dr=0.19, x0r=24., yyr=7.77; //rails
39 Double_t rs=50., ds=0.001; // something belonging to the ITS (screen ?)
41 if (t->GetX() > riw) {
42 if (!t->PropagateTo(riw,diw,x0iw)) return 1;
43 if (TMath::Abs(t->GetY())>yr) t->CorrectForMaterial(dr);
44 if (TMath::Abs(t->GetZ())<zm) t->CorrectForMaterial(dm);
45 if (!t->PropagateTo(rcd,dcd,x0cd)) return 1;
46 //Double_t x,y,z; t->GetGlobalXYZat(rr,x,y,z);
47 //if (TMath::Abs(y)<yyr) t->PropagateTo(rr,dr,x0r);
48 if (!t->PropagateTo(rs,ds)) return 1;
49 } else if (t->GetX() < rs) {
50 if (!t->PropagateTo(rs,-ds)) return 1;
51 //Double_t x,y,z; t->GetGlobalXYZat(rr,x,y,z);
52 //if (TMath::Abs(y)<yyr) t->PropagateTo(rr,-dr,x0r);
53 if (!t->PropagateTo(rcd,-dcd,x0cd)) return 1;
54 if (!t->PropagateTo(riw+0.001,-diw,x0iw)) return 1;
56 // ::Error("CorrectForDeadZoneMaterial","track is already in the dead zone !");
63 Int_t AliHLTITStracker::Clusters2Tracks(AliESD *event) {
64 //--------------------------------------------------------------------
65 // This functions reconstructs HLT ITS tracks
66 //--------------------------------------------------------------------
67 TObjArray itsTracks(15000);
69 {/* Read HLT ESD tracks */
71 nentr=event->GetNumberOfTracks();
72 Info("Clusters2Tracks", "Number of ESD HLT tracks: %d\n", nentr);
75 AliESDtrack *esd=event->GetTrack(nentr);
79 t=new AliHLTITStrack(*esd);
80 } catch (const Char_t *msg) {
81 Warning("Clusters2Tracks",msg);
85 if (TMath::Abs(t->GetD(GetX(),GetY()))>5) {
90 if (CorrectForDeadZoneMaterial(t)!=0) {
91 Warning("Clusters2Tracks",
92 "failed to correct for the material in the dead zone !\n");
98 } /* End Read HLT ESD tracks */
101 Int_t nentr=itsTracks.GetEntriesFast();
102 Info("Clusters2Tracks", "Number of Selected for tracking HLT ESD tracks: %d\n", nentr);
105 for (fPass=0; fPass<2; fPass++) {
106 Int_t &constraint=fConstraint[fPass]; if (constraint<0) continue;
107 for (Int_t i=0; i<nentr; i++) {
108 AliHLTITStrack *t=(AliHLTITStrack*)itsTracks.UncheckedAt(i);
109 if (t==0) continue; //this track has been already tracked
110 Int_t tpcLabel=t->GetLabel(); //save the TPC track label
111 ResetTrackToFollow(*t);
114 for (FollowProlongation(); fI<kMaxLayer; fI++) {
115 while (TakeNextProlongation()) FollowProlongation();
118 if (fBestTrack.GetNumberOfClusters() == 0) continue;
120 if (fConstraint[fPass]) {
121 ResetTrackToFollow(*t);
122 if (!RefitAt(3.7, &fTrackToFollow, &fBestTrack)) continue;
126 if (!fBestTrack.PropagateTo(3.,0.0028,65.19)) continue;
127 if (!fBestTrack.PropagateToVertex(event->GetVertex())) continue;
128 fBestTrack.SetLabel(tpcLabel);
129 fBestTrack.CookdEdx();
130 CookLabel(&fBestTrack,0.); //For comparison only
131 fBestTrack.UpdateESDtrack(AliESDtrack::kITSin);
133 AliESDtrack *esdTrack =fBestTrack.GetESDtrack();
134 Float_t r[3]={0.,0.,0.};
136 esdTrack->RelateToVertex(event->GetVertex(),GetBz(r),maxD);
138 UseClusters(&fBestTrack);
139 delete itsTracks.RemoveAt(i);
146 Info("Clusters2Tracks","Number of prolonged tracks: %d\n",ntrk);
151 Int_t AliHLTITStracker::PropagateBack(AliESD *event) {
152 //--------------------------------------------------------------------
153 // This functions propagates reconstructed ITS tracks back
154 //--------------------------------------------------------------------
155 Int_t nentr=event->GetNumberOfTracks();
156 Info("PropagateBack", "The method is not yet implemented! %d\n", nentr);
160 Int_t AliHLTITStracker::RefitInward(AliESD *event) {
161 //--------------------------------------------------------------------
162 // This functions refits ITS tracks using the
163 // "inward propagated" TPC tracks
164 //--------------------------------------------------------------------
166 Int_t nentr=event->GetNumberOfTracks();
167 Info("RefitInward", "The method is not yet implemented! %d",nentr);