]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/ITS/AliL3ITStracker.cxx
Major update of the HLT Hough transform and ITS tracking code. Hough transform tracks...
[u/mrichter/AliRoot.git] / HLT / ITS / AliL3ITStracker.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, 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 //               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 //-------------------------------------------------------------------------
22
23 #include "AliESD.h"
24 #include "AliL3ITStrack.h"
25 #include "AliL3ITStracker.h"
26
27 ClassImp(AliL3ITStracker)
28
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 ?)
40
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;
55   } else {
56     //    ::Error("CorrectForDeadZoneMaterial","track is already in the dead zone !");
57     return 1;
58   }
59   
60   return 0;
61 }
62
63 Int_t AliL3ITStracker::Clusters2Tracks(AliESD *event) {
64   //--------------------------------------------------------------------
65   // This functions reconstructs HLT ITS tracks
66   //--------------------------------------------------------------------
67   TObjArray itsTracks(15000);
68
69   {/* Read HLT ESD tracks */
70     Int_t nentr;
71     nentr=event->GetNumberOfTracks();
72     Info("Clusters2Tracks", "Number of ESD HLT tracks: %d\n", nentr);
73     while (nentr--) {
74
75       AliESDtrack *esd=event->GetTrack(nentr);
76
77       AliL3ITStrack *t=0;
78       try {
79         t=new AliL3ITStrack(*esd);
80       } catch (const Char_t *msg) {
81         Warning("Clusters2Tracks",msg);
82         delete t;
83         continue;
84       }
85       if (TMath::Abs(t->GetD())>5) {
86         delete t;
87         continue;
88       }
89
90       if (CorrectForDeadZoneMaterial(t)!=0) {
91         Warning("Clusters2Tracks",
92                 "failed to correct for the material in the dead zone !\n");
93         delete t;
94         continue;
95       }
96       itsTracks.AddLast(t);
97     }
98   } /* End Read HLT ESD tracks */
99
100   itsTracks.Sort();
101   Int_t nentr=itsTracks.GetEntriesFast();
102   Info("Clusters2Tracks", "Number of Selected for tracking HLT ESD tracks: %d\n", nentr);
103
104   Int_t ntrk=0;
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        AliL3ITStrack *t=(AliL3ITStrack*)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);
112        ResetBestTrack();
113
114        for (FollowProlongation(); fI<kMaxLayer; fI++) {
115          while (TakeNextProlongation()) FollowProlongation();
116        }
117
118        if (fBestTrack.GetNumberOfClusters() == 0) continue;
119        
120        if (fConstraint[fPass]) {
121          ResetTrackToFollow(*t);
122          if (!RefitAt(3.7, &fTrackToFollow, &fBestTrack)) continue;
123          ResetBestTrack();
124        }
125        
126        if (!fBestTrack.PropagateTo(3.,0.0028,65.19)) continue;
127        if (!fBestTrack.PropagateToVertex()) continue;
128        fBestTrack.SetLabel(tpcLabel);
129        fBestTrack.CookdEdx();
130        CookLabel(&fBestTrack,0.); //For comparison only
131        fBestTrack.UpdateESDtrack(AliESDtrack::kITSin);
132        UseClusters(&fBestTrack);
133        delete itsTracks.RemoveAt(i);
134        ntrk++;
135      }
136   }
137
138   itsTracks.Delete();
139
140   Info("Clusters2Tracks","Number of prolonged tracks: %d\n",ntrk);
141
142   return 0;
143 }
144
145 Int_t AliL3ITStracker::PropagateBack(AliESD *event) {
146   //--------------------------------------------------------------------
147   // This functions propagates reconstructed ITS tracks back
148   //--------------------------------------------------------------------
149   Int_t nentr=event->GetNumberOfTracks();
150   Info("PropagateBack", "The method is not yet implemented! %d\n", nentr);
151   return 0;
152 }
153
154 Int_t AliL3ITStracker::RefitInward(AliESD *event) {
155   //--------------------------------------------------------------------
156   // This functions refits ITS tracks using the 
157   // "inward propagated" TPC tracks
158   //--------------------------------------------------------------------
159
160   Int_t nentr=event->GetNumberOfTracks();
161   Info("RefitInward", "The method is not yet implemented! %d",nentr);
162   return 0;
163 }